download the original source code.
  1 /*
  2    Example 2
  3 
  4    Interface:    Structured interface (Struct)
  5 
  6    Compile with: make ex2
  7 
  8    Sample run:   mpirun -np 2 ex2
  9 
 10    Description:  This is a two processor example and is similar to the previous
 11                  structured interface example (Example 1). However, in
 12                  this case the grid boxes are exactly those in the example
 13                  diagram in the struct interface chapter of the User's Manual.
 14                  (Processor 0 owns two boxes and processor 1 owns one box.)
 15                  The solver is PCG with SMG preconditioner.
 16 
 17                  We recommend viewing example 1 before viewing this
 18                  example.
 19 */
 20 
 21 #include <stdio.h>
 22 
 23 /* Struct linear solvers header */
 24 #include "HYPRE_struct_ls.h"
 25 
 26 #include "vis.c"
 27 
 28 int main (int argc, char *argv[])
 29 {
 30    int i, j;
 31 
 32    int myid, num_procs;
 33 
 34    int vis = 0;
 35 
 36    HYPRE_StructGrid     grid;
 37    HYPRE_StructStencil  stencil;
 38    HYPRE_StructMatrix   A;
 39    HYPRE_StructVector   b;
 40    HYPRE_StructVector   x;
 41    HYPRE_StructSolver   solver;
 42    HYPRE_StructSolver   precond;
 43 
 44    /* Initialize MPI */
 45    MPI_Init(&argc, &argv);
 46    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 47    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 48 
 49    if (num_procs != 2)
 50    {
 51       if (myid == 0) printf("Must run with 2 processors!\n");
 52       MPI_Finalize();
 53 
 54       return(0);
 55    }
 56 
 57    /* Parse command line */
 58    {
 59       int arg_index = 0;
 60       int print_usage = 0;
 61 
 62       while (arg_index < argc)
 63       {
 64          if ( strcmp(argv[arg_index], "-vis") == 0 )
 65          {
 66             arg_index++;
 67             vis = 1;
 68          }
 69          else if ( strcmp(argv[arg_index], "-help") == 0 )
 70          {
 71             print_usage = 1;
 72             break;
 73          }
 74          else
 75          {
 76             arg_index++;
 77          }
 78       }
 79 
 80       if ((print_usage) && (myid == 0))
 81       {
 82          printf("\n");
 83          printf("Usage: %s [<options>]\n", argv[0]);
 84          printf("\n");
 85          printf("  -vis : save the solution for GLVis visualization\n");
 86          printf("\n");
 87       }
 88 
 89       if (print_usage)
 90       {
 91          MPI_Finalize();
 92          return (0);
 93       }
 94    }
 95 
 96    /* 1. Set up a grid */
 97    {
 98       /* Create an empty 2D grid object */
 99       HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
100 
101       /* Processor 0 owns two boxes in the grid. */
102       if (myid == 0)
103       {
104          /* Add a new box to the grid */
105          {
106             int ilower[2] = {-3, 1};
107             int iupper[2] = {-1, 2};
108 
109             HYPRE_StructGridSetExtents(grid, ilower, iupper);
110          }
111 
112          /* Add a new box to the grid */
113          {
114             int ilower[2] = {0, 1};
115             int iupper[2] = {2, 4};
116 
117             HYPRE_StructGridSetExtents(grid, ilower, iupper);
118          }
119       }
120 
121       /* Processor 1 owns one box in the grid. */
122       else if (myid == 1)
123       {
124          /* Add a new box to the grid */
125          {
126             int ilower[2] = {3, 1};
127             int iupper[2] = {6, 4};
128 
129             HYPRE_StructGridSetExtents(grid, ilower, iupper);
130          }
131       }
132 
133       /* This is a collective call finalizing the grid assembly.
134          The grid is now ``ready to be used'' */
135       HYPRE_StructGridAssemble(grid);
136    }
137 
138    /* 2. Define the discretization stencil */
139    {
140       /* Create an empty 2D, 5-pt stencil object */
141       HYPRE_StructStencilCreate(2, 5, &stencil);
142 
143       /* Define the geometry of the stencil. Each represents a
144          relative offset (in the index space). */
145       {
146          int entry;
147          int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
148 
149          /* Assign each of the 5 stencil entries */
150          for (entry = 0; entry < 5; entry++)
151             HYPRE_StructStencilSetElement(stencil, entry, offsets[entry]);
152       }
153    }
154 
155    /* 3. Set up a Struct Matrix */
156    {
157       /* Create an empty matrix object */
158       HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
159 
160       /* Indicate that the matrix coefficients are ready to be set */
161       HYPRE_StructMatrixInitialize(A);
162 
163       if (myid == 0)
164       {
165          /* Set the matrix coefficients for some set of stencil entries
166             over all the gridpoints in my first box (account for boundary
167             grid points later) */
168          {
169             int ilower[2] = {-3, 1};
170             int iupper[2] = {-1, 2};
171 
172             int nentries = 5;
173             int nvalues  = 30; /* 6 grid points, each with 5 stencil entries */
174             double values[30];
175 
176             int stencil_indices[5];
177             for (j = 0; j < nentries; j++) /* label the stencil indices -
178                                               these correspond to the offsets
179                                               defined above */
180                stencil_indices[j] = j;
181 
182             for (i = 0; i < nvalues; i += nentries)
183             {
184                values[i] = 4.0;
185                for (j = 1; j < nentries; j++)
186                   values[i+j] = -1.0;
187             }
188 
189             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
190                                            stencil_indices, values);
191          }
192 
193          /* Set the matrix coefficients for some set of stencil entries
194             over the gridpoints in my second box */
195          {
196             int ilower[2] = {0, 1};
197             int iupper[2] = {2, 4};
198 
199             int nentries = 5;
200             int nvalues  = 60; /* 12 grid points, each with 5 stencil entries */
201             double values[60];
202 
203             int stencil_indices[5];
204             for (j = 0; j < nentries; j++)
205                stencil_indices[j] = j;
206 
207             for (i = 0; i < nvalues; i += nentries)
208             {
209                values[i] = 4.0;
210                for (j = 1; j < nentries; j++)
211                   values[i+j] = -1.0;
212             }
213 
214             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
215                                            stencil_indices, values);
216          }
217       }
218       else if (myid == 1)
219       {
220          /* Set the matrix coefficients for some set of stencil entries
221             over the gridpoints in my box */
222          {
223             int ilower[2] = {3, 1};
224             int iupper[2] = {6, 4};
225 
226             int nentries = 5;
227             int nvalues  = 80; /* 16 grid points, each with 5 stencil entries */
228             double values[80];
229 
230             int stencil_indices[5];
231             for (j = 0; j < nentries; j++)
232                stencil_indices[j] = j;
233 
234             for (i = 0; i < nvalues; i += nentries)
235             {
236                values[i] = 4.0;
237                for (j = 1; j < nentries; j++)
238                   values[i+j] = -1.0;
239             }
240 
241             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
242                                            stencil_indices, values);
243          }
244       }
245 
246       /* For each box, set any coefficients that reach ouside of the
247          boundary to 0 */
248       if (myid == 0)
249       {
250          int maxnvalues = 6;
251          double values[6];
252 
253          for (i = 0; i < maxnvalues; i++)
254             values[i] = 0.0;
255 
256          {
257             /* Values below our first AND second box */
258             int ilower[2] = {-3, 1};
259             int iupper[2] = { 2, 1};
260 
261             int stencil_indices[1] = {3};
262 
263             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
264                                            stencil_indices, values);
265          }
266 
267          {
268             /* Values to the left of our first box */
269             int ilower[2] = {-3, 1};
270             int iupper[2] = {-3, 2};
271 
272             int stencil_indices[1] = {1};
273 
274             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
275                                            stencil_indices, values);
276          }
277 
278          {
279             /* Values above our first box */
280             int ilower[2] = {-3, 2};
281             int iupper[2] = {-1, 2};
282 
283             int stencil_indices[1] = {4};
284 
285             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
286                                            stencil_indices, values);
287          }
288 
289          {
290             /* Values to the left of our second box (that do not border the
291                first box). */
292             int ilower[2] = { 0, 3};
293             int iupper[2] = { 0, 4};
294 
295             int stencil_indices[1] = {1};
296 
297             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
298                                            stencil_indices, values);
299          }
300 
301          {
302             /* Values above our second box */
303             int ilower[2] = { 0, 4};
304             int iupper[2] = { 2, 4};
305 
306             int stencil_indices[1] = {4};
307 
308             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
309                                            stencil_indices, values);
310          }
311       }
312       else if (myid == 1)
313       {
314          int maxnvalues = 4;
315          double values[4];
316          for (i = 0; i < maxnvalues; i++)
317             values[i] = 0.0;
318 
319          {
320             /* Values below our box */
321             int ilower[2] = { 3, 1};
322             int iupper[2] = { 6, 1};
323 
324             int stencil_indices[1] = {3};
325 
326             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
327                                            stencil_indices, values);
328          }
329 
330          {
331             /* Values to the right of our box */
332             int ilower[2] = { 6, 1};
333             int iupper[2] = { 6, 4};
334 
335             int stencil_indices[1] = {2};
336 
337             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
338                                            stencil_indices, values);
339          }
340 
341          {
342             /* Values above our box */
343             int ilower[2] = { 3, 4};
344             int iupper[2] = { 6, 4};
345 
346             int stencil_indices[1] = {4};
347 
348             HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
349                                            stencil_indices, values);
350          }
351       }
352 
353       /* This is a collective call finalizing the matrix assembly.
354          The matrix is now ``ready to be used'' */
355       HYPRE_StructMatrixAssemble(A);
356    }
357 
358    /* 4. Set up Struct Vectors for b and x */
359    {
360       /* Create an empty vector object */
361       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
362       HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
363 
364       /* Indicate that the vector coefficients are ready to be set */
365       HYPRE_StructVectorInitialize(b);
366       HYPRE_StructVectorInitialize(x);
367 
368       if (myid == 0)
369       {
370          /* Set the vector coefficients over the gridpoints in my first box */
371          {
372             int ilower[2] = {-3, 1};
373             int iupper[2] = {-1, 2};
374 
375             int nvalues = 6;  /* 6 grid points */
376             double values[6];
377 
378             for (i = 0; i < nvalues; i ++)
379                values[i] = 1.0;
380             HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
381 
382             for (i = 0; i < nvalues; i ++)
383                values[i] = 0.0;
384             HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
385          }
386 
387          /* Set the vector coefficients over the gridpoints in my second box */
388          {
389             int ilower[2] = { 0, 1};
390             int iupper[2] = { 2, 4};
391 
392             int nvalues = 12; /* 12 grid points */
393             double values[12];
394 
395             for (i = 0; i < nvalues; i ++)
396                values[i] = 1.0;
397             HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
398 
399             for (i = 0; i < nvalues; i ++)
400                values[i] = 0.0;
401             HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
402          }
403       }
404       else if (myid == 1)
405       {
406          /* Set the vector coefficients over the gridpoints in my box */
407          {
408             int ilower[2] = { 3, 1};
409             int iupper[2] = { 6, 4};
410 
411             int nvalues = 16; /* 16 grid points */
412             double values[16];
413 
414             for (i = 0; i < nvalues; i ++)
415                values[i] = 1.0;
416             HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
417 
418             for (i = 0; i < nvalues; i ++)
419                values[i] = 0.0;
420             HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
421          }
422       }
423 
424       /* This is a collective call finalizing the vector assembly.
425          The vectors are now ``ready to be used'' */
426       HYPRE_StructVectorAssemble(b);
427       HYPRE_StructVectorAssemble(x);
428    }
429 
430 
431    /* 5. Set up and use a solver (See the Reference Manual for descriptions
432       of all of the options.) */
433    {
434       /* Create an empty PCG Struct solver */
435       HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
436 
437       /* Set PCG parameters */
438       HYPRE_StructPCGSetTol(solver, 1.0e-06);
439       HYPRE_StructPCGSetPrintLevel(solver, 2);
440       HYPRE_StructPCGSetMaxIter(solver, 50);
441 
442       /* Use symmetric SMG as preconditioner */
443       HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
444       HYPRE_StructSMGSetMaxIter(precond, 1);
445       HYPRE_StructSMGSetTol(precond, 0.0);
446       HYPRE_StructSMGSetZeroGuess(precond);
447       HYPRE_StructSMGSetNumPreRelax(precond, 1);
448       HYPRE_StructSMGSetNumPostRelax(precond, 1);
449 
450       /* Set preconditioner and solve */
451       HYPRE_StructPCGSetPrecond(solver, HYPRE_StructSMGSolve,
452                                 HYPRE_StructSMGSetup, precond);
453       HYPRE_StructPCGSetup(solver, A, b, x);
454       HYPRE_StructPCGSolve(solver, A, b, x);
455    }
456 
457    /* Save the solution for GLVis visualization, see vis/glvis-ex2.sh */
458    if (vis)
459    {
460       GLVis_PrintStructGrid(grid, "vis/ex2.mesh", myid, NULL, NULL);
461       GLVis_PrintStructVector(x, "vis/ex2.sol", myid);
462       GLVis_PrintData("vis/ex2.data", myid, num_procs);
463    }
464 
465    /* Free memory */
466    HYPRE_StructGridDestroy(grid);
467    HYPRE_StructStencilDestroy(stencil);
468    HYPRE_StructMatrixDestroy(A);
469    HYPRE_StructVectorDestroy(b);
470    HYPRE_StructVectorDestroy(x);
471    HYPRE_StructPCGDestroy(solver);
472    HYPRE_StructSMGDestroy(precond);
473 
474    /* Finalize MPI */
475    MPI_Finalize();
476 
477    return (0);
478 }


syntax highlighted by Code2HTML, v. 0.9.1