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 /* See MatrixSetObjectType below */
258 object_type = HYPRE_SSTRUCT;
259 HYPRE_SStructGraphSetObjectType(graph, object_type);
260
261 /* Use the 5-pt stencil on part 0 */
262 part = 0;
263 HYPRE_SStructGraphSetStencil(graph, part, var, stencil_5pt);
264
265 /* Use the 9-pt stencil on part 1 */
266 part = 1;
267 HYPRE_SStructGraphSetStencil(graph, part, var, stencil_9pt);
268
269 /* Use the 5-pt stencil on part 2 */
270 part = 2;
271 HYPRE_SStructGraphSetStencil(graph, part, var, stencil_5pt);
272
273 /* Since we have only stencil connections between parts, we don't need to
274 call HYPRE_SStructGraphAddEntries. */
275
276 /* Assemble the graph */
277 HYPRE_SStructGraphAssemble(graph);
278 }
279
280 /* 4. Set up a SStruct Matrix */
281 {
282 int i,j;
283 int part;
284 int var = 0;
285
286 /* Create the empty matrix object */
287 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
288
289 /* Set the object type (by default HYPRE_SSTRUCT). This determines the
290 data structure used to store the matrix. If you want to use unstructured
291 solvers, e.g. BoomerAMG, the object type should be HYPRE_PARCSR.
292 If the problem is purely structured (with one part), you may want to use
293 HYPRE_STRUCT to access the structured solvers. Since we have two parts
294 with different stencils, we set the object type to HYPRE_SSTRUCT. */
295 object_type = HYPRE_SSTRUCT;
296 HYPRE_SStructMatrixSetObjectType(A, object_type);
297
298 /* Get ready to set values */
299 HYPRE_SStructMatrixInitialize(A);
300
301 /* Each processor must set the stencil values for their boxes on each part.
302 In this example, we only set stencil entries and therefore use
303 HYPRE_SStructMatrixSetBoxValues. If we need to set non-stencil entries,
304 we have to use HYPRE_SStructMatrixSetValues. */
305
306 if (myid == 0)
307 {
308 /* Set the matrix coefficients for some set of stencil entries
309 over all the gridpoints in my first box (account for boundary
310 grid points later) */
311 {
312 int ilower[2] = {-3, 1};
313 int iupper[2] = {-1, 2};
314
315 int nentries = 5;
316 int nvalues = 30; /* 6 grid points, each with 5 stencil entries */
317 double values[30];
318
319 int stencil_indices[5];
320 for (j = 0; j < nentries; j++) /* label the stencil indices -
321 these correspond to the offsets
322 defined above */
323 stencil_indices[j] = j;
324
325 for (i = 0; i < nvalues; i += nentries)
326 {
327 values[i] = 4.0;
328 for (j = 1; j < nentries; j++)
329 values[i+j] = -1.0;
330 }
331
332 part = 0;
333 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
334 var, nentries,
335 stencil_indices, values);
336 }
337
338 /* Set the matrix coefficients for some set of stencil entries
339 over the gridpoints in my second box */
340 {
341 int ilower[2] = {0, 1};
342 int iupper[2] = {2, 4};
343
344 int nentries = 9;
345 int nvalues = 108; /* 12 grid points, each with 5 stencil entries */
346 double values[108];
347
348 int stencil_indices[9];
349 for (j = 0; j < nentries; j++)
350 stencil_indices[j] = j;
351
352 for (i = 0; i < nvalues; i += nentries)
353 {
354 values[i] = 8./3.;
355 for (j = 1; j < nentries; j++)
356 values[i+j] = -1./3.;
357 }
358
359 part = 1;
360 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
361 var, nentries,
362 stencil_indices, values);
363 }
364 }
365 else if (myid == 1)
366 {
367 /* Set the matrix coefficients for some set of stencil entries
368 over the gridpoints in my box */
369 {
370 int ilower[2] = {3, 1};
371 int iupper[2] = {6, 4};
372
373 int nentries = 5;
374 int nvalues = 80; /* 16 grid points, each with 5 stencil entries */
375 double values[80];
376
377 int stencil_indices[5];
378 for (j = 0; j < nentries; j++)
379 stencil_indices[j] = j;
380
381 for (i = 0; i < nvalues; i += nentries)
382 {
383 values[i] = 4.0;
384 for (j = 1; j < nentries; j++)
385 values[i+j] = -1.0;
386 }
387
388 part = 2;
389 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
390 var, nentries,
391 stencil_indices, values);
392 }
393 }
394
395 /* Modify the 9-pt stencil on the boundary between parts to ensure
396 symmetry and good global approximation. */
397 if (myid == 0)
398 {
399 int nentries = 6;
400 int nvalues = 24; /* 4 grid points, each with 6 stencil entries */
401 double values[24];
402
403 part = 1;
404
405 for (i = 0; i < nvalues; i += nentries)
406 {
407 values[i] = 10./3.;
408 values[i+1] = -1.;
409 values[i+2] = -2./3.;
410 values[i+3] = -2./3.;
411 values[i+4] = 0.0;
412 values[i+5] = 0.0;
413 }
414
415 {
416 /* Values to the right of the second box */
417 int ilower[2] = { 2, 1};
418 int iupper[2] = { 2, 4};
419
420 int stencil_indices[6] = {0,2,3,4,6,7};
421
422 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
423 var, nentries,
424 stencil_indices, values);
425 }
426
427 {
428 /* Values to the left of the second box */
429 int ilower[2] = { 0, 1};
430 int iupper[2] = { 0, 4};
431
432 int stencil_indices[6] = {0,1,3,4,5,8};
433
434 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
435 var, nentries,
436 stencil_indices, values);
437 }
438 }
439
440 /* For each box, set any coefficients that reach ouside of the
441 boundary to 0 */
442 if (myid == 0)
443 {
444 int maxnvalues = 9;
445 double values[9];
446
447 for (i = 0; i < maxnvalues; i++)
448 values[i] = 0.0;
449
450 part = 0;
451
452 {
453 /* Values below our first box */
454 int ilower[2] = {-3, 1};
455 int iupper[2] = {-1, 1};
456
457 int stencil_indices[1] = {3};
458
459 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
460 var, 1,
461 stencil_indices, values);
462 }
463
464 {
465 /* Values to the left of our first box */
466 int ilower[2] = {-3, 1};
467 int iupper[2] = {-3, 2};
468
469 int stencil_indices[1] = {1};
470
471 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
472 var, 1,
473 stencil_indices, values);
474 }
475
476 {
477 /* Values above our first box */
478 int ilower[2] = {-3, 2};
479 int iupper[2] = {-1, 2};
480
481 int stencil_indices[1] = {4};
482
483 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
484 var, 1,
485 stencil_indices, values);
486 }
487
488 part = 1;
489
490 {
491 /* Values below our second box */
492 int ilower[2] = { 0, 1};
493 int iupper[2] = { 2, 1};
494
495 int stencil_indices[3] = {3,5,6};
496
497 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
498 var, 3,
499 stencil_indices, values);
500 }
501
502 {
503 /* Values to the left of our second box (that do not border the
504 first box). */
505 int ilower[2] = { 0, 3};
506 int iupper[2] = { 0, 4};
507
508 int stencil_indices[3] = {1,5,8};
509
510 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
511 var, 3,
512 stencil_indices, values);
513 }
514
515 {
516 /* Values above our second box */
517 int ilower[2] = { 0, 4};
518 int iupper[2] = { 2, 4};
519
520 int stencil_indices[3] = {4,7,8};
521
522 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
523 var, 3,
524 stencil_indices, values);
525 }
526 }
527 else if (myid == 1)
528 {
529 int maxnvalues = 4;
530 double values[4];
531
532 for (i = 0; i < maxnvalues; i++)
533 values[i] = 0.0;
534
535 part = 2;
536
537 {
538 /* Values below our box */
539 int ilower[2] = { 3, 1};
540 int iupper[2] = { 6, 1};
541
542 int stencil_indices[1] = {3};
543
544 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
545 var, 1,
546 stencil_indices, values);
547 }
548
549 {
550 /* Values to the right of our box */
551 int ilower[2] = { 6, 1};
552 int iupper[2] = { 6, 4};
553
554 int stencil_indices[1] = {2};
555
556 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
557 var, 1,
558 stencil_indices, values);
559 }
560
561 {
562 /* Values above our box */
563 int ilower[2] = { 3, 4};
564 int iupper[2] = { 6, 4};
565
566 int stencil_indices[1] = {4};
567
568 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
569 var, 1,
570 stencil_indices, values);
571 }
572 }
573
574 /* This is a collective call finalizing the matrix assembly.
575 The matrix is now ``ready to be used'' */
576 HYPRE_SStructMatrixAssemble(A);
577 }
578
579 /* 5. Set up SStruct Vectors for b and x */
580 {
581 int i;
582 int part;
583 int var = 0;
584
585 /* Create an empty vector object */
586 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
587 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
588
589 /* As with the matrix, set the object type for the vectors
590 to be the sstruct type */
591 object_type = HYPRE_SSTRUCT;
592 HYPRE_SStructVectorSetObjectType(b, object_type);
593 HYPRE_SStructVectorSetObjectType(x, object_type);
594
595 /* Indicate that the vector coefficients are ready to be set */
596 HYPRE_SStructVectorInitialize(b);
597 HYPRE_SStructVectorInitialize(x);
598
599 if (myid == 0)
600 {
601 /* Set the vector coefficients over the gridpoints in my first box */
602 {
603 int ilower[2] = {-3, 1};
604 int iupper[2] = {-1, 2};
605
606 int nvalues = 6; /* 6 grid points */
607 double values[6];
608
609 part = 0;
610
611 for (i = 0; i < nvalues; i ++)
612 values[i] = 1.0;
613 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
614
615 for (i = 0; i < nvalues; i ++)
616 values[i] = 0.0;
617 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
618 }
619
620 /* Set the vector coefficients over the gridpoints in my second box */
621 {
622 int ilower[2] = { 0, 1};
623 int iupper[2] = { 2, 4};
624
625 int nvalues = 12; /* 12 grid points */
626 double values[12];
627
628 part = 1;
629
630 for (i = 0; i < nvalues; i ++)
631 values[i] = 1.0;
632 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
633
634 for (i = 0; i < nvalues; i ++)
635 values[i] = 0.0;
636 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
637 }
638 }
639 else if (myid == 1)
640 {
641 /* Set the vector coefficients over the gridpoints in my box */
642 {
643 int ilower[2] = { 3, 1};
644 int iupper[2] = { 6, 4};
645
646 int nvalues = 16; /* 16 grid points */
647 double values[16];
648
649 part = 2;
650
651 for (i = 0; i < nvalues; i ++)
652 values[i] = 1.0;
653 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
654
655 for (i = 0; i < nvalues; i ++)
656 values[i] = 0.0;
657 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
658 }
659 }
660
661 /* This is a collective call finalizing the vector assembly.
662 The vectors are now ``ready to be used'' */
663 HYPRE_SStructVectorAssemble(b);
664 HYPRE_SStructVectorAssemble(x);
665 }
666
667
668 /* 6. Set up and use a solver (See the Reference Manual for descriptions
669 of all of the options.) */
670 {
671 /* Create an empty PCG Struct solver */
672 HYPRE_SStructPCGCreate(MPI_COMM_WORLD, &solver);
673
674 /* Set PCG parameters */
675 HYPRE_SStructPCGSetTol(solver, 1.0e-6 );
676 HYPRE_SStructPCGSetPrintLevel(solver, 2);
677 HYPRE_SStructPCGSetMaxIter(solver, 50);
678
679 /* Create a split SStruct solver for use as a preconditioner */
680 HYPRE_SStructSplitCreate(MPI_COMM_WORLD, &precond);
681 HYPRE_SStructSplitSetMaxIter(precond, 1);
682 HYPRE_SStructSplitSetTol(precond, 0.0);
683 HYPRE_SStructSplitSetZeroGuess(precond);
684
685 /* Set the preconditioner type to split-SMG */
686 HYPRE_SStructSplitSetStructSolver(precond, HYPRE_SMG);
687
688 /* Set preconditioner and solve */
689 HYPRE_SStructPCGSetPrecond(solver, HYPRE_SStructSplitSolve,
690 HYPRE_SStructSplitSetup, precond);
691 HYPRE_SStructPCGSetup(solver, A, b, x);
692 HYPRE_SStructPCGSolve(solver, A, b, x);
693 }
694
695 /* Free memory */
696 HYPRE_SStructGridDestroy(grid);
697 HYPRE_SStructStencilDestroy(stencil_5pt);
698 HYPRE_SStructStencilDestroy(stencil_9pt);
699 HYPRE_SStructGraphDestroy(graph);
700 HYPRE_SStructMatrixDestroy(A);
701 HYPRE_SStructVectorDestroy(b);
702 HYPRE_SStructVectorDestroy(x);
703
704 HYPRE_SStructPCGDestroy(solver);
705 HYPRE_SStructSplitDestroy(precond);
706
707 /* Finalize MPI */
708 MPI_Finalize();
709
710 return (0);
711 }
syntax highlighted by Code2HTML, v. 0.9.1