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 /* See MatrixSetObjectType below */
139 HYPRE_SStructGraphSetObjectType(graph, object_type);
140
141 /* Now we need to tell the graph which stencil to use for each variable on
142 each part (we only have one variable and one part) */
143 HYPRE_SStructGraphSetStencil(graph, part, var, stencil);
144
145 /* Here we could establish connections between parts if we had more than
146 one part using the graph. For example, we could use
147 HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborPart() */
148
149 /* Assemble the graph */
150 HYPRE_SStructGraphAssemble(graph);
151 }
152
153 /* 4. Set up a SStruct Matrix */
154 {
155 /* Create an empty matrix object */
156 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
157
158 /* Set the object type (by default HYPRE_SSTRUCT). This determines the
159 data structure used to store the matrix. For PFMG we need to use
160 HYPRE_STRUCT, and for BoomerAMG we need HYPRE_PARCSR (set above). */
161 HYPRE_SStructMatrixSetObjectType(A, object_type);
162
163 /* Get ready to set values */
164 HYPRE_SStructMatrixInitialize(A);
165
166 /* Set the matrix coefficients. Each processor assigns coefficients for
167 the boxes in the grid that it owns. Note that the coefficients
168 associated with each stencil entry may vary from grid point to grid
169 point if desired. Here, we first set the same stencil entries for each
170 grid point. Then we make modifications to grid points near the
171 boundary. Note that the ilower values are different from those used in
172 ex1 because of the way nodal variables are referenced. Also note that
173 some of the stencil values are set on both processor 0 and processor 1.
174 See the User and Reference manuals for more details. */
175 if (myid == 0)
176 {
177 int ilower[2]={-4,0}, iupper[2]={-1,2};
178 int stencil_indices[5] = {0,1,2,3,4}; /* labels for the stencil entries -
179 these correspond to the offsets
180 defined above */
181 int nentries = 5;
182 int nvalues = 60; /* 12 grid points, each with 5 stencil entries */
183 double values[60];
184
185 for (i = 0; i < nvalues; i += nentries)
186 {
187 values[i] = 4.0;
188 for (j = 1; j < nentries; j++)
189 values[i+j] = -1.0;
190 }
191
192 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, nentries,
193 stencil_indices, values);
194 }
195 else if (myid == 1)
196 {
197 int ilower[2]={-1,0}, iupper[2]={2,4};
198 int stencil_indices[5] = {0,1,2,3,4};
199 int nentries = 5;
200 int nvalues = 100; /* 20 grid points, each with 5 stencil entries */
201 double values[100];
202
203 for (i = 0; i < nvalues; i += nentries)
204 {
205 values[i] = 4.0;
206 for (j = 1; j < nentries; j++)
207 values[i+j] = -1.0;
208 }
209
210 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, nentries,
211 stencil_indices, values);
212 }
213
214 /* Set the coefficients reaching outside of the boundary to 0. Note that
215 * both ilower *and* iupper may be different from those in ex1. */
216 if (myid == 0)
217 {
218 double values[4];
219 for (i = 0; i < 4; i++)
220 values[i] = 0.0;
221 {
222 /* values below our box */
223 int ilower[2]={-4,0}, iupper[2]={-1,0};
224 int stencil_indices[1] = {3};
225 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
226 stencil_indices, values);
227 }
228 {
229 /* values to the left of our box */
230 int ilower[2]={-4,0}, iupper[2]={-4,2};
231 int stencil_indices[1] = {1};
232 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
233 stencil_indices, values);
234 }
235 {
236 /* values above our box */
237 int ilower[2]={-4,2}, iupper[2]={-2,2};
238 int stencil_indices[1] = {4};
239 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
240 stencil_indices, values);
241 }
242 }
243 else if (myid == 1)
244 {
245 double values[5];
246 for (i = 0; i < 5; i++)
247 values[i] = 0.0;
248 {
249 /* values below our box */
250 int ilower[2]={-1,0}, iupper[2]={2,0};
251 int stencil_indices[1] = {3};
252 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
253 stencil_indices, values);
254 }
255 {
256 /* values to the right of our box */
257 int ilower[2]={2,0}, iupper[2]={2,4};
258 int stencil_indices[1] = {2};
259 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
260 stencil_indices, values);
261 }
262 {
263 /* values above our box */
264 int ilower[2]={-1,4}, iupper[2]={2,4};
265 int stencil_indices[1] = {4};
266 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
267 stencil_indices, values);
268 }
269 {
270 /* values to the left of our box
271 (that do not border the other box on proc. 0) */
272 int ilower[2]={-1,3}, iupper[2]={-1,4};
273 int stencil_indices[1] = {1};
274 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var, 1,
275 stencil_indices, values);
276 }
277 }
278
279 /* This is a collective call finalizing the matrix assembly.
280 The matrix is now ``ready to be used'' */
281 HYPRE_SStructMatrixAssemble(A);
282 }
283
284 /* 5. Set up SStruct Vectors for b and x. */
285 {
286 /* Create an empty vector object */
287 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
288 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
289
290 /* As with the matrix, set the appropriate object type for the vectors */
291 HYPRE_SStructVectorSetObjectType(b, object_type);
292 HYPRE_SStructVectorSetObjectType(x, object_type);
293
294 /* Indicate that the vector coefficients are ready to be set */
295 HYPRE_SStructVectorInitialize(b);
296 HYPRE_SStructVectorInitialize(x);
297
298 /* Set the vector coefficients. Again, note that the ilower values are
299 different from those used in ex1, and some of the values are set on
300 both processors. */
301 if (myid == 0)
302 {
303 int ilower[2]={-4,0}, iupper[2]={-1,2};
304 double values[12]; /* 12 grid points */
305
306 for (i = 0; i < 12; i ++)
307 values[i] = 1.0;
308 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
309
310 for (i = 0; i < 12; i ++)
311 values[i] = 0.0;
312 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
313 }
314 else if (myid == 1)
315 {
316 int ilower[2]={0,1}, iupper[2]={2,4};
317 double values[20]; /* 20 grid points */
318
319 for (i = 0; i < 20; i ++)
320 values[i] = 1.0;
321 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var, values);
322
323 for (i = 0; i < 20; i ++)
324 values[i] = 0.0;
325 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
326 }
327
328 /* This is a collective call finalizing the vector assembly.
329 The vectors are now ``ready to be used'' */
330 HYPRE_SStructVectorAssemble(b);
331 HYPRE_SStructVectorAssemble(x);
332 }
333
334 /* 6. Set up and use a solver (See the Reference Manual for descriptions
335 of all of the options.) */
336 if (precond_id == 1) /* PFMG */
337 {
338 HYPRE_StructMatrix sA;
339 HYPRE_StructVector sb;
340 HYPRE_StructVector sx;
341
342 HYPRE_StructSolver solver;
343 HYPRE_StructSolver precond;
344
345 /* Because we are using a struct solver, we need to get the
346 object of the matrix and vectors to pass in to the struct solvers */
347 HYPRE_SStructMatrixGetObject(A, (void **) &sA);
348 HYPRE_SStructVectorGetObject(b, (void **) &sb);
349 HYPRE_SStructVectorGetObject(x, (void **) &sx);
350
351 /* Create an empty PCG Struct solver */
352 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
353
354 /* Set PCG parameters */
355 HYPRE_StructPCGSetTol(solver, 1.0e-06);
356 HYPRE_StructPCGSetPrintLevel(solver, 2);
357 HYPRE_StructPCGSetMaxIter(solver, 50);
358
359 /* Create the Struct PFMG solver for use as a preconditioner */
360 HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
361
362 /* Set PFMG parameters */
363 HYPRE_StructPFMGSetMaxIter(precond, 1);
364 HYPRE_StructPFMGSetTol(precond, 0.0);
365 HYPRE_StructPFMGSetZeroGuess(precond);
366 HYPRE_StructPFMGSetNumPreRelax(precond, 2);
367 HYPRE_StructPFMGSetNumPostRelax(precond, 2);
368 /* non-Galerkin coarse grid (more efficient for this problem) */
369 HYPRE_StructPFMGSetRAPType(precond, 1);
370 /* R/B Gauss-Seidel */
371 HYPRE_StructPFMGSetRelaxType(precond, 2);
372 /* skip relaxation on some levels (more efficient for this problem) */
373 HYPRE_StructPFMGSetSkipRelax(precond, 1);
374
375
376 /* Set preconditioner and solve */
377 HYPRE_StructPCGSetPrecond(solver, HYPRE_StructPFMGSolve,
378 HYPRE_StructPFMGSetup, precond);
379 HYPRE_StructPCGSetup(solver, sA, sb, sx);
380 HYPRE_StructPCGSolve(solver, sA, sb, sx);
381
382 /* Free memory */
383 HYPRE_StructPCGDestroy(solver);
384 HYPRE_StructPFMGDestroy(precond);
385 }
386 else if (precond_id == 2) /* BoomerAMG */
387 {
388 HYPRE_ParCSRMatrix parA;
389 HYPRE_ParVector parb;
390 HYPRE_ParVector parx;
391
392 HYPRE_Solver solver;
393 HYPRE_Solver precond;
394
395 /* Because we are using a struct solver, we need to get the
396 object of the matrix and vectors to pass in to the struct solvers */
397 HYPRE_SStructMatrixGetObject(A, (void **) &parA);
398 HYPRE_SStructVectorGetObject(b, (void **) &parb);
399 HYPRE_SStructVectorGetObject(x, (void **) &parx);
400
401 /* Create an empty PCG Struct solver */
402 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
403
404 /* Set PCG parameters */
405 HYPRE_ParCSRPCGSetTol(solver, 1.0e-06);
406 HYPRE_ParCSRPCGSetPrintLevel(solver, 2);
407 HYPRE_ParCSRPCGSetMaxIter(solver, 50);
408
409 /* Create the BoomerAMG solver for use as a preconditioner */
410 HYPRE_BoomerAMGCreate(&precond);
411
412 /* Set BoomerAMG parameters */
413 HYPRE_BoomerAMGSetMaxIter(precond, 1);
414 HYPRE_BoomerAMGSetTol(precond, 0.0);
415 HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
416 HYPRE_BoomerAMGSetCoarsenType(precond, 6);
417 HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
418 HYPRE_BoomerAMGSetNumSweeps(precond, 1);
419
420 /* Set preconditioner and solve */
421 HYPRE_ParCSRPCGSetPrecond(solver, HYPRE_BoomerAMGSolve,
422 HYPRE_BoomerAMGSetup, precond);
423 HYPRE_ParCSRPCGSetup(solver, parA, parb, parx);
424 HYPRE_ParCSRPCGSolve(solver, parA, parb, parx);
425
426 /* Free memory */
427 HYPRE_ParCSRPCGDestroy(solver);
428 HYPRE_BoomerAMGDestroy(precond);
429 }
430
431 /* Free memory */
432 HYPRE_SStructGridDestroy(grid);
433 HYPRE_SStructStencilDestroy(stencil);
434 HYPRE_SStructGraphDestroy(graph);
435 HYPRE_SStructMatrixDestroy(A);
436 HYPRE_SStructVectorDestroy(b);
437 HYPRE_SStructVectorDestroy(x);
438
439 /* Finalize MPI */
440 MPI_Finalize();
441
442 return (0);
443 }
syntax highlighted by Code2HTML, v. 0.9.1