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


syntax highlighted by Code2HTML, v. 0.9.1