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 HYPRE_SStructGrid edge_grid;
197 HYPRE_SStructGraph A_graph;
198 HYPRE_SStructMatrix A;
199 HYPRE_SStructVector b;
200 HYPRE_SStructVector x;
201 HYPRE_SStructGrid node_grid;
202 HYPRE_SStructGraph G_graph;
203 HYPRE_SStructStencil G_stencil[3];
204 HYPRE_SStructMatrix G;
205 HYPRE_SStructVector xcoord, ycoord, zcoord;
206
207 HYPRE_Solver solver, precond;
208
209 /* Initialize MPI */
210 MPI_Init(&argc, &argv);
211 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
212 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
213
214 /* Set default parameters */
215 n = 10;
216 print_solution = 0;
217 optionAlpha = 0;
218 optionBeta = 0;
219 maxit = 20;
220 tol = 1e-6;
221 cycle_type = 13;
222 rlx_type = 2;
223 rlx_sweeps = 1;
224 rlx_weight = 1.0;
225 rlx_omega = 1.0;
226 amg_coarsen_type = 10;
227 amg_agg_levels = 1;
228 amg_rlx_type = 6;
229 theta = 0.25;
230 amg_interp_type = 6;
231 amg_Pmax = 4;
232 singular_problem = 0;
233
234 /* Parse command line */
235 {
236 int arg_index = 0;
237 int print_usage = 0;
238
239 while (arg_index < argc)
240 {
241 if ( strcmp(argv[arg_index], "-n") == 0 )
242 {
243 arg_index++;
244 n = atoi(argv[arg_index++]);
245 }
246 else if ( strcmp(argv[arg_index], "-a") == 0 )
247 {
248 arg_index++;
249 optionAlpha = atoi(argv[arg_index++]);
250 }
251 else if ( strcmp(argv[arg_index], "-b") == 0 )
252 {
253 arg_index++;
254 optionBeta = atoi(argv[arg_index++]);
255 }
256 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
257 {
258 arg_index++;
259 print_solution = 1;
260 }
261 else if ( strcmp(argv[arg_index], "-maxit") == 0 )
262 {
263 arg_index++;
264 maxit = atoi(argv[arg_index++]);
265 }
266 else if ( strcmp(argv[arg_index], "-tol") == 0 )
267 {
268 arg_index++;
269 tol = atof(argv[arg_index++]);
270 }
271 else if ( strcmp(argv[arg_index], "-type") == 0 )
272 {
273 arg_index++;
274 cycle_type = atoi(argv[arg_index++]);
275 }
276 else if ( strcmp(argv[arg_index], "-rlx") == 0 )
277 {
278 arg_index++;
279 rlx_type = atoi(argv[arg_index++]);
280 }
281 else if ( strcmp(argv[arg_index], "-rlxn") == 0 )
282 {
283 arg_index++;
284 rlx_sweeps = atoi(argv[arg_index++]);
285 }
286 else if ( strcmp(argv[arg_index], "-rlxw") == 0 )
287 {
288 arg_index++;
289 rlx_weight = atof(argv[arg_index++]);
290 }
291 else if ( strcmp(argv[arg_index], "-rlxo") == 0 )
292 {
293 arg_index++;
294 rlx_omega = atof(argv[arg_index++]);
295 }
296 else if ( strcmp(argv[arg_index], "-ctype") == 0 )
297 {
298 arg_index++;
299 amg_coarsen_type = atoi(argv[arg_index++]);
300 }
301 else if ( strcmp(argv[arg_index], "-amgrlx") == 0 )
302 {
303 arg_index++;
304 amg_rlx_type = atoi(argv[arg_index++]);
305 }
306 else if ( strcmp(argv[arg_index], "-agg") == 0 )
307 {
308 arg_index++;
309 amg_agg_levels = atoi(argv[arg_index++]);
310 }
311 else if ( strcmp(argv[arg_index], "-itype") == 0 )
312 {
313 arg_index++;
314 amg_interp_type = atoi(argv[arg_index++]);
315 }
316 else if ( strcmp(argv[arg_index], "-pmax") == 0 )
317 {
318 arg_index++;
319 amg_Pmax = atoi(argv[arg_index++]);
320 }
321 else if ( strcmp(argv[arg_index], "-sing") == 0 )
322 {
323 arg_index++;
324 singular_problem = 1;
325 }
326 else if ( strcmp(argv[arg_index], "-theta") == 0 )
327 {
328 arg_index++;
329 theta = atof(argv[arg_index++]);
330 }
331
332 else if ( strcmp(argv[arg_index], "-help") == 0 )
333 {
334 print_usage = 1;
335 break;
336 }
337 else
338 {
339 arg_index++;
340 }
341 }
342
343 if ((print_usage) && (myid == 0))
344 {
345 printf("\n");
346 printf("Usage: %s [<options>]\n", argv[0]);
347 printf("\n");
348 printf(" -n <n> : problem size per processor (default: 10)\n");
349 printf(" -a <alpha_opt> : choice for the curl-curl coefficient (default: 1.0)\n");
350 printf(" -b <beta_opt> : choice for the mass coefficient (default: 1.0)\n");
351 printf(" -print_solution : print the solution vector\n");
352 printf("\n");
353 printf("PCG-AMS solver options: \n");
354 printf(" -maxit <num> : maximum number of iterations (100) \n");
355 printf(" -tol <num> : convergence tolerance (1e-6) \n");
356 printf(" -type <num> : 3-level cycle type (0-8, 11-14) \n");
357 printf(" -theta <num> : BoomerAMG threshold (0.25) \n");
358 printf(" -ctype <num> : BoomerAMG coarsening type \n");
359 printf(" -agg <num> : Levels of BoomerAMG agg. coarsening \n");
360 printf(" -amgrlx <num> : BoomerAMG relaxation type \n");
361 printf(" -itype <num> : BoomerAMG interpolation type \n");
362 printf(" -pmax <num> : BoomerAMG interpolation truncation \n");
363 printf(" -rlx <num> : relaxation type \n");
364 printf(" -rlxn <num> : number of relaxation sweeps \n");
365 printf(" -rlxw <num> : damping parameter (usually <=1) \n");
366 printf(" -rlxo <num> : SOR parameter (usually in (0,2)) \n");
367 printf(" -sing : curl-curl only (singular) problem \n");
368 printf("\n");
369 printf("\n");
370 }
371
372 if (print_usage)
373 {
374 MPI_Finalize();
375 return (0);
376 }
377 }
378
379 /* Figure out the processor grid (N x N x N). The local problem size is n^3,
380 while pi, pj and pk indicate the position in the processor grid. */
381 N = pow(num_procs,1.0/3.0) + 0.5;
382 if (num_procs != N*N*N)
383 {
384 if (myid == 0) printf("Can't run on %d processors, try %d.\n",
385 num_procs, N*N*N);
386 MPI_Finalize();
387 exit(1);
388 }
389 h = 1.0 / (N*n);
390 pk = myid / (N*N);
391 pj = myid/N - pk*N;
392 pi = myid - pj*N - pk*N*N;
393
394 /* 1. Set up the edge and nodal grids. Note that we do this simultaneously
395 to make sure that they have the same extends. For simplicity we use
396 only one part to represent the unit cube. */
397 {
398 int ndim = 3;
399 int nparts = 1;
400
401 /* Create empty 2D grid objects */
402 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &node_grid);
403 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &edge_grid);
404
405 /* Set the extents of the grid - each processor sets its grid boxes. */
406 {
407 int part = 0;
408 int ilower[3] = {1 + pi*n, 1 + pj*n, 1 + pk*n};
409 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
410
411 HYPRE_SStructGridSetExtents(node_grid, part, ilower, iupper);
412 HYPRE_SStructGridSetExtents(edge_grid, part, ilower, iupper);
413 }
414
415 /* Set the variable type and number of variables on each grid. */
416 {
417 int i;
418 int nnodevars = 1;
419 int nedgevars = 3;
420
421 HYPRE_SStructVariable nodevars[1] = {HYPRE_SSTRUCT_VARIABLE_NODE};
422 HYPRE_SStructVariable edgevars[3] = {HYPRE_SSTRUCT_VARIABLE_XEDGE,
423 HYPRE_SSTRUCT_VARIABLE_YEDGE,
424 HYPRE_SSTRUCT_VARIABLE_ZEDGE};
425 for (i = 0; i < nparts; i++)
426 {
427 HYPRE_SStructGridSetVariables(node_grid, i, nnodevars, nodevars);
428 HYPRE_SStructGridSetVariables(edge_grid, i, nedgevars, edgevars);
429 }
430 }
431
432 /* Since there is only one part, there is no need to call the
433 SetNeighborPart or SetSharedPart functions, which determine the spatial
434 relation between the parts. See Examples 12, 13 and 14 for
435 illustrations of these calls. */
436
437 /* Now the grids are ready to be used */
438 HYPRE_SStructGridAssemble(node_grid);
439 HYPRE_SStructGridAssemble(edge_grid);
440 }
441
442 /* 2. Create the finite element stiffness matrix A and load vector b. */
443 {
444 int part = 0; /* this problem has only one part */
445
446 /* Set the ordering of the variables in the finite element problem. This
447 is done by listing the variable offset directions relative to the
448 element's center. See the Reference Manual for more details. */
449 {
450 int ordering[48] = { 0, 0, -1, -1, /* x-edge [0]-[1] */
451 1, +1, 0, -1, /* y-edge [1]-[2] */
452 /* [7]------[6] */ 0, 0, +1, -1, /* x-edge [3]-[2] */
453 /* /| /| */ 1, -1, 0, -1, /* y-edge [0]-[3] */
454 /* / | / | */ 0, 0, -1, +1, /* x-edge [4]-[5] */
455 /* [4]------[5] | */ 1, +1, 0, +1, /* y-edge [5]-[6] */
456 /* | [3]----|-[2] */ 0, 0, +1, +1, /* x-edge [7]-[6] */
457 /* | / | / */ 1, -1, 0, +1, /* y-edge [4]-[7] */
458 /* |/ |/ */ 2, -1, -1, 0, /* z-edge [0]-[4] */
459 /* [0]------[1] */ 2, +1, -1, 0, /* z-edge [1]-[5] */
460 2, +1, +1, 0, /* z-edge [2]-[6] */
461 2, -1, +1, 0 }; /* z-edge [3]-[7] */
462
463 HYPRE_SStructGridSetFEMOrdering(edge_grid, part, ordering);
464 }
465
466 /* Set up the Graph - this determines the non-zero structure of the
467 matrix. */
468 {
469 int part = 0;
470
471 /* Create the graph object */
472 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, edge_grid, &A_graph);
473
474 /* Indicate that this problem uses finite element stiffness matrices and
475 load vectors, instead of stencils. */
476 HYPRE_SStructGraphSetFEM(A_graph, part);
477
478 /* The edge finite element matrix is full, so there is no need to call the
479 HYPRE_SStructGraphSetFEMSparsity() function. */
480
481 /* Assemble the graph */
482 HYPRE_SStructGraphAssemble(A_graph);
483 }
484
485 /* Set up the SStruct Matrix and right-hand side vector */
486 {
487 /* Create the matrix object */
488 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, A_graph, &A);
489 /* Use a ParCSR storage */
490 HYPRE_SStructMatrixSetObjectType(A, HYPRE_PARCSR);
491 /* Indicate that the matrix coefficients are ready to be set */
492 HYPRE_SStructMatrixInitialize(A);
493
494 /* Create an empty vector object */
495 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, edge_grid, &b);
496 /* Use a ParCSR storage */
497 HYPRE_SStructVectorSetObjectType(b, HYPRE_PARCSR);
498 /* Indicate that the vector coefficients are ready to be set */
499 HYPRE_SStructVectorInitialize(b);
500 }
501
502 /* Set the matrix and vector entries by finite element assembly */
503 {
504 /* local stiffness matrix and load vector */
505 double S[12][12], F[12];
506
507 int i, j, k;
508 int index[3];
509
510 for (i = 1; i <= n; i++)
511 for (j = 1; j <= n; j++)
512 for (k = 1; k <= n; k++)
513 {
514 /* Compute the FEM matrix and r.h.s. for cell (i,j,k) with
515 coefficients evaluated at the cell center. */
516 index[0] = i + pi*n; index[1] = j + pj*n; index[2] = k + pk*n;
517 ComputeFEMND1(S,F,(pi*n+i)*h-h/2,(pj*n+j)*h-h/2,(pk*n+k)*h-h/2,h);
518
519 /* Eliminate boundary conditions on x = 0 */
520 if (index[0] == 1)
521 {
522 int ii, jj, bc_edges[4] = { 3, 11, 7, 8 };
523 for (ii = 0; ii < 4; ii++)
524 {
525 for (jj = 0; jj < 12; jj++)
526 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
527 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
528 F[bc_edges[ii]] = 0.0;
529 }
530 }
531 /* Eliminate boundary conditions on y = 0 */
532 if (index[1] == 1)
533 {
534 int ii, jj, bc_edges[4] = { 0, 9, 4, 8 };
535 for (ii = 0; ii < 4; ii++)
536 {
537 for (jj = 0; jj < 12; jj++)
538 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
539 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
540 F[bc_edges[ii]] = 0.0;
541 }
542 }
543 /* Eliminate boundary conditions on z = 0 */
544 if (index[2] == 1)
545 {
546 int ii, jj, bc_edges[4] = { 0, 1, 2, 3 };
547 for (ii = 0; ii < 4; ii++)
548 {
549 for (jj = 0; jj < 12; jj++)
550 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
551 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
552 F[bc_edges[ii]] = 0.0;
553 }
554 }
555 /* Eliminate boundary conditions on x = 1 */
556 if (index[0] == N*n)
557 {
558 int ii, jj, bc_edges[4] = { 1, 10, 5, 9 };
559 for (ii = 0; ii < 4; ii++)
560 {
561 for (jj = 0; jj < 12; jj++)
562 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
563 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
564 F[bc_edges[ii]] = 0.0;
565 }
566 }
567 /* Eliminate boundary conditions on y = 1 */
568 if (index[1] == N*n)
569 {
570 int ii, jj, bc_edges[4] = { 2, 10, 6, 11 };
571 for (ii = 0; ii < 4; ii++)
572 {
573 for (jj = 0; jj < 12; jj++)
574 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
575 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
576 F[bc_edges[ii]] = 0.0;
577 }
578 }
579 /* Eliminate boundary conditions on z = 1 */
580 if (index[2] == N*n)
581 {
582 int ii, jj, bc_edges[4] = { 4, 5, 6, 7 };
583 for (ii = 0; ii < 4; ii++)
584 {
585 for (jj = 0; jj < 12; jj++)
586 S[bc_edges[ii]][jj] = S[jj][bc_edges[ii]] = 0.0;
587 S[bc_edges[ii]][bc_edges[ii]] = 1.0;
588 F[bc_edges[ii]] = 0.0;
589 }
590 }
591
592 /* Assemble the matrix */
593 HYPRE_SStructMatrixAddFEMValues(A, part, index, &S[0][0]);
594
595 /* Assemble the vector */
596 HYPRE_SStructVectorAddFEMValues(b, part, index, F);
597 }
598 }
599
600 /* Collective calls finalizing the matrix and vector assembly */
601 HYPRE_SStructMatrixAssemble(A);
602 HYPRE_SStructVectorAssemble(b);
603 }
604
605 /* 3. Create the discrete gradient matrix G, which is needed in AMS. */
606 {
607 int part = 0;
608 int stencil_size = 2;
609
610 /* Define the discretization stencil relating the edges and nodes of the
611 grid. */
612 {
613 int ndim = 3;
614 int entry;
615 int var = 0; /* the node variable */
616
617 /* The discrete gradient stencils connect edge to node variables. */
618 int Gx_offsets[2][3] = {{-1,0,0},{0,0,0}}; /* x-edge [7]-[6] */
619 int Gy_offsets[2][3] = {{0,-1,0},{0,0,0}}; /* y-edge [5]-[6] */
620 int Gz_offsets[2][3] = {{0,0,-1},{0,0,0}}; /* z-edge [2]-[6] */
621
622 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[0]);
623 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[1]);
624 HYPRE_SStructStencilCreate(ndim, stencil_size, &G_stencil[2]);
625
626 for (entry = 0; entry < stencil_size; entry++)
627 {
628 HYPRE_SStructStencilSetEntry(G_stencil[0], entry, Gx_offsets[entry], var);
629 HYPRE_SStructStencilSetEntry(G_stencil[1], entry, Gy_offsets[entry], var);
630 HYPRE_SStructStencilSetEntry(G_stencil[2], entry, Gz_offsets[entry], var);
631 }
632 }
633
634 /* Set up the Graph - this determines the non-zero structure of the
635 matrix. */
636 {
637 int nvars = 3;
638 int var; /* the edge variables */
639
640 /* Create the discrete gradient graph object */
641 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, edge_grid, &G_graph);
642
643 /* Since the discrete gradient relates edge and nodal variables (it is a
644 rectangular matrix), we have to specify the domain (column) grid. */
645 HYPRE_SStructGraphSetDomainGrid(G_graph, node_grid);
646
647 /* Tell the graph which stencil to use for each edge variable on each
648 part (we only have one part). */
649 for (var = 0; var < nvars; var++)
650 HYPRE_SStructGraphSetStencil(G_graph, part, var, G_stencil[var]);
651
652 /* Assemble the graph */
653 HYPRE_SStructGraphAssemble(G_graph);
654 }
655
656 /* Set up the SStruct Matrix */
657 {
658 /* Create the matrix object */
659 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, G_graph, &G);
660 /* Use a ParCSR storage */
661 HYPRE_SStructMatrixSetObjectType(G, HYPRE_PARCSR);
662 /* Indicate that the matrix coefficients are ready to be set */
663 HYPRE_SStructMatrixInitialize(G);
664 }
665
666 /* Set the discrete gradient values, assuming a "natural" orientation of
667 the edges (i.e. one in agreement with the coordinate directions). */
668 {
669 int i;
670 int nedges = n*(n+1)*(n+1);
671 double *values;
672 int stencil_indices[2] = {0,1}; /* the nodes of each edge */
673
674 values = calloc(2*nedges, sizeof(double));
675
676 /* The edge orientation is fixed: from first to second node */
677 for (i = 0; i < nedges; i++)
678 {
679 values[2*i] = -1.0;
680 values[2*i+1] = 1.0;
681 }
682
683 /* Set the values in the discrete gradient x-edges */
684 {
685 int var = 0;
686 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
687 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
688 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
689 stencil_size, stencil_indices,
690 values);
691 }
692 /* Set the values in the discrete gradient y-edges */
693 {
694 int var = 1;
695 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
696 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
697 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
698 stencil_size, stencil_indices,
699 values);
700 }
701 /* Set the values in the discrete gradient z-edges */
702 {
703 int var = 2;
704 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
705 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
706 HYPRE_SStructMatrixSetBoxValues(G, part, ilower, iupper, var,
707 stencil_size, stencil_indices,
708 values);
709 }
710
711 free(values);
712 }
713
714 /* Finalize the matrix assembly */
715 HYPRE_SStructMatrixAssemble(G);
716 }
717
718 /* 4. Create the vectors of nodal coordinates xcoord, ycoord and zcoord,
719 which are needed in AMS. */
720 {
721 int i, j, k;
722 int part = 0;
723 int var = 0; /* the node variable */
724 int index[3];
725 double xval, yval, zval;
726
727 /* Create empty vector objects */
728 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &xcoord);
729 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &ycoord);
730 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, node_grid, &zcoord);
731 /* Set the object type to ParCSR */
732 HYPRE_SStructVectorSetObjectType(xcoord, HYPRE_PARCSR);
733 HYPRE_SStructVectorSetObjectType(ycoord, HYPRE_PARCSR);
734 HYPRE_SStructVectorSetObjectType(zcoord, HYPRE_PARCSR);
735 /* Indicate that the vector coefficients are ready to be set */
736 HYPRE_SStructVectorInitialize(xcoord);
737 HYPRE_SStructVectorInitialize(ycoord);
738 HYPRE_SStructVectorInitialize(zcoord);
739
740 /* Compute and set the coordinates of the nodes */
741 for (i = 0; i <= n; i++)
742 for (j = 0; j <= n; j++)
743 for (k = 0; k <= n; k++)
744 {
745 index[0] = i + pi*n; index[1] = j + pj*n; index[2] = k + pk*n;
746
747 xval = index[0]*h;
748 yval = index[1]*h;
749 zval = index[2]*h;
750
751 HYPRE_SStructVectorSetValues(xcoord, part, index, var, &xval);
752 HYPRE_SStructVectorSetValues(ycoord, part, index, var, &yval);
753 HYPRE_SStructVectorSetValues(zcoord, part, index, var, &zval);
754 }
755
756 /* Finalize the vector assembly */
757 HYPRE_SStructVectorAssemble(xcoord);
758 HYPRE_SStructVectorAssemble(ycoord);
759 HYPRE_SStructVectorAssemble(zcoord);
760 }
761
762 /* 5. Set up a SStruct Vector for the solution vector x */
763 {
764 int part = 0;
765 int nvalues = n*(n+1)*(n+1);
766 double *values;
767
768 values = calloc(nvalues, sizeof(double));
769
770 /* Create an empty vector object */
771 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, edge_grid, &x);
772 /* Set the object type to ParCSR */
773 HYPRE_SStructVectorSetObjectType(x, HYPRE_PARCSR);
774 /* Indicate that the vector coefficients are ready to be set */
775 HYPRE_SStructVectorInitialize(x);
776
777 /* Set the values for the initial guess x-edge */
778 {
779 int var = 0;
780 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
781 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
782 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
783 }
784 /* Set the values for the initial guess y-edge */
785 {
786 int var = 1;
787 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
788 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
789 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
790 }
791 /* Set the values for the initial guess z-edge */
792 {
793 int var = 2;
794 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
795 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
796 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
797 }
798
799 free(values);
800
801 /* Finalize the vector assembly */
802 HYPRE_SStructVectorAssemble(x);
803 }
804
805 /* 6. Set up and call the PCG-AMS solver (Solver options can be found in the
806 Reference Manual.) */
807 {
808 double final_res_norm;
809 int its;
810
811 HYPRE_ParCSRMatrix par_A;
812 HYPRE_ParVector par_b;
813 HYPRE_ParVector par_x;
814
815 HYPRE_ParCSRMatrix par_G;
816 HYPRE_ParVector par_xcoord;
817 HYPRE_ParVector par_ycoord;
818 HYPRE_ParVector par_zcoord;
819
820 /* Extract the ParCSR objects needed in the solver */
821 HYPRE_SStructMatrixGetObject(A, (void **) &par_A);
822 HYPRE_SStructVectorGetObject(b, (void **) &par_b);
823 HYPRE_SStructVectorGetObject(x, (void **) &par_x);
824 HYPRE_SStructMatrixGetObject(G, (void **) &par_G);
825 HYPRE_SStructVectorGetObject(xcoord, (void **) &par_xcoord);
826 HYPRE_SStructVectorGetObject(ycoord, (void **) &par_ycoord);
827 HYPRE_SStructVectorGetObject(zcoord, (void **) &par_zcoord);
828
829 /* Create solver */
830 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
831
832 /* Set some parameters (See Reference Manual for more parameters) */
833 HYPRE_PCGSetMaxIter(solver, maxit); /* max iterations */
834 HYPRE_PCGSetTol(solver, tol); /* conv. tolerance */
835 HYPRE_PCGSetTwoNorm(solver, 0); /* use the two norm as the stopping criteria */
836 HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
837 HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
838
839 /* Create AMS preconditioner */
840 HYPRE_AMSCreate(&precond);
841
842 /* Set AMS parameters */
843 HYPRE_AMSSetMaxIter(precond, 1);
844 HYPRE_AMSSetTol(precond, 0.0);
845 HYPRE_AMSSetCycleType(precond, cycle_type);
846 HYPRE_AMSSetPrintLevel(precond, 1);
847
848 /* Set discrete gradient */
849 HYPRE_AMSSetDiscreteGradient(precond, par_G);
850
851 /* Set vertex coordinates */
852 HYPRE_AMSSetCoordinateVectors(precond,
853 par_xcoord, par_ycoord, par_zcoord);
854
855 if (singular_problem)
856 HYPRE_AMSSetBetaPoissonMatrix(precond, NULL);
857
858 /* Smoothing and AMG options */
859 HYPRE_AMSSetSmoothingOptions(precond,
860 rlx_type, rlx_sweeps,
861 rlx_weight, rlx_omega);
862 HYPRE_AMSSetAlphaAMGOptions(precond,
863 amg_coarsen_type, amg_agg_levels,
864 amg_rlx_type, theta, amg_interp_type,
865 amg_Pmax);
866 HYPRE_AMSSetBetaAMGOptions(precond,
867 amg_coarsen_type, amg_agg_levels,
868 amg_rlx_type, theta, amg_interp_type,
869 amg_Pmax);
870
871 /* Set the PCG preconditioner */
872 HYPRE_PCGSetPrecond(solver,
873 (HYPRE_PtrToSolverFcn) HYPRE_AMSSolve,
874 (HYPRE_PtrToSolverFcn) HYPRE_AMSSetup,
875 precond);
876
877 /* Call the setup */
878 HYPRE_ParCSRPCGSetup(solver, par_A, par_b, par_x);
879
880 /* Call the solve */
881 HYPRE_ParCSRPCGSolve(solver, par_A, par_b, par_x);
882
883 /* Get some info */
884 HYPRE_PCGGetNumIterations(solver, &its);
885 HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
886
887 /* Clean up */
888 HYPRE_AMSDestroy(precond);
889 HYPRE_ParCSRPCGDestroy(solver);
890
891 /* Gather the solution vector */
892 HYPRE_SStructVectorGather(x);
893
894 /* Print the solution with replicated shared data */
895 if (print_solution)
896 {
897 FILE *file;
898 char filename[255];
899
900 int part = 0;
901 int nvalues = n*(n+1)*(n+1);
902 double *xvalues, *yvalues, *zvalues;
903
904 xvalues = calloc(nvalues, sizeof(double));
905 yvalues = calloc(nvalues, sizeof(double));
906 zvalues = calloc(nvalues, sizeof(double));
907
908 /* Get local solution in the x-edges */
909 {
910 int var = 0;
911 int ilower[3] = {1 + pi*n, 0 + pj*n, 0 + pk*n};
912 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
913 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
914 var, xvalues);
915 }
916 /* Get local solution in the y-edges */
917 {
918 int var = 1;
919 int ilower[3] = {0 + pi*n, 1 + pj*n, 0 + pk*n};
920 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
921 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
922 var, yvalues);
923 }
924 /* Get local solution in the z-edges */
925 {
926 int var = 2;
927 int ilower[3] = {0 + pi*n, 0 + pj*n, 1 + pk*n};
928 int iupper[3] = {n + pi*n, n + pj*n, n + pk*n};
929 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
930 var, zvalues);
931 }
932
933 sprintf(filename, "sstruct.out.x.00.00.%05d", myid);
934 if ((file = fopen(filename, "w")) == NULL)
935 {
936 printf("Error: can't open output file %s\n", filename);
937 MPI_Finalize();
938 exit(1);
939 }
940
941 /* Save the edge values, element by element, using the same numbering
942 as the local finite element degrees of freedom. */
943 {
944 int i, j, k, s;
945
946 /* Initial x-, y- and z-edge indices in the values arrays */
947 int oi[4] = { 0, n, n*(n+1), n*(n+1)+n }; /* e_0, e_2, e_4, e_6 */
948 int oj[4] = { 0, 1, n*(n+1), n*(n+1)+1 }; /* e_3, e_1, e_7, e_5 */
949 int ok[4] = { 0, 1, n+1, n+2 }; /* e_8, e_9, e_11, e_10 */
950 /* Loop over the cells while updating the above offsets */
951 for (k = 0; k < n; k++)
952 {
953 for (j = 0; j < n; j++)
954 {
955 for (i = 0; i < n; i++)
956 {
957 fprintf(file,
958 "%.14e\n%.14e\n%.14e\n%.14e\n"
959 "%.14e\n%.14e\n%.14e\n%.14e\n"
960 "%.14e\n%.14e\n%.14e\n%.14e\n",
961 xvalues[oi[0]], yvalues[oj[1]], xvalues[oi[1]], yvalues[oj[0]],
962 xvalues[oi[2]], yvalues[oj[3]], xvalues[oi[3]], yvalues[oj[2]],
963 zvalues[ok[0]], zvalues[ok[1]], zvalues[ok[3]], zvalues[ok[2]]);
964
965 for (s=0; s<4; s++) oi[s]++, oj[s]++, ok[s]++;
966 }
967 for (s=0; s<4; s++) oj[s]++, ok[s]++;
968 }
969 for (s=0; s<4; s++) oi[s]+=n, ok[s]+=n+1;
970 }
971 }
972
973 fflush(file);
974 fclose(file);
975
976 free(xvalues);
977 free(yvalues);
978 free(zvalues);
979 }
980
981 if (myid == 0)
982 {
983 printf("\n");
984 printf("Iterations = %d\n", its);
985 printf("Final Relative Residual Norm = %g\n", final_res_norm);
986 printf("\n");
987 }
988 }
989
990 /* Free memory */
991 HYPRE_SStructGridDestroy(edge_grid);
992 HYPRE_SStructGraphDestroy(A_graph);
993 HYPRE_SStructMatrixDestroy(A);
994 HYPRE_SStructVectorDestroy(b);
995 HYPRE_SStructVectorDestroy(x);
996 HYPRE_SStructGridDestroy(node_grid);
997 HYPRE_SStructGraphDestroy(G_graph);
998 HYPRE_SStructStencilDestroy(G_stencil[0]);
999 HYPRE_SStructStencilDestroy(G_stencil[1]);
1000 HYPRE_SStructStencilDestroy(G_stencil[2]);
1001 HYPRE_SStructMatrixDestroy(G);
1002 HYPRE_SStructVectorDestroy(xcoord);
1003 HYPRE_SStructVectorDestroy(ycoord);
1004 HYPRE_SStructVectorDestroy(zcoord);
1005
1006 /* Finalize MPI */
1007 MPI_Finalize();
1008
1009 return 0;
1010 }
syntax highlighted by Code2HTML, v. 0.9.1