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


syntax highlighted by Code2HTML, v. 0.9.1