download the original source code.
1 /*
2 Example 6
3
4 Interface: Semi-Structured interface (SStruct)
5
6 Compile with: make ex6
7
8 Sample run: mpirun -np 2 ex6
9
10 Description: This is a two processor example and is the same problem
11 as is solved with the structured interface in Example 2.
12 (The grid boxes are exactly those in the example
13 diagram in the struct interface chapter of the User's Manual.
14 Processor 0 owns two boxes and processor 1 owns one box.)
15
16 This is the simplest sstruct example, and it demonstrates how
17 the semi-structured interface can be used for structured problems.
18 There is one part and one variable. The solver is PCG with SMG
19 preconditioner. We use a structured solver for this example.
20 */
21
22 #include <stdio.h>
23
24 /* SStruct linear solvers headers */
25 #include "HYPRE_sstruct_ls.h"
26
27 #include "vis.c"
28
29 int main (int argc, char *argv[])
30 {
31 int myid, num_procs;
32
33 int vis = 0;
34
35 HYPRE_SStructGrid grid;
36 HYPRE_SStructGraph graph;
37 HYPRE_SStructStencil stencil;
38 HYPRE_SStructMatrix A;
39 HYPRE_SStructVector b;
40 HYPRE_SStructVector x;
41
42 /* We are using struct solvers for this example */
43 HYPRE_StructSolver solver;
44 HYPRE_StructSolver precond;
45
46 int object_type;
47
48 /* Initialize MPI */
49 MPI_Init(&argc, &argv);
50 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
51 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
52
53 if (num_procs != 2)
54 {
55 if (myid ==0) printf("Must run with 2 processors!\n");
56 MPI_Finalize();
57
58 return(0);
59 }
60
61 /* Parse command line */
62 {
63 int arg_index = 0;
64 int print_usage = 0;
65
66 while (arg_index < argc)
67 {
68 if ( strcmp(argv[arg_index], "-vis") == 0 )
69 {
70 arg_index++;
71 vis = 1;
72 }
73 else if ( strcmp(argv[arg_index], "-help") == 0 )
74 {
75 print_usage = 1;
76 break;
77 }
78 else
79 {
80 arg_index++;
81 }
82 }
83
84 if ((print_usage) && (myid == 0))
85 {
86 printf("\n");
87 printf("Usage: %s [<options>]\n", argv[0]);
88 printf("\n");
89 printf(" -vis : save the solution for GLVis visualization\n");
90 printf("\n");
91 }
92
93 if (print_usage)
94 {
95 MPI_Finalize();
96 return (0);
97 }
98 }
99
100 /* 1. Set up the 2D grid. This gives the index space in each part.
101 Here we only use one part and one variable. (So the part id is 0
102 and the variable id is 0) */
103 {
104 int ndim = 2;
105 int nparts = 1;
106 int part = 0;
107
108 /* Create an empty 2D grid object */
109 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &grid);
110
111 /* Set the extents of the grid - each processor sets its grid
112 boxes. Each part has its own relative index space numbering,
113 but in this example all boxes belong to the same part. */
114
115 /* Processor 0 owns two boxes in the grid. */
116 if (myid == 0)
117 {
118 /* Add a new box to the grid */
119 {
120 int ilower[2] = {-3, 1};
121 int iupper[2] = {-1, 2};
122
123 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
124 }
125
126 /* Add a new box to the grid */
127 {
128 int ilower[2] = {0, 1};
129 int iupper[2] = {2, 4};
130
131 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
132 }
133 }
134
135 /* Processor 1 owns one box in the grid. */
136 else if (myid == 1)
137 {
138 /* Add a new box to the grid */
139 {
140 int ilower[2] = {3, 1};
141 int iupper[2] = {6, 4};
142
143 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
144 }
145 }
146
147 /* Set the variable type and number of variables on each part. */
148 {
149 int i;
150 int nvars = 1;
151 HYPRE_SStructVariable vartypes[1] = {HYPRE_SSTRUCT_VARIABLE_CELL};
152
153 for (i = 0; i< nparts; i++)
154 HYPRE_SStructGridSetVariables(grid, i, nvars, vartypes);
155 }
156
157 /* Now the grid is ready to use */
158 HYPRE_SStructGridAssemble(grid);
159 }
160
161 /* 2. Define the discretization stencil(s) */
162 {
163 /* Create an empty 2D, 5-pt stencil object */
164 HYPRE_SStructStencilCreate(2, 5, &stencil);
165
166 /* Define the geometry of the stencil. Each represents a
167 relative offset (in the index space). */
168 {
169 int entry;
170 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
171 int var = 0;
172
173 /* Assign numerical values to the offsets so that we can
174 easily refer to them - the last argument indicates the
175 variable for which we are assigning this stencil - we are
176 just using one variable in this example so it is the first one (0) */
177 for (entry = 0; entry < 5; entry++)
178 HYPRE_SStructStencilSetEntry(stencil, entry, offsets[entry], var);
179 }
180 }
181
182 /* 3. Set up the Graph - this determines the non-zero structure
183 of the matrix and allows non-stencil relationships between the parts */
184 {
185 int var = 0;
186 int part = 0;
187
188 /* Create the graph object */
189 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
190
191 /* See MatrixSetObjectType below */
192 object_type = HYPRE_STRUCT;
193 HYPRE_SStructGraphSetObjectType(graph, object_type);
194
195 /* Now we need to tell the graph which stencil to use for each
196 variable on each part (we only have one variable and one part) */
197 HYPRE_SStructGraphSetStencil(graph, part, var, stencil);
198
199 /* Here we could establish connections between parts if we
200 had more than one part using the graph. For example, we could
201 use HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborBox() */
202
203 /* Assemble the graph */
204 HYPRE_SStructGraphAssemble(graph);
205 }
206
207 /* 4. Set up a SStruct Matrix */
208 {
209 int i,j;
210 int part = 0;
211 int var = 0;
212
213 /* Create the empty matrix object */
214 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
215
216 /* Set the object type (by default HYPRE_SSTRUCT). This determines the
217 data structure used to store the matrix. If you want to use unstructured
218 solvers, e.g. BoomerAMG, the object type should be HYPRE_PARCSR.
219 If the problem is purely structured (with one part), you may want to use
220 HYPRE_STRUCT to access the structured solvers. Here we have a purely
221 structured example. */
222 object_type = HYPRE_STRUCT;
223 HYPRE_SStructMatrixSetObjectType(A, object_type);
224
225 /* Get ready to set values */
226 HYPRE_SStructMatrixInitialize(A);
227
228 /* Each processor must set the stencil values for their boxes on each part.
229 In this example, we only set stencil entries and therefore use
230 HYPRE_SStructMatrixSetBoxValues. If we need to set non-stencil entries,
231 we have to use HYPRE_SStructMatrixSetValues (shown in a later example). */
232
233 if (myid == 0)
234 {
235 /* Set the matrix coefficients for some set of stencil entries
236 over all the gridpoints in my first box (account for boundary
237 grid points later) */
238 {
239 int ilower[2] = {-3, 1};
240 int iupper[2] = {-1, 2};
241
242 int nentries = 5;
243 int nvalues = 30; /* 6 grid points, each with 5 stencil entries */
244 double values[30];
245
246 int stencil_indices[5];
247 for (j = 0; j < nentries; j++) /* label the stencil indices -
248 these correspond to the offsets
249 defined above */
250 stencil_indices[j] = j;
251
252 for (i = 0; i < nvalues; i += nentries)
253 {
254 values[i] = 4.0;
255 for (j = 1; j < nentries; j++)
256 values[i+j] = -1.0;
257 }
258
259 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
260 var, nentries,
261 stencil_indices, values);
262 }
263
264 /* Set the matrix coefficients for some set of stencil entries
265 over the gridpoints in my second box */
266 {
267 int ilower[2] = {0, 1};
268 int iupper[2] = {2, 4};
269
270 int nentries = 5;
271 int nvalues = 60; /* 12 grid points, each with 5 stencil entries */
272 double values[60];
273
274 int stencil_indices[5];
275 for (j = 0; j < nentries; j++)
276 stencil_indices[j] = j;
277
278 for (i = 0; i < nvalues; i += nentries)
279 {
280 values[i] = 4.0;
281 for (j = 1; j < nentries; j++)
282 values[i+j] = -1.0;
283 }
284
285 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
286 var, nentries,
287 stencil_indices, values);
288 }
289 }
290 else if (myid == 1)
291 {
292 /* Set the matrix coefficients for some set of stencil entries
293 over the gridpoints in my box */
294 {
295 int ilower[2] = {3, 1};
296 int iupper[2] = {6, 4};
297
298 int nentries = 5;
299 int nvalues = 80; /* 16 grid points, each with 5 stencil entries */
300 double values[80];
301
302 int stencil_indices[5];
303 for (j = 0; j < nentries; j++)
304 stencil_indices[j] = j;
305
306 for (i = 0; i < nvalues; i += nentries)
307 {
308 values[i] = 4.0;
309 for (j = 1; j < nentries; j++)
310 values[i+j] = -1.0;
311 }
312
313 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
314 var, nentries,
315 stencil_indices, values);
316 }
317 }
318
319 /* For each box, set any coefficients that reach ouside of the
320 boundary to 0 */
321 if (myid == 0)
322 {
323 int maxnvalues = 6;
324 double values[6];
325
326 for (i = 0; i < maxnvalues; i++)
327 values[i] = 0.0;
328
329 {
330 /* Values below our first AND second box */
331 int ilower[2] = {-3, 1};
332 int iupper[2] = { 2, 1};
333
334 int stencil_indices[1] = {3};
335
336 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
337 var, 1,
338 stencil_indices, values);
339 }
340
341 {
342 /* Values to the left of our first box */
343 int ilower[2] = {-3, 1};
344 int iupper[2] = {-3, 2};
345
346 int stencil_indices[1] = {1};
347
348 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
349 var, 1,
350 stencil_indices, values);
351 }
352
353 {
354 /* Values above our first box */
355 int ilower[2] = {-3, 2};
356 int iupper[2] = {-1, 2};
357
358 int stencil_indices[1] = {4};
359
360 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
361 var, 1,
362 stencil_indices, values);
363 }
364
365 {
366 /* Values to the left of our second box (that do not border the
367 first box). */
368 int ilower[2] = { 0, 3};
369 int iupper[2] = { 0, 4};
370
371 int stencil_indices[1] = {1};
372
373 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
374 var, 1,
375 stencil_indices, values);
376 }
377
378 {
379 /* Values above our second box */
380 int ilower[2] = { 0, 4};
381 int iupper[2] = { 2, 4};
382
383 int stencil_indices[1] = {4};
384
385 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
386 var, 1,
387 stencil_indices, values);
388 }
389 }
390 else if (myid == 1)
391 {
392 int maxnvalues = 4;
393 double values[4];
394 for (i = 0; i < maxnvalues; i++)
395 values[i] = 0.0;
396
397 {
398 /* Values below our box */
399 int ilower[2] = { 3, 1};
400 int iupper[2] = { 6, 1};
401
402 int stencil_indices[1] = {3};
403
404 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
405 var, 1,
406 stencil_indices, values);
407 }
408
409 {
410 /* Values to the right of our box */
411 int ilower[2] = { 6, 1};
412 int iupper[2] = { 6, 4};
413
414 int stencil_indices[1] = {2};
415
416 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
417 var, 1,
418 stencil_indices, values);
419 }
420
421 {
422 /* Values above our box */
423 int ilower[2] = { 3, 4};
424 int iupper[2] = { 6, 4};
425
426 int stencil_indices[1] = {4};
427
428 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
429 var, 1,
430 stencil_indices, values);
431 }
432 }
433
434 /* This is a collective call finalizing the matrix assembly.
435 The matrix is now ``ready to be used'' */
436 HYPRE_SStructMatrixAssemble(A);
437 }
438
439
440 /* 5. Set up SStruct Vectors for b and x */
441 {
442 int i;
443
444 /* We have one part and one variable. */
445 int part = 0;
446 int var = 0;
447
448 /* Create an empty vector object */
449 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
450 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
451
452 /* As with the matrix, set the object type for the vectors
453 to be the struct type */
454 object_type = HYPRE_STRUCT;
455 HYPRE_SStructVectorSetObjectType(b, object_type);
456 HYPRE_SStructVectorSetObjectType(x, object_type);
457
458 /* Indicate that the vector coefficients are ready to be set */
459 HYPRE_SStructVectorInitialize(b);
460 HYPRE_SStructVectorInitialize(x);
461
462 if (myid == 0)
463 {
464 /* Set the vector coefficients over the gridpoints in my first box */
465 {
466 int ilower[2] = {-3, 1};
467 int iupper[2] = {-1, 2};
468
469 int nvalues = 6; /* 6 grid points */
470 double values[6];
471
472 for (i = 0; i < nvalues; i ++)
473 values[i] = 1.0;
474 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
475
476 for (i = 0; i < nvalues; i ++)
477 values[i] = 0.0;
478 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
479 }
480
481 /* Set the vector coefficients over the gridpoints in my second box */
482 {
483 int ilower[2] = { 0, 1};
484 int iupper[2] = { 2, 4};
485
486 int nvalues = 12; /* 12 grid points */
487 double values[12];
488
489 for (i = 0; i < nvalues; i ++)
490 values[i] = 1.0;
491 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
492
493 for (i = 0; i < nvalues; i ++)
494 values[i] = 0.0;
495 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
496 }
497 }
498 else if (myid == 1)
499 {
500 /* Set the vector coefficients over the gridpoints in my box */
501 {
502 int ilower[2] = { 3, 1};
503 int iupper[2] = { 6, 4};
504
505 int nvalues = 16; /* 16 grid points */
506 double values[16];
507
508 for (i = 0; i < nvalues; i ++)
509 values[i] = 1.0;
510 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
511
512 for (i = 0; i < nvalues; i ++)
513 values[i] = 0.0;
514 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
515 }
516 }
517
518 /* This is a collective call finalizing the vector assembly.
519 The vectors are now ``ready to be used'' */
520 HYPRE_SStructVectorAssemble(b);
521 HYPRE_SStructVectorAssemble(x);
522 }
523
524 /* 6. Set up and use a solver (See the Reference Manual for descriptions
525 of all of the options.) */
526 {
527 HYPRE_StructMatrix sA;
528 HYPRE_StructVector sb;
529 HYPRE_StructVector sx;
530
531 /* Because we are using a struct solver, we need to get the
532 object of the matrix and vectors to pass in to the struct solvers */
533 HYPRE_SStructMatrixGetObject(A, (void **) &sA);
534 HYPRE_SStructVectorGetObject(b, (void **) &sb);
535 HYPRE_SStructVectorGetObject(x, (void **) &sx);
536
537 /* Create an empty PCG Struct solver */
538 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
539
540 /* Set PCG parameters */
541 HYPRE_StructPCGSetTol(solver, 1.0e-06);
542 HYPRE_StructPCGSetPrintLevel(solver, 2);
543 HYPRE_StructPCGSetMaxIter(solver, 50);
544
545 /* Create the Struct SMG solver for use as a preconditioner */
546 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
547
548 /* Set SMG parameters */
549 HYPRE_StructSMGSetMaxIter(precond, 1);
550 HYPRE_StructSMGSetTol(precond, 0.0);
551 HYPRE_StructSMGSetZeroGuess(precond);
552 HYPRE_StructSMGSetNumPreRelax(precond, 1);
553 HYPRE_StructSMGSetNumPostRelax(precond, 1);
554
555 /* Set preconditioner and solve */
556 HYPRE_StructPCGSetPrecond(solver, HYPRE_StructSMGSolve,
557 HYPRE_StructSMGSetup, precond);
558 HYPRE_StructPCGSetup(solver, sA, sb, sx);
559 HYPRE_StructPCGSolve(solver, sA, sb, sx);
560 }
561
562 /* Save the solution for GLVis visualization, see vis/glvis-ex6.sh */
563 if (vis)
564 {
565 GLVis_PrintSStructGrid(grid, "vis/ex6.mesh", myid, NULL, NULL);
566 GLVis_PrintSStructVector(x, 0, "vis/ex6.sol", myid);
567 GLVis_PrintData("vis/ex6.data", myid, num_procs);
568 }
569
570 /* Free memory */
571 HYPRE_SStructGridDestroy(grid);
572 HYPRE_SStructStencilDestroy(stencil);
573 HYPRE_SStructGraphDestroy(graph);
574 HYPRE_SStructMatrixDestroy(A);
575 HYPRE_SStructVectorDestroy(b);
576 HYPRE_SStructVectorDestroy(x);
577
578 HYPRE_StructPCGDestroy(solver);
579 HYPRE_StructSMGDestroy(precond);
580
581 /* Finalize MPI */
582 MPI_Finalize();
583
584 return (0);
585 }
syntax highlighted by Code2HTML, v. 0.9.1