download the original source code.
  1 /*
  2    Example 6
  3 
  4    Interface:    Semi-Structured interface (SStruct).  Babel version.
  5 
  6    Compile with: make ex6
  7 
  8    Sample run:   mpirun -np 2 ex6
  9 
 10    Description:  This is a two processor example and is the same problem
 11                  as is solved with the structured interface in Example 2.
 12                  (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 
 16                  This is the simplest sstruct example, and it demonstrates how
 17                  the semi-structured interface can be used for structured problems.
 18                  There is one part and one variable.  The solver is PCG with SMG
 19                  preconditioner. We use a structured solver for this example.
 20 */
 21 
 22 #include <stdio.h>
 23 
 24 #include <mpi.h>
 25 #include "bHYPRE.h"
 26 #include "HYPRE.h"
 27 
 28 int main (int argc, char *argv[])
 29 {
 30    int myid, num_procs;
 31 
 32    bHYPRE_SStructGrid     grid;
 33    bHYPRE_SStructGraph    graph;
 34    bHYPRE_SStructStencil  stencil;
 35    bHYPRE_SStructMatrix   A;
 36    bHYPRE_SStructVector   b;
 37    bHYPRE_SStructVector   x;
 38 
 39    /* We are using struct solvers for this example */
 40    bHYPRE_PCG PCGsolver;
 41    bHYPRE_StructSMG SMGprecond;
 42    bHYPRE_Solver precond;
 43 
 44    sidl_BaseInterface _ex = NULL;
 45    MPI_Comm mpicommworld = MPI_COMM_WORLD;
 46    bHYPRE_MPICommunicator mpi_comm;
 47    int object_type;
 48    int ndim = 2;
 49 
 50    /* Initialize MPI */
 51    MPI_Init(&argc, &argv);
 52    MPI_Comm_rank(MPI_COMM_WORLD, &myid);
 53    MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
 54    mpi_comm = bHYPRE_MPICommunicator_CreateC( &mpicommworld, &_ex );
 55 
 56    if (num_procs != 2)
 57    {
 58       if (myid ==0) printf("Must run with 2 processors!\n");
 59       MPI_Finalize();
 60 
 61       return(0);
 62    }
 63 
 64    /* 1. Set up the 2D grid.  This gives the index space in each part.
 65       Here we only use one part and one variable. (So the part id is 0
 66       and the variable id is 0) */
 67    {
 68       int nparts = 1;
 69       int part = 0;
 70 
 71       /* Create an empty 2D grid object */
 72       grid = bHYPRE_SStructGrid_Create( mpi_comm, ndim, nparts, &_ex );
 73 
 74       /* Set the extents of the grid - each processor sets its grid
 75          boxes.  Each part has its own relative index space numbering,
 76          but in this example all boxes belong to the same part. */
 77 
 78       /* Processor 0 owns two boxes in the grid. */
 79       if (myid == 0)
 80       {
 81          /* Add a new box to the grid */
 82          {
 83             int ilower[2] = {-3, 1};
 84             int iupper[2] = {-1, 2};
 85 
 86             bHYPRE_SStructGrid_SetExtents( grid, part, ilower, iupper, ndim, &_ex );
 87          }
 88 
 89          /* Add a new box to the grid */
 90          {
 91             int ilower[2] = {0, 1};
 92             int iupper[2] = {2, 4};
 93 
 94             bHYPRE_SStructGrid_SetExtents( grid, part, ilower, iupper, ndim, &_ex );
 95          }
 96       }
 97 
 98       /* Processor 1 owns one box in the grid. */
 99       else if (myid == 1)
100       {
101          /* Add a new box to the grid */
102          {
103             int ilower[2] = {3, 1};
104             int iupper[2] = {6, 4};
105 
106             bHYPRE_SStructGrid_SetExtents( grid, part, ilower, iupper, ndim, &_ex );
107          }
108       }
109 
110       /* Set the variable type and number of variables on each part. */
111       {
112          int i;
113          int nvars = 1;
114          int var = 0;
115          enum bHYPRE_SStructVariable__enum vartypes[1] = {bHYPRE_SStructVariable_CELL};
116 
117          for (i = 0; i< nparts; i++)
118             bHYPRE_SStructGrid_SetVariable( grid, i, var, nvars, vartypes[1], &_ex );
119          /* ... if nvars>1 we would need a number of AddVariable calls */
120       }
121 
122       /* Now the grid is ready to use */
123       bHYPRE_SStructGrid_Assemble( grid, &_ex );
124    }
125 
126    /* 2. Define the discretization stencil(s) */
127    {
128       /* Create an empty 2D, 5-pt stencil object */
129       stencil = bHYPRE_SStructStencil_Create( 2, 5, &_ex );
130 
131       /* Define the geometry of the stencil. Each represents a
132          relative offset (in the index space). */
133       {
134          int entry;
135          int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
136          int var = 0;
137 
138          /* Assign numerical values to the offsets so that we can
139             easily refer to them  - the last argument indicates the
140             variable for which we are assigning this stencil - we are
141             just using one variable in this example so it is the first one (0) */
142          for (entry = 0; entry < 5; entry++)
143             bHYPRE_SStructStencil_SetEntry( stencil, entry, offsets[entry], ndim, var, &_ex );
144       }
145    }
146 
147    /* 3. Set up the Graph  - this determines the non-zero structure
148       of the matrix and allows non-stencil relationships between the parts */
149    {
150       int var = 0;
151       int part = 0;
152 
153       /* Create the graph object */
154       graph = bHYPRE_SStructGraph_Create( mpi_comm, grid, &_ex );
155 
156       /* Now we need to tell the graph which stencil to use for each
157          variable on each part (we only have one variable and one part) */
158       bHYPRE_SStructGraph_SetStencil( graph, part, var, stencil, &_ex );
159 
160       /* Here we could establish connections between parts if we
161          had more than one part using the graph. For example, we could
162          use HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborBox() */
163 
164       /* Assemble the graph */
165       bHYPRE_SStructGraph_Assemble( graph, &_ex );
166    }
167 
168    /* 4. Set up a SStruct Matrix */
169    {
170       int i,j;
171       int part = 0;
172       int var = 0;
173 
174       /* Create the empty matrix object */
175       A = bHYPRE_SStructMatrix_Create( mpi_comm, graph, &_ex );
176 
177       /* Set the object type (by default HYPRE_SSTRUCT). This determines the
178          data structure used to store the matrix.  If you want to use unstructured
179          solvers, e.g. BoomerAMG, the object type should be HYPRE_PARCSR.
180          If the problem is purely structured (with one part), you may want to use
181          HYPRE_STRUCT to access the structured solvers. Here we have a purely
182          structured example. */
183       object_type = HYPRE_STRUCT;
184       bHYPRE_SStructMatrix_SetObjectType( A, object_type, &_ex );
185 
186       /* Get ready to set values */
187       bHYPRE_SStructMatrix_Initialize( A, &_ex );
188 
189       /* Each processor must set the stencil values for their boxes on each part.
190          In this example, we only set stencil entries and therefore use
191          HYPRE_SStructMatrixSetBoxValues.  If we need to set non-stencil entries,
192          we have to use HYPRE_SStructMatrixSetValues (shown in a later example). */
193 
194       if (myid == 0)
195       {
196          /* Set the matrix coefficients for some set of stencil entries
197             over all the gridpoints in my first box (account for boundary
198             grid points later) */
199          {
200             int ilower[2] = {-3, 1};
201             int iupper[2] = {-1, 2};
202 
203             int nentries = 5;
204             int nvalues  = 30; /* 6 grid points, each with 5 stencil entries */
205             double values[30];
206 
207             int stencil_indices[5];
208             for (j = 0; j < nentries; j++) /* label the stencil indices -
209                                               these correspond to the offsets
210                                               defined above */
211                stencil_indices[j] = j;
212 
213             for (i = 0; i < nvalues; i += nentries)
214             {
215                values[i] = 4.0;
216                for (j = 1; j < nentries; j++)
217                   values[i+j] = -1.0;
218             }
219 
220             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
221                                                var, nentries, stencil_indices,
222                                                values, nvalues, &_ex );
223          }
224 
225          /* Set the matrix coefficients for some set of stencil entries
226             over the gridpoints in my second box */
227          {
228             int ilower[2] = {0, 1};
229             int iupper[2] = {2, 4};
230 
231             int nentries = 5;
232             int nvalues  = 60; /* 12 grid points, each with 5 stencil entries */
233             double values[60];
234 
235             int stencil_indices[5];
236             for (j = 0; j < nentries; j++)
237                stencil_indices[j] = j;
238 
239             for (i = 0; i < nvalues; i += nentries)
240             {
241                values[i] = 4.0;
242                for (j = 1; j < nentries; j++)
243                   values[i+j] = -1.0;
244             }
245 
246             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
247                                                var, nentries, stencil_indices,
248                                                values, nvalues, &_ex );
249          }
250       }
251       else if (myid == 1)
252       {
253          /* Set the matrix coefficients for some set of stencil entries
254             over the gridpoints in my box */
255          {
256             int ilower[2] = {3, 1};
257             int iupper[2] = {6, 4};
258 
259             int nentries = 5;
260             int nvalues  = 80; /* 16 grid points, each with 5 stencil entries */
261             double values[80];
262 
263             int stencil_indices[5];
264             for (j = 0; j < nentries; j++)
265                stencil_indices[j] = j;
266 
267             for (i = 0; i < nvalues; i += nentries)
268             {
269                values[i] = 4.0;
270                for (j = 1; j < nentries; j++)
271                   values[i+j] = -1.0;
272             }
273 
274             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
275                                                var, nentries, stencil_indices,
276                                                values, nvalues, &_ex );
277          }
278       }
279 
280       /* For each box, set any coefficients that reach ouside of the
281          boundary to 0 */
282       if (myid == 0)
283       {
284          int maxnvalues = 6;
285          double values[6];
286 
287          for (i = 0; i < maxnvalues; i++)
288             values[i] = 0.0;
289 
290          {
291             /* Values below our first AND second box */
292             int ilower[2] = {-3, 1};
293             int iupper[2] = { 2, 1};
294             int nvalues = 6;
295             int stencil_indices[1] = {3};
296 
297             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
298                                                var, 1, stencil_indices,
299                                                values, nvalues, &_ex );
300          }
301 
302          {
303             /* Values to the left of our first box */
304             int ilower[2] = {-3, 1};
305             int iupper[2] = {-3, 2};
306             int nvalues = 2;
307             int stencil_indices[1] = {1};
308 
309             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
310                                                var, 1, stencil_indices,
311                                                values, nvalues, &_ex );
312          }
313 
314          {
315             /* Values above our first box */
316             int ilower[2] = {-3, 2};
317             int iupper[2] = {-1, 2};
318             int nvalues = 3;
319             int stencil_indices[1] = {4};
320 
321             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
322                                                var, 1, stencil_indices,
323                                                values, nvalues, &_ex );
324          }
325 
326          {
327             /* Values to the left of our second box (that do not border the
328                first box). */
329             int ilower[2] = { 0, 3};
330             int iupper[2] = { 0, 4};
331             int nvalues = 2;
332             int stencil_indices[1] = {1};
333 
334             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
335                                                var, 1, stencil_indices,
336                                                values, nvalues, &_ex );
337          }
338 
339          {
340             /* Values above our second box */
341             int ilower[2] = { 0, 4};
342             int iupper[2] = { 2, 4};
343             int nvalues = 3;
344             int stencil_indices[1] = {4};
345 
346             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
347                                                var, 1, stencil_indices,
348                                                values, nvalues, &_ex );
349          }
350       }
351       else if (myid == 1)
352       {
353          int maxnvalues = 4;
354          double values[4];
355          for (i = 0; i < maxnvalues; i++)
356             values[i] = 0.0;
357 
358          {
359             /* Values below our box */
360             int ilower[2] = { 3, 1};
361             int iupper[2] = { 6, 1};
362             int nvalues = 4;
363             int stencil_indices[1] = {3};
364 
365             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
366                                                var, 1, stencil_indices,
367                                                values, nvalues, &_ex );
368          }
369 
370          {
371             /* Values to the right of our box */
372             int ilower[2] = { 6, 1};
373             int iupper[2] = { 6, 4};
374             int nvalues = 4;
375             int stencil_indices[1] = {2};
376 
377             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
378                                                var, 1, stencil_indices,
379                                                values, nvalues, &_ex );
380          }
381 
382          {
383             /* Values above our box */
384             int ilower[2] = { 3, 4};
385             int iupper[2] = { 6, 4};
386             int nvalues = 4;
387             int stencil_indices[1] = {4};
388 
389             bHYPRE_SStructMatrix_SetBoxValues( A, part, ilower, iupper, ndim,
390                                                var, 1, stencil_indices,
391                                                values, nvalues, &_ex );
392          }
393       }
394 
395       /* This is a collective call finalizing the matrix assembly.
396          The matrix is now ``ready to be used'' */
397       bHYPRE_SStructMatrix_Assemble( A, &_ex );
398    }
399 
400 
401    /* 5. Set up SStruct Vectors for b and x */
402    {
403       int i;
404 
405       /* We have one part and one variable. */
406       int part = 0;
407       int var = 0;
408 
409       /* Create an empty vector object */
410       b = bHYPRE_SStructVector_Create( mpi_comm, grid, &_ex );
411       x = bHYPRE_SStructVector_Create( mpi_comm, grid, &_ex );
412 
413       /* As with the matrix,  set the object type for the vectors
414          to be the struct type */
415       object_type = HYPRE_STRUCT;
416       bHYPRE_SStructVector_SetObjectType( b, object_type, &_ex );
417       bHYPRE_SStructVector_SetObjectType( x, object_type, &_ex );
418 
419       /* Indicate that the vector coefficients are ready to be set */
420       bHYPRE_SStructVector_Initialize( b, &_ex );
421       bHYPRE_SStructVector_Initialize( x, &_ex );
422 
423       if (myid == 0)
424       {
425          /* Set the vector coefficients over the gridpoints in my first box */
426          {
427             int ilower[2] = {-3, 1};
428             int iupper[2] = {-1, 2};
429 
430             int nvalues = 6;  /* 6 grid points */
431             double values[6];
432 
433             for (i = 0; i < nvalues; i ++)
434                values[i] = 1.0;
435             bHYPRE_SStructVector_SetBoxValues( b, part, ilower, iupper, ndim,
436                                                var, values, nvalues, &_ex );
437 
438             for (i = 0; i < nvalues; i ++)
439                values[i] = 0.0;
440             bHYPRE_SStructVector_SetBoxValues( x, part, ilower, iupper, ndim,
441                                                var, values, nvalues, &_ex );
442          }
443 
444          /* Set the vector coefficients over the gridpoints in my second box */
445          {
446             int ilower[2] = { 0, 1};
447             int iupper[2] = { 2, 4};
448 
449             int nvalues = 12; /* 12 grid points */
450             double values[12];
451 
452             for (i = 0; i < nvalues; i ++)
453                values[i] = 1.0;
454             bHYPRE_SStructVector_SetBoxValues( b, part, ilower, iupper, ndim,
455                                                var, values, nvalues, &_ex );
456 
457             for (i = 0; i < nvalues; i ++)
458                values[i] = 0.0;
459             bHYPRE_SStructVector_SetBoxValues( x, part, ilower, iupper, ndim,
460                                                var, values, nvalues, &_ex );
461          }
462       }
463       else if (myid == 1)
464       {
465          /* Set the vector coefficients over the gridpoints in my box */
466          {
467             int ilower[2] = { 3, 1};
468             int iupper[2] = { 6, 4};
469 
470             int nvalues = 16; /* 16 grid points */
471             double values[16];
472 
473             for (i = 0; i < nvalues; i ++)
474                values[i] = 1.0;
475             bHYPRE_SStructVector_SetBoxValues( b, part, ilower, iupper, ndim,
476                                                var, values, nvalues, &_ex );
477 
478             for (i = 0; i < nvalues; i ++)
479                values[i] = 0.0;
480             bHYPRE_SStructVector_SetBoxValues( x, part, ilower, iupper, ndim,
481                                                var, values, nvalues, &_ex );
482          }
483       }
484 
485       /* This is a collective call finalizing the vector assembly.
486          The vectors are now ``ready to be used'' */
487       bHYPRE_SStructVector_Assemble( b, &_ex );
488       bHYPRE_SStructVector_Assemble( x, &_ex );
489    }
490 
491 
492    /* 6. Set up and use a solver (See the Reference Manual for descriptions
493       of all of the options.) */
494    {
495       bHYPRE_StructMatrix sA;
496       bHYPRE_Vector vb;
497       bHYPRE_Vector vx;
498       sidl_BaseInterface dummy;
499       bHYPRE_Operator opA;
500 
501       /* Because we are using a struct solver, we need to get the
502          object of the matrix and vectors to pass in to the struct solvers */
503       bHYPRE_SStructMatrix_GetObject( A, &dummy, &_ex );
504       sA = bHYPRE_StructMatrix__cast( dummy, &_ex );
505       sidl_BaseInterface_deleteRef( dummy, &_ex );
506       bHYPRE_SStructVector_GetObject( b, &dummy, &_ex );
507       vb = bHYPRE_Vector__cast( dummy, &_ex );
508       sidl_BaseInterface_deleteRef( dummy, &_ex );
509       bHYPRE_SStructVector_GetObject( x, &dummy, &_ex );
510       vx = bHYPRE_Vector__cast( dummy, &_ex );
511       sidl_BaseInterface_deleteRef( dummy, &_ex );
512 
513       /* Create an empty PCG Struct solver */
514       opA = bHYPRE_Operator__cast( sA, &_ex );
515       PCGsolver = bHYPRE_PCG_Create( mpi_comm, opA, &_ex );
516 
517       /* Set PCG parameters */
518       bHYPRE_PCG_SetDoubleParameter( PCGsolver, "Tolerance", 1.0e-6, &_ex );
519       bHYPRE_PCG_SetIntParameter( PCGsolver, "PrintLevel", 2, &_ex );
520       bHYPRE_PCG_SetIntParameter( PCGsolver, "MaxIter", 50, &_ex );
521 
522       /* Create the Struct SMG solver for use as a preconditioner */
523       SMGprecond = bHYPRE_StructSMG_Create( mpi_comm, sA, &_ex );
524 
525       /* Set SMG parameters */
526       bHYPRE_StructSMG_SetIntParameter( SMGprecond, "MaxIter", 1, &_ex );
527       bHYPRE_StructSMG_SetDoubleParameter( SMGprecond, "Tolerance", 0.0, &_ex );
528       bHYPRE_StructSMG_SetIntParameter( SMGprecond, "ZeroGuess", 1, &_ex );
529       bHYPRE_StructSMG_SetIntParameter( SMGprecond, "NumPreRelax", 1, &_ex );
530       bHYPRE_StructSMG_SetIntParameter( SMGprecond, "NumPostRelax", 1, &_ex );
531 
532       /* Set preconditioner and solve */
533       precond = bHYPRE_Solver__cast( SMGprecond, &_ex );
534       bHYPRE_PCG_SetPreconditioner( PCGsolver, precond, &_ex );
535 
536       bHYPRE_PCG_Setup( PCGsolver, vb, vx, &_ex );
537       bHYPRE_PCG_Apply( PCGsolver, vb, &vx, &_ex );
538 
539       bHYPRE_Operator_deleteRef( opA, &_ex );
540       bHYPRE_Vector_deleteRef( vx, &_ex );
541       bHYPRE_Vector_deleteRef( vb, &_ex );
542       bHYPRE_StructMatrix_deleteRef( sA, &_ex );
543    }
544 
545    /* Free memory */
546    bHYPRE_Solver_deleteRef( precond, &_ex );
547    bHYPRE_StructSMG_deleteRef( SMGprecond, &_ex );
548    bHYPRE_PCG_deleteRef( PCGsolver, &_ex );
549    bHYPRE_SStructVector_deleteRef( x, &_ex );
550    bHYPRE_SStructVector_deleteRef( b, &_ex );
551    bHYPRE_SStructMatrix_deleteRef( A, &_ex );
552    bHYPRE_SStructGraph_deleteRef( graph, &_ex );
553    bHYPRE_SStructStencil_deleteRef( stencil, &_ex );
554    bHYPRE_SStructGrid_deleteRef( grid, &_ex );
555    bHYPRE_MPICommunicator_deleteRef( mpi_comm, &_ex );
556 
557    /* Finalize MPI */
558    MPI_Finalize();
559 
560    return (0);
561 }


syntax highlighted by Code2HTML, v. 0.9.1