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