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


syntax highlighted by Code2HTML, v. 0.9.1