download the original source code.
1 /*
2 Example 3
3
4 Interface: Structured interface (Struct)
5
6 Compile with: make ex3
7
8 Sample run: mpirun -np 16 ex3 -n 33 -solver 0 -v 1 1
9
10 To see options: ex3 -help
11
12 Description: This code solves a system corresponding to a discretization
13 of the Laplace equation with zero boundary conditions on the
14 unit square. The domain is split into an N x N processor grid.
15 Thus, the given number of processors should be a perfect square.
16 Each processor's piece of the grid has n x n cells with n x n
17 nodes connected by the standard 5-point stencil. Note that the
18 struct interface assumes a cell-centered grid, and, therefore,
19 the nodes are not shared. This example demonstrates more
20 features than the previous two struct examples (Example 1 and
21 Example 2). Two solvers are available.
22
23 To incorporate the boundary conditions, we do the following:
24 Let x_i and x_b be the interior and boundary parts of the
25 solution vector x. We can split the matrix A as
26 A = [A_ii A_ib; A_bi A_bb].
27 Let u_0 be the Dirichlet B.C. We can simply say that x_b = u_0.
28 If b_i is the right-hand side, then we just need to solve in
29 the interior:
30 A_ii x_i = b_i - A_ib u_0.
31 For this partitcular example, u_0 = 0, so we are just solving
32 A_ii x_i = b_i.
33
34 We recommend viewing examples 1 and 2 before viewing this
35 example.
36 */
37
38 #include <math.h>
39 #include "_hypre_utilities.h"
40 #include "HYPRE_struct_ls.h"
41
42 int main (int argc, char *argv[])
43 {
44 int i, j;
45
46 int myid, num_procs;
47
48 int n, N, pi, pj;
49 double h, h2;
50 int ilower[2], iupper[2];
51
52 int solver_id;
53 int n_pre, n_post;
54
55 HYPRE_StructGrid grid;
56 HYPRE_StructStencil stencil;
57 HYPRE_StructMatrix A;
58 HYPRE_StructVector b;
59 HYPRE_StructVector x;
60 HYPRE_StructSolver solver;
61 HYPRE_StructSolver precond;
62
63 int num_iterations;
64 double final_res_norm;
65
66 int print_solution;
67
68 /* Initialize MPI */
69 MPI_Init(&argc, &argv);
70 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
71 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
72
73 /* Set defaults */
74 n = 33;
75 solver_id = 0;
76 n_pre = 1;
77 n_post = 1;
78 print_solution = 0;
79
80 /* Parse command line */
81 {
82 int arg_index = 0;
83 int print_usage = 0;
84
85 while (arg_index < argc)
86 {
87 if ( strcmp(argv[arg_index], "-n") == 0 )
88 {
89 arg_index++;
90 n = atoi(argv[arg_index++]);
91 }
92 else if ( strcmp(argv[arg_index], "-solver") == 0 )
93 {
94 arg_index++;
95 solver_id = atoi(argv[arg_index++]);
96 }
97 else if ( strcmp(argv[arg_index], "-v") == 0 )
98 {
99 arg_index++;
100 n_pre = atoi(argv[arg_index++]);
101 n_post = atoi(argv[arg_index++]);
102 }
103 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
104 {
105 arg_index++;
106 print_solution = 1;
107 }
108 else if ( strcmp(argv[arg_index], "-help") == 0 )
109 {
110 print_usage = 1;
111 break;
112 }
113 else
114 {
115 arg_index++;
116 }
117 }
118
119 if ((print_usage) && (myid == 0))
120 {
121 printf("\n");
122 printf("Usage: %s [<options>]\n", argv[0]);
123 printf("\n");
124 printf(" -n <n> : problem size per processor (default: 33)\n");
125 printf(" -solver <ID> : solver ID\n");
126 printf(" 0 - PCG with SMG precond (default)\n");
127 printf(" 1 - SMG\n");
128 printf(" -v <n_pre> <n_post> : number of pre and post relaxations (default: 1 1)\n");
129 printf(" -print_solution : print the solution vector\n");
130 printf("\n");
131 }
132
133 if (print_usage)
134 {
135 MPI_Finalize();
136 return (0);
137 }
138 }
139
140 /* Figure out the processor grid (N x N). The local problem
141 size for the interior nodes is indicated by n (n x n).
142 pi and pj indicate position in the processor grid. */
143 N = sqrt(num_procs);
144 h = 1.0 / (N*n+1); /* note that when calculating h we must
145 remember to count the boundary nodes */
146 h2 = h*h;
147 pj = myid / N;
148 pi = myid - pj*N;
149
150 /* Figure out the extents of each processor's piece of the grid. */
151 ilower[0] = pi*n;
152 ilower[1] = pj*n;
153
154 iupper[0] = ilower[0] + n-1;
155 iupper[1] = ilower[1] + n-1;
156
157 /* 1. Set up a grid */
158 {
159 /* Create an empty 2D grid object */
160 HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
161
162 /* Add a new box to the grid */
163 HYPRE_StructGridSetExtents(grid, ilower, iupper);
164
165 /* This is a collective call finalizing the grid assembly.
166 The grid is now ``ready to be used'' */
167 HYPRE_StructGridAssemble(grid);
168 }
169
170 /* 2. Define the discretization stencil */
171 {
172 /* Create an empty 2D, 5-pt stencil object */
173 HYPRE_StructStencilCreate(2, 5, &stencil);
174
175 /* Define the geometry of the stencil */
176 {
177 int entry;
178 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
179
180 for (entry = 0; entry < 5; entry++)
181 HYPRE_StructStencilSetElement(stencil, entry, offsets[entry]);
182 }
183 }
184
185 /* 3. Set up a Struct Matrix */
186 {
187 int nentries = 5;
188 int nvalues = nentries*n*n;
189 double *values;
190 int stencil_indices[5];
191
192 /* Create an empty matrix object */
193 HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
194
195 /* Indicate that the matrix coefficients are ready to be set */
196 HYPRE_StructMatrixInitialize(A);
197
198 values = calloc(nvalues, sizeof(double));
199
200 for (j = 0; j < nentries; j++)
201 stencil_indices[j] = j;
202
203 /* Set the standard stencil at each grid point,
204 we will fix the boundaries later */
205 for (i = 0; i < nvalues; i += nentries)
206 {
207 values[i] = 4.0;
208 for (j = 1; j < nentries; j++)
209 values[i+j] = -1.0;
210 }
211
212 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
213 stencil_indices, values);
214
215 free(values);
216 }
217
218 /* 4. Incorporate the zero boundary conditions: go along each edge of
219 the domain and set the stencil entry that reaches to the boundary to
220 zero.*/
221 {
222 int bc_ilower[2];
223 int bc_iupper[2];
224 int nentries = 1;
225 int nvalues = nentries*n; /* number of stencil entries times the length
226 of one side of my grid box */
227 double *values;
228 int stencil_indices[1];
229
230 values = calloc(nvalues, sizeof(double));
231 for (j = 0; j < nvalues; j++)
232 values[j] = 0.0;
233
234 /* Recall: pi and pj describe position in the processor grid */
235 if (pj == 0)
236 {
237 /* Bottom row of grid points */
238 bc_ilower[0] = pi*n;
239 bc_ilower[1] = pj*n;
240
241 bc_iupper[0] = bc_ilower[0] + n-1;
242 bc_iupper[1] = bc_ilower[1];
243
244 stencil_indices[0] = 3;
245
246 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
247 stencil_indices, values);
248 }
249
250 if (pj == N-1)
251 {
252 /* upper row of grid points */
253 bc_ilower[0] = pi*n;
254 bc_ilower[1] = pj*n + n-1;
255
256 bc_iupper[0] = bc_ilower[0] + n-1;
257 bc_iupper[1] = bc_ilower[1];
258
259 stencil_indices[0] = 4;
260
261 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
262 stencil_indices, values);
263 }
264
265 if (pi == 0)
266 {
267 /* Left row of grid points */
268 bc_ilower[0] = pi*n;
269 bc_ilower[1] = pj*n;
270
271 bc_iupper[0] = bc_ilower[0];
272 bc_iupper[1] = bc_ilower[1] + n-1;
273
274 stencil_indices[0] = 1;
275
276 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
277 stencil_indices, values);
278 }
279
280 if (pi == N-1)
281 {
282 /* Right row of grid points */
283 bc_ilower[0] = pi*n + n-1;
284 bc_ilower[1] = pj*n;
285
286 bc_iupper[0] = bc_ilower[0];
287 bc_iupper[1] = bc_ilower[1] + n-1;
288
289 stencil_indices[0] = 2;
290
291 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
292 stencil_indices, values);
293 }
294
295 free(values);
296 }
297
298 /* This is a collective call finalizing the matrix assembly.
299 The matrix is now ``ready to be used'' */
300 HYPRE_StructMatrixAssemble(A);
301
302 /* 5. Set up Struct Vectors for b and x */
303 {
304 int nvalues = n*n;
305 double *values;
306
307 values = calloc(nvalues, sizeof(double));
308
309 /* Create an empty vector object */
310 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
311 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
312
313 /* Indicate that the vector coefficients are ready to be set */
314 HYPRE_StructVectorInitialize(b);
315 HYPRE_StructVectorInitialize(x);
316
317 /* Set the values */
318 for (i = 0; i < nvalues; i ++)
319 values[i] = h2;
320 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
321
322 for (i = 0; i < nvalues; i ++)
323 values[i] = 0.0;
324 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
325
326 free(values);
327
328 /* This is a collective call finalizing the vector assembly.
329 The vector is now ``ready to be used'' */
330 HYPRE_StructVectorAssemble(b);
331 HYPRE_StructVectorAssemble(x);
332 }
333
334 /* 6. Set up and use a struct solver
335 (Solver options can be found in the Reference Manual.) */
336 if (solver_id == 0)
337 {
338 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
339 HYPRE_StructPCGSetMaxIter(solver, 50 );
340 HYPRE_StructPCGSetTol(solver, 1.0e-06 );
341 HYPRE_StructPCGSetTwoNorm(solver, 1 );
342 HYPRE_StructPCGSetRelChange(solver, 0 );
343 HYPRE_StructPCGSetPrintLevel(solver, 2 ); /* print each CG iteration */
344 HYPRE_StructPCGSetLogging(solver, 1);
345
346 /* Use symmetric SMG as preconditioner */
347 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
348 HYPRE_StructSMGSetMemoryUse(precond, 0);
349 HYPRE_StructSMGSetMaxIter(precond, 1);
350 HYPRE_StructSMGSetTol(precond, 0.0);
351 HYPRE_StructSMGSetZeroGuess(precond);
352 HYPRE_StructSMGSetNumPreRelax(precond, 1);
353 HYPRE_StructSMGSetNumPostRelax(precond, 1);
354
355 /* Set the preconditioner and solve */
356 HYPRE_StructPCGSetPrecond(solver, HYPRE_StructSMGSolve,
357 HYPRE_StructSMGSetup, precond);
358 HYPRE_StructPCGSetup(solver, A, b, x);
359 HYPRE_StructPCGSolve(solver, A, b, x);
360
361 /* Get some info on the run */
362 HYPRE_StructPCGGetNumIterations(solver, &num_iterations);
363 HYPRE_StructPCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
364
365 /* Clean up */
366 HYPRE_StructPCGDestroy(solver);
367 }
368
369 if (solver_id == 1)
370 {
371 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &solver);
372 HYPRE_StructSMGSetMemoryUse(solver, 0);
373 HYPRE_StructSMGSetMaxIter(solver, 50);
374 HYPRE_StructSMGSetTol(solver, 1.0e-06);
375 HYPRE_StructSMGSetRelChange(solver, 0);
376 HYPRE_StructSMGSetNumPreRelax(solver, n_pre);
377 HYPRE_StructSMGSetNumPostRelax(solver, n_post);
378 /* Logging must be on to get iterations and residual norm info below */
379 HYPRE_StructSMGSetLogging(solver, 1);
380
381 /* Setup and solve */
382 HYPRE_StructSMGSetup(solver, A, b, x);
383 HYPRE_StructSMGSolve(solver, A, b, x);
384
385 /* Get some info on the run */
386 HYPRE_StructSMGGetNumIterations(solver, &num_iterations);
387 HYPRE_StructSMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
388
389 /* Clean up */
390 HYPRE_StructSMGDestroy(solver);
391 }
392
393 /* Print the solution and other info */
394 if (print_solution)
395 HYPRE_StructVectorPrint("struct.out.x", x, 0);
396
397 if (myid == 0)
398 {
399 printf("\n");
400 printf("Iterations = %d\n", num_iterations);
401 printf("Final Relative Residual Norm = %g\n", final_res_norm);
402 printf("\n");
403 }
404
405 /* Free memory */
406 HYPRE_StructGridDestroy(grid);
407 HYPRE_StructStencilDestroy(stencil);
408 HYPRE_StructMatrixDestroy(A);
409 HYPRE_StructVectorDestroy(b);
410 HYPRE_StructVectorDestroy(x);
411
412 /* Finalize MPI */
413 MPI_Finalize();
414
415 return (0);
416 }
syntax highlighted by Code2HTML, v. 0.9.1