download the original source code.
1 /*
2 Example 13
3
4 Interface: Semi-Structured interface (SStruct)
5
6 Compile with: make ex13
7
8 Sample run: mpirun -np 6 ex13 -n 10
9
10 To see options: ex13 -help
11
12 Description: This code solves the 2D Laplace equation using bilinear
13 finite element discretization on a mesh with an "enhanced
14 connectivity" point. Specifically, we solve -Delta u = 1
15 with zero boundary conditions on a star-shaped domain
16 consisting of identical rhombic parts each meshed with a
17 uniform n x n grid. Every part is assigned to a different
18 processor and all parts meet at the origin, equally
19 subdividing the 2*pi angle there. The case of six processors
20 (parts) looks as follows:
21
22 +
23 / \
24 / \
25 / \
26 +--------+ 1 +---------+
27 \ \ / /
28 \ 2 \ / 0 /
29 \ \ / /
30 +--------+---------+
31 / / \ \
32 / 3 / \ 5 \
33 / / \ \
34 +--------+ 4 +---------+
35 \ /
36 \ /
37 \ /
38 +
39
40 Note that in this problem we use nodal variables, which are
41 shared between the different parts. The node at the origin,
42 for example, belongs to all parts as illustrated below:
43
44 .
45 / \
46 . .
47 / \ / \
48 o . *
49 .---.---o \ / \ / *---.---.
50 \ \ \ o * / / /
51 .---.---o \ / *---.---.
52 \ \ \ x / / /
53 @---@---x x---z---z
54 @---@---x x---z---z
55 / / / x \ \ \
56 .---.---a / \ #---.---.
57 / / / a # \ \ \
58 .---.---a / \ / \ #---.---.
59 a . #
60 \ / \ /
61 . .
62 \ /
63 .
64
65 We recommend viewing the Struct examples before viewing this
66 and the other SStruct examples. The primary role of this
67 particular SStruct example is to demonstrate a stencil-based
68 way to set up finite element problems in SStruct, and
69 specifically to show how to handle problems with an "enhanced
70 connectivity" point.
71 */
72
73 #include <math.h>
74 #include "_hypre_utilities.h"
75 #include "HYPRE_sstruct_mv.h"
76 #include "HYPRE_sstruct_ls.h"
77 #include "HYPRE.h"
78
79 #ifndef M_PI
80 #define M_PI 3.14159265358979
81 #endif
82
83
84 /*
85 This routine computes the bilinear finite element stiffness matrix and
86 load vector on a rhombus with angle gamma. Specifically, let R be the
87 rhombus
88 [3]------[2]
89 / /
90 / /
91 [0]------[1]
92
93 with sides of length h. The finite element stiffness matrix
94
95 S_ij = (grad phi_i,grad phi_j)_R
96
97 with bilinear finite element functions {phi_i} has the form
98
99 / 4-k -1 -2+k -1 \
100 alpha . | -1 4+k -1 -2-k |
101 | -2+k -1 4-k -1 |
102 \ -1 -2-k -1 4+k /
103
104 where alpha = 1/(6*sin(gamma)) and k = 3*cos(gamma). The load vector
105 corresponding to a right-hand side of 1 is
106
107 F_j = (1,phi_j)_R = h^2/4 * sin(gamma)
108 */
109 void ComputeFEMRhombus (double S[4][4], double F[4], double gamma, double h)
110 {
111 int i, j;
112
113 double h2_4 = h*h/4;
114 double sing = sin(gamma);
115 double alpha = 1/(6*sing);
116 double k = 3*cos(gamma);
117
118 S[0][0] = alpha * (4-k);
119 S[0][1] = alpha * (-1);
120 S[0][2] = alpha * (-2+k);
121 S[0][3] = alpha * (-1);
122 S[1][1] = alpha * (4+k);
123 S[1][2] = alpha * (-1);
124 S[1][3] = alpha * (-2-k);
125 S[2][2] = alpha * (4-k);
126 S[2][3] = alpha * (-1);
127 S[3][3] = alpha * (4+k);
128
129 /* The stiffness matrix is symmetric */
130 for (i = 1; i < 4; i++)
131 for (j = 0; j < i; j++)
132 S[i][j] = S[j][i];
133
134 for (i = 0; i < 4; i++)
135 F[i] = h2_4*sing;
136 }
137
138
139 int main (int argc, char *argv[])
140 {
141 int myid, num_procs;
142 int n;
143 double gamma, h;
144 int print_solution;
145
146 HYPRE_SStructGrid grid;
147 HYPRE_SStructGraph graph;
148 HYPRE_SStructStencil stencil;
149 HYPRE_SStructMatrix A;
150 HYPRE_SStructVector b;
151 HYPRE_SStructVector x;
152
153 HYPRE_Solver solver;
154
155 /* Initialize MPI */
156 MPI_Init(&argc, &argv);
157 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
158 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
159
160 /* Set default parameters */
161 n = 10;
162 print_solution = 0;
163
164 /* Parse command line */
165 {
166 int arg_index = 0;
167 int print_usage = 0;
168
169 while (arg_index < argc)
170 {
171 if ( strcmp(argv[arg_index], "-n") == 0 )
172 {
173 arg_index++;
174 n = atoi(argv[arg_index++]);
175 }
176 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
177 {
178 arg_index++;
179 print_solution = 1;
180 }
181 else if ( strcmp(argv[arg_index], "-help") == 0 )
182 {
183 print_usage = 1;
184 break;
185 }
186 else
187 {
188 arg_index++;
189 }
190 }
191
192 if ((print_usage) && (myid == 0))
193 {
194 printf("\n");
195 printf("Usage: %s [<options>]\n", argv[0]);
196 printf("\n");
197 printf(" -n <n> : problem size per processor (default: 10)\n");
198 printf(" -print_solution : print the solution vector\n");
199 printf("\n");
200 }
201
202 if (print_usage)
203 {
204 MPI_Finalize();
205 return (0);
206 }
207 }
208
209 /* Set the rhombus angle, gamma, and the mesh size, h, depending on the
210 number of processors np and the given n */
211 if (num_procs < 3)
212 {
213 if (myid ==0) printf("Must run with at least 3 processors!\n");
214 MPI_Finalize();
215 exit(1);
216 }
217 gamma = 2*M_PI/num_procs;
218 h = 1.0/n;
219
220 /* 1. Set up the grid. We will set up the grid so that processor X owns
221 part X. Note that each part has its own index space numbering. Later
222 we relate the parts to each other. */
223 {
224 int ndim = 2;
225 int nparts = num_procs;
226
227 /* Create an empty 2D grid object */
228 HYPRE_SStructGridCreate(MPI_COMM_WORLD, ndim, nparts, &grid);
229
230 /* Set the extents of the grid - each processor sets its grid boxes. Each
231 part has its own relative index space numbering */
232 {
233 int part = myid;
234 int ilower[2] = {1,1}; /* lower-left cell touching the origin */
235 int iupper[2] = {n,n}; /* upper-right cell */
236
237 HYPRE_SStructGridSetExtents(grid, part, ilower, iupper);
238 }
239
240 /* Set the variable type and number of variables on each part. These need
241 to be set in each part which is neighboring or contains boxes owned by
242 the processor. */
243 {
244 int i;
245 int nvars = 1;
246
247 HYPRE_SStructVariable vartypes[1] = {HYPRE_SSTRUCT_VARIABLE_NODE};
248 for (i = 0; i < nparts; i++)
249 HYPRE_SStructGridSetVariables(grid, i, nvars, vartypes);
250 }
251
252 /* Now we need to set the spatial relation between each of the parts.
253 Since we are using nodal variables, we have to use SetSharedPart to
254 establish the connection at the origin. */
255 {
256 /* Relation to the clockwise-previous neighbor part, e.g. 0 and 1 for
257 the case of 6 parts. Note that we could have used SetNeighborPart
258 here instead of SetSharedPart. */
259 {
260 int part = myid;
261 /* the box of cells intersecting the boundary in the current part */
262 int ilower[2] = {1,1}, iupper[2] = {1,n};
263 /* share all data on the left side of the box */
264 int offset[2] = {-1,0};
265
266 int shared_part = (myid+1) % num_procs;
267 /* the box of cells intersecting the boundary in the neighbor */
268 int shared_ilower[2] = {1,1}, shared_iupper[2] = {n,1};
269 /* share all data on the bottom of the box */
270 int shared_offset[2] = {0,-1};
271
272 /* x/y-direction on the current part is -y/x on the neighbor */
273 int index_map[2] = {1,0};
274 int index_dir[2] = {-1,1};
275
276 HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
277 shared_part, shared_ilower,
278 shared_iupper, shared_offset,
279 index_map, index_dir);
280 }
281
282 /* Relation to the clockwise-following neighbor part, e.g. 0 and 5 for
283 the case of 6 parts. Note that we could have used SetNeighborPart
284 here instead of SetSharedPart. */
285 {
286 int part = myid;
287 /* the box of cells intersecting the boundary in the current part */
288 int ilower[2] = {1,1}, iupper[2] = {n,1};
289 /* share all data on the bottom of the box */
290 int offset[2] = {0,-1};
291
292 int shared_part = (myid+num_procs-1) % num_procs;
293 /* the box of cells intersecting the boundary in the neighbor */
294 int shared_ilower[2] = {1,1}, shared_iupper[2] = {1,n};
295 /* share all data on the left side of the box */
296 int shared_offset[2] = {-1,0};
297
298 /* x/y-direction on the current part is y/-x on the neighbor */
299 int index_map[2] = {1,0};
300 int index_dir[2] = {1,-1};
301
302 HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
303 shared_part, shared_ilower,
304 shared_iupper, shared_offset,
305 index_map, index_dir);
306 }
307
308 /* Relation to all other parts, e.g. 0 and 2,3,4. This can be
309 described only by SetSharedPart. */
310 {
311 int part = myid;
312 /* the (one cell) box that touches the origin */
313 int ilower[2] = {1,1}, iupper[2] = {1,1};
314 /* share all data in the bottom left corner (i.e. the origin) */
315 int offset[2] = {-1,-1};
316
317 int shared_part;
318 /* the box of one cell that touches the origin */
319 int shared_ilower[2] = {1,1}, shared_iupper[2] = {1,1};
320 /* share all data in the bottom left corner (i.e. the origin) */
321 int shared_offset[2] = {-1,-1};
322
323 /* x/y-direction on the current part is -x/-y on the neighbor, but
324 in this case the arguments are not really important since we are
325 only sharing a point */
326 int index_map[2] = {0,1};
327 int index_dir[2] = {-1,-1};
328
329 for (shared_part = 0; shared_part < myid-1; shared_part++)
330 HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
331 shared_part, shared_ilower,
332 shared_iupper, shared_offset,
333 index_map, index_dir);
334
335 for (shared_part = myid+2; shared_part < num_procs; shared_part++)
336 HYPRE_SStructGridSetSharedPart(grid, part, ilower, iupper, offset,
337 shared_part, shared_ilower,
338 shared_iupper, shared_offset,
339 index_map, index_dir);
340 }
341 }
342
343 /* Now the grid is ready to be used */
344 HYPRE_SStructGridAssemble(grid);
345 }
346
347 /* 2. Define the discretization stencils. Since this is a finite element
348 discretization we define here a full 9-point stencil. We will later
349 use four sub-stencils for the rows of the local stiffness matrix. */
350 {
351 int ndim = 2;
352 int var = 0;
353 int entry;
354
355 /* Define the geometry of the 9-point stencil */
356 int stencil_size = 9;
357 int offsets[9][2] = {{0,0}, /* [8] [4] [7] */
358 {-1,0}, {1,0}, /* \ | / */
359 {0,-1}, {0,1}, /* [1]-[0]-[2] */
360 {-1,-1}, {1,-1}, /* / | \ */
361 {1,1}, {-1,1}}; /* [5] [3] [6] */
362
363 HYPRE_SStructStencilCreate(ndim, stencil_size, &stencil);
364
365 for (entry = 0; entry < stencil_size; entry++)
366 HYPRE_SStructStencilSetEntry(stencil, entry, offsets[entry], var);
367 }
368
369 /* 3. Set up the Graph - this determines the non-zero structure of the
370 matrix. */
371 {
372 int part;
373 int var = 0;
374
375 /* Create the graph object */
376 HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, &graph);
377
378 /* Now we need to tell the graph which stencil to use for each
379 variable on each part (we only have one variable) */
380 for (part = 0; part < num_procs; part++)
381 HYPRE_SStructGraphSetStencil(graph, part, var, stencil);
382
383 /* Assemble the graph */
384 HYPRE_SStructGraphAssemble(graph);
385 }
386
387 /* 4. Set up the SStruct Matrix and right-hand side vector */
388 {
389 int part = myid;
390 int var = 0;
391
392 /* Create the matrix object */
393 HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, &A);
394 /* Use a ParCSR storage */
395 HYPRE_SStructMatrixSetObjectType(A, HYPRE_PARCSR);
396 /* Indicate that the matrix coefficients are ready to be set */
397 HYPRE_SStructMatrixInitialize(A);
398
399 /* Create an empty vector object */
400 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &b);
401 /* Use a ParCSR storage */
402 HYPRE_SStructVectorSetObjectType(b, HYPRE_PARCSR);
403 /* Indicate that the vector coefficients are ready to be set */
404 HYPRE_SStructVectorInitialize(b);
405
406 /* Set the matrix and vector entries by finite element assembly */
407 {
408 /* local stifness matrix and load vector */
409 double S[4][4], F[4];
410
411 /* The index of the local nodes 0-3 relative to the cell index,
412 i.e. node k in cell (i,j) is in the upper-right corner of the
413 cell (i,j) + node_index_offset[k]. */
414 int node_index_offset[4][2] = {{-1,-1},{0,-1},{0,0},{-1,0}};
415
416 /* The cell sub-stencils of nodes 0-3 indexed from the full stencil,
417 i.e. we take the full stencil in each node of a fixed cell, and
418 restrict it to that as is done in the finite element stiffness
419 matrix:
420 [4] [7] [8] [4] [1]-[0] [0]-[2]
421 | / \ | / | | \
422 [0]-[2] , [1]-[0] , [5] [3] , [3] [6]
423
424 Note that the ordering of the local nodes remains fixed, and
425 therefore the above sub-stencil at node k corresponds to the kth row
426 of the local stiffness matrix and the kth entry of the local load
427 vector. */
428 int node_stencil[4][4] = {{0,2,7,4},{1,0,4,8},{5,3,0,1},{3,6,2,0}};
429
430 int i, j, k;
431 int index[2];
432 int nentries = 4;
433
434 /* set the values in the interior cells */
435 {
436 ComputeFEMRhombus(S, F, gamma, h);
437
438 for (i = 1; i <= n; i++)
439 for (j = 1; j <= n; j++)
440 for (k = 0; k < 4; k++) /* node k in cell (i,j) */
441 {
442 index[0] = i + node_index_offset[k][0];
443 index[1] = j + node_index_offset[k][1];
444 HYPRE_SStructMatrixAddToValues(A, part, index, var,
445 nentries, node_stencil[k],
446 &S[k][0]);
447 HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
448 }
449 }
450
451 /* cells having nodes 1,2 on the domain boundary */
452 {
453 ComputeFEMRhombus(S, F, gamma, h);
454
455 /* eliminate nodes 1,2 from S and F */
456 for (k = 0; k < 4; k++)
457 {
458 S[1][k] = S[k][1] = 0.0;
459 S[2][k] = S[k][2] = 0.0;
460 }
461 S[1][1] = 1.0;
462 S[2][2] = 1.0;
463 F[1] = 0.0;
464 F[2] = 0.0;
465
466 for (i = n; i <= n; i++)
467 for (j = 1; j <= n; j++)
468 for (k = 0; k < 4; k++) /* node k in cell (n,j) */
469 {
470 index[0] = i + node_index_offset[k][0];
471 index[1] = j + node_index_offset[k][1];
472 HYPRE_SStructMatrixAddToValues(A, part, index, var,
473 nentries, node_stencil[k],
474 &S[k][0]);
475 HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
476 }
477 }
478
479 /* cells having nodes 2,3 on the domain boundary */
480 {
481 ComputeFEMRhombus(S, F, gamma, h);
482
483 /* eliminate nodes 2,3 from S and F */
484 for (k = 0; k < 4; k++)
485 {
486 S[2][k] = S[k][2] = 0.0;
487 S[3][k] = S[k][3] = 0.0;
488 }
489 S[2][2] = 1.0;
490 S[3][3] = 1.0;
491 F[2] = 0.0;
492 F[3] = 0.0;
493
494 for (i = 1; i <= n; i++)
495 for (j = n; j <= n; j++)
496 for (k = 0; k < 4; k++) /* node k in cell (i,n) */
497 {
498 index[0] = i + node_index_offset[k][0];
499 index[1] = j + node_index_offset[k][1];
500 HYPRE_SStructMatrixAddToValues(A, part, index, var,
501 nentries, node_stencil[k],
502 &S[k][0]);
503 HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
504 }
505 }
506
507 /* cells having nodes 1,2,3 on the domain boundary */
508 {
509 ComputeFEMRhombus(S, F, gamma, h);
510
511 /* eliminate nodes 2,3 from S and F */
512 for (k = 0; k < 4; k++)
513 {
514 S[1][k] = S[k][1] = 0.0;
515 S[2][k] = S[k][2] = 0.0;
516 S[3][k] = S[k][3] = 0.0;
517 }
518 S[1][1] = 1.0;
519 S[2][2] = 1.0;
520 S[3][3] = 1.0;
521 F[1] = 0.0;
522 F[2] = 0.0;
523 F[3] = 0.0;
524
525 for (i = n; i <= n; i++)
526 for (j = n; j <= n; j++)
527 for (k = 0; k < 4; k++) /* node k in cell (n,n) */
528 {
529 index[0] = i + node_index_offset[k][0];
530 index[1] = j + node_index_offset[k][1];
531 HYPRE_SStructMatrixAddToValues(A, part, index, var,
532 nentries, node_stencil[k],
533 &S[k][0]);
534 HYPRE_SStructVectorAddToValues(b, part, index, var, &F[k]);
535 }
536 }
537 }
538 }
539
540 /* Collective calls finalizing the matrix and vector assembly */
541 HYPRE_SStructMatrixAssemble(A);
542 HYPRE_SStructVectorAssemble(b);
543
544 /* 5. Set up SStruct Vector for the solution vector x */
545 {
546 int part = myid;
547 int var = 0;
548 int nvalues = (n+1)*(n+1);
549 double *values;
550
551 /* Since the SetBoxValues() calls below set the values of the nodes in
552 the upper-right corners of the cells, the nodal box should start
553 from (0,0) instead of (1,1). */
554 int ilower[2] = {0,0};
555 int iupper[2] = {n,n};
556
557 values = calloc(nvalues, sizeof(double));
558
559 /* Create an empty vector object */
560 HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, &x);
561 /* Set the object type to ParCSR */
562 HYPRE_SStructVectorSetObjectType(x, HYPRE_PARCSR);
563 /* Indicate that the vector coefficients are ready to be set */
564 HYPRE_SStructVectorInitialize(x);
565 /* Set the values for the initial guess */
566 HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper, var, values);
567
568 free(values);
569
570 /* Finalize the vector assembly */
571 HYPRE_SStructVectorAssemble(x);
572 }
573
574 /* 6. Set up and call the solver (Solver options can be found in the
575 Reference Manual.) */
576 {
577 double final_res_norm;
578 int its;
579
580 HYPRE_ParCSRMatrix par_A;
581 HYPRE_ParVector par_b;
582 HYPRE_ParVector par_x;
583
584 /* Extract the ParCSR objects needed in the solver */
585 HYPRE_SStructMatrixGetObject(A, (void **) &par_A);
586 HYPRE_SStructVectorGetObject(b, (void **) &par_b);
587 HYPRE_SStructVectorGetObject(x, (void **) &par_x);
588
589 /* Here we construct a BoomerAMG solver. See the other SStruct examples
590 as well as the Reference manual for additional solver choices. */
591 HYPRE_BoomerAMGCreate(&solver);
592 HYPRE_BoomerAMGSetCoarsenType(solver, 6);
593 HYPRE_BoomerAMGSetStrongThreshold(solver, 0.25);
594 HYPRE_BoomerAMGSetTol(solver, 1e-6);
595 HYPRE_BoomerAMGSetPrintLevel(solver, 2);
596 HYPRE_BoomerAMGSetMaxIter(solver, 50);
597
598 /* call the setup */
599 HYPRE_BoomerAMGSetup(solver, par_A, par_b, par_x);
600
601 /* call the solve */
602 HYPRE_BoomerAMGSolve(solver, par_A, par_b, par_x);
603
604 /* get some info */
605 HYPRE_BoomerAMGGetNumIterations(solver, &its);
606 HYPRE_BoomerAMGGetFinalRelativeResidualNorm(solver,
607 &final_res_norm);
608 /* clean up */
609 HYPRE_BoomerAMGDestroy(solver);
610
611 /* Gather the solution vector */
612 HYPRE_SStructVectorGather(x);
613
614 /* Print the solution with replicated shared data */
615 if (print_solution)
616 {
617 FILE *file;
618 char filename[255];
619
620 int i, part = myid, var = 0;
621 int nvalues = (n+1)*(n+1);
622 double *values;
623 int ilower[2] = {0,0};
624 int iupper[2] = {n,n};
625
626 values = calloc(nvalues, sizeof(double));
627
628 /* get all local data (including a local copy of the shared values) */
629 HYPRE_SStructVectorGetBoxValues(x, part, ilower, iupper,
630 var, values);
631
632 sprintf(filename, "sstruct.out.x.00.00.%05d", part);
633 if ((file = fopen(filename, "w")) == NULL)
634 {
635 printf("Error: can't open output file %s\n", filename);
636 MPI_Finalize();
637 exit(1);
638 }
639 for (i = 0; i < nvalues; i++)
640 fprintf(file, "%.14e\n", values[i]);
641 fflush(file);
642 fclose(file);
643
644 free(values);
645 }
646
647 if (myid == 0)
648 {
649 printf("\n");
650 printf("Iterations = %d\n", its);
651 printf("Final Relative Residual Norm = %g\n", final_res_norm);
652 printf("\n");
653 }
654 }
655
656 /* Free memory */
657 HYPRE_SStructGridDestroy(grid);
658 HYPRE_SStructStencilDestroy(stencil);
659 HYPRE_SStructGraphDestroy(graph);
660 HYPRE_SStructMatrixDestroy(A);
661 HYPRE_SStructVectorDestroy(b);
662 HYPRE_SStructVectorDestroy(x);
663
664 /* Finalize MPI */
665 MPI_Finalize();
666
667 return 0;
668 }
syntax highlighted by Code2HTML, v. 0.9.1