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