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