download the original source code.
  1 /*
  2    Example 5
  3 
  4    Interface:    Linear-Algebraic (IJ), Babel-based version
  5 
  6    Compile with: make ex5b
  7 
  8    Sample run:   mpirun -np 4 ex5b
  9 
 10    Description:  This example solves the 2-D
 11                  Laplacian problem with zero boundary conditions
 12                  on an nxn grid.  The number of unknowns is N=n^2.
 13                  The standard 5-point stencil is used, and we solve
 14                  for the interior nodes only.
 15 
 16                  This example solves the same problem as Example 3.
 17                  Available solvers are AMG, PCG, and PCG with AMG or
 18                  Parasails preconditioners.
 19 */
 20 
 21 #include <math.h>
 22 #include <assert.h>
 23 #include "utilities.h"
 24 #include "krylov.h"
 25 #include "HYPRE.h"
 26 #include "HYPRE_parcsr_ls.h"
 27 
 28 /* Babel interface headers */
 29 #include "bHYPRE.h"
 30 #include "bHYPRE_Vector.h"
 31 #include "bHYPRE_IJParCSRMatrix.h"
 32 #include "bHYPRE_IJParCSRVector.h"
 33 #include "bHYPRE_ParCSRDiagScale.h"
 34 #include "bHYPRE_BoomerAMG.h"
 35 
 36 int main (int argc, char *argv[])
 37 {
 38    int i;
 39    int myid, num_procs;
 40    int N, n;
 41 
 42    int ilower, iupper;
 43    int local_size, extra;
 44 
 45    int solver_id;
 46    int print_solution;
 47 
 48    double h, h2;
 49    MPI_Comm mpicommworld = MPI_COMM_WORLD;
 50 
 51    int ierr = 0;
 52    /* If this gets set to anything else, it's an error.
 53     Most functions return error flags, 0 unless there's an error.
 54     For clarity, they aren't checked much in this file. */
 55 
 56    bHYPRE_MPICommunicator mpi_comm;
 57    bHYPRE_IJParCSRMatrix parcsr_A;
 58    bHYPRE_Operator           op_A;
 59    bHYPRE_IJParCSRVector par_b;
 60    bHYPRE_IJParCSRVector par_x;
 61    bHYPRE_Vector vec_b;
 62    bHYPRE_Vector vec_x;
 63 
 64    bHYPRE_Solver precond;
 65    bHYPRE_BoomerAMG amg_solver;
 66    bHYPRE_ParaSails ps_solver;
 67    bHYPRE_PCG pcg_solver;
 68    bHYPRE_IdentitySolver identity;
 69 
 70    /* Initialize MPI */
 71    MPI_Init(&argc, &argv);
 72    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 73    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 74    mpi_comm = bHYPRE_MPICommunicator_CreateC( &mpicommworld );
 75 
 76    /* Default problem parameters */
 77    n = 33;
 78    solver_id = 0;
 79    print_solution  = 0;
 80 
 81    /* Parse command line */
 82    {
 83       int arg_index = 0;
 84       int print_usage = 0;
 85 
 86       while (arg_index < argc)
 87       {
 88          if ( strcmp(argv[arg_index], "-n") == 0 )
 89          {
 90             arg_index++;
 91             n = atoi(argv[arg_index++]);
 92          }
 93          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 94          {
 95             arg_index++;
 96             solver_id = atoi(argv[arg_index++]);
 97          }
 98          else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
 99          {
100             arg_index++;
101             print_solution = 1;
102          }
103          else if ( strcmp(argv[arg_index], "-help") == 0 )
104          {
105             print_usage = 1;
106             break;
107          }
108          else
109          {
110             arg_index++;
111          }
112       }
113 
114       if ((print_usage) && (myid == 0))
115       {
116          printf("\n");
117          printf("Usage: %s [<options>]\n", argv[0]);
118          printf("\n");
119          printf("  -n <n>              : problem size in each direction (default: 33)\n");
120          printf("  -solver <ID>        : solver ID\n");
121          printf("                        0  - AMG (default) \n");
122          printf("                        1  - AMG-PCG\n");
123          printf("                        8  - ParaSails-PCG\n");
124          printf("                        50 - PCG\n");
125          printf("  -print_solution     : print the solution vector\n");
126          printf("\n");
127       }
128 
129       if (print_usage)
130       {
131          MPI_Finalize();
132          return (0);
133       }
134    }
135 
136    /* Preliminaries: want at least one processor per row */
137    if (n*n < num_procs) n = sqrt(num_procs) + 1;
138    N = n*n; /* global number of rows */
139    h = 1.0/(n+1); /* mesh size*/
140    h2 = h*h;
141 
142    /* Each processor knows only of its own rows - the range is denoted by ilower
143       and upper.  Here we partition the rows. We account for the fact that
144       N may not divide evenly by the number of processors. */
145    local_size = N/num_procs;
146    extra = N - local_size*num_procs;
147 
148    ilower = local_size*myid;
149    ilower += hypre_min(myid, extra);
150 
151    iupper = local_size*(myid+1);
152    iupper += hypre_min(myid+1, extra);
153    iupper = iupper - 1;
154 
155    /* How many rows do I have? */
156    local_size = iupper - ilower + 1;
157 
158    /* Create the matrix.
159       Note that this is a square matrix, so we indicate the row partition
160       size twice (since number of rows = number of cols) */
161    parcsr_A = bHYPRE_IJParCSRMatrix_Create( mpi_comm,
162                                             ilower, iupper, ilower, iupper );
163 
164    op_A = bHYPRE_Operator__cast( parcsr_A ); /* needed later as a function argument */
165 
166    /* Choose a parallel csr format storage (see the User's Manual) */
167    /* Note: Here the HYPRE interface requires a SetObjectType call.
168       I am using the bHYPRE interface in a way which does not because
169       the object type is already specified through the class name. */
170 
171    /* Initialize before setting coefficients */
172    bHYPRE_IJParCSRMatrix_Initialize( parcsr_A );
173 
174    /* Now go through my local rows and set the matrix entries.
175       Each row has at most 5 entries. For example, if n=3:
176 
177       A = [M -I 0; -I M -I; 0 -I M]
178       M = [4 -1 0; -1 4 -1; 0 -1 4]
179 
180       Note that here we are setting one row at a time, though
181       one could set all the rows together (see the User's Manual).
182    */
183    {
184       int nnz;
185       double values[5];
186       int cols[5];
187 
188       for (i = ilower; i <= iupper; i++)
189       {
190          nnz = 0;
191 
192          /* The left identity block:position i-n */
193          if ((i-n)>=0)
194          {
195 	    cols[nnz] = i-n;
196 	    values[nnz] = -1.0;
197 	    nnz++;
198          }
199 
200          /* The left -1: position i-1 */
201          if (i%n)
202          {
203             cols[nnz] = i-1;
204             values[nnz] = -1.0;
205             nnz++;
206          }
207 
208          /* Set the diagonal: position i */
209          cols[nnz] = i;
210          values[nnz] = 4.0;
211          nnz++;
212 
213          /* The right -1: position i+1 */
214          if ((i+1)%n)
215          {
216             cols[nnz] = i+1;
217             values[nnz] = -1.0;
218             nnz++;
219          }
220 
221          /* The right identity block:position i+n */
222          if ((i+n)< N)
223          {
224             cols[nnz] = i+n;
225             values[nnz] = -1.0;
226             nnz++;
227          }
228 
229          /* Set the values for row i */
230          bHYPRE_IJParCSRMatrix_SetValues( parcsr_A, 1, &nnz, &i, cols, values, 5 );
231       }
232    }
233 
234    /* Assemble after setting the coefficients */
235    bHYPRE_IJParCSRMatrix_Assemble( parcsr_A );
236 
237    /* Create the rhs and solution */
238    par_b = bHYPRE_IJParCSRVector_Create( mpi_comm, ilower, iupper );
239 
240    vec_b = bHYPRE_Vector__cast( par_b ); /* needed later for function arguments */
241 
242    bHYPRE_IJParCSRVector_Initialize( par_b );
243 
244    par_x = bHYPRE_IJParCSRVector_Create( mpi_comm, ilower, iupper );
245 
246    vec_x = bHYPRE_Vector__cast( par_x ); /* needed later for function arguments */
247 
248    bHYPRE_IJParCSRVector_Initialize( par_x );
249 
250    /* Set the rhs values to h^2 and the solution to zero */
251    {
252       double *rhs_values, *x_values;
253       int    *rows;
254 
255       rhs_values = calloc(local_size, sizeof(double));
256       x_values = calloc(local_size, sizeof(double));
257       rows = calloc(local_size, sizeof(int));
258 
259       for (i=0; i<local_size; i++)
260       {
261          rhs_values[i] = h2;
262          x_values[i] = 0.0;
263          rows[i] = ilower + i;
264       }
265 
266       bHYPRE_IJParCSRVector_SetValues( par_b, local_size, rows, rhs_values );
267       bHYPRE_IJParCSRVector_SetValues( par_x, local_size, rows, x_values );
268 
269       free(x_values);
270       free(rhs_values);
271       free(rows);
272    }
273 
274    bHYPRE_IJParCSRVector_Assemble( par_b );
275    bHYPRE_IJParCSRVector_Assemble( par_x );
276 
277    /* Choose a solver and solve the system */
278 
279    /* AMG */
280    if (solver_id == 0)
281    {
282       int num_iterations;
283       double final_res_norm;
284 
285       /* Create solver */
286       amg_solver = bHYPRE_BoomerAMG_Create( mpi_comm, parcsr_A );
287 
288       /* Set some parameters (See Reference Manual for more parameters) */
289       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "PrintLevel", 3 );  /* print solve info + parameters */
290       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "CoarsenType", 6); /* Falgout coarsening */
291       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "RelaxType", 3);   /* G-S/Jacobi hybrid relaxation */
292       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "NumSweeps", 1);   /* Sweeeps on each level */
293       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "MaxLevels", 20);  /* maximum number of levels */
294       bHYPRE_BoomerAMG_SetDoubleParameter( amg_solver, "Tolerance", 1e-7);      /* conv. tolerance */
295 
296       /* Now setup and solve! */
297       bHYPRE_BoomerAMG_Setup( amg_solver, vec_b, vec_x );
298       bHYPRE_BoomerAMG_Apply( amg_solver, vec_b, &vec_x );
299 
300       /* Run info - needed logging turned on */
301       ierr += bHYPRE_BoomerAMG_GetIntValue( amg_solver, "NumIterations", &num_iterations );
302       bHYPRE_BoomerAMG_GetDoubleValue( amg_solver, "RelResidualNorm", &final_res_norm );
303 
304       if (myid == 0)
305       {
306          printf("\n");
307          printf("Iterations = %d\n", num_iterations);
308          printf("Final Relative Residual Norm = %e\n", final_res_norm);
309          printf("\n");
310       }
311 
312       /* Destroy solver */
313       bHYPRE_BoomerAMG_deleteRef( amg_solver );
314    }
315 
316    /* PCG */
317    else if (solver_id == 50)
318    {
319       int num_iterations;
320       double final_res_norm;
321 
322       /* Create solver */
323       pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A );
324 
325       /* Set some parameters (See Reference Manual for more parameters) */
326       bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000 ); /* max iterations */
327       bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7 ); /* conv. tolerance */
328       bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
329       bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2 ); /* prints out the iteration info */
330       bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1 ); /* needed to get run info later */
331 
332       identity = bHYPRE_IdentitySolver_Create( mpi_comm );
333       precond = bHYPRE_Solver__cast( identity );
334       bHYPRE_PCG_SetPreconditioner( pcg_solver, precond );
335 
336       /* Now setup and solve! */
337       bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x );
338       bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x );
339 
340       /* Run info - needed logging turned on */
341       bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations );
342       bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm );
343       if (myid == 0)
344       {
345          printf("\n");
346          printf("Iterations = %d\n", num_iterations);
347          printf("Final Relative Residual Norm = %e\n", final_res_norm);
348          printf("\n");
349       }
350 
351       /* Destroy solvers */
352       bHYPRE_PCG_deleteRef( pcg_solver );
353       bHYPRE_Solver_deleteRef( precond );
354    }
355    /* PCG with AMG preconditioner */
356    else if (solver_id == 1)
357    {
358       int num_iterations;
359       double final_res_norm;
360 
361       /* Create solver */
362       pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A );
363 
364       /* Set some parameters (See Reference Manual for more parameters) */
365       bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000 ); /* max iterations */
366       bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7 ); /* conv. tolerance */
367       bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
368       bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2 ); /* prints out the iteration info */
369       bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1 ); /* needed to get run info later */
370 
371       /* Now set up the AMG preconditioner and specify any parameters */
372       amg_solver = bHYPRE_BoomerAMG_Create( mpi_comm, parcsr_A );
373       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "PrintLevel", 1 ); /* print amg solution info*/
374       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "CoarsenType", 6); /* Falgout coarsening */
375       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "RelaxType", 6);   /* Sym G-S/Jacobi hybrid relaxation */
376       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "NumSweeps", 1);   /* Sweeeps on each level */
377       bHYPRE_BoomerAMG_SetDoubleParameter( amg_solver, "Tolerance", 1e-3);      /* conv. tolerance */
378       bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "MaxIter", 1 ); /* do only one iteration! */
379 
380       /* Set the PCG preconditioner */
381       precond = bHYPRE_Solver__cast( amg_solver );
382       bHYPRE_PCG_SetPreconditioner( pcg_solver, precond );
383 
384       /* Now setup and solve! */
385       bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x );
386       bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x );
387 
388       /* Run info - needed logging turned on */
389       bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations );
390       bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm );
391       if (myid == 0)
392       {
393          printf("\n");
394          printf("Iterations = %d\n", num_iterations);
395          printf("Final Relative Residual Norm = %e\n", final_res_norm);
396          printf("\n");
397       }
398 
399       /* Destroy solver and preconditioner */
400       bHYPRE_PCG_deleteRef( pcg_solver );
401       bHYPRE_Solver_deleteRef( precond );
402    }
403 
404    /* PCG with Parasails Preconditioner */
405    else if (solver_id == 8)
406    {
407       int    num_iterations;
408       double final_res_norm;
409 
410       int      sai_max_levels = 1;
411       double   sai_threshold = 0.1;
412       double   sai_filter = 0.05;
413       int      sai_sym = 1;
414 
415       /* Create solver */
416       pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A );
417 
418       /* Set some parameters (See Reference Manual for more parameters) */
419       bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000 ); /* max iterations */
420       bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7 ); /* conv. tolerance */
421       bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
422       bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2 ); /* prints out the iteration info */
423       bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1 ); /* needed to get run info later */
424 
425       /* Now set up the ParaSails preconditioner and specify any parameters */
426       ps_solver = bHYPRE_ParaSails_Create( mpi_comm, parcsr_A );
427 
428       /* Set some parameters (See Reference Manual for more parameters) */
429       bHYPRE_ParaSails_SetDoubleParameter( ps_solver, "Thresh", sai_threshold );
430       bHYPRE_ParaSails_SetIntParameter( ps_solver, "Nlevels", sai_max_levels );
431       bHYPRE_ParaSails_SetDoubleParameter( ps_solver, "Filter", sai_filter );
432       bHYPRE_ParaSails_SetIntParameter( ps_solver, "Sym", sai_sym );
433       bHYPRE_ParaSails_SetIntParameter( ps_solver, "Logging", 3 );
434 
435       /* Set the PCG preconditioner */
436       precond = bHYPRE_Solver__cast( ps_solver );
437       bHYPRE_PCG_SetPreconditioner( pcg_solver, precond );
438 
439       /* Now setup and solve! */
440       bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x);
441       bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x);
442 
443 
444       /* Run info - needed logging turned on */
445       bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations );
446       bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm );
447       if (myid == 0)
448       {
449          printf("\n");
450          printf("Iterations = %d\n", num_iterations);
451          printf("Final Relative Residual Norm = %e\n", final_res_norm);
452          printf("\n");
453       }
454 
455       /* Destory solver and preconditioner */
456       bHYPRE_PCG_deleteRef( pcg_solver );
457       bHYPRE_Solver_deleteRef( precond );
458    }
459    else
460    {
461       if (myid ==0) printf("Invalid solver id specified.\n");
462    }
463 
464    /* Print the solution */
465    if (print_solution)
466       bHYPRE_IJParCSRVector_Print( par_x, "ij.out.x" );
467 
468    /* Clean up */
469    bHYPRE_IJParCSRMatrix_deleteRef( parcsr_A );
470    bHYPRE_IJParCSRVector_deleteRef( par_b );
471    bHYPRE_IJParCSRVector_deleteRef( par_x );
472    bHYPRE_MPICommunicator_deleteRef( mpi_comm );
473 
474    hypre_assert( ierr == 0 );
475 
476    /* Finalize MPI*/
477    MPI_Finalize();
478 
479    return(0);
480 }


syntax highlighted by Code2HTML, v. 0.9.1