download the original source code.
1 /*
2 Example 17
3
4 Interface: Structured interface (Struct)
5
6 Compile with: make ex17
7
8 Sample run: mpirun -np 16 ex17 -n 10
9
10 To see options: ex17 -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_struct_ls.h"
18
19 #define NDIM 4
20 #define NSTENC (2*NDIM+1)
21
22 int main (int argc, char *argv[])
23 {
24 int d, i, j;
25 int myid, num_procs;
26 int n, N, nvol, div, rem;
27 int p[NDIM], ilower[NDIM], iupper[NDIM];
28
29 int solver_id;
30
31 HYPRE_StructGrid grid;
32 HYPRE_StructStencil stencil;
33 HYPRE_StructMatrix A;
34 HYPRE_StructVector b;
35 HYPRE_StructVector x;
36 HYPRE_StructSolver solver;
37
38 int num_iterations;
39 double final_res_norm;
40
41 /* Initialize MPI */
42 MPI_Init(&argc, &argv);
43 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
44 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
45
46 /* Set defaults */
47 n = 10;
48 solver_id = 0;
49
50 /* Parse command line */
51 {
52 int arg_index = 0;
53 int print_usage = 0;
54
55 while (arg_index < argc)
56 {
57 if ( strcmp(argv[arg_index], "-n") == 0 )
58 {
59 arg_index++;
60 n = atoi(argv[arg_index++]);
61 }
62 else if ( strcmp(argv[arg_index], "-solver") == 0 )
63 {
64 arg_index++;
65 solver_id = atoi(argv[arg_index++]);
66 }
67 else if ( strcmp(argv[arg_index], "-help") == 0 )
68 {
69 print_usage = 1;
70 break;
71 }
72 else
73 {
74 arg_index++;
75 }
76 }
77
78 if ((print_usage) && (myid == 0))
79 {
80 printf("\n");
81 printf("Usage: %s [<options>]\n", argv[0]);
82 printf("\n");
83 printf(" -n <n> : problem size per processor (default: 33)\n");
84 printf(" -solver <ID> : solver ID\n");
85 printf(" 0 - CG (default)\n");
86 printf(" 1 - GMRES\n");
87 printf("\n");
88 }
89
90 if (print_usage)
91 {
92 MPI_Finalize();
93 return (0);
94 }
95 }
96
97 nvol = pow(n, NDIM);
98
99 /* Figure out the processor grid (N x N x N x N). The local problem size for
100 the interior nodes is indicated by n (n x n x n x n). p indicates the
101 position in the processor grid. */
102 N = pow(num_procs, 1.0/NDIM) + 1.0e-6;
103 div = pow(N, NDIM);
104 rem = myid;
105 if (num_procs != div)
106 {
107 printf("Num procs is not a perfect NDIM-th root!\n");
108 MPI_Finalize();
109 exit(1);
110 }
111 for (d = NDIM-1; d >= 0; d--)
112 {
113 div /= N;
114 p[d] = rem / div;
115 rem %= div;
116 }
117
118 /* Figure out the extents of each processor's piece of the grid. */
119 for (d = 0; d < NDIM; d++)
120 {
121 ilower[d] = p[d]*n;
122 iupper[d] = ilower[d] + n-1;
123 }
124
125 /* 1. Set up a grid */
126 {
127 /* Create an empty 2D grid object */
128 HYPRE_StructGridCreate(MPI_COMM_WORLD, NDIM, &grid);
129
130 /* Add a new box to the grid */
131 HYPRE_StructGridSetExtents(grid, ilower, iupper);
132
133 /* This is a collective call finalizing the grid assembly.
134 The grid is now ``ready to be used'' */
135 HYPRE_StructGridAssemble(grid);
136 }
137
138 /* 2. Define the discretization stencil */
139 {
140 /* Create an empty NDIM-D, NSTENC-pt stencil object */
141 HYPRE_StructStencilCreate(NDIM, NSTENC, &stencil);
142
143 /* Define the geometry of the stencil */
144 {
145 int entry;
146 int offset[NDIM];
147
148 entry = 0;
149 for (d = 0; d < NDIM; d++)
150 {
151 offset[d] = 0;
152 }
153 HYPRE_StructStencilSetElement(stencil, entry++, offset);
154 for (d = 0; d < NDIM; d++)
155 {
156 offset[d] = -1;
157 HYPRE_StructStencilSetElement(stencil, entry++, offset);
158 offset[d] = 1;
159 HYPRE_StructStencilSetElement(stencil, entry++, offset);
160 offset[d] = 0;
161 }
162 }
163 }
164
165 /* 3. Set up a Struct Matrix */
166 {
167 int nentries = NSTENC;
168 int nvalues = nentries*nvol;
169 double *values;
170 int stencil_indices[NSTENC];
171
172 /* Create an empty matrix object */
173 HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
174
175 /* Indicate that the matrix coefficients are ready to be set */
176 HYPRE_StructMatrixInitialize(A);
177
178 values = (double*) calloc(nvalues, sizeof(double));
179
180 for (j = 0; j < nentries; j++)
181 {
182 stencil_indices[j] = j;
183 }
184
185 /* Set the standard stencil at each grid point; fix boundaries later */
186 for (i = 0; i < nvalues; i += nentries)
187 {
188 values[i] = NSTENC; /* Use absolute row sum */
189 for (j = 1; j < nentries; j++)
190 {
191 values[i+j] = -1.0;
192 }
193 }
194
195 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
196 stencil_indices, values);
197
198 free(values);
199 }
200
201 /* 4. Incorporate zero boundary conditions: go along each edge of the domain
202 and set the stencil entry that reaches to the boundary to zero.*/
203 {
204 int bc_ilower[NDIM];
205 int bc_iupper[NDIM];
206 int nentries = 1;
207 int nvalues = nentries*nvol/n; /* number of stencil entries times the
208 length of one side of my grid box */
209 double *values;
210 int stencil_indices[1];
211
212 values = (double*) calloc(nvalues, sizeof(double));
213 for (j = 0; j < nvalues; j++)
214 {
215 values[j] = 0.0;
216 }
217
218 for (d = 0; d < NDIM; d++)
219 {
220 bc_ilower[d] = ilower[d];
221 bc_iupper[d] = iupper[d];
222 }
223 stencil_indices[0] = 1;
224 for (d = 0; d < NDIM; d++)
225 {
226 /* lower boundary in dimension d */
227 if (p[d] == 0)
228 {
229 bc_iupper[d] = ilower[d];
230 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
231 stencil_indices, values);
232 bc_iupper[d] = iupper[d];
233 }
234 stencil_indices[0]++;
235
236 /* upper boundary in dimension d */
237 if (p[d] == N-1)
238 {
239 bc_ilower[d] = iupper[d];
240 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
241 stencil_indices, values);
242 bc_ilower[d] = ilower[d];
243 }
244 stencil_indices[0]++;
245 }
246
247 free(values);
248 }
249
250 /* This is a collective call finalizing the matrix assembly.
251 The matrix is now ``ready to be used'' */
252 HYPRE_StructMatrixAssemble(A);
253
254 /* 5. Set up Struct Vectors for b and x */
255 {
256 int nvalues = nvol;
257 double *values;
258
259 values = (double*) calloc(nvalues, sizeof(double));
260
261 /* Create an empty vector object */
262 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
263 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
264
265 /* Indicate that the vector coefficients are ready to be set */
266 HYPRE_StructVectorInitialize(b);
267 HYPRE_StructVectorInitialize(x);
268
269 /* Set the values */
270 for (i = 0; i < nvalues; i ++)
271 {
272 values[i] = 1.0;
273 }
274 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
275
276 for (i = 0; i < nvalues; i ++)
277 {
278 values[i] = 0.0;
279 }
280 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
281
282 free(values);
283
284 /* This is a collective call finalizing the vector assembly.
285 The vector is now ``ready to be used'' */
286 HYPRE_StructVectorAssemble(b);
287 HYPRE_StructVectorAssemble(x);
288 }
289
290 #if 0
291 HYPRE_StructMatrixPrint("ex17.out.A", A, 0);
292 HYPRE_StructVectorPrint("ex17.out.b", b, 0);
293 HYPRE_StructVectorPrint("ex17.out.x0", x, 0);
294 #endif
295
296 /* 6. Set up and use a struct solver
297 (Solver options can be found in the Reference Manual.) */
298 if (solver_id == 0)
299 {
300 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
301 HYPRE_StructPCGSetMaxIter(solver, 100);
302 HYPRE_StructPCGSetTol(solver, 1.0e-06);
303 HYPRE_StructPCGSetTwoNorm(solver, 1);
304 HYPRE_StructPCGSetRelChange(solver, 0);
305 HYPRE_StructPCGSetPrintLevel(solver, 2); /* print each CG iteration */
306 HYPRE_StructPCGSetLogging(solver, 1);
307
308 /* No preconditioner */
309
310 HYPRE_StructPCGSetup(solver, A, b, x);
311 HYPRE_StructPCGSolve(solver, A, b, x);
312
313 /* Get some info on the run */
314 HYPRE_StructPCGGetNumIterations(solver, &num_iterations);
315 HYPRE_StructPCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
316
317 /* Clean up */
318 HYPRE_StructPCGDestroy(solver);
319 }
320
321 if (myid == 0)
322 {
323 printf("\n");
324 printf("Iterations = %d\n", num_iterations);
325 printf("Final Relative Residual Norm = %g\n", final_res_norm);
326 printf("\n");
327 }
328
329 /* Free memory */
330 HYPRE_StructGridDestroy(grid);
331 HYPRE_StructStencilDestroy(stencil);
332 HYPRE_StructMatrixDestroy(A);
333 HYPRE_StructVectorDestroy(b);
334 HYPRE_StructVectorDestroy(x);
335
336 /* Finalize MPI */
337 MPI_Finalize();
338
339 return (0);
340 }
syntax highlighted by Code2HTML, v. 0.9.1