download the original source code.
1 /*
2 Example 1
3
4 Interface: Structured interface (Struct)
5
6 Compile with: make ex1 (may need to edit HYPRE_DIR in Makefile)
7
8 Sample run: mpirun -np 2 ex1
9
10 Description: This is a two processor example. Each processor owns one
11 box in the grid. For reference, the two grid boxes are those
12 in the example diagram in the struct interface chapter
13 of the User's Manual. Note that in this example code, we have
14 used the two boxes shown in the diagram as belonging
15 to processor 0 (and given one box to each processor). The
16 solver is PCG with no preconditioner.
17
18 We recommend viewing examples 1-4 sequentially for
19 a nice overview/tutorial of the struct interface.
20 */
21
22 #include <stdio.h>
23
24 /* Struct linear solvers header */
25 #include "HYPRE_struct_ls.h"
26
27 int main (int argc, char *argv[])
28 {
29 int i, j, myid, num_procs;
30
31 HYPRE_StructGrid grid;
32 HYPRE_StructStencil stencil;
33 HYPRE_StructMatrix A;
34 HYPRE_StructVector b;
35 HYPRE_StructVector x;
36 HYPRE_StructSolver solver;
37
38 /* Initialize MPI */
39 MPI_Init(&argc, &argv);
40 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
41 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
42
43 if (num_procs != 2)
44 {
45 if (myid ==0) printf("Must run with 2 processors!\n");
46 MPI_Finalize();
47
48 return(0);
49 }
50
51 /* 1. Set up a grid. Each processor describes the piece
52 of the grid that it owns. */
53 {
54 /* Create an empty 2D grid object */
55 HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
56
57 /* Add boxes to the grid */
58 if (myid == 0)
59 {
60 int ilower[2]={-3,1}, iupper[2]={-1,2};
61 HYPRE_StructGridSetExtents(grid, ilower, iupper);
62 }
63 else if (myid == 1)
64 {
65 int ilower[2]={0,1}, iupper[2]={2,4};
66 HYPRE_StructGridSetExtents(grid, ilower, iupper);
67 }
68
69 /* This is a collective call finalizing the grid assembly.
70 The grid is now ``ready to be used'' */
71 HYPRE_StructGridAssemble(grid);
72 }
73
74 /* 2. Define the discretization stencil */
75 {
76 /* Create an empty 2D, 5-pt stencil object */
77 HYPRE_StructStencilCreate(2, 5, &stencil);
78
79 /* Define the geometry of the stencil. Each represents a
80 relative offset (in the index space). */
81 {
82 int entry;
83 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
84
85 /* Assign each of the 5 stencil entries */
86 for (entry = 0; entry < 5; entry++)
87 HYPRE_StructStencilSetElement(stencil, entry, offsets[entry]);
88 }
89 }
90
91 /* 3. Set up a Struct Matrix */
92 {
93 /* Create an empty matrix object */
94 HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
95
96 /* Indicate that the matrix coefficients are ready to be set */
97 HYPRE_StructMatrixInitialize(A);
98
99 /* Set the matrix coefficients. Each processor assigns coefficients
100 for the boxes in the grid that it owns. Note that the coefficients
101 associated with each stencil entry may vary from grid point to grid
102 point if desired. Here, we first set the same stencil entries for
103 each grid point. Then we make modifications to grid points near
104 the boundary. */
105 if (myid == 0)
106 {
107 int ilower[2]={-3,1}, iupper[2]={-1,2};
108 int stencil_indices[5] = {0,1,2,3,4}; /* labels for the stencil entries -
109 these correspond to the offsets
110 defined above */
111 int nentries = 5;
112 int nvalues = 30; /* 6 grid points, each with 5 stencil entries */
113 double values[30];
114
115 /* We have 6 grid points, each with 5 stencil entries */
116 for (i = 0; i < nvalues; i += nentries)
117 {
118 values[i] = 4.0;
119 for (j = 1; j < nentries; j++)
120 values[i+j] = -1.0;
121 }
122
123 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
124 stencil_indices, values);
125 }
126 else if (myid == 1)
127 {
128 int ilower[2]={0,1}, iupper[2]={2,4};
129 int stencil_indices[5] = {0,1,2,3,4};
130 int nentries = 5;
131 int nvalues = 60; /* 12 grid points, each with 5 stencil entries */
132 double values[60];
133
134 for (i = 0; i < nvalues; i += nentries)
135 {
136 values[i] = 4.0;
137 for (j = 1; j < nentries; j++)
138 values[i+j] = -1.0;
139 }
140
141 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
142 stencil_indices, values);
143 }
144
145 /* Set the coefficients reaching outside of the boundary to 0 */
146 if (myid == 0)
147 {
148 double values[3];
149 for (i = 0; i < 3; i++)
150 values[i] = 0.0;
151 {
152 /* values below our box */
153 int ilower[2]={-3,1}, iupper[2]={-1,1};
154 int stencil_indices[1] = {3};
155 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
156 stencil_indices, values);
157 }
158 {
159 /* values to the left of our box */
160 int ilower[2]={-3,1}, iupper[2]={-3,2};
161 int stencil_indices[1] = {1};
162 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
163 stencil_indices, values);
164 }
165 {
166 /* values above our box */
167 int ilower[2]={-3,2}, iupper[2]={-1,2};
168 int stencil_indices[1] = {4};
169 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
170 stencil_indices, values);
171 }
172 }
173 else if (myid == 1)
174 {
175 double values[4];
176 for (i = 0; i < 4; i++)
177 values[i] = 0.0;
178 {
179 /* values below our box */
180 int ilower[2]={0,1}, iupper[2]={2,1};
181 int stencil_indices[1] = {3};
182 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
183 stencil_indices, values);
184 }
185 {
186 /* values to the right of our box */
187 int ilower[2]={2,1}, iupper[2]={2,4};
188 int stencil_indices[1] = {2};
189 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
190 stencil_indices, values);
191 }
192 {
193 /* values above our box */
194 int ilower[2]={0,4}, iupper[2]={2,4};
195 int stencil_indices[1] = {4};
196 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
197 stencil_indices, values);
198 }
199 {
200 /* values to the left of our box
201 (that do not border the other box on proc. 0) */
202 int ilower[2]={0,3}, iupper[2]={0,4};
203 int stencil_indices[1] = {1};
204 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
205 stencil_indices, values);
206 }
207 }
208
209 /* This is a collective call finalizing the matrix assembly.
210 The matrix is now ``ready to be used'' */
211 HYPRE_StructMatrixAssemble(A);
212 }
213
214 /* 4. Set up Struct Vectors for b and x. Each processor sets the vectors
215 corresponding to its boxes. */
216 {
217 /* Create an empty vector object */
218 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
219 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
220
221 /* Indicate that the vector coefficients are ready to be set */
222 HYPRE_StructVectorInitialize(b);
223 HYPRE_StructVectorInitialize(x);
224
225 /* Set the vector coefficients */
226 if (myid == 0)
227 {
228 int ilower[2]={-3,1}, iupper[2]={-1,2};
229 double values[6]; /* 6 grid points */
230
231 for (i = 0; i < 6; i ++)
232 values[i] = 1.0;
233 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
234
235 for (i = 0; i < 6; i ++)
236 values[i] = 0.0;
237 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
238 }
239 else if (myid == 1)
240 {
241 int ilower[2]={0,1}, iupper[2]={2,4};
242 double values[12]; /* 12 grid points */
243
244 for (i = 0; i < 12; i ++)
245 values[i] = 1.0;
246 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
247
248 for (i = 0; i < 12; i ++)
249 values[i] = 0.0;
250 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
251 }
252
253 /* This is a collective call finalizing the vector assembly.
254 The vectors are now ``ready to be used'' */
255 HYPRE_StructVectorAssemble(b);
256 HYPRE_StructVectorAssemble(x);
257 }
258
259 /* 5. Set up and use a solver (See the Reference Manual for descriptions
260 of all of the options.) */
261 {
262 /* Create an empty PCG Struct solver */
263 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
264
265 /* Set some parameters */
266 HYPRE_StructPCGSetTol(solver, 1.0e-06); /* convergence tolerance */
267 HYPRE_StructPCGSetPrintLevel(solver, 2); /* amount of info. printed */
268
269 /* Setup and solve */
270 HYPRE_StructPCGSetup(solver, A, b, x);
271 HYPRE_StructPCGSolve(solver, A, b, x);
272 }
273
274 /* Free memory */
275 HYPRE_StructGridDestroy(grid);
276 HYPRE_StructStencilDestroy(stencil);
277 HYPRE_StructMatrixDestroy(A);
278 HYPRE_StructVectorDestroy(b);
279 HYPRE_StructVectorDestroy(x);
280 HYPRE_StructPCGDestroy(solver);
281
282 /* Finalize MPI */
283 MPI_Finalize();
284
285 return (0);
286 }
syntax highlighted by Code2HTML, v. 0.9.1