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 Laplacian problem with zero boundary
 11                  conditions on an n x n grid.  The number of unknowns is N=n^2.
 12                  The standard 5-point stencil is used, and we solve for the
 13                  interior nodes only.
 14 
 15                  This example solves the same problem as Example 3.  Available
 16                  solvers are AMG, PCG, and PCG with AMG or Parasails
 17                  preconditioners.  */
 18 
 19 #include <math.h>
 20 #include "_hypre_utilities.h"
 21 #include "HYPRE_krylov.h"
 22 #include "HYPRE.h"
 23 #include "HYPRE_parcsr_ls.h"
 24 
 25 #include "vis.c"
 26 
 27 int hypre_FlexGMRESModifyPCAMGExample(void *precond_data, int iterations,
 28                                       double rel_residual_norm);
 29 
 30 
 31 int main (int argc, char *argv[])
 32 {
 33    int i;
 34    int myid, num_procs;
 35    int N, n;
 36 
 37    int ilower, iupper;
 38    int local_size, extra;
 39 
 40    int solver_id;
 41    int vis, print_system;
 42 
 43    double h, h2;
 44 
 45    HYPRE_IJMatrix A;
 46    HYPRE_ParCSRMatrix parcsr_A;
 47    HYPRE_IJVector b;
 48    HYPRE_ParVector par_b;
 49    HYPRE_IJVector x;
 50    HYPRE_ParVector par_x;
 51 
 52    HYPRE_Solver solver, precond;
 53 
 54    /* Initialize MPI */
 55    MPI_Init(&argc, &argv);
 56    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 57    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 58 
 59    /* Default problem parameters */
 60    n = 33;
 61    solver_id = 0;
 62    vis = 0;
 63    print_system = 0;
 64 
 65 
 66    /* Parse command line */
 67    {
 68       int arg_index = 0;
 69       int print_usage = 0;
 70 
 71       while (arg_index < argc)
 72       {
 73          if ( strcmp(argv[arg_index], "-n") == 0 )
 74          {
 75             arg_index++;
 76             n = atoi(argv[arg_index++]);
 77          }
 78          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 79          {
 80             arg_index++;
 81             solver_id = atoi(argv[arg_index++]);
 82          }
 83          else if ( strcmp(argv[arg_index], "-vis") == 0 )
 84          {
 85             arg_index++;
 86             vis = 1;
 87          }
 88          else if ( strcmp(argv[arg_index], "-print_system") == 0 )
 89          {
 90             arg_index++;
 91             print_system = 1;
 92          }
 93          else if ( strcmp(argv[arg_index], "-help") == 0 )
 94          {
 95             print_usage = 1;
 96             break;
 97          }
 98          else
 99          {
100             arg_index++;
101          }
102       }
103 
104       if ((print_usage) && (myid == 0))
105       {
106          printf("\n");
107          printf("Usage: %s [<options>]\n", argv[0]);
108          printf("\n");
109          printf("  -n <n>              : problem size in each direction (default: 33)\n");
110          printf("  -solver <ID>        : solver ID\n");
111          printf("                        0  - AMG (default) \n");
112          printf("                        1  - AMG-PCG\n");
113          printf("                        8  - ParaSails-PCG\n");
114          printf("                        50 - PCG\n");
115          printf("                        61 - AMG-FlexGMRES\n");
116          printf("  -vis                : save the solution for GLVis visualization\n");
117          printf("  -print_system       : print the matrix and rhs\n");
118          printf("\n");
119       }
120 
121       if (print_usage)
122       {
123          MPI_Finalize();
124          return (0);
125       }
126    }
127 
128    /* Preliminaries: want at least one processor per row */
129    if (n*n < num_procs) n = sqrt(num_procs) + 1;
130    N = n*n; /* global number of rows */
131    h = 1.0/(n+1); /* mesh size*/
132    h2 = h*h;
133 
134    /* Each processor knows only of its own rows - the range is denoted by ilower
135       and upper.  Here we partition the rows. We account for the fact that
136       N may not divide evenly by the number of processors. */
137    local_size = N/num_procs;
138    extra = N - local_size*num_procs;
139 
140    ilower = local_size*myid;
141    ilower += hypre_min(myid, extra);
142 
143    iupper = local_size*(myid+1);
144    iupper += hypre_min(myid+1, extra);
145    iupper = iupper - 1;
146 
147    /* How many rows do I have? */
148    local_size = iupper - ilower + 1;
149 
150    /* Create the matrix.
151       Note that this is a square matrix, so we indicate the row partition
152       size twice (since number of rows = number of cols) */
153    HYPRE_IJMatrixCreate(MPI_COMM_WORLD, ilower, iupper, ilower, iupper, &A);
154 
155    /* Choose a parallel csr format storage (see the User's Manual) */
156    HYPRE_IJMatrixSetObjectType(A, HYPRE_PARCSR);
157 
158    /* Initialize before setting coefficients */
159    HYPRE_IJMatrixInitialize(A);
160 
161    /* Now go through my local rows and set the matrix entries.
162       Each row has at most 5 entries. For example, if n=3:
163 
164       A = [M -I 0; -I M -I; 0 -I M]
165       M = [4 -1 0; -1 4 -1; 0 -1 4]
166 
167       Note that here we are setting one row at a time, though
168       one could set all the rows together (see the User's Manual).
169    */
170    {
171       int nnz;
172       double values[5];
173       int cols[5];
174 
175       for (i = ilower; i <= iupper; i++)
176       {
177          nnz = 0;
178 
179          /* The left identity block:position i-n */
180          if ((i-n)>=0)
181          {
182             cols[nnz] = i-n;
183             values[nnz] = -1.0;
184             nnz++;
185          }
186 
187          /* The left -1: position i-1 */
188          if (i%n)
189          {
190             cols[nnz] = i-1;
191             values[nnz] = -1.0;
192             nnz++;
193          }
194 
195          /* Set the diagonal: position i */
196          cols[nnz] = i;
197          values[nnz] = 4.0;
198          nnz++;
199 
200          /* The right -1: position i+1 */
201          if ((i+1)%n)
202          {
203             cols[nnz] = i+1;
204             values[nnz] = -1.0;
205             nnz++;
206          }
207 
208          /* The right identity block:position i+n */
209          if ((i+n)< N)
210          {
211             cols[nnz] = i+n;
212             values[nnz] = -1.0;
213             nnz++;
214          }
215 
216          /* Set the values for row i */
217          HYPRE_IJMatrixSetValues(A, 1, &nnz, &i, cols, values);
218       }
219    }
220 
221    /* Assemble after setting the coefficients */
222    HYPRE_IJMatrixAssemble(A);
223 
224    /* Note: for the testing of small problems, one may wish to read
225       in a matrix in IJ format (for the format, see the output files
226       from the -print_system option).
227       In this case, one would use the following routine:
228       HYPRE_IJMatrixRead( <filename>, MPI_COMM_WORLD,
229                           HYPRE_PARCSR, &A );
230       <filename>  = IJ.A.out to read in what has been printed out
231       by -print_system (processor numbers are omitted).
232       A call to HYPRE_IJMatrixRead is an *alternative* to the
233       following sequence of HYPRE_IJMatrix calls:
234       Create, SetObjectType, Initialize, SetValues, and Assemble
235    */
236 
237 
238    /* Get the parcsr matrix object to use */
239    HYPRE_IJMatrixGetObject(A, (void**) &parcsr_A);
240 
241 
242    /* Create the rhs and solution */
243    HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&b);
244    HYPRE_IJVectorSetObjectType(b, HYPRE_PARCSR);
245    HYPRE_IJVectorInitialize(b);
246 
247    HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&x);
248    HYPRE_IJVectorSetObjectType(x, HYPRE_PARCSR);
249    HYPRE_IJVectorInitialize(x);
250 
251    /* Set the rhs values to h^2 and the solution to zero */
252    {
253       double *rhs_values, *x_values;
254       int    *rows;
255 
256       rhs_values = calloc(local_size, sizeof(double));
257       x_values = calloc(local_size, sizeof(double));
258       rows = calloc(local_size, sizeof(int));
259 
260       for (i=0; i<local_size; i++)
261       {
262          rhs_values[i] = h2;
263          x_values[i] = 0.0;
264          rows[i] = ilower + i;
265       }
266 
267       HYPRE_IJVectorSetValues(b, local_size, rows, rhs_values);
268       HYPRE_IJVectorSetValues(x, local_size, rows, x_values);
269 
270       free(x_values);
271       free(rhs_values);
272       free(rows);
273    }
274 
275 
276    HYPRE_IJVectorAssemble(b);
277    /*  As with the matrix, for testing purposes, one may wish to read in a rhs:
278        HYPRE_IJVectorRead( <filename>, MPI_COMM_WORLD,
279                                  HYPRE_PARCSR, &b );
280        as an alternative to the
281        following sequence of HYPRE_IJVectors calls:
282        Create, SetObjectType, Initialize, SetValues, and Assemble
283    */
284    HYPRE_IJVectorGetObject(b, (void **) &par_b);
285 
286    HYPRE_IJVectorAssemble(x);
287    HYPRE_IJVectorGetObject(x, (void **) &par_x);
288 
289 
290   /*  Print out the system  - files names will be IJ.out.A.XXXXX
291        and IJ.out.b.XXXXX, where XXXXX = processor id */
292    if (print_system)
293    {
294       HYPRE_IJMatrixPrint(A, "IJ.out.A");
295       HYPRE_IJVectorPrint(b, "IJ.out.b");
296    }
297 
298 
299    /* Choose a solver and solve the system */
300 
301    /* AMG */
302    if (solver_id == 0)
303    {
304       int num_iterations;
305       double final_res_norm;
306 
307       /* Create solver */
308       HYPRE_BoomerAMGCreate(&solver);
309 
310       /* Set some parameters (See Reference Manual for more parameters) */
311       HYPRE_BoomerAMGSetPrintLevel(solver, 3);  /* print solve info + parameters */
312       HYPRE_BoomerAMGSetCoarsenType(solver, 6); /* Falgout coarsening */
313       HYPRE_BoomerAMGSetRelaxType(solver, 3);   /* G-S/Jacobi hybrid relaxation */
314       HYPRE_BoomerAMGSetNumSweeps(solver, 1);   /* Sweeeps on each level */
315       HYPRE_BoomerAMGSetMaxLevels(solver, 20);  /* maximum number of levels */
316       HYPRE_BoomerAMGSetTol(solver, 1e-7);      /* conv. tolerance */
317 
318       /* Now setup and solve! */
319       HYPRE_BoomerAMGSetup(solver, parcsr_A, par_b, par_x);
320       HYPRE_BoomerAMGSolve(solver, parcsr_A, par_b, par_x);
321 
322       /* Run info - needed logging turned on */
323       HYPRE_BoomerAMGGetNumIterations(solver, &num_iterations);
324       HYPRE_BoomerAMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
325       if (myid == 0)
326       {
327          printf("\n");
328          printf("Iterations = %d\n", num_iterations);
329          printf("Final Relative Residual Norm = %e\n", final_res_norm);
330          printf("\n");
331       }
332 
333       /* Destroy solver */
334       HYPRE_BoomerAMGDestroy(solver);
335    }
336    /* PCG */
337    else if (solver_id == 50)
338    {
339       int num_iterations;
340       double final_res_norm;
341 
342       /* Create solver */
343       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
344 
345       /* Set some parameters (See Reference Manual for more parameters) */
346       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
347       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
348       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
349       HYPRE_PCGSetPrintLevel(solver, 2); /* prints out the iteration info */
350       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
351 
352       /* Now setup and solve! */
353       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
354       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
355 
356       /* Run info - needed logging turned on */
357       HYPRE_PCGGetNumIterations(solver, &num_iterations);
358       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
359       if (myid == 0)
360       {
361          printf("\n");
362          printf("Iterations = %d\n", num_iterations);
363          printf("Final Relative Residual Norm = %e\n", final_res_norm);
364          printf("\n");
365       }
366 
367       /* Destroy solver */
368       HYPRE_ParCSRPCGDestroy(solver);
369    }
370    /* PCG with AMG preconditioner */
371    else if (solver_id == 1)
372    {
373       int num_iterations;
374       double final_res_norm;
375 
376       /* Create solver */
377       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
378 
379       /* Set some parameters (See Reference Manual for more parameters) */
380       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
381       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
382       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
383       HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
384       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
385 
386       /* Now set up the AMG preconditioner and specify any parameters */
387       HYPRE_BoomerAMGCreate(&precond);
388       HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
389       HYPRE_BoomerAMGSetCoarsenType(precond, 6);
390       HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
391       HYPRE_BoomerAMGSetNumSweeps(precond, 1);
392       HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
393       HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
394 
395       /* Set the PCG preconditioner */
396       HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
397                           (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup, precond);
398 
399       /* Now setup and solve! */
400       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
401       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
402 
403       /* Run info - needed logging turned on */
404       HYPRE_PCGGetNumIterations(solver, &num_iterations);
405       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
406       if (myid == 0)
407       {
408          printf("\n");
409          printf("Iterations = %d\n", num_iterations);
410          printf("Final Relative Residual Norm = %e\n", final_res_norm);
411          printf("\n");
412       }
413 
414       /* Destroy solver and preconditioner */
415       HYPRE_ParCSRPCGDestroy(solver);
416       HYPRE_BoomerAMGDestroy(precond);
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       HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
431 
432       /* Set some parameters (See Reference Manual for more parameters) */
433       HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
434       HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
435       HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
436       HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
437       HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
438 
439       /* Now set up the ParaSails preconditioner and specify any parameters */
440       HYPRE_ParaSailsCreate(MPI_COMM_WORLD, &precond);
441 
442       /* Set some parameters (See Reference Manual for more parameters) */
443       HYPRE_ParaSailsSetParams(precond, sai_threshold, sai_max_levels);
444       HYPRE_ParaSailsSetFilter(precond, sai_filter);
445       HYPRE_ParaSailsSetSym(precond, sai_sym);
446       HYPRE_ParaSailsSetLogging(precond, 3);
447 
448       /* Set the PCG preconditioner */
449       HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSolve,
450                           (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSetup, precond);
451 
452       /* Now setup and solve! */
453       HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
454       HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
455 
456 
457       /* Run info - needed logging turned on */
458       HYPRE_PCGGetNumIterations(solver, &num_iterations);
459       HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
460       if (myid == 0)
461       {
462          printf("\n");
463          printf("Iterations = %d\n", num_iterations);
464          printf("Final Relative Residual Norm = %e\n", final_res_norm);
465          printf("\n");
466       }
467 
468       /* Destory solver and preconditioner */
469       HYPRE_ParCSRPCGDestroy(solver);
470       HYPRE_ParaSailsDestroy(precond);
471    }
472    /* Flexible GMRES with  AMG Preconditioner */
473    else if (solver_id == 61)
474    {
475       int    num_iterations;
476       double final_res_norm;
477       int    restart = 30;
478       int    modify = 1;
479 
480 
481       /* Create solver */
482       HYPRE_ParCSRFlexGMRESCreate(MPI_COMM_WORLD, &solver);
483 
484       /* Set some parameters (See Reference Manual for more parameters) */
485       HYPRE_FlexGMRESSetKDim(solver, restart);
486       HYPRE_FlexGMRESSetMaxIter(solver, 1000); /* max iterations */
487       HYPRE_FlexGMRESSetTol(solver, 1e-7); /* conv. tolerance */
488       HYPRE_FlexGMRESSetPrintLevel(solver, 2); /* print solve info */
489       HYPRE_FlexGMRESSetLogging(solver, 1); /* needed to get run info later */
490 
491 
492       /* Now set up the AMG preconditioner and specify any parameters */
493       HYPRE_BoomerAMGCreate(&precond);
494       HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
495       HYPRE_BoomerAMGSetCoarsenType(precond, 6);
496       HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
497       HYPRE_BoomerAMGSetNumSweeps(precond, 1);
498       HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
499       HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
500 
501       /* Set the FlexGMRES preconditioner */
502       HYPRE_FlexGMRESSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
503                           (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup, precond);
504 
505 
506       if (modify)
507       /* this is an optional call  - if you don't call it, hypre_FlexGMRESModifyPCDefault
508          is used - which does nothing.  Otherwise, you can define your own, similar to
509          the one used here */
510          HYPRE_FlexGMRESSetModifyPC( solver,
511                                      (HYPRE_PtrToModifyPCFcn) hypre_FlexGMRESModifyPCAMGExample);
512 
513 
514       /* Now setup and solve! */
515       HYPRE_ParCSRFlexGMRESSetup(solver, parcsr_A, par_b, par_x);
516       HYPRE_ParCSRFlexGMRESSolve(solver, parcsr_A, par_b, par_x);
517 
518       /* Run info - needed logging turned on */
519       HYPRE_FlexGMRESGetNumIterations(solver, &num_iterations);
520       HYPRE_FlexGMRESGetFinalRelativeResidualNorm(solver, &final_res_norm);
521       if (myid == 0)
522       {
523          printf("\n");
524          printf("Iterations = %d\n", num_iterations);
525          printf("Final Relative Residual Norm = %e\n", final_res_norm);
526          printf("\n");
527       }
528 
529       /* Destory solver and preconditioner */
530       HYPRE_ParCSRFlexGMRESDestroy(solver);
531       HYPRE_BoomerAMGDestroy(precond);
532 
533    }
534    else
535    {
536       if (myid ==0) printf("Invalid solver id specified.\n");
537    }
538 
539    /* Save the solution for GLVis visualization, see vis/glvis-ex5.sh */
540    if (vis)
541    {
542       FILE *file;
543       char filename[255];
544 
545       int nvalues = local_size;
546       int *rows = calloc(nvalues, sizeof(int));
547       double *values = calloc(nvalues, sizeof(double));
548 
549       for (i = 0; i < nvalues; i++)
550          rows[i] = ilower + i;
551 
552       /* get the local solution */
553       HYPRE_IJVectorGetValues(x, nvalues, rows, values);
554 
555       sprintf(filename, "%s.%06d", "vis/ex5.sol", myid);
556       if ((file = fopen(filename, "w")) == NULL)
557       {
558          printf("Error: can't open output file %s\n", filename);
559          MPI_Finalize();
560          exit(1);
561       }
562 
563       /* save solution */
564       for (i = 0; i < nvalues; i++)
565          fprintf(file, "%.14e\n", values[i]);
566 
567       fflush(file);
568       fclose(file);
569 
570       free(rows);
571       free(values);
572 
573       /* save global finite element mesh */
574       if (myid == 0)
575          GLVis_PrintGlobalSquareMesh("vis/ex5.mesh", n-1);
576    }
577 
578    /* Clean up */
579    HYPRE_IJMatrixDestroy(A);
580    HYPRE_IJVectorDestroy(b);
581    HYPRE_IJVectorDestroy(x);
582 
583    /* Finalize MPI*/
584    MPI_Finalize();
585 
586    return(0);
587 }
588 
589 /*--------------------------------------------------------------------------
590    hypre_FlexGMRESModifyPCAMGExample -
591 
592     This is an example (not recommended)
593    of how we can modify things about AMG that
594    affect the solve phase based on how FlexGMRES is doing...For
595    another preconditioner it may make sense to modify the tolerance..
596 
597  *--------------------------------------------------------------------------*/
598 
599 int hypre_FlexGMRESModifyPCAMGExample(void *precond_data, int iterations,
600                                    double rel_residual_norm)
601 {
602 
603 
604    if (rel_residual_norm > .1)
605    {
606       HYPRE_BoomerAMGSetNumSweeps(precond_data, 10);
607    }
608    else
609    {
610       HYPRE_BoomerAMGSetNumSweeps(precond_data, 1);
611    }
612 
613 
614    return 0;
615 }


syntax highlighted by Code2HTML, v. 0.9.1