download the original source code.
  1 /*
  2    Example 5
  3 
  4    Interface:    Linear-Algebraic (IJ)
  5 
  6    Compile with: make ex5
  7 
  8    Sample run:   mpirun -np 4 ex5
  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 "_hypre_utilities.h"
 23 #include "HYPRE_krylov.h"
 24 #include "HYPRE.h"
 25 #include "HYPRE_parcsr_ls.h"
 26 
 27 int main (int argc, char *argv[])
 28 {
 29    int i;
 30    int myid, num_procs;
 31    int N, n;
 32 
 33    int ilower, iupper;
 34    int local_size, extra;
 35 
 36    int solver_id;
 37    int print_solution, print_system;
 38 
 39    double h, h2;
 40 
 41    HYPRE_IJMatrix A;
 42    HYPRE_ParCSRMatrix parcsr_A;
 43    HYPRE_IJVector b;
 44    HYPRE_ParVector par_b;
 45    HYPRE_IJVector x;
 46    HYPRE_ParVector par_x;
 47 
 48    HYPRE_Solver solver, precond;
 49 
 50    /* Initialize MPI */
 51    MPI_Init(&argc, &argv);
 52    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 53    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 54 
 55    /* Default problem parameters */
 56    n = 33;
 57    solver_id = 0;
 58    print_solution  = 0;
 59    print_system = 0;
 60    
 61 
 62    /* Parse command line */
 63    {
 64       int arg_index = 0;
 65       int print_usage = 0;
 66 
 67       while (arg_index < argc)
 68       {
 69          if ( strcmp(argv[arg_index], "-n") == 0 )
 70          {
 71             arg_index++;
 72             n = atoi(argv[arg_index++]);
 73          }
 74          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 75          {
 76             arg_index++;
 77             solver_id = atoi(argv[arg_index++]);
 78          }
 79          else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
 80          {
 81             arg_index++;
 82             print_solution = 1;
 83          }
 84          else if ( strcmp(argv[arg_index], "-print_system") == 0 )
 85          {
 86             arg_index++;
 87             print_system = 1;
 88          }
 89 
 90 
 91          else if ( strcmp(argv[arg_index], "-help") == 0 )
 92          {
 93             print_usage = 1;
 94             break;
 95          }
 96          else
 97          {
 98             arg_index++;
 99          }
100       }
101 
102       if ((print_usage) && (myid == 0))
103       {
104          printf("\n");
105          printf("Usage: %s [<options>]\n", argv[0]);
106          printf("\n");
107          printf("  -n <n>              : problem size in each direction (default: 33)\n");
108          printf("  -solver <ID>        : solver ID\n");
109          printf("                        0  - AMG (default) \n");
110          printf("                        1  - AMG-PCG\n");
111          printf("                        8  - ParaSails-PCG\n");
112          printf("                        50 - PCG\n");
113          printf("  -print_solution     : print the solution vector\n");
114          printf("  -print_system       : print the matrix and rhs\n");
115          printf("\n");
116       }
117 
118       if (print_usage)
119       {
120          MPI_Finalize();
121          return (0);
122       }
123    }
124 
125    /* Preliminaries: want at least one processor per row */
126    if (n*n < num_procs) n = sqrt(num_procs) + 1;
127    N = n*n; /* global number of rows */
128    h = 1.0/(n+1); /* mesh size*/
129    h2 = h*h;
130 
131    /* Each processor knows only of its own rows - the range is denoted by ilower
132       and upper.  Here we partition the rows. We account for the fact that
133       N may not divide evenly by the number of processors. */
134    local_size = N/num_procs;
135    extra = N - local_size*num_procs;
136 
137    ilower = local_size*myid;
138    ilower += hypre_min(myid, extra);
139 
140    iupper = local_size*(myid+1);
141    iupper += hypre_min(myid+1, extra);
142    iupper = iupper - 1;
143 
144    /* How many rows do I have? */
145    local_size = iupper - ilower + 1;
146 
147    /* Create the matrix.
148       Note that this is a square matrix, so we indicate the row partition
149       size twice (since number of rows = number of cols) */
150    HYPRE_IJMatrixCreate(MPI_COMM_WORLD, ilower, iupper, ilower, iupper, &A);
151 
152    /* Choose a parallel csr format storage (see the User's Manual) */
153    HYPRE_IJMatrixSetObjectType(A, HYPRE_PARCSR);
154 
155    /* Initialize before setting coefficients */
156    HYPRE_IJMatrixInitialize(A);
157 
158    /* Now go through my local rows and set the matrix entries.
159       Each row has at most 5 entries. For example, if n=3:
160 
161       A = [M -I 0; -I M -I; 0 -I M]
162       M = [4 -1 0; -1 4 -1; 0 -1 4]
163 
164       Note that here we are setting one row at a time, though
165       one could set all the rows together (see the User's Manual).
166    */
167    {
168       int nnz;
169       double values[5];
170       int cols[5];
171 
172       for (i = ilower; i <= iupper; i++)
173       {
174          nnz = 0;
175 
176          /* The left identity block:position i-n */
177          if ((i-n)>=0)
178          {
179 	    cols[nnz] = i-n;
180 	    values[nnz] = -1.0;
181 	    nnz++;
182          }
183 
184          /* The left -1: position i-1 */
185          if (i%n)
186          {
187             cols[nnz] = i-1;
188             values[nnz] = -1.0;
189             nnz++;
190          }
191 
192          /* Set the diagonal: position i */
193          cols[nnz] = i;
194          values[nnz] = 4.0;
195          nnz++;
196 
197          /* The right -1: position i+1 */
198          if ((i+1)%n)
199          {
200             cols[nnz] = i+1;
201             values[nnz] = -1.0;
202             nnz++;
203          }
204 
205          /* The right identity block:position i+n */
206          if ((i+n)< N)
207          {
208             cols[nnz] = i+n;
209             values[nnz] = -1.0;
210             nnz++;
211          }
212 
213          /* Set the values for row i */
214          HYPRE_IJMatrixSetValues(A, 1, &nnz, &i, cols, values);
215       }
216    }
217 
218    /* Assemble after setting the coefficients */
219    HYPRE_IJMatrixAssemble(A);
220 
221    /* Note: for the testing of small problems, one may wish to read
222       in a matrix in IJ format (for the format, see the output files 
223       from the -print_system option).
224       In this case, one would use the following routine:  
225       HYPRE_IJMatrixRead( <filename>, MPI_COMM_WORLD,
226                           HYPRE_PARCSR, &A );
227       <filename>  = IJ.A.out to read in what has been printed out
228       by -print_system (processor numbers are omitted).
229       A call to HYPRE_IJMatrixRead is an *alternative* to the 
230       following sequence of HYPRE_IJMatrix calls: 
231       Create, SetObjectType, Initialize, SetValues, and Assemble
232    */                     
233 
234 
235    /* Get the parcsr matrix object to use */
236    HYPRE_IJMatrixGetObject(A, (void**) &parcsr_A);
237 
238 
239    /* Create the rhs and solution */
240    HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&b);
241    HYPRE_IJVectorSetObjectType(b, HYPRE_PARCSR);
242    HYPRE_IJVectorInitialize(b);
243 
244    HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&x);
245    HYPRE_IJVectorSetObjectType(x, HYPRE_PARCSR);
246    HYPRE_IJVectorInitialize(x);
247 
248    /* Set the rhs values to h^2 and the solution to zero */
249    {
250       double *rhs_values, *x_values;
251       int    *rows;
252 
253       rhs_values = calloc(local_size, sizeof(double));
254       x_values = calloc(local_size, sizeof(double));
255       rows = calloc(local_size, sizeof(int));
256 
257       for (i=0; i<local_size; i++)
258       {
259          rhs_values[i] = h2;
260          x_values[i] = 0.0;
261          rows[i] = ilower + i;
262       }
263 
264       HYPRE_IJVectorSetValues(b, local_size, rows, rhs_values);
265       HYPRE_IJVectorSetValues(x, local_size, rows, x_values);
266 
267       free(x_values);
268       free(rhs_values);
269       free(rows);
270    }
271 
272 
273    HYPRE_IJVectorAssemble(b);
274    /*  As with the matrix, for testing purposes, one may wish to read in a rhs:
275        HYPRE_IJVectorRead( <filename>, MPI_COMM_WORLD, 
276                                  HYPRE_PARCSR, &b ); 
277        as an alternative to the 
278        following sequence of HYPRE_IJVectors calls: 
279        Create, SetObjectType, Initialize, SetValues, and Assemble
280    */
281    HYPRE_IJVectorGetObject(b, (void **) &par_b);
282 
283    HYPRE_IJVectorAssemble(x);
284    HYPRE_IJVectorGetObject(x, (void **) &par_x);
285 
286  
287   /*  Print out the system  - files names will be IJ.out.A.XXXXX
288        and IJ.out.b.XXXXX, where XXXXX = processor id */
289    if (print_system)
290    {
291       HYPRE_IJMatrixPrint(A, "IJ.out.A");
292       HYPRE_IJVectorPrint(b, "IJ.out.b");
293    }
294 
295 
296    /* Choose a solver and solve the system */
297 
298    /* AMG */
299    if (solver_id == 0)
300    {
301       int num_iterations;
302       double final_res_norm;
303 
304       /* Create solver */
305       HYPRE_BoomerAMGCreate(&solver);
306 
307       /* Set some parameters (See Reference Manual for more parameters) */
308       HYPRE_BoomerAMGSetPrintLevel(solver, 3);  /* print solve info + parameters */
309       HYPRE_BoomerAMGSetCoarsenType(solver, 6); /* Falgout coarsening */
310       HYPRE_BoomerAMGSetRelaxType(solver, 3);   /* G-S/Jacobi hybrid relaxation */
311       HYPRE_BoomerAMGSetNumSweeps(solver, 1);   /* Sweeeps on each level */
312       HYPRE_BoomerAMGSetMaxLevels(solver, 20);  /* maximum number of levels */
313       HYPRE_BoomerAMGSetTol(solver, 1e-7);      /* conv. tolerance */
314 
315       /* Now setup and solve! */
316       HYPRE_BoomerAMGSetup(solver, parcsr_A, par_b, par_x);
317       HYPRE_BoomerAMGSolve(solver, parcsr_A, par_b, par_x);
318 
319       /* Run info - needed logging turned on */
320       HYPRE_BoomerAMGGetNumIterations(solver, &num_iterations);
321       HYPRE_BoomerAMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
322       if (myid == 0)
323       {
324          printf("\n");
325          printf("Iterations = %d\n", num_iterations);
326          printf("Final Relative Residual Norm = %e\n", final_res_norm);
327          printf("\n");
328       }
329 
330       /* Destroy solver */
331       HYPRE_BoomerAMGDestroy(solver);
332    }
333    /* PCG */
334    else if (solver_id == 50)
335    {
336       int num_iterations;
337       double final_res_norm;
338 
339       /* Create solver */
340       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
341 
342       /* Set some parameters (See Reference Manual for more parameters) */
343       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
344       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
345       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
346       HYPRE_PCGSetPrintLevel(solver, 2); /* prints out the iteration info */
347       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
348 
349       /* Now setup and solve! */
350       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
351       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
352 
353       /* Run info - needed logging turned on */
354       HYPRE_PCGGetNumIterations(solver, &num_iterations);
355       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &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 solver */
365       HYPRE_ParCSRPCGDestroy(solver);
366    }
367    /* PCG with AMG preconditioner */
368    else if (solver_id == 1)
369    {
370       int num_iterations;
371       double final_res_norm;
372 
373       /* Create solver */
374       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
375 
376       /* Set some parameters (See Reference Manual for more parameters) */
377       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
378       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
379       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
380       HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
381       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
382 
383       /* Now set up the AMG preconditioner and specify any parameters */
384       HYPRE_BoomerAMGCreate(&precond);
385       HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
386       HYPRE_BoomerAMGSetCoarsenType(precond, 6);
387       HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */ 
388       HYPRE_BoomerAMGSetNumSweeps(precond, 1);
389       HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
390       HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
391 
392       /* Set the PCG preconditioner */
393       HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
394                           (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup, precond);
395 
396       /* Now setup and solve! */
397       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
398       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
399 
400       /* Run info - needed logging turned on */
401       HYPRE_PCGGetNumIterations(solver, &num_iterations);
402       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
403       if (myid == 0)
404       {
405          printf("\n");
406          printf("Iterations = %d\n", num_iterations);
407          printf("Final Relative Residual Norm = %e\n", final_res_norm);
408          printf("\n");
409       }
410 
411       /* Destroy solver and preconditioner */
412       HYPRE_ParCSRPCGDestroy(solver);
413       HYPRE_BoomerAMGDestroy(precond);
414    }
415    /* PCG with Parasails Preconditioner */
416    else if (solver_id == 8)
417    {
418       int    num_iterations;
419       double final_res_norm;
420 
421       int      sai_max_levels = 1;
422       double   sai_threshold = 0.1;
423       double   sai_filter = 0.05;
424       int      sai_sym = 1;
425 
426       /* Create solver */
427       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
428 
429       /* Set some parameters (See Reference Manual for more parameters) */
430       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
431       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
432       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
433       HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
434       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
435 
436       /* Now set up the ParaSails preconditioner and specify any parameters */
437       HYPRE_ParaSailsCreate(MPI_COMM_WORLD, &precond);
438 
439       /* Set some parameters (See Reference Manual for more parameters) */
440       HYPRE_ParaSailsSetParams(precond, sai_threshold, sai_max_levels);
441       HYPRE_ParaSailsSetFilter(precond, sai_filter);
442       HYPRE_ParaSailsSetSym(precond, sai_sym);
443       HYPRE_ParaSailsSetLogging(precond, 3);
444 
445       /* Set the PCG preconditioner */
446       HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSolve,
447                           (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSetup, precond);
448 
449       /* Now setup and solve! */
450       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
451       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
452 
453 
454       /* Run info - needed logging turned on */
455       HYPRE_PCGGetNumIterations(solver, &num_iterations);
456       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
457       if (myid == 0)
458       {
459          printf("\n");
460          printf("Iterations = %d\n", num_iterations);
461          printf("Final Relative Residual Norm = %e\n", final_res_norm);
462          printf("\n");
463       }
464 
465       /* Destory solver and preconditioner */
466       HYPRE_ParCSRPCGDestroy(solver);
467       HYPRE_ParaSailsDestroy(precond);
468    }
469    else
470    {
471       if (myid ==0) printf("Invalid solver id specified.\n");
472    }
473 
474    /* Print the solution */
475    if (print_solution)
476       HYPRE_IJVectorPrint(x, "ij.out.x");
477 
478    /* Clean up */
479    HYPRE_IJMatrixDestroy(A);
480    HYPRE_IJVectorDestroy(b);
481    HYPRE_IJVectorDestroy(x);
482 
483    /* Finalize MPI*/
484    MPI_Finalize();
485 
486    return(0);
487 }


syntax highlighted by Code2HTML, v. 0.9.1