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
52 /* 1. Set up a grid. Each processor describes the piece
53 of the grid that it owns. */
54 {
55 /* Create an empty 2D grid object */
56 HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
57
58 /* Add boxes to the grid */
59 if (myid == 0)
60 {
61 int ilower[2]={-3,1}, iupper[2]={-1,2};
62 HYPRE_StructGridSetExtents(grid, ilower, iupper);
63 }
64 else if (myid == 1)
65 {
66 int ilower[2]={0,1}, iupper[2]={2,4};
67 HYPRE_StructGridSetExtents(grid, ilower, iupper);
68 }
69
70 /* This is a collective call finalizing the grid assembly.
71 The grid is now ``ready to be used'' */
72 HYPRE_StructGridAssemble(grid);
73 }
74
75 /* 2. Define the discretization stencil */
76 {
77 /* Create an empty 2D, 5-pt stencil object */
78 HYPRE_StructStencilCreate(2, 5, &stencil);
79
80 /* Define the geometry of the stencil. Each represents a
81 relative offset (in the index space). */
82 {
83 int entry;
84 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
85
86 /* Assign each of the 5 stencil entries */
87 for (entry = 0; entry < 5; entry++)
88 HYPRE_StructStencilSetElement(stencil, entry, offsets[entry]);
89 }
90 }
91
92 /* 3. Set up a Struct Matrix */
93 {
94 /* Create an empty matrix object */
95 HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
96
97 /* Indicate that the matrix coefficients are ready to be set */
98 HYPRE_StructMatrixInitialize(A);
99
100 /* Set the matrix coefficients. Each processor assigns coefficients
101 for the boxes in the grid that it owns. Note that the coefficients
102 associated with each stencil entry may vary from grid point to grid
103 point if desired. Here, we first set the same stencil entries for
104 each grid point. Then we make modifications to grid points near
105 the boundary. */
106 if (myid == 0)
107 {
108 int ilower[2]={-3,1}, iupper[2]={-1,2};
109 int stencil_indices[5] = {0,1,2,3,4}; /* labels for the stencil entries -
110 these correspond to the offsets
111 defined above */
112 int nentries = 5;
113 int nvalues = 30; /* 6 grid points, each with 5 stencil entries */
114 double values[30];
115
116 /* We have 6 grid points, each with 5 stencil entries */
117 for (i = 0; i < nvalues; i += nentries)
118 {
119 values[i] = 4.0;
120 for (j = 1; j < nentries; j++)
121 values[i+j] = -1.0;
122 }
123
124 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
125 stencil_indices, values);
126 }
127 else if (myid == 1)
128 {
129 int ilower[2]={0,1}, iupper[2]={2,4};
130 int stencil_indices[5] = {0,1,2,3,4};
131 int nentries = 5;
132 int nvalues = 60; /* 12 grid points, each with 5 stencil entries */
133 double values[60];
134
135 for (i = 0; i < nvalues; i += nentries)
136 {
137 values[i] = 4.0;
138 for (j = 1; j < nentries; j++)
139 values[i+j] = -1.0;
140 }
141
142 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, nentries,
143 stencil_indices, values);
144 }
145
146 /* Set the coefficients reaching outside of the boundary to 0 */
147 if (myid == 0)
148 {
149 double values[3];
150 for (i = 0; i < 3; i++)
151 values[i] = 0.0;
152 {
153 /* values below our box */
154 int ilower[2]={-3,1}, iupper[2]={-1,1};
155 int stencil_indices[1] = {3};
156 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
157 stencil_indices, values);
158 }
159 {
160 /* values to the left of our box */
161 int ilower[2]={-3,1}, iupper[2]={-3,2};
162 int stencil_indices[1] = {1};
163 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
164 stencil_indices, values);
165 }
166 {
167 /* values above our box */
168 int ilower[2]={-3,2}, iupper[2]={-1,2};
169 int stencil_indices[1] = {4};
170 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
171 stencil_indices, values);
172 }
173 }
174 else if (myid == 1)
175 {
176 double values[4];
177 for (i = 0; i < 4; i++)
178 values[i] = 0.0;
179 {
180 /* values below our box */
181 int ilower[2]={0,1}, iupper[2]={2,1};
182 int stencil_indices[1] = {3};
183 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
184 stencil_indices, values);
185 }
186 {
187 /* values to the right of our box
188 (that do not border the other box on proc. 0) */
189 int ilower[2]={2,1}, iupper[2]={2,4};
190 int stencil_indices[1] = {2};
191 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
192 stencil_indices, values);
193 }
194 {
195 /* values above our box */
196 int ilower[2]={0,4}, iupper[2]={2,4};
197 int stencil_indices[1] = {4};
198 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
199 stencil_indices, values);
200 }
201 {
202 /* values to the left of our box */
203 int ilower[2]={0,3}, iupper[2]={0,4};
204 int stencil_indices[1] = {1};
205 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 1,
206 stencil_indices, values);
207 }
208 }
209
210 /* This is a collective call finalizing the matrix assembly.
211 The matrix is now ``ready to be used'' */
212 HYPRE_StructMatrixAssemble(A);
213 }
214
215 /* 4. Set up Struct Vectors for b and x. Each processor sets the vectors
216 corresponding to its boxes. */
217 {
218 /* Create an empty vector object */
219 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
220 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
221
222 /* Indicate that the vector coefficients are ready to be set */
223 HYPRE_StructVectorInitialize(b);
224 HYPRE_StructVectorInitialize(x);
225
226 /* Set the vector coefficients */
227 if (myid == 0)
228 {
229 int ilower[2]={-3,1}, iupper[2]={-1,2};
230 double values[6]; /* 6 grid points */
231
232 for (i = 0; i < 6; i ++)
233 values[i] = 1.0;
234 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
235
236 for (i = 0; i < 6; i ++)
237 values[i] = 0.0;
238 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
239 }
240 else if (myid == 1)
241 {
242 int ilower[2]={0,1}, iupper[2]={2,4};
243 double values[12]; /* 12 grid points */
244
245 for (i = 0; i < 12; i ++)
246 values[i] = 1.0;
247 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
248
249 for (i = 0; i < 12; i ++)
250 values[i] = 0.0;
251 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
252 }
253
254 /* This is a collective call finalizing the vector assembly.
255 The vectors are now ``ready to be used'' */
256 HYPRE_StructVectorAssemble(b);
257 HYPRE_StructVectorAssemble(x);
258 }
259
260 /* 5. Set up and use a solver (See the Reference Manual for descriptions
261 of all of the options.) */
262 {
263 /* Create an empty PCG Struct solver */
264 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
265
266 /* Set some parameters */
267 HYPRE_StructPCGSetTol(solver, 1.0e-06); /* convergence tolerance */
268 HYPRE_StructPCGSetPrintLevel(solver, 2); /* amount of info. printed */
269
270 /* Setup and solve */
271 HYPRE_StructPCGSetup(solver, A, b, x);
272 HYPRE_StructPCGSolve(solver, A, b, x);
273 }
274
275 /* Free memory */
276 HYPRE_StructGridDestroy(grid);
277 HYPRE_StructStencilDestroy(stencil);
278 HYPRE_StructMatrixDestroy(A);
279 HYPRE_StructVectorDestroy(b);
280 HYPRE_StructVectorDestroy(x);
281 HYPRE_StructPCGDestroy(solver);
282
283 /* Finalize MPI */
284 MPI_Finalize();
285
286 return (0);
287 }
syntax highlighted by Code2HTML, v. 0.9.1