download the original source code.
1 /*
2 Example 12
3
4 Interface: Semi-Structured interface (SStruct)
5
6 Compile with: make ex12 (may need to edit HYPRE_DIR in Makefile)
7
8 Sample runs: mpirun -np 2 ex12 -pfmg
9 mpirun -np 2 ex12 -boomeramg
10
11 Description: The grid layout is the same as ex1, but with nodal unknowns. The
12 solver is PCG preconditioned with either PFMG or BoomerAMG,
13 selected on the command line.
14
15 We recommend viewing the Struct examples before viewing this
16 and the other SStruct examples. This is one of the simplest
17 SStruct examples, used primarily to demonstrate how to set up
18 non-cell-centered problems, and to demonstrate how easy it is
19 to switch between structured solvers (PFMG) and solvers
20 designed for more general settings (AMG).
21 */
22
23 #include <stdio.h>
24 #include <stdlib.h>
25 #include <string.h>
26
27 #include "HYPRE_sstruct_ls.h"
28 #include "HYPRE_parcsr_ls.h"
29 #include "HYPRE_krylov.h"
30
31 int main (int argc, char *argv[])
32 {
33 int i, j, myid, num_procs;
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 only have one part and one variable */
43 int nparts = 1;
44 int nvars = 1;
45 int part = 0;
46 int var = 0;
47
48 int precond_id, object_type;
49
50 /* Initialize MPI */
51 MPI_Init(&argc, &argv);
52 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
53 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
54
55 if (num_procs != 2)
56 {
57 if (myid == 0) printf("Must run with 2 processors!\n");
58 exit(1);
59 }
60
61 /* Parse the command line to determine the solver */
62 if (argc != 2)
63 {
64 if (myid == 0) printf("Must specify a solver!\n");
65 exit(1);
66 }
67 if ( strcmp(argv[1], "-pfmg") == 0 )
68 {
69 precond_id = 1;
70 object_type = HYPRE_STRUCT;
71 }
72 else if ( strcmp(argv[1], "-boomeramg") == 0 )
73 {
74 precond_id = 2;
75 object_type = HYPRE_PARCSR;
76 }
77 else
78 {
79 if (myid == 0) printf("Invalid solver!\n");
80 exit(1);
81 }
82
83 /* 1. Set up the grid. Here we use only one part. Each processor describes
84 the piece of the grid that it owns. */
85 {
86 /* Create an empty 2D grid object */
87 HYPRE_SStructGridCreate(MPI_COMM_WORLD, 2, nparts, &grid);
88
89 /* Add boxes to the grid */
90 if (myid == 0)
91 {
92 int ilower[2]={-3,1}, iupper[2]={-1,2};
93 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
94 }
95 else if (myid == 1)
96 {
97 int ilower[2]={0,1}, iupper[2]={2,4};
98 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
99 }
100
101 /* Set the variable type and number of variables on each part. */
102 {
103 HYPRE_SStructVariable vartypes[1] = {HYPRE_SSTRUCT_VARIABLE_NODE};
104
105 HYPRE_SStructGridSetVariables(grid, part, nvars, vartypes);
106 }
107
108 /* This is a collective call finalizing the grid assembly.
109 The grid is now ``ready to be used'' */
110 HYPRE_SStructGridAssemble(grid);
111 }
112
113 /* 2. Define the discretization stencil */
114 {
115 /* Create an empty 2D, 5-pt stencil object */
116 HYPRE_SStructStencilCreate(2, 5, &stencil);
117
118 /* Define the geometry of the stencil. Each represents a relative offset
119 (in the index space). */
120 {
121 int entry;
122 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
123
124 /* Assign numerical values to the offsets so that we can easily refer
125 to them - the last argument indicates the variable for which we are
126 assigning this stencil */
127 for (entry = 0; entry < 5; entry++)
128 HYPRE_SStructStencilSetEntry(stencil, entry, offsets[entry], var);
129 }
130 }
131
132 /* 3. Set up the Graph - this determines the non-zero structure of the matrix
133 and allows non-stencil relationships between the parts */
134 {
135 /* Create the graph object */
136 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
137
138 /* Now we need to tell the graph which stencil to use for each variable on
139 each part (we only have one variable and one part) */
140 HYPRE_SStructGraphSetStencil(graph, part, var, stencil);
141
142 /* Here we could establish connections between parts if we had more than
143 one part using the graph. For example, we could use
144 HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborPart() */
145
146 /* Assemble the graph */
147 HYPRE_SStructGraphAssemble(graph);
148 }
149
150 /* 4. Set up a SStruct Matrix */
151 {
152 /* Create an empty matrix object */
153 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
154
155 /* Set the object type (by default HYPRE_SSTRUCT). This determines the
156 data structure used to store the matrix. For PFMG we need to use
157 HYPRE_STRUCT, and for BoomerAMG we need HYPRE_PARCSR (set above). */
158 HYPRE_SStructMatrixSetObjectType(A, object_type);
159
160 /* Get ready to set values */
161 HYPRE_SStructMatrixInitialize(A);
162
163 /* Set the matrix coefficients. Each processor assigns coefficients for
164 the boxes in the grid that it owns. Note that the coefficients
165 associated with each stencil entry may vary from grid point to grid
166 point if desired. Here, we first set the same stencil entries for each
167 grid point. Then we make modifications to grid points near the
168 boundary. Note that the ilower values are different from those used in
169 ex1 because of the way nodal variables are referenced. Also note that
170 some of the stencil values are set on both processor 0 and processor 1.
171 See the User and Reference manuals for more details. */
172 if (myid == 0)
173 {
174 int ilower[2]={-4,0}, iupper[2]={-1,2};
175 int stencil_indices[5] = {0,1,2,3,4}; /* labels for the stencil entries -
176 these correspond to the offsets
177 defined above */
178 int nentries = 5;
179 int nvalues = 60; /* 12 grid points, each with 5 stencil entries */
180 double values[60];
181
182 for (i = 0; i < nvalues; i += nentries)
183 {
184 values[i] = 4.0;
185 for (j = 1; j < nentries; j++)
186 values[i+j] = -1.0;
187 }
188
189 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, nentries,
190 stencil_indices, values);
191 }
192 else if (myid == 1)
193 {
194 int ilower[2]={-1,0}, iupper[2]={2,4};
195 int stencil_indices[5] = {0,1,2,3,4};
196 int nentries = 5;
197 int nvalues = 100; /* 20 grid points, each with 5 stencil entries */
198 double values[100];
199
200 for (i = 0; i < nvalues; i += nentries)
201 {
202 values[i] = 4.0;
203 for (j = 1; j < nentries; j++)
204 values[i+j] = -1.0;
205 }
206
207 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, nentries,
208 stencil_indices, values);
209 }
210
211 /* Set the coefficients reaching outside of the boundary to 0. Note that
212 * both ilower *and* iupper may be different from those in ex1. */
213 if (myid == 0)
214 {
215 double values[4];
216 for (i = 0; i < 4; i++)
217 values[i] = 0.0;
218 {
219 /* values below our box */
220 int ilower[2]={-4,0}, iupper[2]={-1,0};
221 int stencil_indices[1] = {3};
222 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
223 stencil_indices, values);
224 }
225 {
226 /* values to the left of our box */
227 int ilower[2]={-4,0}, iupper[2]={-4,2};
228 int stencil_indices[1] = {1};
229 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
230 stencil_indices, values);
231 }
232 {
233 /* values above our box */
234 int ilower[2]={-4,2}, iupper[2]={-2,2};
235 int stencil_indices[1] = {4};
236 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
237 stencil_indices, values);
238 }
239 }
240 else if (myid == 1)
241 {
242 double values[5];
243 for (i = 0; i < 5; i++)
244 values[i] = 0.0;
245 {
246 /* values below our box */
247 int ilower[2]={-1,0}, iupper[2]={2,0};
248 int stencil_indices[1] = {3};
249 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
250 stencil_indices, values);
251 }
252 {
253 /* values to the right of our box */
254 int ilower[2]={2,0}, iupper[2]={2,4};
255 int stencil_indices[1] = {2};
256 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
257 stencil_indices, values);
258 }
259 {
260 /* values above our box */
261 int ilower[2]={-1,4}, iupper[2]={2,4};
262 int stencil_indices[1] = {4};
263 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
264 stencil_indices, values);
265 }
266 {
267 /* values to the left of our box
268 (that do not border the other box on proc. 0) */
269 int ilower[2]={-1,3}, iupper[2]={-1,4};
270 int stencil_indices[1] = {1};
271 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
272 stencil_indices, values);
273 }
274 }
275
276 /* This is a collective call finalizing the matrix assembly.
277 The matrix is now ``ready to be used'' */
278 HYPRE_SStructMatrixAssemble(A);
279 }
280
281 /* 5. Set up SStruct Vectors for b and x. */
282 {
283 /* Create an empty vector object */
284 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
285 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
286
287 /* As with the matrix, set the appropriate object type for the vectors */
288 HYPRE_SStructVectorSetObjectType(b, object_type);
289 HYPRE_SStructVectorSetObjectType(x, object_type);
290
291 /* Indicate that the vector coefficients are ready to be set */
292 HYPRE_SStructVectorInitialize(b);
293 HYPRE_SStructVectorInitialize(x);
294
295 /* Set the vector coefficients. Again, note that the ilower values are
296 different from those used in ex1, and some of the values are set on
297 both processors. */
298 if (myid == 0)
299 {
300 int ilower[2]={-4,0}, iupper[2]={-1,2};
301 double values[12]; /* 12 grid points */
302
303 for (i = 0; i < 12; i ++)
304 values[i] = 1.0;
305 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
306
307 for (i = 0; i < 12; i ++)
308 values[i] = 0.0;
309 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
310 }
311 else if (myid == 1)
312 {
313 int ilower[2]={0,1}, iupper[2]={2,4};
314 double values[20]; /* 20 grid points */
315
316 for (i = 0; i < 20; i ++)
317 values[i] = 1.0;
318 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
319
320 for (i = 0; i < 20; i ++)
321 values[i] = 0.0;
322 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
323 }
324
325 /* This is a collective call finalizing the vector assembly.
326 The vectors are now ``ready to be used'' */
327 HYPRE_SStructVectorAssemble(b);
328 HYPRE_SStructVectorAssemble(x);
329 }
330
331 /* 6. Set up and use a solver (See the Reference Manual for descriptions
332 of all of the options.) */
333 if (precond_id == 1) /* PFMG */
334 {
335 HYPRE_StructMatrix sA;
336 HYPRE_StructVector sb;
337 HYPRE_StructVector sx;
338
339 HYPRE_StructSolver solver;
340 HYPRE_StructSolver precond;
341
342 /* Because we are using a struct solver, we need to get the
343 object of the matrix and vectors to pass in to the struct solvers */
344 HYPRE_SStructMatrixGetObject(A, (void **) &sA);
345 HYPRE_SStructVectorGetObject(b, (void **) &sb);
346 HYPRE_SStructVectorGetObject(x, (void **) &sx);
347
348 /* Create an empty PCG Struct solver */
349 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
350
351 /* Set PCG parameters */
352 HYPRE_StructPCGSetTol(solver, 1.0e-06);
353 HYPRE_StructPCGSetPrintLevel(solver, 2);
354 HYPRE_StructPCGSetMaxIter(solver, 50);
355
356 /* Create the Struct PFMG solver for use as a preconditioner */
357 HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
358
359 /* Set PFMG parameters */
360 HYPRE_StructPFMGSetMaxIter(precond, 1);
361 HYPRE_StructPFMGSetTol(precond, 0.0);
362 HYPRE_StructPFMGSetZeroGuess(precond);
363 HYPRE_StructPFMGSetNumPreRelax(precond, 2);
364 HYPRE_StructPFMGSetNumPostRelax(precond, 2);
365 /* non-Galerkin coarse grid (more efficient for this problem) */
366 HYPRE_StructPFMGSetRAPType(precond, 1);
367 /* R/B Gauss-Seidel */
368 HYPRE_StructPFMGSetRelaxType(precond, 2);
369 /* skip relaxation on some levels (more efficient for this problem) */
370 HYPRE_StructPFMGSetSkipRelax(precond, 1);
371
372
373 /* Set preconditioner and solve */
374 HYPRE_StructPCGSetPrecond(solver, HYPRE_StructPFMGSolve,
375 HYPRE_StructPFMGSetup, precond);
376 HYPRE_StructPCGSetup(solver, sA, sb, sx);
377 HYPRE_StructPCGSolve(solver, sA, sb, sx);
378
379 /* Free memory */
380 HYPRE_StructPCGDestroy(solver);
381 HYPRE_StructPFMGDestroy(precond);
382 }
383 else if (precond_id == 2) /* BoomerAMG */
384 {
385 HYPRE_ParCSRMatrix parA;
386 HYPRE_ParVector parb;
387 HYPRE_ParVector parx;
388
389 HYPRE_Solver solver;
390 HYPRE_Solver precond;
391
392 /* Because we are using a struct solver, we need to get the
393 object of the matrix and vectors to pass in to the struct solvers */
394 HYPRE_SStructMatrixGetObject(A, (void **) &parA);
395 HYPRE_SStructVectorGetObject(b, (void **) &parb);
396 HYPRE_SStructVectorGetObject(x, (void **) &parx);
397
398 /* Create an empty PCG Struct solver */
399 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
400
401 /* Set PCG parameters */
402 HYPRE_ParCSRPCGSetTol(solver, 1.0e-06);
403 HYPRE_ParCSRPCGSetPrintLevel(solver, 2);
404 HYPRE_ParCSRPCGSetMaxIter(solver, 50);
405
406 /* Create the BoomerAMG solver for use as a preconditioner */
407 HYPRE_BoomerAMGCreate(&precond);
408
409 /* Set BoomerAMG parameters */
410 HYPRE_BoomerAMGSetMaxIter(precond, 1);
411 HYPRE_BoomerAMGSetTol(precond, 0.0);
412 HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
413 HYPRE_BoomerAMGSetCoarsenType(precond, 6);
414 HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
415 HYPRE_BoomerAMGSetNumSweeps(precond, 1);
416
417 /* Set preconditioner and solve */
418 HYPRE_ParCSRPCGSetPrecond(solver, HYPRE_BoomerAMGSolve,
419 HYPRE_BoomerAMGSetup, precond);
420 HYPRE_ParCSRPCGSetup(solver, parA, parb, parx);
421 HYPRE_ParCSRPCGSolve(solver, parA, parb, parx);
422
423 /* Free memory */
424 HYPRE_ParCSRPCGDestroy(solver);
425 HYPRE_BoomerAMGDestroy(precond);
426 }
427
428 /* Free memory */
429 HYPRE_SStructGridDestroy(grid);
430 HYPRE_SStructStencilDestroy(stencil);
431 HYPRE_SStructGraphDestroy(graph);
432 HYPRE_SStructMatrixDestroy(A);
433 HYPRE_SStructVectorDestroy(b);
434 HYPRE_SStructVectorDestroy(x);
435
436 /* Finalize MPI */
437 MPI_Finalize();
438
439 return (0);
440 }
syntax highlighted by Code2HTML, v. 0.9.1