download the original source code.
  1 /*
  2    Example 3
  3 
  4    Interface:      Structured interface (Struct)
  5 
  6    Compile with:   make ex3
  7 
  8    Sample run:     mpirun -np 16 ex3 -n 33 -solver 0 -v 1 1
  9 
 10    To see options: ex3 -help
 11 
 12    Description:    This code solves a system corresponding to a discretization
 13                    of the Laplace equation -Delta u = 1 with zero boundary
 14                    conditions on the unit square.  The domain is split into
 15                    an N x N processor grid.  Thus, the given number of processors
 16                    should be a perfect square.  Each processor's piece of the
 17                    grid has n x n cells with n x n nodes connected by the
 18                    standard 5-point stencil. Note that the struct interface
 19                    assumes a cell-centered grid, and, therefore, the nodes are
 20                    not shared.  This example demonstrates more features than the
 21                    previous two struct examples (Example 1 and Example 2).  Two
 22                    solvers are available.
 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. We can split the matrix A as
 27                                    A = [A_ii A_ib; A_bi A_bb].
 28                    Let u_0 be the Dirichlet B.C.  We can simply say that x_b = u_0.
 29                    If b_i is the right-hand side, then we just need to solve in
 30                    the interior:
 31                                     A_ii x_i = b_i - A_ib u_0.
 32                    For this partitcular example, u_0 = 0, so we are just solving
 33                    A_ii x_i = b_i.
 34 
 35                    We recommend viewing examples 1 and 2 before viewing this
 36                    example.
 37 */
 38 
 39 #include <math.h>
 40 #include "_hypre_utilities.h"
 41 #include "HYPRE_struct_ls.h"
 42 
 43 #include "vis.c"
 44 
 45 int main (int argc, char *argv[])
 46 {
 47    int i, j, k;
 48 
 49    int myid, num_procs;
 50 
 51    int n, N, pi, pj;
 52    double h, h2;
 53    int ilower[2], iupper[2];
 54 
 55    int solver_id;
 56    int n_pre, n_post;
 57 
 58    HYPRE_StructGrid     grid;
 59    HYPRE_StructStencil  stencil;
 60    HYPRE_StructMatrix   A;
 61    HYPRE_StructVector   b;
 62    HYPRE_StructVector   x;
 63    HYPRE_StructSolver   solver;
 64    HYPRE_StructSolver   precond;
 65 
 66    int num_iterations;
 67    double final_res_norm;
 68 
 69    int vis;
 70 
 71    /* Initialize MPI */
 72    MPI_Init(&argc, &argv);
 73    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 74    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 75 
 76    /* Set defaults */
 77    n = 33;
 78    solver_id = 0;
 79    n_pre  = 1;
 80    n_post = 1;
 81    vis = 0;
 82 
 83    /* Parse command line */
 84    {
 85       int arg_index = 0;
 86       int print_usage = 0;
 87 
 88       while (arg_index < argc)
 89       {
 90          if ( strcmp(argv[arg_index], "-n") == 0 )
 91          {
 92             arg_index++;
 93             n = atoi(argv[arg_index++]);
 94          }
 95          else if ( strcmp(argv[arg_index], "-solver") == 0 )
 96          {
 97             arg_index++;
 98             solver_id = atoi(argv[arg_index++]);
 99          }
100          else if ( strcmp(argv[arg_index], "-v") == 0 )
101          {
102             arg_index++;
103             n_pre = atoi(argv[arg_index++]);
104             n_post = atoi(argv[arg_index++]);
105          }
106          else if ( strcmp(argv[arg_index], "-vis") == 0 )
107          {
108             arg_index++;
109             vis = 1;
110          }
111          else if ( strcmp(argv[arg_index], "-help") == 0 )
112          {
113             print_usage = 1;
114             break;
115          }
116          else
117          {
118             arg_index++;
119          }
120       }
121 
122       if ((print_usage) && (myid == 0))
123       {
124          printf("\n");
125          printf("Usage: %s [<options>]\n", argv[0]);
126          printf("\n");
127          printf("  -n <n>              : problem size per processor (default: 33)\n");
128          printf("  -solver <ID>        : solver ID\n");
129          printf("                        0  - PCG with SMG precond (default)\n");
130          printf("                        1  - SMG\n");
131          printf("  -v <n_pre> <n_post> : number of pre and post relaxations (default: 1 1)\n");
132          printf("  -vis                : save the solution for GLVis visualization\n");
133          printf("\n");
134       }
135 
136       if (print_usage)
137       {
138          MPI_Finalize();
139          return (0);
140       }
141    }
142 
143    /* Figure out the processor grid (N x N).  The local problem
144       size for the interior nodes is indicated by n (n x n).
145       pi and pj indicate position in the processor grid. */
146    N  = sqrt(num_procs);
147    h  = 1.0 / (N*n+1); /* note that when calculating h we must
148                           remember to count the boundary nodes */
149    h2 = h*h;
150    pj = myid / N;
151    pi = myid - pj*N;
152 
153    /* Figure out the extents of each processor's piece of the grid. */
154    ilower[0] = pi*n;
155    ilower[1] = pj*n;
156 
157    iupper[0] = ilower[0] + n-1;
158    iupper[1] = ilower[1] + n-1;
159 
160    /* 1. Set up a grid */
161    {
162       /* Create an empty 2D grid object */
163       HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
164 
165       /* Add a new box to the grid */
166       HYPRE_StructGridSetExtents(grid, ilower, iupper);
167 
168       /* This is a collective call finalizing the grid assembly.
169          The grid is now ``ready to be used'' */
170       HYPRE_StructGridAssemble(grid);
171    }
172 
173    /* 2. Define the discretization stencil */
174    {
175       /* Create an empty 2D, 5-pt stencil object */
176       HYPRE_StructStencilCreate(2, 5, &stencil);
177 
178       /* Define the geometry of the stencil */
179       {
180          int entry;
181          int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
182 
183          for (entry = 0; entry < 5; entry++)
184             HYPRE_StructStencilSetElement(stencil, entry, offsets[entry]);
185       }
186    }
187 
188    /* 3. Set up a Struct Matrix */
189    {
190       int nentries = 5;
191       int nvalues = nentries*n*n;
192       double *values;
193       int stencil_indices[5];
194 
195       /* Create an empty matrix object */
196       HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
197 
198       /* Indicate that the matrix coefficients are ready to be set */
199       HYPRE_StructMatrixInitialize(A);
200 
201       values = calloc(nvalues, sizeof(double));
202 
203       for (j = 0; j < nentries; j++)
204          stencil_indices[j] = j;
205 
206       /* Set the standard stencil at each grid point,
207          we will fix the boundaries later */
208       for (i = 0; i < nvalues; i += nentries)
209       {
210          values[i] = 4.0;
211          for (j = 1; j < nentries; j++)
212             values[i+j] = -1.0;
213       }
214 
215       HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
216                                      stencil_indices, values);
217 
218       free(values);
219    }
220 
221    /* 4. Incorporate the zero boundary conditions: go along each edge of
222          the domain and set the stencil entry that reaches to the boundary to
223          zero.*/
224    {
225       int bc_ilower[2];
226       int bc_iupper[2];
227       int nentries = 1;
228       int nvalues  = nentries*n; /*  number of stencil entries times the length
229                                      of one side of my grid box */
230       double *values;
231       int stencil_indices[1];
232 
233       values = calloc(nvalues, sizeof(double));
234       for (j = 0; j < nvalues; j++)
235          values[j] = 0.0;
236 
237       /* Recall: pi and pj describe position in the processor grid */
238       if (pj == 0)
239       {
240          /* Bottom row of grid points */
241          bc_ilower[0] = pi*n;
242          bc_ilower[1] = pj*n;
243 
244          bc_iupper[0] = bc_ilower[0] + n-1;
245          bc_iupper[1] = bc_ilower[1];
246 
247          stencil_indices[0] = 3;
248 
249          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
250                                         stencil_indices, values);
251       }
252 
253       if (pj == N-1)
254       {
255          /* upper row of grid points */
256          bc_ilower[0] = pi*n;
257          bc_ilower[1] = pj*n + n-1;
258 
259          bc_iupper[0] = bc_ilower[0] + n-1;
260          bc_iupper[1] = bc_ilower[1];
261 
262          stencil_indices[0] = 4;
263 
264          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
265                                         stencil_indices, values);
266       }
267 
268       if (pi == 0)
269       {
270          /* Left row of grid points */
271          bc_ilower[0] = pi*n;
272          bc_ilower[1] = pj*n;
273 
274          bc_iupper[0] = bc_ilower[0];
275          bc_iupper[1] = bc_ilower[1] + n-1;
276 
277          stencil_indices[0] = 1;
278 
279          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
280                                         stencil_indices, values);
281       }
282 
283       if (pi == N-1)
284       {
285          /* Right row of grid points */
286          bc_ilower[0] = pi*n + n-1;
287          bc_ilower[1] = pj*n;
288 
289          bc_iupper[0] = bc_ilower[0];
290          bc_iupper[1] = bc_ilower[1] + n-1;
291 
292          stencil_indices[0] = 2;
293 
294          HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
295                                         stencil_indices, values);
296       }
297 
298       free(values);
299    }
300 
301    /* This is a collective call finalizing the matrix assembly.
302       The matrix is now ``ready to be used'' */
303    HYPRE_StructMatrixAssemble(A);
304 
305    /* 5. Set up Struct Vectors for b and x */
306    {
307       int    nvalues = n*n;
308       double *values;
309 
310       values = calloc(nvalues, sizeof(double));
311 
312       /* Create an empty vector object */
313       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
314       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
315 
316       /* Indicate that the vector coefficients are ready to be set */
317       HYPRE_StructVectorInitialize(b);
318       HYPRE_StructVectorInitialize(x);
319 
320      /* Set the values */
321       for (i = 0; i < nvalues; i ++)
322          values[i] = h2;
323       HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
324 
325       for (i = 0; i < nvalues; i ++)
326          values[i] = 0.0;
327       HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
328 
329       free(values);
330 
331       /* This is a collective call finalizing the vector assembly.
332          The vector is now ``ready to be used'' */
333       HYPRE_StructVectorAssemble(b);
334       HYPRE_StructVectorAssemble(x);
335    }
336 
337    /* 6. Set up and use a struct solver
338       (Solver options can be found in the Reference Manual.) */
339    if (solver_id == 0)
340    {
341       HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
342       HYPRE_StructPCGSetMaxIter(solver, 50 );
343       HYPRE_StructPCGSetTol(solver, 1.0e-06 );
344       HYPRE_StructPCGSetTwoNorm(solver, 1 );
345       HYPRE_StructPCGSetRelChange(solver, 0 );
346       HYPRE_StructPCGSetPrintLevel(solver, 2 ); /* print each CG iteration */
347       HYPRE_StructPCGSetLogging(solver, 1);
348 
349       /* Use symmetric SMG as preconditioner */
350       HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
351       HYPRE_StructSMGSetMemoryUse(precond, 0);
352       HYPRE_StructSMGSetMaxIter(precond, 1);
353       HYPRE_StructSMGSetTol(precond, 0.0);
354       HYPRE_StructSMGSetZeroGuess(precond);
355       HYPRE_StructSMGSetNumPreRelax(precond, 1);
356       HYPRE_StructSMGSetNumPostRelax(precond, 1);
357 
358       /* Set the preconditioner and solve */
359       HYPRE_StructPCGSetPrecond(solver, HYPRE_StructSMGSolve,
360                                   HYPRE_StructSMGSetup, precond);
361       HYPRE_StructPCGSetup(solver, A, b, x);
362       HYPRE_StructPCGSolve(solver, A, b, x);
363 
364       /* Get some info on the run */
365       HYPRE_StructPCGGetNumIterations(solver, &num_iterations);
366       HYPRE_StructPCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
367 
368       /* Clean up */
369       HYPRE_StructPCGDestroy(solver);
370    }
371 
372    if (solver_id == 1)
373    {
374       HYPRE_StructSMGCreate(MPI_COMM_WORLD, &solver);
375       HYPRE_StructSMGSetMemoryUse(solver, 0);
376       HYPRE_StructSMGSetMaxIter(solver, 50);
377       HYPRE_StructSMGSetTol(solver, 1.0e-06);
378       HYPRE_StructSMGSetRelChange(solver, 0);
379       HYPRE_StructSMGSetNumPreRelax(solver, n_pre);
380       HYPRE_StructSMGSetNumPostRelax(solver, n_post);
381       /* Logging must be on to get iterations and residual norm info below */
382       HYPRE_StructSMGSetLogging(solver, 1);
383 
384       /* Setup and solve */
385       HYPRE_StructSMGSetup(solver, A, b, x);
386       HYPRE_StructSMGSolve(solver, A, b, x);
387 
388       /* Get some info on the run */
389       HYPRE_StructSMGGetNumIterations(solver, &num_iterations);
390       HYPRE_StructSMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
391 
392       /* Clean up */
393       HYPRE_StructSMGDestroy(solver);
394    }
395 
396    /* Save the solution for GLVis visualization, see vis/glvis-ex3.sh */
397    if (vis)
398    {
399       FILE *file;
400       char filename[255];
401 
402       int nvalues = n*n;
403       double *values = calloc(nvalues, sizeof(double));
404 
405       /* get the local solution */
406       HYPRE_StructVectorGetBoxValues(x, ilower, iupper, values);
407 
408       sprintf(filename, "%s.%06d", "vis/ex3.sol", myid);
409       if ((file = fopen(filename, "w")) == NULL)
410       {
411          printf("Error: can't open output file %s\n", filename);
412          MPI_Finalize();
413          exit(1);
414       }
415 
416       /* save solution with global unknown numbers */
417       k = 0;
418       for (j = 0; j < n; j++)
419          for (i = 0; i < n; i++)
420             fprintf(file, "%06d %.14e\n", pj*N*n*n+pi*n+j*N*n+i, values[k++]);
421 
422       fflush(file);
423       fclose(file);
424       free(values);
425 
426       /* save global finite element mesh */
427       if (myid == 0)
428          GLVis_PrintGlobalSquareMesh("vis/ex3.mesh", N*n-1);
429    }
430 
431    if (myid == 0)
432    {
433       printf("\n");
434       printf("Iterations = %d\n", num_iterations);
435       printf("Final Relative Residual Norm = %g\n", final_res_norm);
436       printf("\n");
437    }
438 
439    /* Free memory */
440    HYPRE_StructGridDestroy(grid);
441    HYPRE_StructStencilDestroy(stencil);
442    HYPRE_StructMatrixDestroy(A);
443    HYPRE_StructVectorDestroy(b);
444    HYPRE_StructVectorDestroy(x);
445 
446    /* Finalize MPI */
447    MPI_Finalize();
448 
449    return (0);
450 }


syntax highlighted by Code2HTML, v. 0.9.1