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


syntax highlighted by Code2HTML, v. 0.9.1