download the original source code.
  1 /*
  2    Example 5
  3 
  4    Interface:    Linear-Algebraic (IJ), Babel-based version
  5 
  6    Compile with: make ex5bxx
  7 
  8    Sample run:   mpirun -np 4 ex5bxx
  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 
 22 #include <math.h>
 23 #include <assert.h>
 24 #include "_hypre_utilities.h"
 25 #include "HYPRE_krylov.h"
 26 #include "HYPRE.h"
 27 #include "HYPRE_parcsr_ls.h"
 28 
 29 /* Babel interface headers */
 30 #include "bHYPRE.hxx"
 31 #include "bHYPRE_Vector.h"
 32 #include "bHYPRE_IJParCSRMatrix.h"
 33 #include "bHYPRE_IJParCSRVector.h"
 34 #include "bHYPRE_ParCSRDiagScale.h"
 35 #include "bHYPRE_BoomerAMG.h"
 36 
 37 int main (int argc, char *argv[])
 38 {
 39    using namespace ::bHYPRE;
 40    int i;
 41    int myid, num_procs;
 42    int N, n;
 43 
 44    int ilower, iupper;
 45    int local_size, extra;
 46 
 47    int solver_id;
 48    int print_solution;
 49 
 50    double h, h2;
 51    MPI_Comm mpicommworld = MPI_COMM_WORLD;
 52 
 53    int ierr = 0;
 54    /* If this gets set to anything else, it's an error.
 55     Most functions return error flags, 0 unless there's an error.
 56     For clarity, they aren't checked much in this file. */
 57 
 58    MPICommunicator mpi_comm;
 59    IJParCSRMatrix parcsr_A;
 60    Operator  op_A;
 61    IJParCSRVector par_b;
 62    IJParCSRVector par_x;
 63    Vector vec_x;
 64 
 65    BoomerAMG amg_solver;
 66    ParaSails ps_solver;
 67    PCG pcg_solver;
 68    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 = 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 = int(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 = IJParCSRMatrix::Create( mpi_comm,
162                                         ilower, iupper, ilower, iupper );
163 
164    op_A = parcsr_A; /* we can eliminate op_A later, it's really needed only in C */
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    parcsr_A.Initialize();
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          parcsr_A.SetValues( 1, &nnz, &i, cols, values, 5 );
231       }
232    }
233 
234    /* Assemble after setting the coefficients */
235    parcsr_A.Assemble();
236 
237    /* Create the rhs and solution */
238    par_b = IJParCSRVector::Create( mpi_comm, ilower, iupper );
239 
240    par_b.Initialize();
241 
242    par_x = IJParCSRVector::Create( mpi_comm, ilower, iupper );
243 
244    par_x.Initialize();
245 
246    /* Set the rhs values to h^2 and the solution to zero */
247    {
248       double *rhs_values, *x_values;
249       int    *rows;
250 
251       rhs_values = new double[local_size];
252       x_values = new double[local_size];
253       rows = new int[local_size];
254 
255       for (i=0; i<local_size; i++)
256       {
257          rhs_values[i] = h2;
258          x_values[i] = 0.0;
259          rows[i] = ilower + i;
260       }
261 
262       par_b.SetValues( local_size, rows, rhs_values );
263       par_x.SetValues( local_size, rows, x_values );
264 
265       delete[] x_values;
266       delete[] rhs_values;
267       delete[] rows;
268    }
269 
270    par_b.Assemble();
271    par_x.Assemble();
272 
273    /* Choose a solver and solve the system */
274 
275    /* AMG */
276    if (solver_id == 0)
277    {
278       int num_iterations;
279       double final_res_norm;
280 
281       /* Create solver */
282       amg_solver = BoomerAMG::Create( mpi_comm, parcsr_A );
283 
284       /* Set some parameters (See Reference Manual for more parameters) */
285       amg_solver.SetIntParameter( "PrintLevel", 3 );  /* print solve info + parameters */
286       amg_solver.SetIntParameter( "CoarsenType", 6); /* Falgout coarsening */
287       amg_solver.SetIntParameter( "RelaxType", 3);   /* G-S/Jacobi hybrid relaxation */
288       amg_solver.SetIntParameter( "NumSweeps", 1);   /* Sweeeps on each level */
289       amg_solver.SetIntParameter( "MaxLevels", 20);  /* maximum number of levels */
290       amg_solver.SetDoubleParameter( "Tolerance", 1e-7);      /* conv. tolerance */
291 
292       /* Now setup and solve! */
293       // In C++, use of the Apply function is slightly more complicated than in other
294       // languages.
295       //   The second argument, the solution (and initial guess), is declared inout in
296       // Interfaces.idl.  For C++, Babel implements an inout (or out) argument as a
297       // reference argument, and obviously not a const one.
298       //   So it is important not to implicitly ask the C++ environment for a type
299       // conversion when supplying the second argument - results will be unpredictable!
300       //   Our solution vector is par_x, declared type IJParCSRVector, derived from
301       // Vector. In other languages we would use par_x directly for the second argument.
302       // In C++ we must declare a Vector, vec_x, for use as the second argument.
303       //   Fortunately, an assignment is all you have to do to set up vec_x.
304       // Babel casts and code other inside the Apply function will recover the original
305       // IJParCSRVector, par_x, to make sure it gets computed correctly.
306       amg_solver.Setup( par_b, par_x );
307       vec_x = par_x;
308       amg_solver.Apply( par_b, vec_x );
309 
310       /* Run info - needed logging turned on */
311 
312       ierr += amg_solver.GetIntValue( "NumIterations", num_iterations );
313       amg_solver.GetDoubleValue( "RelResidualNorm", final_res_norm );
314 
315       if (myid == 0)
316       {
317          printf("\n");
318          printf("Iterations = %d\n", num_iterations);
319          printf("Final Relative Residual Norm = %e\n", final_res_norm);
320          printf("\n");
321       }
322 
323       /* Destroy solver */
324       /* In C++, unlike C, deleteRef of amg_solver is not needed here -
325          it happens automatically. */
326    }
327 
328    /* PCG */
329    else if (solver_id == 50)
330    {
331       int num_iterations;
332       double final_res_norm;
333 
334       /* Create solver */
335       pcg_solver = PCG::Create( mpi_comm, op_A );
336 
337       /* Set some parameters (See Reference Manual for more parameters) */
338       pcg_solver.SetIntParameter( "MaxIter", 1000 ); /* max iterations */
339       pcg_solver.SetDoubleParameter( "Tolerance", 1e-7 ); /* conv. tolerance */
340       pcg_solver.SetIntParameter( "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
341       pcg_solver.SetIntParameter( "PrintLevel", 2 ); /* prints out the iteration info */
342       pcg_solver.SetIntParameter( "Logging", 1 ); /* needed to get run info later */
343 
344       identity = IdentitySolver::Create( mpi_comm );
345       pcg_solver.SetPreconditioner( identity );
346 
347       /* Now setup and solve! */
348       // See comments for solver 0 on Apply and vec_x.
349       pcg_solver.Setup( par_b, par_x );
350       vec_x = par_x;
351       pcg_solver.Apply( par_b, vec_x );
352 
353       /* Run info - needed logging turned on */
354       pcg_solver.GetIntValue( "NumIterations", num_iterations );
355       pcg_solver.GetDoubleValue( "RelResidualNorm", final_res_norm );
356       if (myid == 0)
357       {
358          printf("\n");
359          printf("Iterations = %d\n", num_iterations);
360          printf("Final Relative Residual Norm = %e\n", final_res_norm);
361          printf("\n");
362       }
363 
364       /* Destroy solvers */
365       /* In C++, unlike C, deleteRef's of solvers are not needed here -
366          it happens automatically. */
367    }
368    /* PCG with AMG preconditioner */
369    else if (solver_id == 1)
370    {
371       int num_iterations;
372       double final_res_norm;
373 
374       /* Create solver */
375       pcg_solver = PCG::Create( mpi_comm, op_A );
376 
377       /* Set some parameters (See Reference Manual for more parameters) */
378       pcg_solver.SetIntParameter( "MaxIter", 1000 ); /* max iterations */
379       pcg_solver.SetDoubleParameter( "Tolerance", 1e-7 ); /* conv. tolerance */
380       pcg_solver.SetIntParameter( "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
381       pcg_solver.SetIntParameter( "PrintLevel", 2 ); /* prints out the iteration info */
382       pcg_solver.SetIntParameter( "Logging", 1 ); /* needed to get run info later */
383 
384       /* Now set up the AMG preconditioner and specify any parameters */
385       amg_solver = BoomerAMG::Create( mpi_comm, parcsr_A );
386       amg_solver.SetIntParameter( "PrintLevel", 1 ); /* print amg solution info*/
387       amg_solver.SetIntParameter( "CoarsenType", 6); /* Falgout coarsening */
388       amg_solver.SetIntParameter( "RelaxType", 6);   /* Sym G-S/Jacobi hybrid relaxation */
389       amg_solver.SetIntParameter( "NumSweeps", 1);   /* Sweeeps on each level */
390       amg_solver.SetDoubleParameter( "Tolerance", 0);      /* conv. tolerance */
391       amg_solver.SetIntParameter( "MaxIter", 1 ); /* do only one iteration! */
392 
393       /* Set the PCG preconditioner */
394       pcg_solver.SetPreconditioner( amg_solver );
395 
396       /* Now setup and solve! */
397       // See comments for solver 0 on Apply and vec_x.
398       pcg_solver.Setup( par_b, par_x );
399       vec_x = par_x;
400       pcg_solver.Apply( par_b, vec_x );
401 
402       /* Run info - needed logging turned on */
403       pcg_solver.GetIntValue( "NumIterations", num_iterations );
404       pcg_solver.GetDoubleValue( "RelResidualNorm", final_res_norm );
405       if (myid == 0)
406       {
407          printf("\n");
408          printf("Iterations = %d\n", num_iterations);
409          printf("Final Relative Residual Norm = %e\n", final_res_norm);
410          printf("\n");
411       }
412 
413       /* Destroy solver and preconditioner */
414       /* In C++, unlike C, deleteRef's of solvers are not needed here -
415          it happens automatically. */
416    }
417 
418    /* PCG with Parasails Preconditioner */
419    else if (solver_id == 8)
420    {
421       int    num_iterations;
422       double final_res_norm;
423 
424       int      sai_max_levels = 1;
425       double   sai_threshold = 0.1;
426       double   sai_filter = 0.05;
427       int      sai_sym = 1;
428 
429       /* Create solver */
430       pcg_solver = PCG::Create( mpi_comm, op_A );
431 
432       /* Set some parameters (See Reference Manual for more parameters) */
433       pcg_solver.SetIntParameter( "MaxIter", 1000 ); /* max iterations */
434       pcg_solver.SetDoubleParameter( "Tolerance", 1e-7 ); /* conv. tolerance */
435       pcg_solver.SetIntParameter( "TwoNorm", 1 ); /* use the two norm as the stopping criteria */
436       pcg_solver.SetIntParameter( "PrintLevel", 2 ); /* prints out the iteration info */
437       pcg_solver.SetIntParameter( "Logging", 1 ); /* needed to get run info later */
438 
439       /* Now set up the ParaSails preconditioner and specify any parameters */
440       ps_solver = ParaSails::Create( mpi_comm, parcsr_A );
441 
442       /* Set some parameters (See Reference Manual for more parameters) */
443       ps_solver.SetDoubleParameter( "Thresh", sai_threshold );
444       ps_solver.SetIntParameter( "Nlevels", sai_max_levels );
445       ps_solver.SetDoubleParameter( "Filter", sai_filter );
446       ps_solver.SetIntParameter( "Sym", sai_sym );
447       ps_solver.SetIntParameter( "Logging", 3 );
448 
449       /* Set the PCG preconditioner */
450       pcg_solver.SetPreconditioner( ps_solver );
451 
452       /* Now setup and solve! */
453       // See comments for solver 0 on Apply and vec_x.
454       pcg_solver.Setup( par_b, par_x);
455       vec_x = par_x;
456       pcg_solver.Apply( par_b, vec_x);
457 
458 
459       /* Run info - needed logging turned on */
460       pcg_solver.GetIntValue( "NumIterations", num_iterations );
461       pcg_solver.GetDoubleValue( "RelResidualNorm", final_res_norm );
462       if (myid == 0)
463       {
464          printf("\n");
465          printf("Iterations = %d\n", num_iterations);
466          printf("Final Relative Residual Norm = %e\n", final_res_norm);
467          printf("\n");
468       }
469 
470       /* Destroy solver and preconditioner */
471       /* In C++, unlike C, deleteRef's of solvers are not needed here -
472          it happens automatically. */
473    }
474    else
475    {
476       if (myid ==0) printf("Invalid solver id specified.\n");
477    }
478 
479    /* Print the solution */
480    if (print_solution)
481       par_x.Print( "ij.out.x" );
482 
483    /* Clean up */
484    /* In C++, unlike C, deleteRef gets called automatically, so we don't
485       need any explicit cleanup. */
486 
487    hypre_assert( ierr == 0 );
488 
489    /* Finalize MPI*/
490    MPI_Finalize();
491 
492    return(0);
493 }
syntax highlighted by Code2HTML, v. 0.9.1