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


syntax highlighted by Code2HTML, v. 0.9.1