download the original source code.
  1 /*
  2    Example 9
  3 
  4    Interface:      Semi-Structured interface (SStruct)
  5 
  6    Compile with:   make ex9
  7 
  8    Sample run:     mpirun -np 16 ex9 -n 33 -solver 0 -v 1 1
  9 
 10    To see options: ex9 -help
 11 
 12    Description:    This code solves a system corresponding to a discretization
 13                    of the biharmonic problem treated as a system of equations
 14                    on the unit square.  Specifically, instead of solving
 15                    Delta^2(u) = f with zero boundary conditions for u and
 16                    Delta(u), we solve the system A x = b, where
 17 
 18                    A = [ Delta -I ; 0 Delta], x = [ u ; v] and b = [ 0 ; f]
 19 
 20                    The corresponding boundary conditions are u = 0 and v = 0.
 21 
 22                    The domain is split into an N x N processor grid.  Thus, the
 23                    given number of processors should be a perfect square.
 24                    Each processor's piece of the grid has n x n cells with n x n
 25                    nodes. We use cell-centered variables, and, therefore, the
 26                    nodes are not shared. Note that we have two variables, u and
 27                    v, and need only one part to describe the domain. We use the
 28                    standard 5-point stencil to discretize the Laplace operators.
 29                    The boundary conditions are incorporated as in Example 3.
 30 
 31                    We recommend viewing Examples 3, 6 and 7 before this example.
 32 */
 33 
 34 #include <math.h>
 35 #include "_hypre_utilities.h"
 36 #include "HYPRE_sstruct_ls.h"
 37 #include "HYPRE_krylov.h"
 38 
 39 #include "vis.c"
 40 
 41 int main (int argc, char *argv[])
 42 {
 43    int i, j;
 44 
 45    int myid, num_procs;
 46 
 47    int n, N, pi, pj;
 48    double h, h2;
 49    int ilower[2], iupper[2];
 50 
 51    int solver_id;
 52    int n_pre, n_post;
 53 
 54    int vis;
 55    int object_type;
 56 
 57    HYPRE_SStructGrid     grid;
 58    HYPRE_SStructGraph    graph;
 59    HYPRE_SStructStencil  stencil_v;
 60    HYPRE_SStructStencil  stencil_u;
 61    HYPRE_SStructMatrix   A;
 62    HYPRE_SStructVector   b;
 63    HYPRE_SStructVector   x;
 64 
 65    /* sstruct solvers */
 66    HYPRE_SStructSolver   solver;
 67    HYPRE_SStructSolver   precond;
 68 
 69    /* parcsr solvers */
 70    HYPRE_Solver          par_solver;
 71    HYPRE_Solver          par_precond;
 72 
 73    /* Initialize MPI */
 74    MPI_Init(&argc, &argv);
 75    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 76    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 77 
 78    /* Set defaults */
 79    n = 33;
 80    solver_id = 0;
 81    n_pre  = 1;
 82    n_post = 1;
 83    vis = 0;
 84 
 85    /* Parse command line */
 86    {
 87       int arg_index = 0;
 88       int print_usage = 0;
 89 
 90       while (arg_index < argc)
 91       {
 92          if ( strcmp(argv[arg_index], "-n") == 0 )
 93          {
 94             arg_index++;
 95             n = atoi(argv[arg_index++]);
 96          }
 97          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 98          {
 99             arg_index++;
100             solver_id = atoi(argv[arg_index++]);
101          }
102          else if ( strcmp(argv[arg_index], "-v") == 0 )
103          {
104             arg_index++;
105             n_pre = atoi(argv[arg_index++]);
106             n_post = atoi(argv[arg_index++]);
107          }
108          else if ( strcmp(argv[arg_index], "-vis") == 0 )
109          {
110             arg_index++;
111             vis = 1;
112          }
113          else if ( strcmp(argv[arg_index], "-help") == 0 )
114          {
115             print_usage = 1;
116             break;
117          }
118          else
119          {
120             arg_index++;
121          }
122       }
123 
124       if ((print_usage) && (myid == 0))
125       {
126          printf("\n");
127          printf("Usage: %s [<options>]\n", argv[0]);
128          printf("\n");
129          printf("  -n <n>              : problem size per processor (default: 33)\n");
130          printf("  -solver <ID>        : solver ID\n");
131          printf("                        0  - GMRES with sysPFMG precond (default)\n");
132          printf("                        1  - sysPFMG\n");
133          printf("                        2  - GMRES with AMG precond\n");
134          printf("                        3  - AMG\n");
135          printf("  -v <n_pre> <n_post> : number of pre and post relaxations for SysPFMG (default: 1 1)\n");
136          printf("  -vis                : save the solution for GLVis visualization\n");
137          printf("\n");
138       }
139 
140       if (print_usage)
141       {
142          MPI_Finalize();
143          return (0);
144       }
145    }
146 
147    /* Figure out the processor grid (N x N).  The local problem
148       size for the interior nodes is indicated by n (n x n).
149       pi and pj indicate position in the processor grid. */
150    N  = sqrt(num_procs);
151    h  = 1.0 / (N*n+1); /* note that when calculating h we must
152                           remember to count the boundary nodes */
153    h2 = h*h;
154    pj = myid / N;
155    pi = myid - pj*N;
156 
157    /* Figure out the extents of each processor's piece of the grid. */
158    ilower[0] = pi*n;
159    ilower[1] = pj*n;
160 
161    iupper[0] = ilower[0] + n-1;
162    iupper[1] = ilower[1] + n-1;
163 
164    /* 1. Set up a grid - we have one part and two variables */
165    {
166       int nparts = 1;
167       int part = 0;
168       int ndim = 2;
169 
170       /* Create an empty 2D grid object */
171       HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &grid);
172 
173       /* Add a new box to the grid */
174       HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
175 
176       /* Set the variable type and number of variables on each part.*/
177       {
178          int i;
179          int nvars = 2;
180          HYPRE_SStructVariable vartypes[2] = {HYPRE_SSTRUCT_VARIABLE_CELL,
181                                               HYPRE_SSTRUCT_VARIABLE_CELL };
182 
183          for (i = 0; i< nparts; i++)
184             HYPRE_SStructGridSetVariables(grid, i, nvars, vartypes);
185       }
186 
187       /* This is a collective call finalizing the grid assembly.
188          The grid is now ``ready to be used'' */
189       HYPRE_SStructGridAssemble(grid);
190    }
191 
192    /* 2. Define the discretization stencils */
193    {
194       int entry;
195       int stencil_size;
196       int var;
197       int ndim = 2;
198 
199       /* Stencil object for variable u (labeled as variable 0) */
200       {
201          int offsets[6][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}, {0,0}};
202          stencil_size = 6;
203 
204          HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil_u);
205 
206          /* The first 5 entries are for the u-u connections */
207          var = 0; /* connect to variable 0 */
208          for (entry = 0; entry < stencil_size-1 ; entry++)
209             HYPRE_SStructStencilSetEntry(stencil_u, entry, offsets[entry], var);
210 
211          /* The last entry is for the u-v connection */
212          var = 1;  /* connect to variable 1 */
213          entry = 5;
214          HYPRE_SStructStencilSetEntry(stencil_u, entry, offsets[entry], var);
215       }
216 
217       /* Stencil object for variable v  (variable 1) */
218       {
219          int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
220          stencil_size = 5;
221 
222          HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil_v);
223 
224          /* These are all v-v connections */
225          var = 1; /* Connect to variable 1 */
226          for (entry = 0; entry < stencil_size; entry++)
227             HYPRE_SStructStencilSetEntry(stencil_v, entry, offsets[entry], var);
228       }
229    }
230 
231    /* 3. Set up the Graph  - this determines the non-zero structure
232       of the matrix and allows non-stencil relationships between the parts. */
233    {
234       int var;
235       int part = 0;
236 
237       /* Create the graph object */
238       HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
239 
240       /* See MatrixSetObjectType below */
241       if (solver_id > 1 && solver_id < 4)
242       {
243          object_type = HYPRE_PARCSR;
244       }
245       else
246       {
247          object_type = HYPRE_SSTRUCT;
248       }
249       HYPRE_SStructGraphSetObjectType(graph, object_type);
250 
251       /* Assign the u-stencil we created to variable u (variable 0) */
252       var = 0;
253       HYPRE_SStructGraphSetStencil(graph, part, var, stencil_u);
254 
255       /* Assign the v-stencil we created to variable v (variable 1) */
256       var = 1;
257       HYPRE_SStructGraphSetStencil(graph, part, var, stencil_v);
258 
259       /* Assemble the graph */
260       HYPRE_SStructGraphAssemble(graph);
261    }
262 
263    /* 4. Set up the SStruct Matrix */
264    {
265       int nentries;
266       int nvalues;
267       int var;
268       int part = 0;
269 
270       /* Create an empty matrix object */
271       HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
272 
273       /* Set the object type (by default HYPRE_SSTRUCT). This determines the
274          data structure used to store the matrix.  If you want to use
275          unstructured solvers, e.g. BoomerAMG, the object type should be
276          HYPRE_PARCSR. If the problem is purely structured (with one part), you
277          may want to use HYPRE_STRUCT to access the structured solvers.  */
278       HYPRE_SStructMatrixSetObjectType(A, object_type);
279 
280       /* Indicate that the matrix coefficients are ready to be set */
281       HYPRE_SStructMatrixInitialize(A);
282 
283       /* Each processor must set the stencil values for their boxes on each part.
284          In this example, we only set stencil entries and therefore use
285          HYPRE_SStructMatrixSetBoxValues.  If we need to set non-stencil entries,
286          we have to use HYPRE_SStructMatrixSetValues. */
287 
288       /* First set the u-stencil entries.  Note that
289          HYPRE_SStructMatrixSetBoxValues can only set values corresponding
290          to stencil entries for the same variable. Therefore, we must set the
291          entries for each variable within a stencil with separate function calls.
292          For example, below the u-u connections and u-v connections are handled
293          in separate calls.  */
294       {
295          int     i, j;
296          double *u_values;
297          int     u_v_indices[1] = {5};
298          int     u_u_indices[5] = {0, 1, 2, 3, 4};
299 
300          var = 0; /* Set values for the u connections */
301 
302          /*  First the u-u connections */
303          nentries = 5;
304          nvalues = nentries*n*n;
305          u_values = (double*) calloc(nvalues, sizeof(double));
306 
307          for (i = 0; i < nvalues; i += nentries)
308          {
309             u_values[i] = 4.0;
310             for (j = 1; j < nentries; j++)
311                u_values[i+j] = -1.0;
312          }
313 
314          HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
315                                          var, nentries,
316                                          u_u_indices, u_values);
317          free(u_values);
318 
319          /* Next the u-v connections */
320          nentries = 1;
321          nvalues = nentries*n*n;
322          u_values = (double*) calloc(nvalues, sizeof(double));
323 
324          for (i = 0; i < nvalues; i++)
325          {
326             u_values[i] = -h2;
327          }
328 
329          HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
330                                          var, nentries,
331                                          u_v_indices, u_values);
332 
333          free(u_values);
334       }
335 
336       /*  Now set the v-stencil entries */
337       {
338          int     i, j;
339          double *v_values;
340          int     v_v_indices[5] = {0, 1, 2, 3, 4};
341 
342          var = 1; /* the v connections */
343 
344          /* the v-v connections */
345          nentries = 5;
346          nvalues = nentries*n*n;
347          v_values = (double*) calloc(nvalues, sizeof(double));
348 
349          for (i = 0; i < nvalues; i += nentries)
350          {
351             v_values[i] = 4.0;
352             for (j = 1; j < nentries; j++)
353                v_values[i+j] = -1.0;
354          }
355 
356          HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
357                                          var, nentries,
358                                          v_v_indices, v_values);
359 
360          free(v_values);
361 
362          /* There are no v-u connections to set */
363       }
364    }
365 
366    /* 5. Incorporate the zero boundary conditions: go along each edge of
367          the domain and set the stencil entry that reaches to the boundary
368          to zero.*/
369    {
370       int bc_ilower[2];
371       int bc_iupper[2];
372       int nentries = 1;
373       int nvalues  = nentries*n; /*  number of stencil entries times the length
374                                      of one side of my grid box */
375       int var;
376       double *values;
377       int stencil_indices[1];
378 
379       int part = 0;
380 
381       values = (double*) calloc(nvalues, sizeof(double));
382       for (j = 0; j < nvalues; j++)
383             values[j] = 0.0;
384 
385       /* Recall: pi and pj describe position in the processor grid */
386       if (pj == 0)
387       {
388          /* Bottom row of grid points */
389          bc_ilower[0] = pi*n;
390          bc_ilower[1] = pj*n;
391 
392          bc_iupper[0] = bc_ilower[0] + n-1;
393          bc_iupper[1] = bc_ilower[1];
394 
395          stencil_indices[0] = 3;
396 
397          /* Need to do this for u and for v */
398          var = 0;
399          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
400                                          var, nentries,
401                                          stencil_indices, values);
402 
403          var = 1;
404          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
405                                          var, nentries,
406                                          stencil_indices, values);
407       }
408 
409       if (pj == N-1)
410       {
411          /* upper row of grid points */
412          bc_ilower[0] = pi*n;
413          bc_ilower[1] = pj*n + n-1;
414 
415          bc_iupper[0] = bc_ilower[0] + n-1;
416          bc_iupper[1] = bc_ilower[1];
417 
418          stencil_indices[0] = 4;
419 
420          /* Need to do this for u and for v */
421          var = 0;
422          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
423                                          var, nentries,
424                                          stencil_indices, values);
425 
426          var = 1;
427          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
428                                          var, nentries,
429                                          stencil_indices, values);
430 
431       }
432 
433       if (pi == 0)
434       {
435          /* Left row of grid points */
436          bc_ilower[0] = pi*n;
437          bc_ilower[1] = pj*n;
438 
439          bc_iupper[0] = bc_ilower[0];
440          bc_iupper[1] = bc_ilower[1] + n-1;
441 
442          stencil_indices[0] = 1;
443 
444          /* Need to do this for u and for v */
445          var = 0;
446          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
447                                          var, nentries,
448                                          stencil_indices, values);
449 
450          var = 1;
451          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
452                                          var, nentries,
453                                          stencil_indices, values);
454       }
455 
456       if (pi == N-1)
457       {
458          /* Right row of grid points */
459          bc_ilower[0] = pi*n + n-1;
460          bc_ilower[1] = pj*n;
461 
462          bc_iupper[0] = bc_ilower[0];
463          bc_iupper[1] = bc_ilower[1] + n-1;
464 
465          stencil_indices[0] = 2;
466 
467          /* Need to do this for u and for v */
468          var = 0;
469          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
470                                          var, nentries,
471                                          stencil_indices, values);
472 
473          var = 1;
474          HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper,
475                                          var, nentries,
476                                          stencil_indices, values);
477       }
478 
479       free(values);
480    }
481 
482    /* This is a collective call finalizing the matrix assembly.
483       The matrix is now ``ready to be used'' */
484    HYPRE_SStructMatrixAssemble(A);
485 
486    /* 5. Set up SStruct Vectors for b and x */
487    {
488       int    nvalues = n*n;
489       double *values;
490       int part = 0;
491       int var;
492 
493       values = (double*) calloc(nvalues, sizeof(double));
494 
495       /* Create an empty vector object */
496       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
497       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
498 
499       /* Set the object type for the vectors
500          to be the same as was already set for the matrix */
501       HYPRE_SStructVectorSetObjectType(b, object_type);
502       HYPRE_SStructVectorSetObjectType(x, object_type);
503 
504       /* Indicate that the vector coefficients are ready to be set */
505       HYPRE_SStructVectorInitialize(b);
506       HYPRE_SStructVectorInitialize(x);
507 
508       /* Set the values for b */
509       for (i = 0; i < nvalues; i ++)
510          values[i] = h2;
511       var = 1;
512       HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
513 
514       for (i = 0; i < nvalues; i ++)
515          values[i] = 0.0;
516       var = 0;
517       HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
518 
519       /* Set the values for the initial guess */
520       var = 0;
521       HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
522 
523       var = 1;
524       HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
525 
526       free(values);
527 
528       /* This is a collective call finalizing the vector assembly.
529          The vector is now ``ready to be used'' */
530       HYPRE_SStructVectorAssemble(b);
531       HYPRE_SStructVectorAssemble(x);
532    }
533 
534    /* 6. Set up and use a solver
535       (Solver options can be found in the Reference Manual.) */
536    {
537       double final_res_norm;
538       int its;
539 
540       HYPRE_ParCSRMatrix    par_A;
541       HYPRE_ParVector       par_b;
542       HYPRE_ParVector       par_x;
543 
544       /* If we are using a parcsr solver, we need to get the object for the
545          matrix and vectors. */
546       if (object_type == HYPRE_PARCSR)
547       {
548          HYPRE_SStructMatrixGetObject(A, (void **) &par_A);
549          HYPRE_SStructVectorGetObject(b, (void **) &par_b);
550          HYPRE_SStructVectorGetObject(x, (void **) &par_x);
551       }
552 
553       if (solver_id ==0 ) /* GMRES with SysPFMG - the default*/
554       {
555          HYPRE_SStructGMRESCreate(MPI_COMM_WORLD, &solver);
556 
557          /* GMRES parameters */
558          HYPRE_SStructGMRESSetMaxIter(solver, 50 );
559          HYPRE_SStructGMRESSetTol(solver, 1.0e-06 );
560          HYPRE_SStructGMRESSetPrintLevel(solver, 2 ); /* print each GMRES
561                                                          iteration */
562          HYPRE_SStructGMRESSetLogging(solver, 1);
563 
564          /* use SysPFMG as precondititioner */
565          HYPRE_SStructSysPFMGCreate(MPI_COMM_WORLD, &precond);
566 
567          /* Set sysPFMG parameters */
568          HYPRE_SStructSysPFMGSetTol(precond, 0.0);
569          HYPRE_SStructSysPFMGSetMaxIter(precond, 1);
570          HYPRE_SStructSysPFMGSetNumPreRelax(precond, n_pre);
571          HYPRE_SStructSysPFMGSetNumPostRelax(precond, n_post);
572          HYPRE_SStructSysPFMGSetPrintLevel(precond, 0);
573          HYPRE_SStructSysPFMGSetZeroGuess(precond);
574 
575          /* Set the preconditioner*/
576          HYPRE_SStructGMRESSetPrecond(solver, HYPRE_SStructSysPFMGSolve,
577                                       HYPRE_SStructSysPFMGSetup, precond);
578          /* do the setup */
579          HYPRE_SStructGMRESSetup(solver, A, b, x);
580 
581          /* do the solve */
582          HYPRE_SStructGMRESSolve(solver, A, b, x);
583 
584          /* get some info */
585          HYPRE_SStructGMRESGetFinalRelativeResidualNorm(solver,
586                                                         &final_res_norm);
587          HYPRE_SStructGMRESGetNumIterations(solver, &its);
588 
589          /* clean up */
590          HYPRE_SStructGMRESDestroy(solver);
591       }
592       else if (solver_id == 1) /* SysPFMG */
593       {
594          HYPRE_SStructSysPFMGCreate(MPI_COMM_WORLD, &solver);
595 
596          /* Set sysPFMG parameters */
597          HYPRE_SStructSysPFMGSetTol(solver, 1.0e-6);
598          HYPRE_SStructSysPFMGSetMaxIter(solver, 50);
599          HYPRE_SStructSysPFMGSetNumPreRelax(solver, n_pre);
600          HYPRE_SStructSysPFMGSetNumPostRelax(solver, n_post);
601          HYPRE_SStructSysPFMGSetPrintLevel(solver, 0);
602          HYPRE_SStructSysPFMGSetLogging(solver, 1);
603 
604          /* do the setup */
605          HYPRE_SStructSysPFMGSetup(solver, A, b, x);
606 
607          /* do the solve */
608          HYPRE_SStructSysPFMGSolve(solver, A, b, x);
609 
610          /* get some info */
611          HYPRE_SStructSysPFMGGetFinalRelativeResidualNorm(solver,
612                                                           &final_res_norm);
613          HYPRE_SStructSysPFMGGetNumIterations(solver, &its);
614 
615          /* clean up */
616          HYPRE_SStructSysPFMGDestroy(solver);
617       }
618       else if (solver_id == 2) /* GMRES with AMG */
619       {
620          HYPRE_ParCSRGMRESCreate(MPI_COMM_WORLD, &par_solver);
621 
622          /* set the GMRES paramaters */
623          HYPRE_GMRESSetKDim(par_solver, 5);
624          HYPRE_GMRESSetMaxIter(par_solver, 100);
625          HYPRE_GMRESSetTol(par_solver, 1.0e-06);
626          HYPRE_GMRESSetPrintLevel(par_solver, 2);
627          HYPRE_GMRESSetLogging(par_solver, 1);
628 
629          /* use BoomerAMG as preconditioner */
630          HYPRE_BoomerAMGCreate(&par_precond);
631          HYPRE_BoomerAMGSetCoarsenType(par_precond, 6);
632          HYPRE_BoomerAMGSetStrongThreshold(par_precond, 0.25);
633          HYPRE_BoomerAMGSetTol(par_precond, 0.0);
634          HYPRE_BoomerAMGSetPrintLevel(par_precond, 1);
635          HYPRE_BoomerAMGSetPrintFileName(par_precond, "ex9.out.log");
636          HYPRE_BoomerAMGSetMaxIter(par_precond, 1);
637 
638          /* set the preconditioner */
639          HYPRE_ParCSRGMRESSetPrecond(par_solver,
640                                      HYPRE_BoomerAMGSolve,
641                                      HYPRE_BoomerAMGSetup,
642                                      par_precond);
643 
644          /* do the setup */
645          HYPRE_ParCSRGMRESSetup(par_solver, par_A, par_b, par_x);
646 
647          /* do the solve */
648          HYPRE_ParCSRGMRESSolve(par_solver, par_A, par_b, par_x);
649 
650          /* get some info */
651          HYPRE_GMRESGetNumIterations(par_solver, &its);
652          HYPRE_GMRESGetFinalRelativeResidualNorm(par_solver,
653                                                  &final_res_norm);
654          /* clean up */
655          HYPRE_ParCSRGMRESDestroy(par_solver);
656          HYPRE_BoomerAMGDestroy(par_precond);
657       }
658       else if (solver_id == 3) /* AMG */
659       {
660          HYPRE_BoomerAMGCreate(&par_solver);
661          HYPRE_BoomerAMGSetCoarsenType(par_solver, 6);
662          HYPRE_BoomerAMGSetStrongThreshold(par_solver, 0.25);
663          HYPRE_BoomerAMGSetTol(par_solver, 1.9e-6);
664          HYPRE_BoomerAMGSetPrintLevel(par_solver, 1);
665          HYPRE_BoomerAMGSetPrintFileName(par_solver, "ex9.out.log");
666          HYPRE_BoomerAMGSetMaxIter(par_solver, 50);
667 
668          /* do the setup */
669          HYPRE_BoomerAMGSetup(par_solver, par_A, par_b, par_x);
670 
671          /* do the solve */
672          HYPRE_BoomerAMGSolve(par_solver, par_A, par_b, par_x);
673 
674          /* get some info */
675          HYPRE_BoomerAMGGetNumIterations(par_solver, &its);
676          HYPRE_BoomerAMGGetFinalRelativeResidualNorm(par_solver,
677                                                      &final_res_norm);
678          /* clean up */
679          HYPRE_BoomerAMGDestroy(par_solver);
680       }
681       else
682       {
683          if (myid ==0) printf("\n ERROR: Invalid solver id specified.\n");
684       }
685 
686       /* Gather the solution vector.  This needs to be done if:
687          (1) the  object  type is parcsr OR
688          (2) any one of the variables is NOT cell-centered */
689       if (object_type == HYPRE_PARCSR)
690       {
691          HYPRE_SStructVectorGather(x);
692       }
693 
694       /* Save the solution for GLVis visualization, see vis/glvis-ex7.sh */
695       if (vis)
696       {
697          FILE *file;
698          char filename[255];
699 
700          int k, part = 0, var;
701          int nvalues = n*n;
702          double *values = (double*) calloc(nvalues, sizeof(double));
703 
704          /* save local solution for variable u */
705          var = 0;
706          HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
707                                          var, values);
708 
709          sprintf(filename, "%s.%06d", "vis/ex9-u.sol", myid);
710          if ((file = fopen(filename, "w")) == NULL)
711          {
712             printf("Error: can't open output file %s\n", filename);
713             MPI_Finalize();
714             exit(1);
715          }
716 
717          /* save solution with global unknown numbers */
718          k = 0;
719          for (j = 0; j < n; j++)
720             for (i = 0; i < n; i++)
721                fprintf(file, "%06d %.14e\n", pj*N*n*n+pi*n+j*N*n+i, values[k++]);
722 
723          fflush(file);
724          fclose(file);
725 
726          /* save local solution for variable v */
727          var = 1;
728          HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
729                                          var, values);
730 
731          sprintf(filename, "%s.%06d", "vis/ex9-v.sol", myid);
732          if ((file = fopen(filename, "w")) == NULL)
733          {
734             printf("Error: can't open output file %s\n", filename);
735             MPI_Finalize();
736             exit(1);
737          }
738 
739          /* save solution with global unknown numbers */
740          k = 0;
741          for (j = 0; j < n; j++)
742             for (i = 0; i < n; i++)
743                fprintf(file, "%06d %.14e\n", pj*N*n*n+pi*n+j*N*n+i, values[k++]);
744 
745          fflush(file);
746          fclose(file);
747 
748          free(values);
749 
750          /* save global finite element mesh */
751          if (myid == 0)
752             GLVis_PrintGlobalSquareMesh("vis/ex9.mesh", N*n-1);
753       }
754 
755       if (myid == 0)
756       {
757          printf("\n");
758          printf("Iterations = %d\n", its);
759          printf("Final Relative Residual Norm = %g\n", final_res_norm);
760          printf("\n");
761       }
762    }
763 
764    /* Free memory */
765    HYPRE_SStructGridDestroy(grid);
766    HYPRE_SStructStencilDestroy(stencil_v);
767    HYPRE_SStructStencilDestroy(stencil_u);
768    HYPRE_SStructGraphDestroy(graph);
769    HYPRE_SStructMatrixDestroy(A);
770    HYPRE_SStructVectorDestroy(b);
771    HYPRE_SStructVectorDestroy(x);
772 
773    /* Finalize MPI */
774    MPI_Finalize();
775 
776    return (0);
777 }
syntax highlighted by Code2HTML, v. 0.9.1