download the original source code.
  1 /*
  2    Example 8
  3 
  4    Interface:    Semi-Structured interface (SStruct)
  5 
  6    Compile with: make ex8
  7 
  8    Sample run:   mpirun -np 2 ex8
  9 
 10    Description:  This is a two processor example which solves a similar
 11                  problem to the one in Example 2, and Example 6 (The grid
 12                  boxes are exactly those in the example diagram in the
 13                  struct interface chapter of the User's Manual.)
 14 
 15                  The difference with the previous examples is that we use
 16                  three parts, two with a 5-point and one with a 9-point
 17                  discretization stencil. The solver is PCG with split-SMG
 18                  preconditioner.
 19 */
 20 
 21 #include <stdio.h>
 22 
 23 /* SStruct linear solvers headers */
 24 #include "HYPRE_sstruct_ls.h"
 25 
 26 
 27 int main (int argc, char *argv[])
 28 {
 29    int myid, num_procs;
 30 
 31    HYPRE_SStructGrid     grid;
 32    HYPRE_SStructGraph    graph;
 33    HYPRE_SStructStencil  stencil_5pt;
 34    HYPRE_SStructStencil  stencil_9pt;
 35    HYPRE_SStructMatrix   A;
 36    HYPRE_SStructVector   b;
 37    HYPRE_SStructVector   x;
 38    HYPRE_SStructSolver   solver;
 39    HYPRE_SStructSolver   precond;
 40 
 41    int object_type;
 42 
 43    /* Initialize MPI */
 44    MPI_Init(&argc, &argv);
 45    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 46    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 47 
 48    if (num_procs != 2)
 49    {
 50       if (myid ==0) printf("Must run with 2 processors!\n");
 51       MPI_Finalize();
 52 
 53       return(0);
 54    }
 55 
 56    /* 1. Set up the 2D grid.  This gives the index space in each part.
 57       We have one variable in each part. */
 58    {
 59       int ndim = 2;
 60       int nparts = 3;
 61       int part;
 62 
 63       /* Create an empty 2D grid object */
 64       HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &grid);
 65 
 66       /* Set the extents of the grid - each processor sets its grid
 67          boxes.  Each part has its own relative index space numbering. */
 68 
 69       /* Processor 0 owns two boxes - one in part 0 and one in part 1. */
 70       if (myid == 0)
 71       {
 72          /* Add the first box to the grid in part 0 */
 73          {
 74             int ilower[2] = {-3, 1};
 75             int iupper[2] = {-1, 2};
 76 
 77             part = 0;
 78             HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
 79          }
 80 
 81          /* Add the second box to the grid in part 1 */
 82          {
 83             /* For convenience we use the same index space across all
 84                parts, but this is not a requirement. For example, on this
 85                part we could have used ilower=[23,24] and iupper=[25,27]. */
 86             int ilower[2] = {0, 1};
 87             int iupper[2] = {2, 4};
 88 
 89             part = 1;
 90             HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
 91          }
 92       }
 93 
 94       /* Processor 1 owns one box in part 2. */
 95       else if (myid == 1)
 96       {
 97          /* Add a new box to the grid in part 2 */
 98          {
 99             int ilower[2] = {3, 1};
100             int iupper[2] = {6, 4};
101 
102             part = 2;
103             HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
104          }
105       }
106 
107       /* Set the variable type and number of variables on each part. */
108       {
109          int i;
110          int nvars = 1;
111          HYPRE_SStructVariable vartypes[1] = {HYPRE_SSTRUCT_VARIABLE_CELL};
112 
113          for (i = 0; i< nparts; i++)
114             HYPRE_SStructGridSetVariables(grid, i, nvars, vartypes);
115       }
116 
117       /* Now we need to set the spatial relation between each of the parts.
118          Since we have the same types of variables on both parts, we can
119          use HYPRE_GridSetNeighborPart().  Each processor calls this function
120          for each part on which it owns boxes that border a different part. */
121 
122       if (myid == 0)
123       {
124          /* Relation between part 0 and part 1 on processor 0 */
125          {
126             int part = 0;
127             int nbor_part = 1;
128             /* Cells just outside of the boundary of part 0 in
129                its coordinates */
130             int b_ilower[2] = {0,1}, b_iupper[2] = {0,2};
131             /* The same cells in part 1's coordinates.  Since we use the same
132                index space across all parts, the coordinates coincide. */
133             int nbor_ilower[2] = {0,1}, nbor_iupper[2] = {0,2};
134             /* These parts have the same orientation, so no
135                rotation is necessary */
136             int index_map[2] = {0,1};
137             /* These parts map increasing values to increasing values 
138                for both variables (note: if decreasing maps to increasing, use -1)*/
139             int index_dir[2] = {1,1};
140 
141             HYPRE_SStructGridSetNeighborPart(grid, part, b_ilower, b_iupper,
142                                              nbor_part, nbor_ilower, nbor_iupper,
143                                              index_map, index_dir);
144          }
145 
146          /* Relation between part 1 and part 0 on processor 0 */
147          {
148             int part = 1;
149             int nbor_part = 0;
150             /* Cells just outside of the boundary of part 1 in
151                its coordinates */
152             int b_ilower[2] = {-1,1}, b_iupper[2] = {-1,2};
153             /* The same cells in part 0's coordinates.  Since we use the same
154                index space across all parts, the coordinates coincide. */
155             int nbor_ilower[2] = {-1,1}, nbor_iupper[2] = {-1,2};
156             /* These parts have the same orientation, so no
157                rotation is necessary */
158             int index_map[2] = {0,1};
159             /* These parts map increasing values to increasing values 
160                for both variables (note: if decreasing maps to increasing, use -1)*/
161             int index_dir[2] = {1,1};
162 
163             HYPRE_SStructGridSetNeighborPart(grid, part, b_ilower, b_iupper,
164                                              nbor_part, nbor_ilower, nbor_iupper,
165                                              index_map, index_dir);
166          }
167 
168          /* Relation between part 1 and part 2 on processor 0 */
169          {
170             int part = 1;
171             int nbor_part = 2;
172             /* Cells just outside of the boundary of part 1 in
173                its coordinates */
174             int b_ilower[2] = {3,1}, b_iupper[2] = {3,4};
175             /* The same cells in part 2's coordinates.  Since we use the same
176                index space across all parts, the coordinates coincide. */
177             int nbor_ilower[2] = {3,1}, nbor_iupper[2] = {3,4};
178             /* These parts have the same orientation, so no
179                rotation is necessary */
180             int index_map[2] = {0,1};
181             /* These parts map increasing values to increasing values 
182                for both variables (note: if decreasing maps to increasing, use -1)*/
183             int index_dir[2] = {1,1};
184 
185             HYPRE_SStructGridSetNeighborPart(grid, part, b_ilower, b_iupper,
186                                             nbor_part, nbor_ilower, nbor_iupper,
187                                              index_map, index_dir);
188          }
189       }
190       else if (myid == 1)
191       {
192          /* Relation between part 2 and part 1 on processor 1 */
193          {
194             int part = 2;
195             int nbor_part = 1;
196             /* Cells just outside of the boundary of part 2 in
197                its coordinates */
198             int b_ilower[2] = {2,1}, b_iupper[2] = {2,4};
199             /* The same cells in part 1's coordinates.  Since we use the same
200                index space across all parts, the coordinates coincide. */
201             int nbor_ilower[2] = {2,1}, nbor_iupper[2] = {2,4};
202             /* These parts have the same orientation, so no
203                rotation is necessary */
204             int index_map[2] = {0,1};
205             /* These parts map increasing values to increasing values 
206               for both variables (note: if decreasing maps to increasing, use -1)*/
207             int index_dir[2] = {1,1};
208 
209             HYPRE_SStructGridSetNeighborPart(grid, part, b_ilower, b_iupper,
210                                              nbor_part, nbor_ilower, nbor_iupper,
211                                              index_map, index_dir); 
212          }
213       }
214 
215       /* Now the grid is ready to use */
216       HYPRE_SStructGridAssemble(grid);
217    }
218 
219    /* 2. Define the discretization stencils */
220    {
221       int ndim = 2;
222       int var = 0;
223       int entry;
224 
225       /* the 5-pt stencil in 2D */
226       {
227          int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
228          int stencil_size = 5;
229 
230          HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil_5pt);
231 
232          for (entry = 0; entry < 5; entry++)
233             HYPRE_SStructStencilSetEntry(stencil_5pt, entry, offsets[entry], var);
234       }
235 
236       /* the 9-pt stencil in 2D */
237       {
238          int offsets[9][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1},
239                               {-1,-1}, {1,-1}, {1,1}, {-1,1}};
240          int stencil_size = 9;
241          HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil_9pt);
242 
243          for (entry = 0; entry < stencil_size; entry++)
244             HYPRE_SStructStencilSetEntry(stencil_9pt, entry, offsets[entry], var);
245       }
246    }
247 
248    /* 3. Set up the Graph  - this determines the non-zero structure
249       of the matrix and allows non-stencil relationships between the parts */
250    {
251       int var = 0;
252       int part;
253 
254       /* Create the graph object */
255       HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
256 
257       /* Use the 5-pt stencil on part 0 */
258       part = 0;
259       HYPRE_SStructGraphSetStencil(graph, part, var, stencil_5pt);
260 
261       /* Use the 9-pt stencil on part 1 */
262       part = 1;
263       HYPRE_SStructGraphSetStencil(graph, part, var, stencil_9pt);
264 
265       /* Use the 5-pt stencil on part 2 */
266       part = 2;
267       HYPRE_SStructGraphSetStencil(graph, part, var, stencil_5pt);
268 
269       /*  Since we have only stencil connections between parts, we don't need to
270           call HYPRE_SStructGraphAddEntries. */
271 
272       /* Assemble the graph */
273       HYPRE_SStructGraphAssemble(graph);
274    }
275 
276    /* 4. Set up a SStruct Matrix */
277    {
278       int i,j;
279       int part;
280       int var = 0;
281 
282       /* Create the empty matrix object */
283       HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
284 
285       /* Set the object type (by default HYPRE_SSTRUCT). This determines the
286          data structure used to store the matrix.  If you want to use unstructured
287          solvers, e.g. BoomerAMG, the object type should be HYPRE_PARCSR.
288          If the problem is purely structured (with one part), you may want to use
289          HYPRE_STRUCT to access the structured solvers.  Since we have two parts
290          with different stencils, we set the object type to HYPRE_SSTRUCT. */
291       object_type = HYPRE_SSTRUCT;
292       HYPRE_SStructMatrixSetObjectType(A, object_type);
293 
294       /* Get ready to set values */
295       HYPRE_SStructMatrixInitialize(A);
296 
297       /* Each processor must set the stencil values for their boxes on each part.
298          In this example, we only set stencil entries and therefore use
299          HYPRE_SStructMatrixSetBoxValues.  If we need to set non-stencil entries,
300          we have to use HYPRE_SStructMatrixSetValues. */
301 
302       if (myid == 0)
303       {
304          /* Set the matrix coefficients for some set of stencil entries
305             over all the gridpoints in my first box (account for boundary
306             grid points later) */
307          {
308             int ilower[2] = {-3, 1};
309             int iupper[2] = {-1, 2};
310 
311             int nentries = 5;
312             int nvalues  = 30; /* 6 grid points, each with 5 stencil entries */
313             double values[30];
314 
315             int stencil_indices[5];
316             for (j = 0; j < nentries; j++) /* label the stencil indices -
317                                               these correspond to the offsets
318                                               defined above */
319                stencil_indices[j] = j;
320 
321             for (i = 0; i < nvalues; i += nentries)
322             {
323                values[i] = 4.0;
324                for (j = 1; j < nentries; j++)
325                   values[i+j] = -1.0;
326             }
327 
328             part = 0;
329             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
330                                             var, nentries,
331                                             stencil_indices, values);
332          }
333 
334          /* Set the matrix coefficients for some set of stencil entries
335             over the gridpoints in my second box */
336          {
337             int ilower[2] = {0, 1};
338             int iupper[2] = {2, 4};
339 
340             int nentries = 9;
341             int nvalues  = 108; /* 12 grid points, each with 5 stencil entries */
342             double values[108];
343 
344             int stencil_indices[9];
345             for (j = 0; j < nentries; j++)
346                stencil_indices[j] = j;
347 
348             for (i = 0; i < nvalues; i += nentries)
349             {
350                values[i] = 8./3.;
351                for (j = 1; j < nentries; j++)
352                   values[i+j] = -1./3.;
353             }
354 
355             part = 1;
356             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
357                                             var, nentries,
358                                             stencil_indices, values);
359          }
360       }
361       else if (myid == 1)
362       {
363          /* Set the matrix coefficients for some set of stencil entries
364             over the gridpoints in my box */
365          {
366             int ilower[2] = {3, 1};
367             int iupper[2] = {6, 4};
368 
369             int nentries = 5;
370             int nvalues  = 80; /* 16 grid points, each with 5 stencil entries */
371             double values[80];
372 
373             int stencil_indices[5];
374             for (j = 0; j < nentries; j++)
375                stencil_indices[j] = j;
376 
377             for (i = 0; i < nvalues; i += nentries)
378             {
379                values[i] = 4.0;
380                for (j = 1; j < nentries; j++)
381                   values[i+j] = -1.0;
382             }
383 
384             part = 2;
385             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
386                                             var, nentries,
387                                             stencil_indices, values);
388          }
389       }
390 
391       /* Modify the 9-pt stencil on the boundary between parts to ensure
392          symmetry and good global approximation. */
393       if (myid == 0)
394       {
395          int nentries = 6;
396          int nvalues  = 24; /* 4 grid points, each with 6 stencil entries */
397          double values[24];
398 
399          part = 1;
400 
401          for (i = 0; i < nvalues; i += nentries)
402          {
403             values[i]   = 10./3.;
404             values[i+1] = -1.;
405             values[i+2] = -2./3.;
406             values[i+3] = -2./3.;
407             values[i+4] = 0.0;
408             values[i+5] = 0.0;
409          }
410 
411          {
412             /* Values to the right of the second box */
413             int ilower[2] = { 2, 1};
414             int iupper[2] = { 2, 4};
415 
416             int stencil_indices[6] = {0,2,3,4,6,7};
417 
418             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
419                                             var, nentries,
420                                             stencil_indices, values);
421          }
422 
423          {
424             /* Values to the left of the second box */
425             int ilower[2] = { 0, 1};
426             int iupper[2] = { 0, 4};
427 
428             int stencil_indices[6] = {0,1,3,4,5,8};
429 
430             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
431                                             var, nentries,
432                                             stencil_indices, values);
433          }
434       }
435 
436       /* For each box, set any coefficients that reach ouside of the
437          boundary to 0 */
438       if (myid == 0)
439       {
440          int maxnvalues = 9;
441          double values[9];
442 
443          for (i = 0; i < maxnvalues; i++)
444             values[i] = 0.0;
445 
446          part = 0;
447 
448          {
449             /* Values below our first box */
450             int ilower[2] = {-3, 1};
451             int iupper[2] = {-1, 1};
452 
453             int stencil_indices[1] = {3};
454 
455             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
456                                             var, 1,
457                                             stencil_indices, values);
458          }
459 
460          {
461             /* Values to the left of our first box */
462             int ilower[2] = {-3, 1};
463             int iupper[2] = {-3, 2};
464 
465             int stencil_indices[1] = {1};
466 
467             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
468                                             var, 1,
469                                             stencil_indices, values);
470          }
471 
472          {
473             /* Values above our first box */
474             int ilower[2] = {-3, 2};
475             int iupper[2] = {-1, 2};
476 
477             int stencil_indices[1] = {4};
478 
479             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
480                                             var, 1,
481                                             stencil_indices, values);
482          }
483 
484          part = 1;
485 
486          {
487             /* Values below our second box */
488             int ilower[2] = { 0, 1};
489             int iupper[2] = { 2, 1};
490 
491             int stencil_indices[3] = {3,5,6};
492 
493             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
494                                             var, 3,
495                                             stencil_indices, values);
496          }
497 
498          {
499             /* Values to the left of our second box (that do not border the
500                first box). */
501             int ilower[2] = { 0, 3};
502             int iupper[2] = { 0, 4};
503 
504             int stencil_indices[3] = {1,5,8};
505 
506             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
507                                             var, 3,
508                                             stencil_indices, values);
509          }
510 
511          {
512             /* Values above our second box */
513             int ilower[2] = { 0, 4};
514             int iupper[2] = { 2, 4};
515 
516             int stencil_indices[3] = {4,7,8};
517 
518             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
519                                             var, 3,
520                                             stencil_indices, values);
521          }
522       }
523       else if (myid == 1)
524       {
525          int maxnvalues = 4;
526          double values[4];
527 
528          for (i = 0; i < maxnvalues; i++)
529             values[i] = 0.0;
530 
531          part = 2;
532 
533          {
534             /* Values below our box */
535             int ilower[2] = { 3, 1};
536             int iupper[2] = { 6, 1};
537 
538             int stencil_indices[1] = {3};
539 
540             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
541                                             var, 1,
542                                             stencil_indices, values);
543          }
544 
545          {
546             /* Values to the right of our box */
547             int ilower[2] = { 6, 1};
548             int iupper[2] = { 6, 4};
549 
550             int stencil_indices[1] = {2};
551 
552             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
553                                             var, 1,
554                                             stencil_indices, values);
555          }
556 
557          {
558             /* Values above our box */
559             int ilower[2] = { 3, 4};
560             int iupper[2] = { 6, 4};
561 
562             int stencil_indices[1] = {4};
563 
564             HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
565                                             var, 1,
566                                             stencil_indices, values);
567          }
568       }
569 
570       /* This is a collective call finalizing the matrix assembly.
571          The matrix is now ``ready to be used'' */
572       HYPRE_SStructMatrixAssemble(A);
573    }
574 
575    /* 5. Set up SStruct Vectors for b and x */
576    {
577       int i;
578       int part;
579       int var = 0;
580 
581       /* Create an empty vector object */
582       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
583       HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
584 
585       /* As with the matrix,  set the object type for the vectors
586          to be the sstruct type */
587       object_type = HYPRE_SSTRUCT;
588       HYPRE_SStructVectorSetObjectType(b, object_type);
589       HYPRE_SStructVectorSetObjectType(x, object_type);
590 
591       /* Indicate that the vector coefficients are ready to be set */
592       HYPRE_SStructVectorInitialize(b);
593       HYPRE_SStructVectorInitialize(x);
594 
595       if (myid == 0)
596       {
597          /* Set the vector coefficients over the gridpoints in my first box */
598          {
599             int ilower[2] = {-3, 1};
600             int iupper[2] = {-1, 2};
601 
602             int nvalues = 6;  /* 6 grid points */
603             double values[6];
604 
605             part = 0;
606 
607             for (i = 0; i < nvalues; i ++)
608                values[i] = 1.0;
609             HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
610 
611             for (i = 0; i < nvalues; i ++)
612                values[i] = 0.0;
613             HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
614          }
615 
616          /* Set the vector coefficients over the gridpoints in my second box */
617          {
618             int ilower[2] = { 0, 1};
619             int iupper[2] = { 2, 4};
620 
621             int nvalues = 12; /* 12 grid points */
622             double values[12];
623 
624             part = 1;
625 
626             for (i = 0; i < nvalues; i ++)
627                values[i] = 1.0;
628             HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
629 
630             for (i = 0; i < nvalues; i ++)
631                values[i] = 0.0;
632             HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
633          }
634       }
635       else if (myid == 1)
636       {
637          /* Set the vector coefficients over the gridpoints in my box */
638          {
639             int ilower[2] = { 3, 1};
640             int iupper[2] = { 6, 4};
641 
642             int nvalues = 16; /* 16 grid points */
643             double values[16];
644 
645             part = 2;
646 
647             for (i = 0; i < nvalues; i ++)
648                values[i] = 1.0;
649             HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
650 
651             for (i = 0; i < nvalues; i ++)
652                values[i] = 0.0;
653             HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
654          }
655       }
656 
657       /* This is a collective call finalizing the vector assembly.
658          The vectors are now ``ready to be used'' */
659       HYPRE_SStructVectorAssemble(b);
660       HYPRE_SStructVectorAssemble(x);
661    }
662 
663 
664    /* 6. Set up and use a solver (See the Reference Manual for descriptions
665       of all of the options.) */
666    {
667       /* Create an empty PCG Struct solver */
668       HYPRE_SStructPCGCreate(MPI_COMM_WORLD, &solver);
669 
670       /* Set PCG parameters */
671       HYPRE_SStructPCGSetTol(solver, 1.0e-6 );
672       HYPRE_SStructPCGSetPrintLevel(solver, 2);
673       HYPRE_SStructPCGSetMaxIter(solver, 50);
674 
675       /* Create a split SStruct solver for use as a preconditioner */
676       HYPRE_SStructSplitCreate(MPI_COMM_WORLD, &precond);
677       HYPRE_SStructSplitSetMaxIter(precond, 1);
678       HYPRE_SStructSplitSetTol(precond, 0.0);
679       HYPRE_SStructSplitSetZeroGuess(precond);
680 
681       /* Set the preconditioner type to split-SMG */
682       HYPRE_SStructSplitSetStructSolver(precond, HYPRE_SMG);
683 
684       /* Set preconditioner and solve */
685       HYPRE_SStructPCGSetPrecond(solver, HYPRE_SStructSplitSolve,
686                                  HYPRE_SStructSplitSetup, precond);
687       HYPRE_SStructPCGSetup(solver, A, b, x);
688       HYPRE_SStructPCGSolve(solver, A, b, x);
689    }
690 
691    /* Free memory */
692    HYPRE_SStructGridDestroy(grid);
693    HYPRE_SStructStencilDestroy(stencil_5pt);
694    HYPRE_SStructStencilDestroy(stencil_9pt);
695    HYPRE_SStructGraphDestroy(graph);
696    HYPRE_SStructMatrixDestroy(A);
697    HYPRE_SStructVectorDestroy(b);
698    HYPRE_SStructVectorDestroy(x);
699 
700    HYPRE_SStructPCGDestroy(solver);
701    HYPRE_SStructSplitDestroy(precond);
702 
703    /* Finalize MPI */
704    MPI_Finalize();
705 
706    return (0);
707 }


syntax highlighted by Code2HTML, v. 0.9.1