download the original source code.
  1 /*
  2    Example 13
  3 
  4    Interface:      Semi-Structured interface (SStruct)
  5 
  6    Compile with:   make ex13
  7 
  8    Sample run:     mpirun -np 6 ex13 -n 10
  9 
 10    To see options: ex13 -help
 11 
 12    Description:    This code solves the 2D Laplace equation using bilinear
 13                    finite element discretization on a mesh with an "enhanced
 14                    connectivity" point.  Specifically, we solve -Delta u = 1
 15                    with zero boundary conditions on a star-shaped domain
 16                    consisting of identical rhombic parts each meshed with a
 17                    uniform n x n grid.  Every part is assigned to a different
 18                    processor and all parts meet at the origin, equally
 19                    subdividing the 2*pi angle there. The case of six processors
 20                    (parts) looks as follows:
 21 
 22                                                 +
 23                                                / \
 24                                               /   \
 25                                              /     \
 26                                    +--------+   1   +---------+
 27                                     \        \     /         /
 28                                      \    2   \   /    0    /
 29                                       \        \ /         /
 30                                        +--------+---------+
 31                                       /        / \         \
 32                                      /    3   /   \    5    \
 33                                     /        /     \         \
 34                                    +--------+   4   +---------+
 35                                              \     /
 36                                               \   /
 37                                                \ /
 38                                                 +
 39 
 40                    Note that in this problem we use nodal variables, which are
 41                    shared between the different parts.  The node at the origin,
 42                    for example, belongs to all parts as illustrated below:
 43 
 44                                                 .
 45                                                / \
 46                                               .   .
 47                                              / \ / \
 48                                             o   .   *
 49                                   .---.---o  \ / \ /  *---.---.
 50                                    \   \   \  o   *  /   /   /
 51                                     .---.---o  \ /  *---.---.
 52                                      \   \   \  x  /   /   /
 53                                       @---@---x   x---z---z
 54                                       @---@---x   x---z---z
 55                                      /   /   /  x  \   \   \
 56                                     .---.---a  / \  #---.---.
 57                                    /   /   /  a   #  \   \   \
 58                                   .---.---a  / \ / \  #---.---.
 59                                             a   .   #
 60                                              \ / \ /
 61                                               .   .
 62                                                \ /
 63                                                 .
 64 
 65                    We recommend viewing the Struct examples before viewing this
 66                    and the other SStruct examples.  The primary role of this
 67                    particular SStruct example is to demonstrate a stencil-based
 68                    way to set up finite element problems in SStruct, and
 69                    specifically to show how to handle problems with an "enhanced
 70                    connectivity" point.
 71 */
 72 
 73 #include <math.h>
 74 #include "_hypre_utilities.h"
 75 #include "HYPRE_sstruct_mv.h"
 76 #include "HYPRE_sstruct_ls.h"
 77 #include "HYPRE.h"
 78 
 79 #ifndef M_PI
 80 #define M_PI 3.14159265358979
 81 #endif
 82 
 83 #include "vis.c"
 84 
 85 /*
 86    This routine computes the bilinear finite element stiffness matrix and
 87    load vector on a rhombus with angle gamma. Specifically, let R be the
 88    rhombus
 89                               [3]------[2]
 90                               /        /
 91                              /        /
 92                            [0]------[1]
 93 
 94    with sides of length h. The finite element stiffness matrix
 95 
 96                   S_ij = (grad phi_i,grad phi_j)_R
 97 
 98    with bilinear finite element functions {phi_i} has the form
 99 
100                        /  4-k    -1  -2+k    -1 \
101                alpha . |   -1   4+k    -1  -2-k |
102                        | -2+k    -1   4-k    -1 |
103                        \   -1  -2-k    -1   4+k /
104 
105    where alpha = 1/(6*sin(gamma)) and k = 3*cos(gamma). The load vector
106    corresponding to a right-hand side of 1 is
107 
108                   F_j = (1,phi_j)_R = h^2/4 * sin(gamma)
109 */
110 void ComputeFEMRhombus (double S[4][4], double F[4], double gamma, double h)
111 {
112    int i, j;
113 
114    double h2_4 = h*h/4;
115    double sing = sin(gamma);
116    double alpha = 1/(6*sing);
117    double k = 3*cos(gamma);
118 
119    S[0][0] = alpha * (4-k);
120    S[0][1] = alpha * (-1);
121    S[0][2] = alpha * (-2+k);
122    S[0][3] = alpha * (-1);
123    S[1][1] = alpha * (4+k);
124    S[1][2] = alpha * (-1);
125    S[1][3] = alpha * (-2-k);
126    S[2][2] = alpha * (4-k);
127    S[2][3] = alpha * (-1);
128    S[3][3] = alpha * (4+k);
129 
130    /* The stiffness matrix is symmetric */
131    for (i = 1; i < 4; i++)
132       for (j = 0; j < i; j++)
133          S[i][j] = S[j][i];
134 
135    for (i = 0; i < 4; i++)
136       F[i] = h2_4*sing;
137 }
138 
139 
140 int main (int argc, char *argv[])
141 {
142    int myid, num_procs;
143    int n;
144    double gamma, h;
145    int vis;
146 
147    HYPRE_SStructGrid     grid;
148    HYPRE_SStructGraph    graph;
149    HYPRE_SStructStencil  stencil;
150    HYPRE_SStructMatrix   A;
151    HYPRE_SStructVector   b;
152    HYPRE_SStructVector   x;
153 
154    HYPRE_Solver          solver;
155 
156    /* Initialize MPI */
157    MPI_Init(&argc, &argv);
158    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
159    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
160 
161    /* Set default parameters */
162    n = 10;
163    vis = 0;
164 
165    /* Parse command line */
166    {
167       int arg_index = 0;
168       int print_usage = 0;
169 
170       while (arg_index < argc)
171       {
172          if ( strcmp(argv[arg_index], "-n") == 0 )
173          {
174             arg_index++;
175             n = atoi(argv[arg_index++]);
176          }
177          else if ( strcmp(argv[arg_index], "-vis") == 0 )
178          {
179             arg_index++;
180             vis = 1;
181          }
182          else if ( strcmp(argv[arg_index], "-help") == 0 )
183          {
184             print_usage = 1;
185             break;
186          }
187          else
188          {
189             arg_index++;
190          }
191       }
192 
193       if ((print_usage) && (myid == 0))
194       {
195          printf("\n");
196          printf("Usage: %s [<options>]\n", argv[0]);
197          printf("\n");
198          printf("  -n <n>              : problem size per processor (default: 10)\n");
199          printf("  -vis                : save the solution for GLVis visualization\n");
200          printf("\n");
201       }
202 
203       if (print_usage)
204       {
205          MPI_Finalize();
206          return (0);
207       }
208    }
209 
210    /* Set the rhombus angle, gamma, and the mesh size, h, depending on the
211       number of processors np and the given n */
212    if (num_procs < 3)
213    {
214       if (myid ==0) printf("Must run with at least 3 processors!\n");
215       MPI_Finalize();
216       exit(1);
217    }
218    gamma = 2*M_PI/num_procs;
219    h = 1.0/n;
220 
221    /* 1. Set up the grid.  We will set up the grid so that processor X owns
222          part X.  Note that each part has its own index space numbering. Later
223          we relate the parts to each other. */
224    {
225       int ndim = 2;
226       int nparts = num_procs;
227 
228       /* Create an empty 2D grid object */
229       HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &grid);
230 
231       /* Set the extents of the grid - each processor sets its grid boxes.  Each
232          part has its own relative index space numbering */
233       {
234          int part = myid;
235          int ilower[2] = {1,1}; /* lower-left cell touching the origin */
236          int iupper[2] = {n,n}; /* upper-right cell */
237 
238          HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
239       }
240 
241       /* Set the variable type and number of variables on each part.  These need
242          to be set in each part which is neighboring or contains boxes owned by
243          the processor. */
244       {
245          int i;
246          int nvars = 1;
247 
248          HYPRE_SStructVariable vartypes[1] = {HYPRE_SSTRUCT_VARIABLE_NODE};
249          for (i = 0; i < nparts; i++)
250             HYPRE_SStructGridSetVariables(grid, i, nvars, vartypes);
251       }
252 
253       /* Now we need to set the spatial relation between each of the parts.
254          Since we are using nodal variables, we have to use SetSharedPart to
255          establish the connection at the origin. */
256       {
257          /* Relation to the clockwise-previous neighbor part, e.g. 0 and 1 for
258             the case of 6 parts.  Note that we could have used SetNeighborPart
259             here instead of SetSharedPart. */
260          {
261             int part = myid;
262             /* the box of cells intersecting the boundary in the current part */
263             int ilower[2] = {1,1}, iupper[2] = {1,n};
264             /* share all data on the left side of the box */
265             int offset[2] = {-1,0};
266 
267             int shared_part = (myid+1) % num_procs;
268             /* the box of cells intersecting the boundary in the neighbor */
269             int shared_ilower[2] = {1,1}, shared_iupper[2] = {n,1};
270             /* share all data on the bottom of the box */
271             int shared_offset[2] = {0,-1};
272 
273             /* x/y-direction on the current part is -y/x on the neighbor */
274             int index_map[2] = {1,0};
275             int index_dir[2] = {-1,1};
276 
277             HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
278                                            shared_part, shared_ilower,
279                                            shared_iupper, shared_offset,
280                                            index_map, index_dir);
281          }
282 
283          /* Relation to the clockwise-following neighbor part, e.g. 0 and 5 for
284             the case of 6 parts.  Note that we could have used SetNeighborPart
285             here instead of SetSharedPart. */
286          {
287             int part = myid;
288             /* the box of cells intersecting the boundary in the current part */
289             int ilower[2] = {1,1}, iupper[2] = {n,1};
290             /* share all data on the bottom of the box */
291             int offset[2] = {0,-1};
292 
293             int shared_part = (myid+num_procs-1) % num_procs;
294             /* the box of cells intersecting the boundary in the neighbor */
295             int shared_ilower[2] = {1,1}, shared_iupper[2] = {1,n};
296             /* share all data on the left side of the box */
297             int shared_offset[2] = {-1,0};
298 
299             /* x/y-direction on the current part is y/-x on the neighbor */
300             int index_map[2] = {1,0};
301             int index_dir[2] = {1,-1};
302 
303             HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
304                                            shared_part, shared_ilower,
305                                            shared_iupper, shared_offset,
306                                            index_map, index_dir);
307          }
308 
309          /* Relation to all other parts, e.g. 0 and 2,3,4.  This can be
310             described only by SetSharedPart. */
311          {
312             int part = myid;
313             /* the (one cell) box that touches the origin */
314             int ilower[2] = {1,1}, iupper[2] = {1,1};
315             /* share all data in the bottom left corner (i.e. the origin) */
316             int offset[2] = {-1,-1};
317 
318             int shared_part;
319             /* the box of one cell that touches the origin */
320             int shared_ilower[2] = {1,1}, shared_iupper[2] = {1,1};
321             /* share all data in the bottom left corner (i.e. the origin) */
322             int shared_offset[2] = {-1,-1};
323 
324             /* x/y-direction on the current part is -x/-y on the neighbor, but
325                in this case the arguments are not really important since we are
326                only sharing a point */
327             int index_map[2] = {0,1};
328             int index_dir[2] = {-1,-1};
329 
330             for (shared_part = 0; shared_part < myid-1; shared_part++)
331                HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
332                                               shared_part, shared_ilower,
333                                               shared_iupper, shared_offset,
334                                               index_map, index_dir);
335 
336             for (shared_part = myid+2; shared_part < num_procs; shared_part++)
337                HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
338                                               shared_part, shared_ilower,
339                                               shared_iupper, shared_offset,
340                                               index_map, index_dir);
341          }
342       }
343 
344       /* Now the grid is ready to be used */
345       HYPRE_SStructGridAssemble(grid);
346    }
347 
348    /* 2. Define the discretization stencils.  Since this is a finite element
349          discretization we define here a full 9-point stencil.  We will later
350          use four sub-stencils for the rows of the local stiffness matrix. */
351    {
352       int ndim = 2;
353       int var = 0;
354       int entry;
355 
356       /* Define the geometry of the 9-point stencil */
357       int stencil_size = 9;
358       int offsets[9][2] = {{0,0},           /*  [8] [4] [7]  */
359                            {-1,0}, {1,0},   /*     \ | /     */
360                            {0,-1}, {0,1},   /*  [1]-[0]-[2]  */
361                            {-1,-1}, {1,-1}, /*     / | \     */
362                            {1,1}, {-1,1}};  /*  [5] [3] [6]  */
363 
364       HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil);
365 
366       for (entry = 0; entry < stencil_size; entry++)
367          HYPRE_SStructStencilSetEntry(stencil, entry, offsets[entry], var);
368    }
369 
370    /* 3. Set up the Graph - this determines the non-zero structure of the
371          matrix. */
372    {
373       int part;
374       int var = 0;
375 
376       /* Create the graph object */
377       HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
378 
379       /* See MatrixSetObjectType below */
380       HYPRE_SStructGraphSetObjectType(graph, HYPRE_PARCSR);
381 
382       /* Now we need to tell the graph which stencil to use for each
383          variable on each part (we only have one variable) */
384       for (part = 0; part < num_procs; part++)
385          HYPRE_SStructGraphSetStencil(graph, part, var, stencil);
386 
387       /* Assemble the graph */
388       HYPRE_SStructGraphAssemble(graph);
389    }
390 
391    /* 4. Set up the SStruct Matrix and right-hand side vector */
392    {
393       int part = myid;
394       int var = 0;
395 
396       /* Create the matrix object */
397       HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
398       /* Use a ParCSR storage */
399       HYPRE_SStructMatrixSetObjectType(A, HYPRE_PARCSR);
400       /* Indicate that the matrix coefficients are ready to be set */
401       HYPRE_SStructMatrixInitialize(A);
402 
403       /* Create an empty vector object */
404       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
405       /* Use a ParCSR storage */
406       HYPRE_SStructVectorSetObjectType(b, HYPRE_PARCSR);
407       /* Indicate that the vector coefficients are ready to be set */
408       HYPRE_SStructVectorInitialize(b);
409 
410       /* Set the matrix and vector entries by finite element assembly */
411       {
412          /* local stifness matrix and load vector */
413          double S[4][4], F[4];
414 
415          /* The index of the local nodes 0-3 relative to the cell index,
416             i.e. node k in cell (i,j) is in the upper-right corner of the
417             cell (i,j) + node_index_offset[k]. */
418          int node_index_offset[4][2] = {{-1,-1},{0,-1},{0,0},{-1,0}};
419 
420          /* The cell sub-stencils of nodes 0-3 indexed from the full stencil,
421             i.e. we take the full stencil in each node of a fixed cell, and
422             restrict it to that as is done in the finite element stiffness
423             matrix:
424                          [4] [7]   [8] [4]   [1]-[0]   [0]-[2]
425                           | /         \ |       / |     | \
426                          [0]-[2] , [1]-[0] , [5] [3] , [3] [6]
427 
428             Note that the ordering of the local nodes remains fixed, and
429             therefore the above sub-stencil at node k corresponds to the kth row
430             of the local stiffness matrix and the kth entry of the local load
431             vector. */
432          int node_stencil[4][4] = {{0,2,7,4},{1,0,4,8},{5,3,0,1},{3,6,2,0}};
433 
434          int i, j, k;
435          int index[2];
436          int nentries = 4;
437 
438          /* set the values in the interior cells */
439          {
440             ComputeFEMRhombus(S, F, gamma, h);
441 
442             for (i = 1; i <= n; i++)
443                for (j = 1; j <= n; j++)
444                   for (k = 0; k < 4; k++) /* node k in cell (i,j) */
445                   {
446                      index[0] = i + node_index_offset[k][0];
447                      index[1] = j + node_index_offset[k][1];
448                      HYPRE_SStructMatrixAddToValues(A, part, index, var,
449                                                     nentries, node_stencil[k],
450                                                     &S[k][0]);
451                      HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
452                   }
453          }
454 
455          /* cells having nodes 1,2 on the domain boundary */
456          {
457             ComputeFEMRhombus(S, F, gamma, h);
458 
459             /* eliminate nodes 1,2 from S and F */
460             for (k = 0; k < 4; k++)
461             {
462                S[1][k] = S[k][1] = 0.0;
463                S[2][k] = S[k][2] = 0.0;
464             }
465             S[1][1] = 1.0;
466             S[2][2] = 1.0;
467             F[1] = 0.0;
468             F[2] = 0.0;
469 
470             for (i = n; i <= n; i++)
471                for (j = 1; j <= n; j++)
472                   for (k = 0; k < 4; k++) /* node k in cell (n,j) */
473                   {
474                      index[0] = i + node_index_offset[k][0];
475                      index[1] = j + node_index_offset[k][1];
476                      HYPRE_SStructMatrixAddToValues(A, part, index, var,
477                                                     nentries, node_stencil[k],
478                                                     &S[k][0]);
479                      HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
480                   }
481          }
482 
483          /* cells having nodes 2,3 on the domain boundary */
484          {
485             ComputeFEMRhombus(S, F, gamma, h);
486 
487             /* eliminate nodes 2,3 from S and F */
488             for (k = 0; k < 4; k++)
489             {
490                S[2][k] = S[k][2] = 0.0;
491                S[3][k] = S[k][3] = 0.0;
492             }
493             S[2][2] = 1.0;
494             S[3][3] = 1.0;
495             F[2] = 0.0;
496             F[3] = 0.0;
497 
498             for (i = 1; i <= n; i++)
499                for (j = n; j <= n; j++)
500                   for (k = 0; k < 4; k++) /* node k in cell (i,n) */
501                   {
502                      index[0] = i + node_index_offset[k][0];
503                      index[1] = j + node_index_offset[k][1];
504                      HYPRE_SStructMatrixAddToValues(A, part, index, var,
505                                                     nentries, node_stencil[k],
506                                                     &S[k][0]);
507                      HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
508                   }
509          }
510 
511          /* cells having nodes 1,2,3 on the domain boundary */
512          {
513             ComputeFEMRhombus(S, F, gamma, h);
514 
515             /* eliminate nodes 2,3 from S and F */
516             for (k = 0; k < 4; k++)
517             {
518                S[1][k] = S[k][1] = 0.0;
519                S[2][k] = S[k][2] = 0.0;
520                S[3][k] = S[k][3] = 0.0;
521             }
522             S[1][1] = 1.0;
523             S[2][2] = 1.0;
524             S[3][3] = 1.0;
525             F[1] = 0.0;
526             F[2] = 0.0;
527             F[3] = 0.0;
528 
529             for (i = n; i <= n; i++)
530                for (j = n; j <= n; j++)
531                   for (k = 0; k < 4; k++) /* node k in cell (n,n) */
532                   {
533                      index[0] = i + node_index_offset[k][0];
534                      index[1] = j + node_index_offset[k][1];
535                      HYPRE_SStructMatrixAddToValues(A, part, index, var,
536                                                     nentries, node_stencil[k],
537                                                     &S[k][0]);
538                      HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
539                   }
540          }
541       }
542    }
543 
544    /* Collective calls finalizing the matrix and vector assembly */
545    HYPRE_SStructMatrixAssemble(A);
546    HYPRE_SStructVectorAssemble(b);
547 
548    /* 5. Set up SStruct Vector for the solution vector x */
549    {
550       int part = myid;
551       int var = 0;
552       int nvalues = (n+1)*(n+1);
553       double *values;
554 
555       /* Since the SetBoxValues() calls below set the values of the nodes in
556          the upper-right corners of the cells, the nodal box should start
557          from (0,0) instead of (1,1). */
558       int ilower[2] = {0,0};
559       int iupper[2] = {n,n};
560 
561       values = (double*) calloc(nvalues, sizeof(double));
562 
563       /* Create an empty vector object */
564       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
565       /* Set the object type to ParCSR */
566       HYPRE_SStructVectorSetObjectType(x, HYPRE_PARCSR);
567       /* Indicate that the vector coefficients are ready to be set */
568       HYPRE_SStructVectorInitialize(x);
569       /* Set the values for the initial guess */
570       HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
571 
572       free(values);
573 
574       /* Finalize the vector assembly */
575       HYPRE_SStructVectorAssemble(x);
576    }
577 
578    /* 6. Set up and call the solver (Solver options can be found in the
579          Reference Manual.) */
580    {
581       double final_res_norm;
582       int its;
583 
584       HYPRE_ParCSRMatrix    par_A;
585       HYPRE_ParVector       par_b;
586       HYPRE_ParVector       par_x;
587 
588       /* Extract the ParCSR objects needed in the solver */
589       HYPRE_SStructMatrixGetObject(A, (void **) &par_A);
590       HYPRE_SStructVectorGetObject(b, (void **) &par_b);
591       HYPRE_SStructVectorGetObject(x, (void **) &par_x);
592 
593       /* Here we construct a BoomerAMG solver.  See the other SStruct examples
594          as well as the Reference manual for additional solver choices. */
595       HYPRE_BoomerAMGCreate(&solver);
596       HYPRE_BoomerAMGSetCoarsenType(solver, 6);
597       HYPRE_BoomerAMGSetStrongThreshold(solver, 0.25);
598       HYPRE_BoomerAMGSetTol(solver, 1e-6);
599       HYPRE_BoomerAMGSetPrintLevel(solver, 2);
600       HYPRE_BoomerAMGSetMaxIter(solver, 50);
601 
602       /* call the setup */
603       HYPRE_BoomerAMGSetup(solver, par_A, par_b, par_x);
604 
605       /* call the solve */
606       HYPRE_BoomerAMGSolve(solver, par_A, par_b, par_x);
607 
608       /* get some info */
609       HYPRE_BoomerAMGGetNumIterations(solver, &its);
610       HYPRE_BoomerAMGGetFinalRelativeResidualNorm(solver,
611                                                   &final_res_norm);
612       /* clean up */
613       HYPRE_BoomerAMGDestroy(solver);
614 
615       /* Gather the solution vector */
616       HYPRE_SStructVectorGather(x);
617 
618       /* Save the solution for GLVis visualization, see vis/glvis-ex13.sh */
619       if (vis)
620       {
621          FILE *file;
622          char filename[255];
623 
624          int i, part = myid, var = 0;
625          int nvalues = (n+1)*(n+1);
626          double *values = (double*) calloc(nvalues, sizeof(double));
627          int ilower[2] = {0,0};
628          int iupper[2] = {n,n};
629 
630          /* get all local data (including a local copy of the shared values) */
631          HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
632                                          var, values);
633 
634          sprintf(filename, "%s.%06d", "vis/ex13.sol", myid);
635          if ((file = fopen(filename, "w")) == NULL)
636          {
637             printf("Error: can't open output file %s\n", filename);
638             MPI_Finalize();
639             exit(1);
640          }
641 
642          /* finite element space header */
643          fprintf(file, "FiniteElementSpace\n");
644          fprintf(file, "FiniteElementCollection: H1_2D_P1\n");
645          fprintf(file, "VDim: 1\n");
646          fprintf(file, "Ordering: 0\n\n");
647 
648          /* save solution */
649          for (i = 0; i < nvalues; i++)
650             fprintf(file, "%.14e\n", values[i]);
651 
652          fflush(file);
653          fclose(file);
654          free(values);
655 
656          /* save local finite element mesh */
657          GLVis_PrintLocalRhombusMesh("vis/ex13.mesh", n, myid, gamma);
658 
659          /* additional visualization data */
660          GLVis_PrintData("vis/ex13.data", myid, num_procs);
661       }
662 
663       if (myid == 0)
664       {
665          printf("\n");
666          printf("Iterations = %d\n", its);
667          printf("Final Relative Residual Norm = %g\n", final_res_norm);
668          printf("\n");
669       }
670    }
671 
672    /* Free memory */
673    HYPRE_SStructGridDestroy(grid);
674    HYPRE_SStructStencilDestroy(stencil);
675    HYPRE_SStructGraphDestroy(graph);
676    HYPRE_SStructMatrixDestroy(A);
677    HYPRE_SStructVectorDestroy(b);
678    HYPRE_SStructVectorDestroy(x);
679 
680    /* Finalize MPI */
681    MPI_Finalize();
682 
683    return 0;
684 }
syntax highlighted by Code2HTML, v. 0.9.1