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


syntax highlighted by Code2HTML, v. 0.9.1