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