download the original source code.
1 /*
2 Example 15
3
4 Interface: Semi-Structured interface (SStruct)
5
6 Compile with: make ex15
7
8 Sample run: mpirun -np 8 ex15 -n 10
9
10 To see options: ex15 -help
11
12 Description: This code solves a 3D electromagnetic diffusion (definite
13 curl-curl) problem using the lowest order Nedelec, or "edge"
14 finite element discretization on a uniform hexahedral meshing
15 of the unit cube. The right-hand-side corresponds to a unit
16 vector force and we use uniform zero Dirichlet boundary
17 conditions. The overall problem reads:
18 curl alpha curl E + beta E = 1,
19 with E x n = 0 on the boundary, where alpha and beta are
20 piecewise-constant material coefficients.
21
22 The linear system is split in parallel using the SStruct
23 interface with an n x n x n grid on each processors, and
24 similar N x N x N processor grid. Therefore, the number of
25 processors should be a perfect cube.
26
27 This example code is mainly meant as an illustration of using
28 the Auxiliary-space Maxwell Solver (AMS) through the SStruct
29 interface. It is also an example of setting up a finite
30 element discretization in the SStruct interface, and we
31 recommend viewing Example 13 and Example 14 before viewing
32 this example.
33 */
34
35 #include <math.h>
36 #include "_hypre_utilities.h"
37 #include "HYPRE_sstruct_mv.h"
38 #include "HYPRE_sstruct_ls.h"
39 #include "_hypre_parcsr_ls.h"
40 #include "HYPRE.h"
41
42
43 int optionAlpha, optionBeta;
44
45 /* Curl-curl coefficient alpha = mu^{-1} */
46 double alpha(double x, double y, double z)
47 {
48 switch (optionAlpha)
49 {
50 case 0: /* uniform coefficient */
51 return 1.0;
52 case 1: /* smooth coefficient */
53 return x*x+exp(y)+sin(z);
54 case 2: /* small outside of an interior cube */
55 if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25) && (fabs(z-0.5) < 0.25))
56 return 1.0;
57 else
58 return 1.0e-6;
59 case 3: /* small outside of an interior ball */
60 if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)+(z-0.5)*(z-0.5)) < 0.0625)
61 return 1.0;
62 else
63 return 1.0e-6;
64 case 4: /* random coefficient */
65 return hypre_Rand();
66 default:
67 return 1.0;
68 }
69 }
70
71 /* Mass coefficient beta = sigma */
72 double beta(double x, double y, double z)
73 {
74 switch (optionBeta)
75 {
76 case 0: /* uniform coefficient */
77 return 1.0;
78 case 1: /* smooth coefficient */
79 return x*x+exp(y)+sin(z);
80 case 2:/* small outside of interior cube */
81 if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25) && (fabs(z-0.5) < 0.25))
82 return 1.0;
83 else
84 return 1.0e-6;
85 case 3: /* small outside of an interior ball */
86 if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)+(z-0.5)*(z-0.5)) < 0.0625)
87 return 1.0;
88 else
89 return 1.0e-6;
90 case 4: /* random coefficient */
91 return hypre_Rand();
92 default:
93 return 1.0;
94 }
95 }
96
97 /*
98 This routine computes the lowest order Nedelec, or "edge" finite element
99 stiffness matrix and load vector on a cube of size h. The 12 edges {e_i}
100 are numbered in terms of the vertices as follows:
101
102 [7]------[6]
103 /| /| e_0 = 01, e_1 = 12, e_2 = 32, e_3 = 03,
104 / | / | e_4 = 45, e_5 = 56, e_6 = 76, e_7 = 47,
105 [4]------[5] | e_8 = 04, e_9 = 15, e_10 = 26, e_11 = 37.
106 | [3]----|-[2]
107 | / | / The edges are oriented from first to the
108 |/ |/ second vertex, e.g. e_0 is from [0] to [1].
109 [0]------[1]
110
111 We allow for different scaling of the curl-curl and the mass parts of the
112 matrix with coefficients alpha and beta respectively:
113
114 S_ij = alpha (curl phi_i,curl phi_j) + beta (phi_i, phi_j).
115
116 The load vector corresponding to a right-hand side of {1,1,1} is
117
118 F_j = (1,phi_j) = h^2/4.
119 */
120 void ComputeFEMND1(double S[12][12], double F[12],
121 double x, double y, double z, double h)
122 {
123 int i, j;
124
125 double h2_4 = h*h/4;
126
127 double cS1 = alpha(x,y,z)/(6.0*h), cS2 = 2*cS1, cS4 = 2*cS2;
128 double cM1 = beta(x,y,z)*h/36.0, cM2 = 2*cM1, cM4 = 2*cM2;
129
130 S[ 0][ 0] = cS4 + cM4; S[ 0][ 1] = cS2; S[ 0][ 2] = -cS1 + cM2;
131 S[ 0][ 3] = -cS2; S[ 0][ 4] = -cS1 + cM2; S[ 0][ 5] = cS1;
132 S[ 0][ 6] = -cS2 + cM1; S[ 0][ 7] = -cS1; S[ 0][ 8] = -cS2;
133 S[ 0][ 9] = cS2; S[ 0][10] = cS1; S[ 0][11] = -cS1;
134
135 S[ 1][ 1] = cS4 + cM4; S[ 1][ 2] = -cS2; S[ 1][ 3] = -cS1 + cM2;
136 S[ 1][ 4] = cS1; S[ 1][ 5] = -cS1 + cM2; S[ 1][ 6] = -cS1;
137 S[ 1][ 7] = -cS2 + cM1; S[ 1][ 8] = -cS1; S[ 1][ 9] = -cS2;
138 S[ 1][10] = cS2; S[ 1][11] = cS1;
139
140 S[ 2][ 2] = cS4 + cM4; S[ 2][ 3] = cS2; S[ 2][ 4] = -cS2 + cM1;
141 S[ 2][ 5] = -cS1; S[ 2][ 6] = -cS1 + cM2; S[ 2][ 7] = cS1;
142 S[ 2][ 8] = -cS1; S[ 2][ 9] = cS1; S[ 2][10] = cS2;
143 S[ 2][11] = -cS2;
144
145 S[ 3][ 3] = cS4 + cM4; S[ 3][ 4] = -cS1; S[ 3][ 5] = -cS2 + cM1;
146 S[ 3][ 6] = cS1; S[ 3][ 7] = -cS1 + cM2; S[ 3][ 8] = -cS2;
147 S[ 3][ 9] = -cS1; S[ 3][10] = cS1; S[ 3][11] = cS2;
148
149 S[ 4][ 4] = cS4 + cM4; S[ 4][ 5] = cS2; S[ 4][ 6] = -cS1 + cM2;
150 S[ 4][ 7] = -cS2; S[ 4][ 8] = cS2; S[ 4][ 9] = -cS2;
151 S[ 4][10] = -cS1; S[ 4][11] = cS1;
152
153 S[ 5][ 5] = cS4 + cM4; S[ 5][ 6] = -cS2; S[ 5][ 7] = -cS1 + cM2;
154 S[ 5][ 8] = cS1; S[ 5][ 9] = cS2; S[ 5][10] = -cS2;
155 S[ 5][11] = -cS1;
156
157 S[ 6][ 6] = cS4 + cM4; S[ 6][ 7] = cS2; S[ 6][ 8] = cS1;
158 S[ 6][ 9] = -cS1; S[ 6][10] = -cS2; S[ 6][11] = cS2;
159
160 S[ 7][ 7] = cS4 + cM4; S[ 7][ 8] = cS2; S[ 7][ 9] = cS1;
161 S[ 7][10] = -cS1; S[ 7][11] = -cS2;
162
163 S[ 8][ 8] = cS4 + cM4; S[ 8][ 9] = -cS1 + cM2; S[ 8][10] = -cS2 + cM1;
164 S[ 8][11] = -cS1 + cM2;
165
166 S[ 9][ 9] = cS4 + cM4; S[ 9][10] = -cS1 + cM2; S[ 9][11] = -cS2 + cM1;
167
168 S[10][10] = cS4 + cM4; S[10][11] = -cS1 + cM2;
169
170 S[11][11] = cS4 + cM4;
171
172 /* The stiffness matrix is symmetric */
173 for (i = 1; i < 12; i++)
174 for (j = 0; j < i; j++)
175 S[i][j] = S[j][i];
176
177 for (i = 0; i < 12; i++)
178 F[i] = h2_4;
179 }
180
181
182 int main (int argc, char *argv[])
183 {
184 int myid, num_procs;
185 int n, N, pi, pj, pk;
186 double h;
187 int print_solution;
188
189 double tol, theta;
190 int maxit, cycle_type;
191 int rlx_type, rlx_sweeps, rlx_weight, rlx_omega;
192 int amg_coarsen_type, amg_agg_levels, amg_rlx_type;
193 int amg_interp_type, amg_Pmax;
194 int singular_problem ;
195
196 int time_index;
197
198 HYPRE_SStructGrid edge_grid;
199 HYPRE_SStructGraph A_graph;
200 HYPRE_SStructMatrix A;
201 HYPRE_SStructVector b;
202 HYPRE_SStructVector x;
203 HYPRE_SStructGrid node_grid;
204 HYPRE_SStructGraph G_graph;
205 HYPRE_SStructStencil G_stencil[3];
206 HYPRE_SStructMatrix G;
207 HYPRE_SStructVector xcoord, ycoord, zcoord;
208
209 HYPRE_Solver solver, precond;
210
211 /* Initialize MPI */
212 MPI_Init(&argc, &argv);
213 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
214 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
215
216 /* Set default parameters */
217 n = 10;
218 print_solution = 0;
219 optionAlpha = 0;
220 optionBeta = 0;
221 maxit = 100;
222 tol = 1e-6;
223 cycle_type = 13;
224 rlx_type = 2;
225 rlx_sweeps = 1;
226 rlx_weight = 1.0;
227 rlx_omega = 1.0;
228 amg_coarsen_type = 10;
229 amg_agg_levels = 1;
230 amg_rlx_type = 6;
231 theta = 0.25;
232 amg_interp_type = 6;
233 amg_Pmax = 4;
234 singular_problem = 0;
235
236 /* Parse command line */
237 {
238 int arg_index = 0;
239 int print_usage = 0;
240
241 while (arg_index < argc)
242 {
243 if ( strcmp(argv[arg_index], "-n") == 0 )
244 {
245 arg_index++;
246 n = atoi(argv[arg_index++]);
247 }
248 else if ( strcmp(argv[arg_index], "-a") == 0 )
249 {
250 arg_index++;
251 optionAlpha = atoi(argv[arg_index++]);
252 }
253 else if ( strcmp(argv[arg_index], "-b") == 0 )
254 {
255 arg_index++;
256 optionBeta = atoi(argv[arg_index++]);
257 }
258 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
259 {
260 arg_index++;
261 print_solution = 1;
262 }
263 else if ( strcmp(argv[arg_index], "-maxit") == 0 )
264 {
265 arg_index++;
266 maxit = atoi(argv[arg_index++]);
267 }
268 else if ( strcmp(argv[arg_index], "-tol") == 0 )
269 {
270 arg_index++;
271 tol = atof(argv[arg_index++]);
272 }
273 else if ( strcmp(argv[arg_index], "-type") == 0 )
274 {
275 arg_index++;
276 cycle_type = atoi(argv[arg_index++]);
277 }
278 else if ( strcmp(argv[arg_index], "-rlx") == 0 )
279 {
280 arg_index++;
281 rlx_type = atoi(argv[arg_index++]);
282 }
283 else if ( strcmp(argv[arg_index], "-rlxn") == 0 )
284 {
285 arg_index++;
286 rlx_sweeps = atoi(argv[arg_index++]);
287 }
288 else if ( strcmp(argv[arg_index], "-rlxw") == 0 )
289 {
290 arg_index++;
291 rlx_weight = atof(argv[arg_index++]);
292 }
293 else if ( strcmp(argv[arg_index], "-rlxo") == 0 )
294 {
295 arg_index++;
296 rlx_omega = atof(argv[arg_index++]);
297 }
298 else if ( strcmp(argv[arg_index], "-ctype") == 0 )
299 {
300 arg_index++;
301 amg_coarsen_type = atoi(argv[arg_index++]);
302 }
303 else if ( strcmp(argv[arg_index], "-amgrlx") == 0 )
304 {
305 arg_index++;
306 amg_rlx_type = atoi(argv[arg_index++]);
307 }
308 else if ( strcmp(argv[arg_index], "-agg") == 0 )
309 {
310 arg_index++;
311 amg_agg_levels = atoi(argv[arg_index++]);
312 }
313 else if ( strcmp(argv[arg_index], "-itype") == 0 )
314 {
315 arg_index++;
316 amg_interp_type = atoi(argv[arg_index++]);
317 }
318 else if ( strcmp(argv[arg_index], "-pmax") == 0 )
319 {
320 arg_index++;
321 amg_Pmax = atoi(argv[arg_index++]);
322 }
323 else if ( strcmp(argv[arg_index], "-sing") == 0 )
324 {
325 arg_index++;
326 singular_problem = 1;
327 }
328 else if ( strcmp(argv[arg_index], "-theta") == 0 )
329 {
330 arg_index++;
331 theta = atof(argv[arg_index++]);
332 }
333
334 else if ( strcmp(argv[arg_index], "-help") == 0 )
335 {
336 print_usage = 1;
337 break;
338 }
339 else
340 {
341 arg_index++;
342 }
343 }
344
345 if ((print_usage) && (myid == 0))
346 {
347 printf("\n");
348 printf("Usage: %s [<options>]\n", argv[0]);
349 printf("\n");
350 printf(" -n <n> : problem size per processor (default: 10)\n");
351 printf(" -a <alpha_opt> : choice for the curl-curl coefficient (default: 1.0)\n");
352 printf(" -b <beta_opt> : choice for the mass coefficient (default: 1.0)\n");
353 printf(" -print_solution : print the solution vector\n");
354 printf("\n");
355 printf("PCG-AMS solver options: \n");
356 printf(" -maxit <num> : maximum number of iterations (100) \n");
357 printf(" -tol <num> : convergence tolerance (1e-6) \n");
358 printf(" -type <num> : 3-level cycle type (0-8, 11-14) \n");
359 printf(" -theta <num> : BoomerAMG threshold (0.25) \n");
360 printf(" -ctype <num> : BoomerAMG coarsening type \n");
361 printf(" -agg <num> : Levels of BoomerAMG agg. coarsening \n");
362 printf(" -amgrlx <num> : BoomerAMG relaxation type \n");
363 printf(" -itype <num> : BoomerAMG interpolation type \n");
364 printf(" -pmax <num> : BoomerAMG interpolation truncation \n");
365 printf(" -rlx <num> : relaxation type \n");
366 printf(" -rlxn <num> : number of relaxation sweeps \n");
367 printf(" -rlxw <num> : damping parameter (usually <=1) \n");
368 printf(" -rlxo <num> : SOR parameter (usually in (0,2)) \n");
369 printf(" -sing : curl-curl only (singular) problem \n");
370 printf("\n");
371 printf("\n");
372 }
373
374 if (print_usage)
375 {
376 MPI_Finalize();
377 return (0);
378 }
379 }
380
381 /* Figure out the processor grid (N x N x N). The local problem size is n^3,
382 while pi, pj and pk indicate the position in the processor grid. */
383 N = pow(num_procs,1.0/3.0) + 0.5;
384 if (num_procs != N*N*N)
385 {
386 if (myid == 0) printf("Can't run on %d processors, try %d.\n",
387 num_procs, N*N*N);
388 MPI_Finalize();
389 exit(1);
390 }
391 h = 1.0 / (N*n);
392 pk = myid / (N*N);
393 pj = myid/N - pk*N;
394 pi = myid - pj*N - pk*N*N;
395
396 /* Start timing */
397 time_index = hypre_InitializeTiming("SStruct Setup");
398 hypre_BeginTiming(time_index);
399
400 /* 1. Set up the edge and nodal grids. Note that we do this simultaneously
401 to make sure that they have the same extends. For simplicity we use
402 only one part to represent the unit cube. */
403 {
404 int ndim = 3;
405 int nparts = 1;
406
407 /* Create empty 2D grid objects */
408 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &node_grid);
409 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &edge_grid);
410
411 /* Set the extents of the grid - each processor sets its grid boxes. */
412 {
413 int part = 0;
414 int ilower[3] = {1 + pi*n, 1 + pj*n, 1 + pk*n};
415 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
416
417 HYPRE_SStructGridSetExtents(node_grid, part, ilower, iupper);
418 HYPRE_SStructGridSetExtents(edge_grid, part, ilower, iupper);
419 }
420
421 /* Set the variable type and number of variables on each grid. */
422 {
423 int i;
424 int nnodevars = 1;
425 int nedgevars = 3;
426
427 HYPRE_SStructVariable nodevars[1] = {HYPRE_SSTRUCT_VARIABLE_NODE};
428 HYPRE_SStructVariable edgevars[3] = {HYPRE_SSTRUCT_VARIABLE_XEDGE,
429 HYPRE_SSTRUCT_VARIABLE_YEDGE,
430 HYPRE_SSTRUCT_VARIABLE_ZEDGE};
431 for (i = 0; i < nparts; i++)
432 {
433 HYPRE_SStructGridSetVariables(node_grid, i, nnodevars, nodevars);
434 HYPRE_SStructGridSetVariables(edge_grid, i, nedgevars, edgevars);
435 }
436 }
437
438 /* Since there is only one part, there is no need to call the
439 SetNeighborPart or SetSharedPart functions, which determine the spatial
440 relation between the parts. See Examples 12, 13 and 14 for
441 illustrations of these calls. */
442
443 /* Now the grids are ready to be used */
444 HYPRE_SStructGridAssemble(node_grid);
445 HYPRE_SStructGridAssemble(edge_grid);
446 }
447
448 /* 2. Create the finite element stiffness matrix A and load vector b. */
449 {
450 int part = 0; /* this problem has only one part */
451
452 /* Set the ordering of the variables in the finite element problem. This
453 is done by listing the variable offset directions relative to the
454 element's center. See the Reference Manual for more details. */
455 {
456 int ordering[48] = { 0, 0, -1, -1, /* x-edge [0]-[1] */
457 1, +1, 0, -1, /* y-edge [1]-[2] */
458 /* [7]------[6] */ 0, 0, +1, -1, /* x-edge [3]-[2] */
459 /* /| /| */ 1, -1, 0, -1, /* y-edge [0]-[3] */
460 /* / | / | */ 0, 0, -1, +1, /* x-edge [4]-[5] */
461 /* [4]------[5] | */ 1, +1, 0, +1, /* y-edge [5]-[6] */
462 /* | [3]----|-[2] */ 0, 0, +1, +1, /* x-edge [7]-[6] */
463 /* | / | / */ 1, -1, 0, +1, /* y-edge [4]-[7] */
464 /* |/ |/ */ 2, -1, -1, 0, /* z-edge [0]-[4] */
465 /* [0]------[1] */ 2, +1, -1, 0, /* z-edge [1]-[5] */
466 2, +1, +1, 0, /* z-edge [2]-[6] */
467 2, -1, +1, 0 }; /* z-edge [3]-[7] */
468
469 HYPRE_SStructGridSetFEMOrdering(edge_grid, part, ordering);
470 }
471
472 /* Set up the Graph - this determines the non-zero structure of the
473 matrix. */
474 {
475 int part = 0;
476
477 /* Create the graph object */
478 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, edge_grid, &A_graph);
479
480 /* Indicate that this problem uses finite element stiffness matrices and
481 load vectors, instead of stencils. */
482 HYPRE_SStructGraphSetFEM(A_graph, part);
483
484 /* The edge finite element matrix is full, so there is no need to call the
485 HYPRE_SStructGraphSetFEMSparsity() function. */
486
487 /* Assemble the graph */
488 HYPRE_SStructGraphAssemble(A_graph);
489 }
490
491 /* Set up the SStruct Matrix and right-hand side vector */
492 {
493 /* Create the matrix object */
494 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, A_graph, &A);
495 /* Use a ParCSR storage */
496 HYPRE_SStructMatrixSetObjectType(A, HYPRE_PARCSR);
497 /* Indicate that the matrix coefficients are ready to be set */
498 HYPRE_SStructMatrixInitialize(A);
499
500 /* Create an empty vector object */
501 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, edge_grid, &b);
502 /* Use a ParCSR storage */
503 HYPRE_SStructVectorSetObjectType(b, HYPRE_PARCSR);
504 /* Indicate that the vector coefficients are ready to be set */
505 HYPRE_SStructVectorInitialize(b);
506 }
507
508 /* Set the matrix and vector entries by finite element assembly */
509 {
510 /* local stiffness matrix and load vector */
511 double S[12][12], F[12];
512
513 int i, j, k;
514 int index[3];
515
516 for (i = 1; i <= n; i++)
517 for (j = 1; j <= n; j++)
518 for (k = 1; k <= n; k++)
519 {
520 /* Compute the FEM matrix and r.h.s. for cell (i,j,k) with
521 coefficients evaluated at the cell center. */
522 index[0] = i + pi*n; index[1] = j + pj*n; index[2] = k + pk*n;
523 ComputeFEMND1(S,F,(pi*n+i)*h-h/2,(pj*n+j)*h-h/2,(pk*n+k)*h-h/2,h);
524
525 /* Eliminate boundary conditions on x = 0 */
526 if (index[0] == 1)
527 {
528 int ii, jj, bc_edges[4] = { 3, 11, 7, 8 };
529 for (ii = 0; ii < 4; ii++)
530 {
531 for (jj = 0; jj < 12; jj++)
532 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
533 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
534 F[bc_edges[ii]] = 0.0;
535 }
536 }
537 /* Eliminate boundary conditions on y = 0 */
538 if (index[1] == 1)
539 {
540 int ii, jj, bc_edges[4] = { 0, 9, 4, 8 };
541 for (ii = 0; ii < 4; ii++)
542 {
543 for (jj = 0; jj < 12; jj++)
544 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
545 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
546 F[bc_edges[ii]] = 0.0;
547 }
548 }
549 /* Eliminate boundary conditions on z = 0 */
550 if (index[2] == 1)
551 {
552 int ii, jj, bc_edges[4] = { 0, 1, 2, 3 };
553 for (ii = 0; ii < 4; ii++)
554 {
555 for (jj = 0; jj < 12; jj++)
556 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
557 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
558 F[bc_edges[ii]] = 0.0;
559 }
560 }
561 /* Eliminate boundary conditions on x = 1 */
562 if (index[0] == N*n)
563 {
564 int ii, jj, bc_edges[4] = { 1, 10, 5, 9 };
565 for (ii = 0; ii < 4; ii++)
566 {
567 for (jj = 0; jj < 12; jj++)
568 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
569 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
570 F[bc_edges[ii]] = 0.0;
571 }
572 }
573 /* Eliminate boundary conditions on y = 1 */
574 if (index[1] == N*n)
575 {
576 int ii, jj, bc_edges[4] = { 2, 10, 6, 11 };
577 for (ii = 0; ii < 4; ii++)
578 {
579 for (jj = 0; jj < 12; jj++)
580 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
581 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
582 F[bc_edges[ii]] = 0.0;
583 }
584 }
585 /* Eliminate boundary conditions on z = 1 */
586 if (index[2] == N*n)
587 {
588 int ii, jj, bc_edges[4] = { 4, 5, 6, 7 };
589 for (ii = 0; ii < 4; ii++)
590 {
591 for (jj = 0; jj < 12; jj++)
592 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
593 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
594 F[bc_edges[ii]] = 0.0;
595 }
596 }
597
598 /* Assemble the matrix */
599 HYPRE_SStructMatrixAddFEMValues(A, part, index, &S[0][0]);
600
601 /* Assemble the vector */
602 HYPRE_SStructVectorAddFEMValues(b, part, index, F);
603 }
604 }
605
606 /* Collective calls finalizing the matrix and vector assembly */
607 HYPRE_SStructMatrixAssemble(A);
608 HYPRE_SStructVectorAssemble(b);
609 }
610
611 /* 3. Create the discrete gradient matrix G, which is needed in AMS. */
612 {
613 int part = 0;
614 int stencil_size = 2;
615
616 /* Define the discretization stencil relating the edges and nodes of the
617 grid. */
618 {
619 int ndim = 3;
620 int entry;
621 int var = 0; /* the node variable */
622
623 /* The discrete gradient stencils connect edge to node variables. */
624 int Gx_offsets[2][3] = {{-1,0,0},{0,0,0}}; /* x-edge [7]-[6] */
625 int Gy_offsets[2][3] = {{0,-1,0},{0,0,0}}; /* y-edge [5]-[6] */
626 int Gz_offsets[2][3] = {{0,0,-1},{0,0,0}}; /* z-edge [2]-[6] */
627
628 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[0]);
629 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[1]);
630 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[2]);
631
632 for (entry = 0; entry < stencil_size; entry++)
633 {
634 HYPRE_SStructStencilSetEntry(G_stencil[0], entry, Gx_offsets[entry], var);
635 HYPRE_SStructStencilSetEntry(G_stencil[1], entry, Gy_offsets[entry], var);
636 HYPRE_SStructStencilSetEntry(G_stencil[2], entry, Gz_offsets[entry], var);
637 }
638 }
639
640 /* Set up the Graph - this determines the non-zero structure of the
641 matrix. */
642 {
643 int nvars = 3;
644 int var; /* the edge variables */
645
646 /* Create the discrete gradient graph object */
647 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, edge_grid, &G_graph);
648
649 /* Since the discrete gradient relates edge and nodal variables (it is a
650 rectangular matrix), we have to specify the domain (column) grid. */
651 HYPRE_SStructGraphSetDomainGrid(G_graph, node_grid);
652
653 /* Tell the graph which stencil to use for each edge variable on each
654 part (we only have one part). */
655 for (var = 0; var < nvars; var++)
656 HYPRE_SStructGraphSetStencil(G_graph, part, var, G_stencil[var]);
657
658 /* Assemble the graph */
659 HYPRE_SStructGraphAssemble(G_graph);
660 }
661
662 /* Set up the SStruct Matrix */
663 {
664 /* Create the matrix object */
665 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, G_graph, &G);
666 /* Use a ParCSR storage */
667 HYPRE_SStructMatrixSetObjectType(G, HYPRE_PARCSR);
668 /* Indicate that the matrix coefficients are ready to be set */
669 HYPRE_SStructMatrixInitialize(G);
670 }
671
672 /* Set the discrete gradient values, assuming a "natural" orientation of
673 the edges (i.e. one in agreement with the coordinate directions). */
674 {
675 int i;
676 int nedges = n*(n+1)*(n+1);
677 double *values;
678 int stencil_indices[2] = {0,1}; /* the nodes of each edge */
679
680 values = calloc(2*nedges, sizeof(double));
681
682 /* The edge orientation is fixed: from first to second node */
683 for (i = 0; i < nedges; i++)
684 {
685 values[2*i] = -1.0;
686 values[2*i+1] = 1.0;
687 }
688
689 /* Set the values in the discrete gradient x-edges */
690 {
691 int var = 0;
692 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
693 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
694 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
695 stencil_size, stencil_indices,
696 values);
697 }
698 /* Set the values in the discrete gradient y-edges */
699 {
700 int var = 1;
701 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
702 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
703 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
704 stencil_size, stencil_indices,
705 values);
706 }
707 /* Set the values in the discrete gradient z-edges */
708 {
709 int var = 2;
710 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
711 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
712 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
713 stencil_size, stencil_indices,
714 values);
715 }
716
717 free(values);
718 }
719
720 /* Finalize the matrix assembly */
721 HYPRE_SStructMatrixAssemble(G);
722 }
723
724 /* 4. Create the vectors of nodal coordinates xcoord, ycoord and zcoord,
725 which are needed in AMS. */
726 {
727 int i, j, k;
728 int part = 0;
729 int var = 0; /* the node variable */
730 int index[3];
731 double xval, yval, zval;
732
733 /* Create empty vector objects */
734 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &xcoord);
735 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &ycoord);
736 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &zcoord);
737 /* Set the object type to ParCSR */
738 HYPRE_SStructVectorSetObjectType(xcoord, HYPRE_PARCSR);
739 HYPRE_SStructVectorSetObjectType(ycoord, HYPRE_PARCSR);
740 HYPRE_SStructVectorSetObjectType(zcoord, HYPRE_PARCSR);
741 /* Indicate that the vector coefficients are ready to be set */
742 HYPRE_SStructVectorInitialize(xcoord);
743 HYPRE_SStructVectorInitialize(ycoord);
744 HYPRE_SStructVectorInitialize(zcoord);
745
746 /* Compute and set the coordinates of the nodes */
747 for (i = 0; i <= n; i++)
748 for (j = 0; j <= n; j++)
749 for (k = 0; k <= n; k++)
750 {
751 index[0] = i + pi*n; index[1] = j + pj*n; index[2] = k + pk*n;
752
753 xval = index[0]*h;
754 yval = index[1]*h;
755 zval = index[2]*h;
756
757 HYPRE_SStructVectorSetValues(xcoord, part, index, var, &xval);
758 HYPRE_SStructVectorSetValues(ycoord, part, index, var, &yval);
759 HYPRE_SStructVectorSetValues(zcoord, part, index, var, &zval);
760 }
761
762 /* Finalize the vector assembly */
763 HYPRE_SStructVectorAssemble(xcoord);
764 HYPRE_SStructVectorAssemble(ycoord);
765 HYPRE_SStructVectorAssemble(zcoord);
766 }
767
768 /* 5. Set up a SStruct Vector for the solution vector x */
769 {
770 int part = 0;
771 int nvalues = n*(n+1)*(n+1);
772 double *values;
773
774 values = calloc(nvalues, sizeof(double));
775
776 /* Create an empty vector object */
777 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, edge_grid, &x);
778 /* Set the object type to ParCSR */
779 HYPRE_SStructVectorSetObjectType(x, HYPRE_PARCSR);
780 /* Indicate that the vector coefficients are ready to be set */
781 HYPRE_SStructVectorInitialize(x);
782
783 /* Set the values for the initial guess x-edge */
784 {
785 int var = 0;
786 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
787 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
788 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
789 }
790 /* Set the values for the initial guess y-edge */
791 {
792 int var = 1;
793 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
794 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
795 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
796 }
797 /* Set the values for the initial guess z-edge */
798 {
799 int var = 2;
800 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
801 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
802 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
803 }
804
805 free(values);
806
807 /* Finalize the vector assembly */
808 HYPRE_SStructVectorAssemble(x);
809 }
810
811 /* Finalize current timing */
812 hypre_EndTiming(time_index);
813 hypre_PrintTiming("SStruct phase times", MPI_COMM_WORLD);
814 hypre_FinalizeTiming(time_index);
815 hypre_ClearTiming();
816
817 /* 6. Set up and call the PCG-AMS solver (Solver options can be found in the
818 Reference Manual.) */
819 {
820 double final_res_norm;
821 int its;
822
823 HYPRE_ParCSRMatrix par_A;
824 HYPRE_ParVector par_b;
825 HYPRE_ParVector par_x;
826
827 HYPRE_ParCSRMatrix par_G;
828 HYPRE_ParVector par_xcoord;
829 HYPRE_ParVector par_ycoord;
830 HYPRE_ParVector par_zcoord;
831
832 /* Extract the ParCSR objects needed in the solver */
833 HYPRE_SStructMatrixGetObject(A, (void **) &par_A);
834 HYPRE_SStructVectorGetObject(b, (void **) &par_b);
835 HYPRE_SStructVectorGetObject(x, (void **) &par_x);
836 HYPRE_SStructMatrixGetObject(G, (void **) &par_G);
837 HYPRE_SStructVectorGetObject(xcoord, (void **) &par_xcoord);
838 HYPRE_SStructVectorGetObject(ycoord, (void **) &par_ycoord);
839 HYPRE_SStructVectorGetObject(zcoord, (void **) &par_zcoord);
840
841 if (myid == 0)
842 printf("Problem size: %d\n\n",
843 hypre_ParCSRMatrixGlobalNumRows((hypre_ParCSRMatrix*)par_A));
844
845 /* Start timing */
846 time_index = hypre_InitializeTiming("AMS Setup");
847 hypre_BeginTiming(time_index);
848
849 /* Create solver */
850 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
851
852 /* Set some parameters (See Reference Manual for more parameters) */
853 HYPRE_PCGSetMaxIter(solver, maxit); /* max iterations */
854 HYPRE_PCGSetTol(solver, tol); /* conv. tolerance */
855 HYPRE_PCGSetTwoNorm(solver, 0); /* use the two norm as the stopping criteria */
856 HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
857 HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
858
859 /* Create AMS preconditioner */
860 HYPRE_AMSCreate(&precond);
861
862 /* Set AMS parameters */
863 HYPRE_AMSSetMaxIter(precond, 1);
864 HYPRE_AMSSetTol(precond, 0.0);
865 HYPRE_AMSSetCycleType(precond, cycle_type);
866 HYPRE_AMSSetPrintLevel(precond, 1);
867
868 /* Set discrete gradient */
869 HYPRE_AMSSetDiscreteGradient(precond, par_G);
870
871 /* Set vertex coordinates */
872 HYPRE_AMSSetCoordinateVectors(precond,
873 par_xcoord, par_ycoord, par_zcoord);
874
875 if (singular_problem)
876 HYPRE_AMSSetBetaPoissonMatrix(precond, NULL);
877
878 /* Smoothing and AMG options */
879 HYPRE_AMSSetSmoothingOptions(precond,
880 rlx_type, rlx_sweeps,
881 rlx_weight, rlx_omega);
882 HYPRE_AMSSetAlphaAMGOptions(precond,
883 amg_coarsen_type, amg_agg_levels,
884 amg_rlx_type, theta, amg_interp_type,
885 amg_Pmax);
886 HYPRE_AMSSetBetaAMGOptions(precond,
887 amg_coarsen_type, amg_agg_levels,
888 amg_rlx_type, theta, amg_interp_type,
889 amg_Pmax);
890
891 /* Set the PCG preconditioner */
892 HYPRE_PCGSetPrecond(solver,
893 (HYPRE_PtrToSolverFcn) HYPRE_AMSSolve,
894 (HYPRE_PtrToSolverFcn) HYPRE_AMSSetup,
895 precond);
896
897 /* Call the setup */
898 HYPRE_ParCSRPCGSetup(solver, par_A, par_b, par_x);
899
900 /* Finalize current timing */
901 hypre_EndTiming(time_index);
902 hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
903 hypre_FinalizeTiming(time_index);
904 hypre_ClearTiming();
905
906 /* Start timing again */
907 time_index = hypre_InitializeTiming("AMS Solve");
908 hypre_BeginTiming(time_index);
909
910 /* Call the solve */
911 HYPRE_ParCSRPCGSolve(solver, par_A, par_b, par_x);
912
913 /* Finalize current timing */
914 hypre_EndTiming(time_index);
915 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
916 hypre_FinalizeTiming(time_index);
917 hypre_ClearTiming();
918
919 /* Get some info */
920 HYPRE_PCGGetNumIterations(solver, &its);
921 HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
922
923 /* Clean up */
924 HYPRE_AMSDestroy(precond);
925 HYPRE_ParCSRPCGDestroy(solver);
926
927 /* Gather the solution vector */
928 HYPRE_SStructVectorGather(x);
929
930 /* Print the solution with replicated shared data */
931 if (print_solution)
932 {
933 FILE *file;
934 char filename[255];
935
936 int part = 0;
937 int nvalues = n*(n+1)*(n+1);
938 double *xvalues, *yvalues, *zvalues;
939
940 xvalues = calloc(nvalues, sizeof(double));
941 yvalues = calloc(nvalues, sizeof(double));
942 zvalues = calloc(nvalues, sizeof(double));
943
944 /* Get local solution in the x-edges */
945 {
946 int var = 0;
947 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
948 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
949 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
950 var, xvalues);
951 }
952 /* Get local solution in the y-edges */
953 {
954 int var = 1;
955 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
956 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
957 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
958 var, yvalues);
959 }
960 /* Get local solution in the z-edges */
961 {
962 int var = 2;
963 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
964 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
965 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
966 var, zvalues);
967 }
968
969 sprintf(filename, "sstruct.out.x.00.00.%05d", myid);
970 if ((file = fopen(filename, "w")) == NULL)
971 {
972 printf("Error: can't open output file %s\n", filename);
973 MPI_Finalize();
974 exit(1);
975 }
976
977 /* Save the edge values, element by element, using the same numbering
978 as the local finite element degrees of freedom. */
979 {
980 int i, j, k, s;
981
982 /* Initial x-, y- and z-edge indices in the values arrays */
983 int oi[4] = { 0, n, n*(n+1), n*(n+1)+n }; /* e_0, e_2, e_4, e_6 */
984 int oj[4] = { 0, 1, n*(n+1), n*(n+1)+1 }; /* e_3, e_1, e_7, e_5 */
985 int ok[4] = { 0, 1, n+1, n+2 }; /* e_8, e_9, e_11, e_10 */
986 /* Loop over the cells while updating the above offsets */
987 for (k = 0; k < n; k++)
988 {
989 for (j = 0; j < n; j++)
990 {
991 for (i = 0; i < n; i++)
992 {
993 fprintf(file,
994 "%.14e\n%.14e\n%.14e\n%.14e\n"
995 "%.14e\n%.14e\n%.14e\n%.14e\n"
996 "%.14e\n%.14e\n%.14e\n%.14e\n",
997 xvalues[oi[0]], yvalues[oj[1]], xvalues[oi[1]], yvalues[oj[0]],
998 xvalues[oi[2]], yvalues[oj[3]], xvalues[oi[3]], yvalues[oj[2]],
999 zvalues[ok[0]], zvalues[ok[1]], zvalues[ok[3]], zvalues[ok[2]]);
1000
1001 for (s=0; s<4; s++) oi[s]++, oj[s]++, ok[s]++;
1002 }
1003 for (s=0; s<4; s++) oj[s]++, ok[s]++;
1004 }
1005 for (s=0; s<4; s++) oi[s]+=n, ok[s]+=n+1;
1006 }
1007 }
1008
1009 fflush(file);
1010 fclose(file);
1011
1012 free(xvalues);
1013 free(yvalues);
1014 free(zvalues);
1015 }
1016
1017 if (myid == 0)
1018 {
1019 printf("\n");
1020 printf("Iterations = %d\n", its);
1021 printf("Final Relative Residual Norm = %g\n", final_res_norm);
1022 printf("\n");
1023 }
1024 }
1025
1026 /* Free memory */
1027 HYPRE_SStructGridDestroy(edge_grid);
1028 HYPRE_SStructGraphDestroy(A_graph);
1029 HYPRE_SStructMatrixDestroy(A);
1030 HYPRE_SStructVectorDestroy(b);
1031 HYPRE_SStructVectorDestroy(x);
1032 HYPRE_SStructGridDestroy(node_grid);
1033 HYPRE_SStructGraphDestroy(G_graph);
1034 HYPRE_SStructStencilDestroy(G_stencil[0]);
1035 HYPRE_SStructStencilDestroy(G_stencil[1]);
1036 HYPRE_SStructStencilDestroy(G_stencil[2]);
1037 HYPRE_SStructMatrixDestroy(G);
1038 HYPRE_SStructVectorDestroy(xcoord);
1039 HYPRE_SStructVectorDestroy(ycoord);
1040 HYPRE_SStructVectorDestroy(zcoord);
1041
1042 /* Finalize MPI */
1043 MPI_Finalize();
1044
1045 return 0;
1046 }
syntax highlighted by Code2HTML, v. 0.9.1