download the original source code.
   1 /*
   2    Example 4
   3 
   4    Interface:      Structured interface (Struct)
   5 
   6    Compile with:   make ex4
   7 
   8    Sample run:     mpirun -np 16 ex4 -n 33 -solver 10 -K 3 -B 0 -C 1 -U0 2 -F 4
   9 
  10    To see options: ex4 -help
  11 
  12    Description:    This example differs from the previous structured example
  13                    (Example 3) in that a more sophisticated stencil and
  14                    boundary conditions are implemented. The method illustrated
  15                    here to implement the boundary conditions is much more general
  16                    than that in the previous example.  Also symmetric storage is
  17                    utilized when applicable.
  18 
  19                    This code solves the convection-reaction-diffusion problem
  20                    div (-K grad u + B u) + C u = F in the unit square with
  21                    boundary condition u = U0.  The domain is split into N x N
  22                    processor grid.  Thus, the given number of processors should
  23                    be a perfect square. Each processor has a n x n grid, with
  24                    nodes connected by a 5-point stencil. Note that the struct
  25                    interface assumes a cell-centered grid, and, therefore, the
  26                    nodes are not shared.
  27 
  28                    To incorporate the boundary conditions, we do the following:
  29                    Let x_i and x_b be the interior and boundary parts of the
  30                    solution vector x. If we split the matrix A as
  31                              A = [A_ii A_ib; A_bi A_bb],
  32                    then we solve
  33                              [A_ii 0; 0 I] [x_i ; x_b] = [b_i - A_ib u_0; u_0].
  34                    Note that this differs from the previous example in that we
  35                    are actually solving for the boundary conditions (so they
  36                    may not be exact as in ex3, where we only solved for the
  37                    interior).  This approach is useful for more general types
  38                    of b.c.
  39 
  40                    A number of solvers are available. More information can be
  41                    found in the Solvers and Preconditioners chapter of the
  42                    User's Manual.
  43 
  44                    We recommend viewing examples 1, 2, and 3 before viewing this
  45                    example.
  46 */
  47 
  48 #include <math.h>
  49 #include "_hypre_utilities.h"
  50 #include "HYPRE_krylov.h"
  51 #include "HYPRE_struct_ls.h"
  52 
  53 #ifdef M_PI
  54   #define PI M_PI
  55 #else
  56   #define PI 3.14159265358979
  57 #endif
  58 
  59 #include "vis.c"
  60 
  61 /* Macro to evaluate a function F in the grid point (i,j) */
  62 #define Eval(F,i,j) (F( (ilower[0]+(i))*h, (ilower[1]+(j))*h ))
  63 #define bcEval(F,i,j) (F( (bc_ilower[0]+(i))*h, (bc_ilower[1]+(j))*h ))
  64 
  65 int optionK, optionB, optionC, optionU0, optionF;
  66 
  67 /* Diffusion coefficient */
  68 double K(double x, double y)
  69 {
  70    switch (optionK)
  71    {
  72       case 0:
  73          return 1.0;
  74       case 1:
  75          return x*x+exp(y);
  76       case 2:
  77          if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25))
  78             return 100.0;
  79          else
  80             return 1.0;
  81       case 3:
  82          if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)) < 0.0625)
  83             return 10.0;
  84          else
  85             return 1.0;
  86       default:
  87          return 1.0;
  88    }
  89 }
  90 
  91 /* Convection vector, first component */
  92 double B1(double x, double y)
  93 {
  94    switch (optionB)
  95    {
  96       case 0:
  97          return 0.0;
  98       case 1:
  99          return -0.1;
 100       case 2:
 101          return 0.25;
 102       case 3:
 103          return 1.0;
 104       default:
 105          return 0.0;
 106    }
 107 }
 108 
 109 /* Convection vector, second component */
 110 double B2(double x, double y)
 111 {
 112    switch (optionB)
 113    {
 114       case 0:
 115          return 0.0;
 116       case 1:
 117          return 0.1;
 118       case 2:
 119          return -0.25;
 120       case 3:
 121          return 1.0;
 122       default:
 123          return 0.0;
 124    }
 125 }
 126 
 127 /* Reaction coefficient */
 128 double C(double x, double y)
 129 {
 130    switch (optionC)
 131    {
 132       case 0:
 133          return 0.0;
 134       case 1:
 135          return 10.0;
 136       case 2:
 137          return 100.0;
 138       default:
 139          return 0.0;
 140    }
 141 }
 142 
 143 /* Boundary condition */
 144 double U0(double x, double y)
 145 {
 146    switch (optionU0)
 147    {
 148       case 0:
 149          return 0.0;
 150       case 1:
 151          return (x+y)/100;
 152       case 2:
 153          return (sin(5*PI*x)+sin(5*PI*y))/1000;
 154       default:
 155          return 0.0;
 156    }
 157 }
 158 
 159 /* Right-hand side */
 160 double F(double x, double y)
 161 {
 162    switch (optionF)
 163    {
 164       case 0:
 165          return 1.0;
 166       case 1:
 167          return 0.0;
 168       case 2:
 169          return 2*PI*PI*sin(PI*x)*sin(PI*y);
 170       case 3:
 171          if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25))
 172             return -1.0;
 173          else
 174             return 1.0;
 175       case 4:
 176          if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)) < 0.0625)
 177             return -1.0;
 178          else
 179             return 1.0;
 180       default:
 181          return 1.0;
 182    }
 183 }
 184 
 185 int main (int argc, char *argv[])
 186 {
 187    int i, j, k;
 188 
 189    int myid, num_procs;
 190 
 191    int n, N, pi, pj;
 192    double h, h2;
 193    int ilower[2], iupper[2];
 194 
 195    int solver_id;
 196    int n_pre, n_post;
 197    int rap, relax, skip, sym;
 198    int time_index;
 199 
 200    int num_iterations;
 201    double final_res_norm;
 202 
 203    int vis;
 204 
 205    HYPRE_StructGrid     grid;
 206    HYPRE_StructStencil  stencil;
 207    HYPRE_StructMatrix   A;
 208    HYPRE_StructVector   b;
 209    HYPRE_StructVector   x;
 210    HYPRE_StructSolver   solver;
 211    HYPRE_StructSolver   precond;
 212 
 213    /* Initialize MPI */
 214    MPI_Init(&argc, &argv);
 215    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 216    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 217 
 218    /* Set default parameters */
 219    n         = 33;
 220    optionK   = 0;
 221    optionB   = 0;
 222    optionC   = 0;
 223    optionU0  = 0;
 224    optionF   = 0;
 225    solver_id = 10;
 226    n_pre     = 1;
 227    n_post    = 1;
 228    rap       = 0;
 229    relax     = 1;
 230    skip      = 0;
 231    sym       = 0;
 232 
 233    vis       = 0;
 234 
 235    /* Parse command line */
 236    {
 237       int arg_index = 0;
 238       int print_usage = 0;
 239 
 240       while (arg_index < argc)
 241       {
 242          if ( strcmp(argv[arg_index], "-n") == 0 )
 243          {
 244             arg_index++;
 245             n = atoi(argv[arg_index++]);
 246          }
 247          else if ( strcmp(argv[arg_index], "-K") == 0 )
 248          {
 249             arg_index++;
 250             optionK = atoi(argv[arg_index++]);
 251          }
 252          else if ( strcmp(argv[arg_index], "-B") == 0 )
 253          {
 254             arg_index++;
 255             optionB = atoi(argv[arg_index++]);
 256          }
 257          else if ( strcmp(argv[arg_index], "-C") == 0 )
 258          {
 259             arg_index++;
 260             optionC = atoi(argv[arg_index++]);
 261          }
 262          else if ( strcmp(argv[arg_index], "-U0") == 0 )
 263          {
 264             arg_index++;
 265             optionU0 = atoi(argv[arg_index++]);
 266          }
 267          else if ( strcmp(argv[arg_index], "-F") == 0 )
 268          {
 269             arg_index++;
 270             optionF = atoi(argv[arg_index++]);
 271          }
 272          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 273          {
 274             arg_index++;
 275             solver_id = atoi(argv[arg_index++]);
 276          }
 277          else if ( strcmp(argv[arg_index], "-v") == 0 )
 278          {
 279             arg_index++;
 280             n_pre = atoi(argv[arg_index++]);
 281             n_post = atoi(argv[arg_index++]);
 282          }
 283          else if ( strcmp(argv[arg_index], "-rap") == 0 )
 284          {
 285             arg_index++;
 286             rap = atoi(argv[arg_index++]);
 287          }
 288          else if ( strcmp(argv[arg_index], "-relax") == 0 )
 289          {
 290             arg_index++;
 291             relax = atoi(argv[arg_index++]);
 292          }
 293          else if ( strcmp(argv[arg_index], "-skip") == 0 )
 294          {
 295             arg_index++;
 296             skip = atoi(argv[arg_index++]);
 297          }
 298          else if ( strcmp(argv[arg_index], "-sym") == 0 )
 299          {
 300             arg_index++;
 301             sym = atoi(argv[arg_index++]);
 302          }
 303          else if ( strcmp(argv[arg_index], "-vis") == 0 )
 304          {
 305             arg_index++;
 306             vis = 1;
 307          }
 308          else if ( strcmp(argv[arg_index], "-help") == 0 )
 309          {
 310             print_usage = 1;
 311             break;
 312          }
 313          else
 314          {
 315             arg_index++;
 316          }
 317       }
 318 
 319       if ((print_usage) && (myid == 0))
 320       {
 321          printf("\n");
 322          printf("Usage: %s [<options>]\n", argv[0]);
 323          printf("\n");
 324          printf("  -n  <n>             : problem size per processor (default: 8)\n");
 325          printf("  -K  <K>             : choice for the diffusion coefficient (default: 1)\n");
 326          printf("  -B  <B>             : choice for the convection vector (default: 0)\n");
 327          printf("  -C  <C>             : choice for the reaction coefficient (default: 0)\n");
 328          printf("  -U0 <U0>            : choice for the boundary condition (default: 0)\n");
 329          printf("  -F  <F>             : choice for the right-hand side (default: 1) \n");
 330          printf("  -solver <ID>        : solver ID\n");
 331          printf("                        0  - SMG \n");
 332          printf("                        1  - PFMG\n");
 333          printf("                        10 - CG with SMG precond (default)\n");
 334          printf("                        11 - CG with PFMG precond\n");
 335          printf("                        17 - CG with 2-step Jacobi\n");
 336          printf("                        18 - CG with diagonal scaling\n");
 337          printf("                        19 - CG\n");
 338          printf("                        30 - GMRES with SMG precond\n");
 339          printf("                        31 - GMRES with PFMG precond\n");
 340          printf("                        37 - GMRES with 2-step Jacobi\n");
 341          printf("                        38 - GMRES with diagonal scaling\n");
 342          printf("                        39 - GMRES\n");
 343          printf("  -v <n_pre> <n_post> : number of pre and post relaxations\n");
 344          printf("  -rap <r>            : coarse grid operator type\n");
 345          printf("                        0 - Galerkin (default)\n");
 346          printf("                        1 - non-Galerkin ParFlow operators\n");
 347          printf("                        2 - Galerkin, general operators\n");
 348          printf("  -relax <r>          : relaxation type\n");
 349          printf("                        0 - Jacobi\n");
 350          printf("                        1 - Weighted Jacobi (default)\n");
 351          printf("                        2 - R/B Gauss-Seidel\n");
 352          printf("                        3 - R/B Gauss-Seidel (nonsymmetric)\n");
 353          printf("  -skip <s>           : skip levels in PFMG (0 or 1)\n");
 354          printf("  -sym <s>            : symmetric storage (1) or not (0)\n");
 355          printf("  -vis                : save the solution for GLVis visualization\n");
 356          printf("\n");
 357       }
 358 
 359       if (print_usage)
 360       {
 361          MPI_Finalize();
 362          return (0);
 363       }
 364    }
 365 
 366    /* Convection produces non-symmetric matrices */
 367    if (optionB && sym)
 368       optionB = 0;
 369 
 370    /* Figure out the processor grid (N x N).  The local
 371       problem size is indicated by n (n x n). pi and pj
 372       indicate position in the processor grid. */
 373    N  = sqrt(num_procs);
 374    h  = 1.0 / (N*n-1);
 375    h2 = h*h;
 376    pj = myid / N;
 377    pi = myid - pj*N;
 378 
 379    /* Define the nodes owned by the current processor (each processor's
 380       piece of the global grid) */
 381    ilower[0] = pi*n;
 382    ilower[1] = pj*n;
 383    iupper[0] = ilower[0] + n-1;
 384    iupper[1] = ilower[1] + n-1;
 385 
 386    /* 1. Set up a grid */
 387    {
 388       /* Create an empty 2D grid object */
 389       HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
 390 
 391       /* Add a new box to the grid */
 392       HYPRE_StructGridSetExtents(grid, ilower, iupper);
 393 
 394       /* This is a collective call finalizing the grid assembly.
 395          The grid is now ``ready to be used'' */
 396       HYPRE_StructGridAssemble(grid);
 397    }
 398 
 399    /* 2. Define the discretization stencil */
 400    if (sym == 0)
 401    {
 402       /* Define the geometry of the stencil */
 403       int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
 404 
 405       /* Create an empty 2D, 5-pt stencil object */
 406       HYPRE_StructStencilCreate(2, 5, &stencil);
 407 
 408       /* Assign stencil entries */
 409       for (i = 0; i < 5; i++)
 410          HYPRE_StructStencilSetElement(stencil, i, offsets[i]);
 411    }
 412    else /* Symmetric storage */
 413    {
 414       /* Define the geometry of the stencil */
 415       int offsets[3][2] = {{0,0}, {1,0}, {0,1}};
 416 
 417       /* Create an empty 2D, 3-pt stencil object */
 418       HYPRE_StructStencilCreate(2, 3, &stencil);
 419 
 420       /* Assign stencil entries */
 421       for (i = 0; i < 3; i++)
 422          HYPRE_StructStencilSetElement(stencil, i, offsets[i]);
 423    }
 424 
 425    /* 3. Set up Struct Vectors for b and x */
 426    {
 427       double *values;
 428 
 429       /* Create an empty vector object */
 430       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
 431       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
 432 
 433       /* Indicate that the vector coefficients are ready to be set */
 434       HYPRE_StructVectorInitialize(b);
 435       HYPRE_StructVectorInitialize(x);
 436 
 437       values = calloc((n*n), sizeof(double));
 438 
 439       /* Set the values of b in left-to-right, bottom-to-top order */
 440       for (k = 0, j = 0; j < n; j++)
 441          for (i = 0; i < n; i++, k++)
 442             values[k] = h2 * Eval(F,i,j);
 443       HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
 444 
 445       /* Set x = 0 */
 446       for (i = 0; i < (n*n); i ++)
 447          values[i] = 0.0;
 448       HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
 449 
 450       free(values);
 451 
 452       /* Assembling is postponed since the vectors will be further modified */
 453    }
 454 
 455    /* 4. Set up a Struct Matrix */
 456    {
 457       /* Create an empty matrix object */
 458       HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
 459 
 460       /* Use symmetric storage? */
 461       HYPRE_StructMatrixSetSymmetric(A, sym);
 462 
 463       /* Indicate that the matrix coefficients are ready to be set */
 464       HYPRE_StructMatrixInitialize(A);
 465 
 466       /* Set the stencil values in the interior. Here we set the values
 467          at every node. We will modify the boundary nodes later. */
 468       if (sym == 0)
 469       {
 470          int stencil_indices[5] = {0, 1, 2, 3, 4}; /* labels correspond
 471                                                       to the offsets */
 472          double *values;
 473 
 474          values = calloc(5*(n*n), sizeof(double));
 475 
 476          /* The order is left-to-right, bottom-to-top */
 477          for (k = 0, j = 0; j < n; j++)
 478             for (i = 0; i < n; i++, k+=5)
 479             {
 480                values[k+1] = - Eval(K,i-0.5,j) - Eval(B1,i-0.5,j);
 481 
 482                values[k+2] = - Eval(K,i+0.5,j) + Eval(B1,i+0.5,j);
 483 
 484                values[k+3] = - Eval(K,i,j-0.5) - Eval(B2,i,j-0.5);
 485 
 486                values[k+4] = - Eval(K,i,j+0.5) + Eval(B2,i,j+0.5);
 487 
 488                values[k] = h2 * Eval(C,i,j)
 489                   + Eval(K ,i-0.5,j) + Eval(K ,i+0.5,j)
 490                   + Eval(K ,i,j-0.5) + Eval(K ,i,j+0.5)
 491                   - Eval(B1,i-0.5,j) + Eval(B1,i+0.5,j)
 492                   - Eval(B2,i,j-0.5) + Eval(B2,i,j+0.5);
 493             }
 494 
 495          HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 5,
 496                                         stencil_indices, values);
 497 
 498          free(values);
 499       }
 500       else /* Symmetric storage */
 501       {
 502          int stencil_indices[3] = {0, 1, 2};
 503          double *values;
 504 
 505          values = calloc(3*(n*n), sizeof(double));
 506 
 507          /* The order is left-to-right, bottom-to-top */
 508          for (k = 0, j = 0; j < n; j++)
 509             for (i = 0; i < n; i++, k+=3)
 510             {
 511                values[k+1] = - Eval(K,i+0.5,j);
 512                values[k+2] = - Eval(K,i,j+0.5);
 513                values[k] = h2 * Eval(C,i,j)
 514                   + Eval(K,i+0.5,j) + Eval(K,i,j+0.5)
 515                   + Eval(K,i-0.5,j) + Eval(K,i,j-0.5);
 516             }
 517 
 518          HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 3,
 519                                         stencil_indices, values);
 520 
 521          free(values);
 522       }
 523    }
 524 
 525    /* 5. Set the boundary conditions, while eliminating the coefficients
 526          reaching ouside of the domain boundary. We must modify the matrix
 527          stencil and the corresponding rhs entries. */
 528    {
 529       int bc_ilower[2];
 530       int bc_iupper[2];
 531 
 532       int stencil_indices[5] = {0, 1, 2, 3, 4};
 533       double *values, *bvalues;
 534 
 535       int nentries;
 536       if (sym == 0)
 537          nentries = 5;
 538       else
 539          nentries = 3;
 540 
 541       values  = calloc(nentries*n, sizeof(double));
 542       bvalues = calloc(n, sizeof(double));
 543 
 544       /* The stencil at the boundary nodes is 1-0-0-0-0. Because
 545          we have I x_b = u_0; */
 546       for (i = 0; i < nentries*n; i += nentries)
 547       {
 548          values[i] = 1.0;
 549          for (j = 1; j < nentries; j++)
 550             values[i+j] = 0.0;
 551       }
 552 
 553       /* Processors at y = 0 */
 554       if (pj == 0)
 555       {
 556          bc_ilower[0] = pi*n;
 557          bc_ilower[1] = pj*n;
 558 
 559          bc_iupper[0] = bc_ilower[0] + n-1;
 560          bc_iupper[1] = bc_ilower[1];
 561 
 562          /* Modify the matrix */
 563          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
 564                                         stencil_indices, values);
 565 
 566          /* Put the boundary conditions in b */
 567          for (i = 0; i < n; i++)
 568             bvalues[i] = bcEval(U0,i,0);
 569 
 570          HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
 571       }
 572 
 573       /* Processors at y = 1 */
 574       if (pj == N-1)
 575       {
 576          bc_ilower[0] = pi*n;
 577          bc_ilower[1] = pj*n + n-1;
 578 
 579          bc_iupper[0] = bc_ilower[0] + n-1;
 580          bc_iupper[1] = bc_ilower[1];
 581 
 582          /* Modify the matrix */
 583          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
 584                                         stencil_indices, values);
 585 
 586          /* Put the boundary conditions in b */
 587          for (i = 0; i < n; i++)
 588             bvalues[i] = bcEval(U0,i,0);
 589 
 590          HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
 591       }
 592 
 593       /* Processors at x = 0 */
 594       if (pi == 0)
 595       {
 596          bc_ilower[0] = pi*n;
 597          bc_ilower[1] = pj*n;
 598 
 599          bc_iupper[0] = bc_ilower[0];
 600          bc_iupper[1] = bc_ilower[1] + n-1;
 601 
 602          /* Modify the matrix */
 603          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
 604                                         stencil_indices, values);
 605 
 606          /* Put the boundary conditions in b */
 607          for (j = 0; j < n; j++)
 608             bvalues[j] = bcEval(U0,0,j);
 609 
 610          HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
 611       }
 612 
 613       /* Processors at x = 1 */
 614       if (pi == N-1)
 615       {
 616          bc_ilower[0] = pi*n + n-1;
 617          bc_ilower[1] = pj*n;
 618 
 619          bc_iupper[0] = bc_ilower[0];
 620          bc_iupper[1] = bc_ilower[1] + n-1;
 621 
 622          /* Modify the matrix */
 623          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
 624                                         stencil_indices, values);
 625 
 626          /* Put the boundary conditions in b */
 627          for (j = 0; j < n; j++)
 628             bvalues[j] = bcEval(U0,0,j);
 629 
 630          HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
 631       }
 632 
 633       /* Recall that the system we are solving is:
 634          [A_ii 0; 0 I] [x_i ; x_b] = [b_i - A_ib u_0; u_0].
 635          This requires removing the connections between the interior
 636          and boundary nodes that we have set up when we set the
 637          5pt stencil at each node. We adjust for removing
 638          these connections by appropriately modifying the rhs.
 639          For the symm ordering scheme, just do the top and right
 640          boundary */
 641 
 642       /* Processors at y = 0, neighbors of boundary nodes */
 643       if (pj == 0)
 644       {
 645          bc_ilower[0] = pi*n;
 646          bc_ilower[1] = pj*n + 1;
 647 
 648          bc_iupper[0] = bc_ilower[0] + n-1;
 649          bc_iupper[1] = bc_ilower[1];
 650 
 651          stencil_indices[0] = 3;
 652 
 653          /* Modify the matrix */
 654          for (i = 0; i < n; i++)
 655             bvalues[i] = 0.0;
 656 
 657          if (sym == 0)
 658             HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
 659                                            stencil_indices, bvalues);
 660 
 661          /* Eliminate the boundary conditions in b */
 662          for (i = 0; i < n; i++)
 663             bvalues[i] = bcEval(U0,i,-1) * (bcEval(K,i,-0.5)+bcEval(B2,i,-0.5));
 664 
 665          if (pi == 0)
 666             bvalues[0] = 0.0;
 667 
 668          if (pi == N-1)
 669             bvalues[n-1] = 0.0;
 670 
 671          /* Note the use of AddToBoxValues (because we have already set values
 672             at these nodes) */
 673          HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
 674       }
 675 
 676       /* Processors at x = 0, neighbors of boundary nodes */
 677       if (pi == 0)
 678       {
 679          bc_ilower[0] = pi*n + 1;
 680          bc_ilower[1] = pj*n;
 681 
 682          bc_iupper[0] = bc_ilower[0];
 683          bc_iupper[1] = bc_ilower[1] + n-1;
 684 
 685          stencil_indices[0] = 1;
 686 
 687          /* Modify the matrix */
 688          for (j = 0; j < n; j++)
 689             bvalues[j] = 0.0;
 690 
 691          if (sym == 0)
 692             HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
 693                                            stencil_indices, bvalues);
 694 
 695          /* Eliminate the boundary conditions in b */
 696          for (j = 0; j < n; j++)
 697             bvalues[j] = bcEval(U0,-1,j) * (bcEval(K,-0.5,j)+bcEval(B1,-0.5,j));
 698 
 699          if (pj == 0)
 700             bvalues[0] = 0.0;
 701 
 702          if (pj == N-1)
 703             bvalues[n-1] = 0.0;
 704 
 705          HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
 706       }
 707 
 708       /* Processors at y = 1, neighbors of boundary nodes */
 709       if (pj == N-1)
 710       {
 711          bc_ilower[0] = pi*n;
 712          bc_ilower[1] = pj*n + (n-1) -1;
 713 
 714          bc_iupper[0] = bc_ilower[0] + n-1;
 715          bc_iupper[1] = bc_ilower[1];
 716 
 717          if (sym == 0)
 718             stencil_indices[0] = 4;
 719          else
 720             stencil_indices[0] = 2;
 721 
 722          /* Modify the matrix */
 723          for (i = 0; i < n; i++)
 724             bvalues[i] = 0.0;
 725 
 726          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
 727                                         stencil_indices, bvalues);
 728 
 729          /* Eliminate the boundary conditions in b */
 730          for (i = 0; i < n; i++)
 731             bvalues[i] = bcEval(U0,i,1) * (bcEval(K,i,0.5)+bcEval(B2,i,0.5));
 732 
 733          if (pi == 0)
 734             bvalues[0] = 0.0;
 735 
 736          if (pi == N-1)
 737             bvalues[n-1] = 0.0;
 738 
 739          HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
 740       }
 741 
 742       /* Processors at x = 1, neighbors of boundary nodes */
 743       if (pi == N-1)
 744       {
 745          bc_ilower[0] = pi*n + (n-1) - 1;
 746          bc_ilower[1] = pj*n;
 747 
 748          bc_iupper[0] = bc_ilower[0];
 749          bc_iupper[1] = bc_ilower[1] + n-1;
 750 
 751          if (sym == 0)
 752             stencil_indices[0] = 2;
 753          else
 754             stencil_indices[0] = 1;
 755 
 756          /* Modify the matrix */
 757          for (j = 0; j < n; j++)
 758             bvalues[j] = 0.0;
 759 
 760          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
 761                                         stencil_indices, bvalues);
 762 
 763          /* Eliminate the boundary conditions in b */
 764          for (j = 0; j < n; j++)
 765             bvalues[j] = bcEval(U0,1,j) * (bcEval(K,0.5,j)+bcEval(B1,0.5,j));
 766 
 767          if (pj == 0)
 768             bvalues[0] = 0.0;
 769 
 770          if (pj == N-1)
 771             bvalues[n-1] = 0.0;
 772 
 773          HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
 774       }
 775 
 776       free(values);
 777       free(bvalues);
 778    }
 779 
 780    /* Finalize the vector and matrix assembly */
 781    HYPRE_StructMatrixAssemble(A);
 782    HYPRE_StructVectorAssemble(b);
 783    HYPRE_StructVectorAssemble(x);
 784 
 785    /* 6. Set up and use a solver */
 786    if (solver_id == 0) /* SMG */
 787    {
 788       /* Start timing */
 789       time_index = hypre_InitializeTiming("SMG Setup");
 790       hypre_BeginTiming(time_index);
 791 
 792       /* Options and setup */
 793       HYPRE_StructSMGCreate(MPI_COMM_WORLD, &solver);
 794       HYPRE_StructSMGSetMemoryUse(solver, 0);
 795       HYPRE_StructSMGSetMaxIter(solver, 50);
 796       HYPRE_StructSMGSetTol(solver, 1.0e-06);
 797       HYPRE_StructSMGSetRelChange(solver, 0);
 798       HYPRE_StructSMGSetNumPreRelax(solver, n_pre);
 799       HYPRE_StructSMGSetNumPostRelax(solver, n_post);
 800       HYPRE_StructSMGSetPrintLevel(solver, 1);
 801       HYPRE_StructSMGSetLogging(solver, 1);
 802       HYPRE_StructSMGSetup(solver, A, b, x);
 803 
 804       /* Finalize current timing */
 805       hypre_EndTiming(time_index);
 806       hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
 807       hypre_FinalizeTiming(time_index);
 808       hypre_ClearTiming();
 809 
 810       /* Start timing again */
 811       time_index = hypre_InitializeTiming("SMG Solve");
 812       hypre_BeginTiming(time_index);
 813 
 814       /* Solve */
 815       HYPRE_StructSMGSolve(solver, A, b, x);
 816 
 817       /* Finalize current timing */
 818       hypre_EndTiming(time_index);
 819       hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
 820       hypre_FinalizeTiming(time_index);
 821       hypre_ClearTiming();
 822 
 823       /* Get info and release memory */
 824       HYPRE_StructSMGGetNumIterations(solver, &num_iterations);
 825       HYPRE_StructSMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
 826       HYPRE_StructSMGDestroy(solver);
 827    }
 828 
 829    if (solver_id == 1) /* PFMG */
 830    {
 831       /* Start timing */
 832       time_index = hypre_InitializeTiming("PFMG Setup");
 833       hypre_BeginTiming(time_index);
 834 
 835       /* Options and setup */
 836       HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &solver);
 837       HYPRE_StructPFMGSetMaxIter(solver, 50);
 838       HYPRE_StructPFMGSetTol(solver, 1.0e-06);
 839       HYPRE_StructPFMGSetRelChange(solver, 0);
 840       HYPRE_StructPFMGSetRAPType(solver, rap);
 841       HYPRE_StructPFMGSetRelaxType(solver, relax);
 842       HYPRE_StructPFMGSetNumPreRelax(solver, n_pre);
 843       HYPRE_StructPFMGSetNumPostRelax(solver, n_post);
 844       HYPRE_StructPFMGSetSkipRelax(solver, skip);
 845       HYPRE_StructPFMGSetPrintLevel(solver, 1);
 846       HYPRE_StructPFMGSetLogging(solver, 1);
 847       HYPRE_StructPFMGSetup(solver, A, b, x);
 848 
 849       /* Finalize current timing */
 850       hypre_EndTiming(time_index);
 851       hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
 852       hypre_FinalizeTiming(time_index);
 853       hypre_ClearTiming();
 854 
 855       /* Start timing again */
 856       time_index = hypre_InitializeTiming("PFMG Solve");
 857       hypre_BeginTiming(time_index);
 858 
 859       /* Solve */
 860       HYPRE_StructPFMGSolve(solver, A, b, x);
 861 
 862       /* Finalize current timing */
 863       hypre_EndTiming(time_index);
 864       hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
 865       hypre_FinalizeTiming(time_index);
 866       hypre_ClearTiming();
 867 
 868       /* Get info and release memory */
 869       HYPRE_StructPFMGGetNumIterations(solver, &num_iterations);
 870       HYPRE_StructPFMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
 871       HYPRE_StructPFMGDestroy(solver);
 872    }
 873 
 874    /* Preconditioned CG */
 875    if ((solver_id > 9) && (solver_id < 20))
 876    {
 877       time_index = hypre_InitializeTiming("PCG Setup");
 878       hypre_BeginTiming(time_index);
 879 
 880       HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
 881       HYPRE_StructPCGSetMaxIter(solver, 200 );
 882       HYPRE_StructPCGSetTol(solver, 1.0e-06 );
 883       HYPRE_StructPCGSetTwoNorm(solver, 1 );
 884       HYPRE_StructPCGSetRelChange(solver, 0 );
 885       HYPRE_StructPCGSetPrintLevel(solver, 2 );
 886 
 887       if (solver_id == 10)
 888       {
 889          /* use symmetric SMG as preconditioner */
 890          HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
 891          HYPRE_StructSMGSetMemoryUse(precond, 0);
 892          HYPRE_StructSMGSetMaxIter(precond, 1);
 893          HYPRE_StructSMGSetTol(precond, 0.0);
 894          HYPRE_StructSMGSetZeroGuess(precond);
 895          HYPRE_StructSMGSetNumPreRelax(precond, n_pre);
 896          HYPRE_StructSMGSetNumPostRelax(precond, n_post);
 897          HYPRE_StructSMGSetPrintLevel(precond, 0);
 898          HYPRE_StructSMGSetLogging(precond, 0);
 899          HYPRE_StructPCGSetPrecond(solver,
 900                                    HYPRE_StructSMGSolve,
 901                                    HYPRE_StructSMGSetup,
 902                                    precond);
 903       }
 904 
 905       else if (solver_id == 11)
 906       {
 907          /* use symmetric PFMG as preconditioner */
 908          HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
 909          HYPRE_StructPFMGSetMaxIter(precond, 1);
 910          HYPRE_StructPFMGSetTol(precond, 0.0);
 911          HYPRE_StructPFMGSetZeroGuess(precond);
 912          HYPRE_StructPFMGSetRAPType(precond, rap);
 913          HYPRE_StructPFMGSetRelaxType(precond, relax);
 914          HYPRE_StructPFMGSetNumPreRelax(precond, n_pre);
 915          HYPRE_StructPFMGSetNumPostRelax(precond, n_post);
 916          HYPRE_StructPFMGSetSkipRelax(precond, skip);
 917          HYPRE_StructPFMGSetPrintLevel(precond, 0);
 918          HYPRE_StructPFMGSetLogging(precond, 0);
 919          HYPRE_StructPCGSetPrecond(solver,
 920                                    HYPRE_StructPFMGSolve,
 921                                    HYPRE_StructPFMGSetup,
 922                                    precond);
 923       }
 924 
 925       else if (solver_id == 17)
 926       {
 927          /* use two-step Jacobi as preconditioner */
 928          HYPRE_StructJacobiCreate(MPI_COMM_WORLD, &precond);
 929          HYPRE_StructJacobiSetMaxIter(precond, 2);
 930          HYPRE_StructJacobiSetTol(precond, 0.0);
 931          HYPRE_StructJacobiSetZeroGuess(precond);
 932          HYPRE_StructPCGSetPrecond( solver,
 933                                     HYPRE_StructJacobiSolve,
 934                                     HYPRE_StructJacobiSetup,
 935                                     precond);
 936       }
 937 
 938       else if (solver_id == 18)
 939       {
 940          /* use diagonal scaling as preconditioner */
 941          precond = NULL;
 942          HYPRE_StructPCGSetPrecond(solver,
 943                                    HYPRE_StructDiagScale,
 944                                    HYPRE_StructDiagScaleSetup,
 945                                    precond);
 946       }
 947 
 948       /* PCG Setup */
 949       HYPRE_StructPCGSetup(solver, A, b, x );
 950 
 951       hypre_EndTiming(time_index);
 952       hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
 953       hypre_FinalizeTiming(time_index);
 954       hypre_ClearTiming();
 955 
 956       time_index = hypre_InitializeTiming("PCG Solve");
 957       hypre_BeginTiming(time_index);
 958 
 959       /* PCG Solve */
 960       HYPRE_StructPCGSolve(solver, A, b, x);
 961 
 962       hypre_EndTiming(time_index);
 963       hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
 964       hypre_FinalizeTiming(time_index);
 965       hypre_ClearTiming();
 966 
 967       /* Get info and release memory */
 968       HYPRE_StructPCGGetNumIterations( solver, &num_iterations );
 969       HYPRE_StructPCGGetFinalRelativeResidualNorm( solver, &final_res_norm );
 970       HYPRE_StructPCGDestroy(solver);
 971 
 972       if (solver_id == 10)
 973       {
 974          HYPRE_StructSMGDestroy(precond);
 975       }
 976       else if (solver_id == 11 )
 977       {
 978          HYPRE_StructPFMGDestroy(precond);
 979       }
 980       else if (solver_id == 17)
 981       {
 982          HYPRE_StructJacobiDestroy(precond);
 983       }
 984    }
 985 
 986    /* Preconditioned GMRES */
 987    if ((solver_id > 29) && (solver_id < 40))
 988    {
 989       time_index = hypre_InitializeTiming("GMRES Setup");
 990       hypre_BeginTiming(time_index);
 991 
 992       HYPRE_StructGMRESCreate(MPI_COMM_WORLD, &solver);
 993 
 994       /* Note that GMRES can be used with all the interfaces - not
 995          just the struct.  So here we demonstrate the
 996          more generic GMRES interface functions. Since we have chosen
 997          a struct solver then we must type cast to the more generic
 998          HYPRE_Solver when setting options with these generic functions.
 999          Note that one could declare the solver to be
1000          type HYPRE_Solver, and then the casting would not be necessary.*/
1001 
1002       HYPRE_GMRESSetMaxIter((HYPRE_Solver) solver, 500 );
1003       HYPRE_GMRESSetKDim((HYPRE_Solver) solver,30);
1004       HYPRE_GMRESSetTol((HYPRE_Solver) solver, 1.0e-06 );
1005       HYPRE_GMRESSetPrintLevel((HYPRE_Solver) solver, 2 );
1006       HYPRE_GMRESSetLogging((HYPRE_Solver) solver, 1 );
1007 
1008       if (solver_id == 30)
1009       {
1010          /* use symmetric SMG as preconditioner */
1011          HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
1012          HYPRE_StructSMGSetMemoryUse(precond, 0);
1013          HYPRE_StructSMGSetMaxIter(precond, 1);
1014          HYPRE_StructSMGSetTol(precond, 0.0);
1015          HYPRE_StructSMGSetZeroGuess(precond);
1016          HYPRE_StructSMGSetNumPreRelax(precond, n_pre);
1017          HYPRE_StructSMGSetNumPostRelax(precond, n_post);
1018          HYPRE_StructSMGSetPrintLevel(precond, 0);
1019          HYPRE_StructSMGSetLogging(precond, 0);
1020          HYPRE_StructGMRESSetPrecond(solver,
1021                                      HYPRE_StructSMGSolve,
1022                                      HYPRE_StructSMGSetup,
1023                                      precond);
1024       }
1025 
1026       else if (solver_id == 31)
1027       {
1028          /* use symmetric PFMG as preconditioner */
1029          HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
1030          HYPRE_StructPFMGSetMaxIter(precond, 1);
1031          HYPRE_StructPFMGSetTol(precond, 0.0);
1032          HYPRE_StructPFMGSetZeroGuess(precond);
1033          HYPRE_StructPFMGSetRAPType(precond, rap);
1034          HYPRE_StructPFMGSetRelaxType(precond, relax);
1035          HYPRE_StructPFMGSetNumPreRelax(precond, n_pre);
1036          HYPRE_StructPFMGSetNumPostRelax(precond, n_post);
1037          HYPRE_StructPFMGSetSkipRelax(precond, skip);
1038          HYPRE_StructPFMGSetPrintLevel(precond, 0);
1039          HYPRE_StructPFMGSetLogging(precond, 0);
1040          HYPRE_StructGMRESSetPrecond( solver,
1041                                       HYPRE_StructPFMGSolve,
1042                                       HYPRE_StructPFMGSetup,
1043                                       precond);
1044       }
1045 
1046       else if (solver_id == 37)
1047       {
1048          /* use two-step Jacobi as preconditioner */
1049          HYPRE_StructJacobiCreate(MPI_COMM_WORLD, &precond);
1050          HYPRE_StructJacobiSetMaxIter(precond, 2);
1051          HYPRE_StructJacobiSetTol(precond, 0.0);
1052          HYPRE_StructJacobiSetZeroGuess(precond);
1053          HYPRE_StructGMRESSetPrecond( solver,
1054                                       HYPRE_StructJacobiSolve,
1055                                       HYPRE_StructJacobiSetup,
1056                                       precond);
1057       }
1058 
1059       else if (solver_id == 38)
1060       {
1061          /* use diagonal scaling as preconditioner */
1062          precond = NULL;
1063          HYPRE_StructGMRESSetPrecond( solver,
1064                                      HYPRE_StructDiagScale,
1065                                      HYPRE_StructDiagScaleSetup,
1066                                      precond);
1067       }
1068 
1069       /* GMRES Setup */
1070       HYPRE_StructGMRESSetup(solver, A, b, x );
1071 
1072       hypre_EndTiming(time_index);
1073       hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
1074       hypre_FinalizeTiming(time_index);
1075       hypre_ClearTiming();
1076 
1077       time_index = hypre_InitializeTiming("GMRES Solve");
1078       hypre_BeginTiming(time_index);
1079 
1080       /* GMRES Solve */
1081       HYPRE_StructGMRESSolve(solver, A, b, x);
1082 
1083       hypre_EndTiming(time_index);
1084       hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
1085       hypre_FinalizeTiming(time_index);
1086       hypre_ClearTiming();
1087 
1088       /* Get info and release memory */
1089       HYPRE_StructGMRESGetNumIterations(solver, &num_iterations);
1090       HYPRE_StructGMRESGetFinalRelativeResidualNorm(solver, &final_res_norm);
1091       HYPRE_StructGMRESDestroy(solver);
1092 
1093       if (solver_id == 30)
1094       {
1095          HYPRE_StructSMGDestroy(precond);
1096       }
1097       else if (solver_id == 31)
1098       {
1099          HYPRE_StructPFMGDestroy(precond);
1100       }
1101       else if (solver_id == 37)
1102       {
1103          HYPRE_StructJacobiDestroy(precond);
1104       }
1105    }
1106 
1107    /* Save the solution for GLVis visualization, see vis/glvis-ex4.sh */
1108    if (vis)
1109    {
1110       FILE *file;
1111       char filename[255];
1112 
1113       int nvalues = n*n;
1114       double *values = calloc(nvalues, sizeof(double));
1115 
1116       /* get the local solution */
1117       HYPRE_StructVectorGetBoxValues(x, ilower, iupper, values);
1118 
1119       sprintf(filename, "%s.%06d", "vis/ex4.sol", myid);
1120       if ((file = fopen(filename, "w")) == NULL)
1121       {
1122          printf("Error: can't open output file %s\n", filename);
1123          MPI_Finalize();
1124          exit(1);
1125       }
1126 
1127       /* save solution with global unknown numbers */
1128       k = 0;
1129       for (j = 0; j < n; j++)
1130          for (i = 0; i < n; i++)
1131             fprintf(file, "%06d %.14e\n", pj*N*n*n+pi*n+j*N*n+i, values[k++]);
1132 
1133       fflush(file);
1134       fclose(file);
1135       free(values);
1136 
1137       /* save global finite element mesh */
1138       if (myid == 0)
1139          GLVis_PrintGlobalSquareMesh("vis/ex4.mesh", N*n-1);
1140    }
1141 
1142    if (myid == 0)
1143    {
1144       printf("\n");
1145       printf("Iterations = %d\n", num_iterations);
1146       printf("Final Relative Residual Norm = %e\n", final_res_norm);
1147       printf("\n");
1148    }
1149 
1150    /* Free memory */
1151    HYPRE_StructGridDestroy(grid);
1152    HYPRE_StructStencilDestroy(stencil);
1153    HYPRE_StructMatrixDestroy(A);
1154    HYPRE_StructVectorDestroy(b);
1155    HYPRE_StructVectorDestroy(x);
1156 
1157    /* Finalize MPI */
1158    MPI_Finalize();
1159 
1160    return (0);
1161 }
syntax highlighted by Code2HTML, v. 0.9.1