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


syntax highlighted by Code2HTML, v. 0.9.1