download the original source code.
1 /*
2 Example 18
3
4 Interface: SStructured interface (SStruct)
5
6 Compile with: make ex18
7
8 Sample run: mpirun -np 16 ex18 -n 4
9
10 To see options: ex18 -help
11
12 Description: This code solves an "NDIM-D Laplacian" using CG.
13 */
14
15 #include <math.h>
16 #include "_hypre_utilities.h"
17 #include "HYPRE_sstruct_ls.h"
18
19 #define NDIM 4
20 #define NPARTS 1
21 #define NVARS 2
22 #define NSTENC NVARS*(2*NDIM+1)
23
24 int main (int argc, char *argv[])
25 {
26 int d, i, j;
27 int myid, num_procs;
28 int n, N, nvol, div, rem;
29 int p[NDIM], ilower[NDIM], iupper[NDIM];
30
31 int solver_id, object_type = HYPRE_SSTRUCT;
32
33 HYPRE_SStructGrid grid;
34 HYPRE_SStructStencil stencil0, stencil1;
35 HYPRE_SStructGraph graph;
36 HYPRE_SStructMatrix A;
37 HYPRE_SStructVector b;
38 HYPRE_SStructVector x;
39
40 HYPRE_SStructSolver solver;
41
42 int num_iterations;
43 double final_res_norm;
44
45 /* Initialize MPI */
46 MPI_Init(&argc, &argv);
47 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
48 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
49
50 /* Set defaults */
51 n = 4;
52 solver_id = 0;
53
54 /* Parse command line */
55 {
56 int arg_index = 0;
57 int print_usage = 0;
58
59 while (arg_index < argc)
60 {
61 if ( strcmp(argv[arg_index], "-n") == 0 )
62 {
63 arg_index++;
64 n = atoi(argv[arg_index++]);
65 }
66 else if ( strcmp(argv[arg_index], "-solver") == 0 )
67 {
68 arg_index++;
69 solver_id = atoi(argv[arg_index++]);
70 }
71 else if ( strcmp(argv[arg_index], "-help") == 0 )
72 {
73 print_usage = 1;
74 break;
75 }
76 else
77 {
78 arg_index++;
79 }
80 }
81
82 if ((print_usage) && (myid == 0))
83 {
84 printf("\n");
85 printf("Usage: %s [<options>]\n", argv[0]);
86 printf("\n");
87 printf(" -n <n> : problem size per processor (default: 4)\n");
88 printf(" -solver <ID> : solver ID\n");
89 printf(" 0 - CG (default)\n");
90 printf(" 1 - GMRES\n");
91 printf("\n");
92 }
93
94 if (print_usage)
95 {
96 MPI_Finalize();
97 return (0);
98 }
99 }
100
101 nvol = pow(n, NDIM);
102
103 /* Figure out the processor grid (N x N x N x N). The local problem size for
104 the interior nodes is indicated by n (n x n x n x n). p indicates the
105 position in the processor grid. */
106 N = pow(num_procs, 1.0/NDIM) + 1.0e-6;
107 div = pow(N, NDIM);
108 rem = myid;
109 if (num_procs != div)
110 {
111 printf("Num procs is not a perfect NDIM-th root!\n");
112 MPI_Finalize();
113 exit(1);
114 }
115 for (d = NDIM-1; d >= 0; d--)
116 {
117 div /= N;
118 p[d] = rem / div;
119 rem %= div;
120 }
121
122 /* Figure out the extents of each processor's piece of the grid. */
123 for (d = 0; d < NDIM; d++)
124 {
125 ilower[d] = p[d]*n;
126 iupper[d] = ilower[d] + n-1;
127 }
128
129 /* 1. Set up a grid */
130 {
131 int part = 0;
132 HYPRE_SStructVariable vartypes[NVARS] = {HYPRE_SSTRUCT_VARIABLE_CELL,
133 HYPRE_SSTRUCT_VARIABLE_CELL};
134
135 /* Create an empty 2D grid object */
136 HYPRE_SStructGridCreate(MPI_COMM_WORLD, NDIM, NPARTS, &grid);
137
138 /* Add a new box to the grid */
139 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
140
141 /* Set the variable type and number of variables on each part. */
142 HYPRE_SStructGridSetVariables(grid, part, NVARS, vartypes);
143
144 /* The grid is now ready to use */
145 HYPRE_SStructGridAssemble(grid);
146 }
147
148 /* 2. Define the discretization stencil */
149 {
150 /* Create two empty NDIM-D, NSTENC-pt stencil objects */
151 HYPRE_SStructStencilCreate(NDIM, NSTENC, &stencil0);
152 HYPRE_SStructStencilCreate(NDIM, NSTENC, &stencil1);
153
154 /* Define the geometry of the stencil */
155 {
156 int entry, var0 = 0, var1 = 1;
157 int offset[NDIM];
158
159 entry = 0;
160 for (d = 0; d < NDIM; d++)
161 {
162 offset[d] = 0;
163 }
164 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var0);
165 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var1);
166 entry++;
167 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var1);
168 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var0);
169 entry++;
170 for (d = 0; d < NDIM; d++)
171 {
172 offset[d] = -1;
173 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var0);
174 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var1);
175 entry++;
176 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var1);
177 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var0);
178 entry++;
179 offset[d] = 1;
180 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var0);
181 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var1);
182 entry++;
183 HYPRE_SStructStencilSetEntry(stencil0, entry, offset, var1);
184 HYPRE_SStructStencilSetEntry(stencil1, entry, offset, var0);
185 entry++;
186 offset[d] = 0;
187 }
188 }
189 }
190
191 /* 3. Set up the Graph */
192 {
193 int part = 0;
194 int var0 = 0, var1 = 1;
195
196 /* Create the graph object */
197 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
198
199 /* Set up the object type (see Matrix and VectorSetObjectType below) */
200 HYPRE_SStructGraphSetObjectType(graph, object_type);
201
202 /* Set the stencil */
203 HYPRE_SStructGraphSetStencil(graph, part, var0, stencil0);
204 HYPRE_SStructGraphSetStencil(graph, part, var1, stencil1);
205
206 /* Assemble the graph */
207 HYPRE_SStructGraphAssemble(graph);
208 }
209
210 /* 4. Set up the Matrix */
211 {
212 int part = 0;
213 int var0 = 0, var1 = 1;
214 int nentries = NSTENC/NVARS;
215 int nvalues = nentries*nvol;
216 double *values;
217 int stencil_indices[NSTENC];
218
219 /* Create an empty matrix object */
220 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
221
222 /* Set up the object type */
223 HYPRE_SStructMatrixSetObjectType(A, object_type);
224
225 /* Get ready to set values */
226 HYPRE_SStructMatrixInitialize(A);
227
228 values = (double*) calloc(nvalues, sizeof(double));
229
230 /* Set intra-variable values; fix boundaries later */
231 for (j = 0; j < nentries; j++)
232 {
233 stencil_indices[j] = 2*j;
234 }
235 for (i = 0; i < nvalues; i += nentries)
236 {
237 values[i] = 1.1*(NSTENC/NVARS); /* Diagonal: Use absolute row sum */
238 for (j = 1; j < nentries; j++)
239 {
240 values[i+j] = -1.0;
241 }
242 }
243 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var0,
244 nentries, stencil_indices, values);
245 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var1,
246 nentries, stencil_indices, values);
247
248 /* Set inter-variable values; fix boundaries later */
249 for (j = 0; j < nentries; j++)
250 {
251 stencil_indices[j] = 2*j+1;
252 }
253 for (i = 0; i < nvalues; i += nentries)
254 {
255 values[i] = -0.1;
256 for (j = 1; j < nentries; j++)
257 {
258 values[i+j] = -0.1;
259 }
260 }
261 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var0,
262 nentries, stencil_indices, values);
263 HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper, var1,
264 nentries, stencil_indices, values);
265
266 free(values);
267 }
268
269 /* 5. Incorporate zero boundary conditions: go along each edge of the domain
270 and set the stencil entry that reaches to the boundary to zero.*/
271 {
272 int part = 0;
273 int var0 = 0, var1 = 1;
274 int bc_ilower[NDIM];
275 int bc_iupper[NDIM];
276 int nentries = 1;
277 int nvalues = nentries*nvol/n; /* number of stencil entries times the
278 length of one side of my grid box */
279 double *values;
280 int stencil_indices[1];
281
282 values = (double*) calloc(nvalues, sizeof(double));
283 for (j = 0; j < nvalues; j++)
284 {
285 values[j] = 0.0;
286 }
287
288 for (d = 0; d < NDIM; d++)
289 {
290 bc_ilower[d] = ilower[d];
291 bc_iupper[d] = iupper[d];
292 }
293 stencil_indices[0] = NVARS;
294 for (d = 0; d < NDIM; d++)
295 {
296 /* lower boundary in dimension d */
297 if (p[d] == 0)
298 {
299 bc_iupper[d] = ilower[d];
300 for (i = 0; i < NVARS; i++)
301 {
302 HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper, var0,
303 nentries, stencil_indices, values);
304 HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper, var1,
305 nentries, stencil_indices, values);
306 stencil_indices[0]++;
307 }
308 bc_iupper[d] = iupper[d];
309 }
310 else
311 {
312 stencil_indices[0] += NVARS;
313 }
314
315 /* upper boundary in dimension d */
316 if (p[d] == N-1)
317 {
318 bc_ilower[d] = iupper[d];
319 for (i = 0; i < NVARS; i++)
320 {
321 HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper, var0,
322 nentries, stencil_indices, values);
323 HYPRE_SStructMatrixSetBoxValues(A, part, bc_ilower, bc_iupper, var1,
324 nentries, stencil_indices, values);
325 stencil_indices[0]++;
326 }
327 bc_ilower[d] = ilower[d];
328 }
329 else
330 {
331 stencil_indices[0] += NVARS;
332 }
333 }
334
335 free(values);
336 }
337
338 /* The matrix is now ready to use */
339 HYPRE_SStructMatrixAssemble(A);
340
341 /* 6. Set up Vectors for b and x */
342 {
343 int part = 0;
344 int var0 = 0, var1 = 1;
345 int nvalues = NVARS*nvol;
346 double *values;
347
348 values = (double*) calloc(nvalues, sizeof(double));
349
350 /* Create an empty vector object */
351 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
352 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
353
354 /* Set up the object type */
355 HYPRE_SStructVectorSetObjectType(b, object_type);
356 HYPRE_SStructVectorSetObjectType(x, object_type);
357
358 /* Indicate that the vector coefficients are ready to be set */
359 HYPRE_SStructVectorInitialize(b);
360 HYPRE_SStructVectorInitialize(x);
361
362 /* Set the values */
363 for (i = 0; i < nvalues; i ++)
364 {
365 values[i] = 1.0;
366 }
367 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var0, values);
368 HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper, var1, values);
369
370 for (i = 0; i < nvalues; i ++)
371 {
372 values[i] = 0.0;
373 }
374 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var0, values);
375 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var1, values);
376
377 free(values);
378
379 /* The vector is now ready to use */
380 HYPRE_SStructVectorAssemble(b);
381 HYPRE_SStructVectorAssemble(x);
382 }
383
384 #if 0
385 HYPRE_SStructMatrixPrint("ex18.out.A", A, 0);
386 HYPRE_SStructVectorPrint("ex18.out.b", b, 0);
387 HYPRE_SStructVectorPrint("ex18.out.x0", x, 0);
388 #endif
389
390 /* 7. Set up and use a struct solver */
391 if (solver_id == 0)
392 {
393 HYPRE_SStructPCGCreate(MPI_COMM_WORLD, &solver);
394 HYPRE_SStructPCGSetMaxIter(solver, 100);
395 HYPRE_SStructPCGSetTol(solver, 1.0e-06);
396 HYPRE_SStructPCGSetTwoNorm(solver, 1);
397 HYPRE_SStructPCGSetRelChange(solver, 0);
398 HYPRE_SStructPCGSetPrintLevel(solver, 2); /* print each CG iteration */
399 HYPRE_SStructPCGSetLogging(solver, 1);
400
401 /* No preconditioner */
402
403 HYPRE_SStructPCGSetup(solver, A, b, x);
404 HYPRE_SStructPCGSolve(solver, A, b, x);
405
406 /* Get some info on the run */
407 HYPRE_SStructPCGGetNumIterations(solver, &num_iterations);
408 HYPRE_SStructPCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
409
410 /* Clean up */
411 HYPRE_SStructPCGDestroy(solver);
412 }
413
414 if (myid == 0)
415 {
416 printf("\n");
417 printf("Iterations = %d\n", num_iterations);
418 printf("Final Relative Residual Norm = %g\n", final_res_norm);
419 printf("\n");
420 }
421
422 /* Free memory */
423 HYPRE_SStructGridDestroy(grid);
424 HYPRE_SStructGraphDestroy(graph);
425 HYPRE_SStructStencilDestroy(stencil0);
426 HYPRE_SStructStencilDestroy(stencil1);
427 HYPRE_SStructMatrixDestroy(A);
428 HYPRE_SStructVectorDestroy(b);
429 HYPRE_SStructVectorDestroy(x);
430
431 /* Finalize MPI */
432 MPI_Finalize();
433
434 return (0);
435 }
syntax highlighted by Code2HTML, v. 0.9.1