download the original source code.
1 /*
2 Example 5
3
4 Interface: Linear-Algebraic (IJ), Babel-based version
5
6 Compile with: make ex5b
7
8 Sample run: mpirun -np 4 ex5b
9
10 Description: This example solves the 2-D
11 Laplacian problem with zero boundary conditions
12 on an nxn grid. The number of unknowns is N=n^2.
13 The standard 5-point stencil is used, and we solve
14 for the interior nodes only.
15
16 This example solves the same problem as Example 3.
17 Available solvers are AMG, PCG, and PCG with AMG or
18 Parasails preconditioners.
19 */
20
21 #include <math.h>
22 #include <assert.h>
23 #include "_hypre_utilities.h"
24 #include "HYPRE_krylov.h"
25 #include "HYPRE.h"
26 #include "HYPRE_parcsr_ls.h"
27
28 /* Babel interface headers */
29 #include "bHYPRE.h"
30 #include "bHYPRE_Vector.h"
31 #include "bHYPRE_IJParCSRMatrix.h"
32 #include "bHYPRE_IJParCSRVector.h"
33 #include "bHYPRE_ParCSRDiagScale.h"
34 #include "bHYPRE_BoomerAMG.h"
35 #include "sidl_Exception.h"
36
37 int main (int argc, char *argv[])
38 {
39 int i;
40 int myid, num_procs;
41 int N, n;
42
43 int ilower, iupper;
44 int local_size, extra;
45
46 int solver_id;
47 int print_solution;
48
49 double h, h2;
50 MPI_Comm mpicommworld = MPI_COMM_WORLD;
51
52 int ierr = 0;
53 /* If this gets set to anything else, it's an error.
54 Most functions return error flags, 0 unless there's an error.
55 For clarity, they aren't checked much in this file. */
56
57 bHYPRE_MPICommunicator mpi_comm;
58 bHYPRE_IJParCSRMatrix parcsr_A;
59 bHYPRE_Operator op_A;
60 bHYPRE_IJParCSRVector par_b;
61 bHYPRE_IJParCSRVector par_x;
62 bHYPRE_Vector vec_b;
63 bHYPRE_Vector vec_x;
64
65 bHYPRE_Solver precond;
66 bHYPRE_BoomerAMG amg_solver;
67 bHYPRE_ParaSails ps_solver;
68 bHYPRE_PCG pcg_solver;
69 bHYPRE_IdentitySolver identity;
70 sidl_BaseInterface _ex = NULL;
71
72 /* Initialize MPI */
73 MPI_Init(&argc, &argv);
74 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
75 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
76 mpi_comm = bHYPRE_MPICommunicator_CreateC( &mpicommworld, &_ex );
77
78 /* Default problem parameters */
79 n = 33;
80 solver_id = 0;
81 print_solution = 0;
82
83 /* Parse command line */
84 {
85 int arg_index = 0;
86 int print_usage = 0;
87
88 while (arg_index < argc)
89 {
90 if ( strcmp(argv[arg_index], "-n") == 0 )
91 {
92 arg_index++;
93 n = atoi(argv[arg_index++]);
94 }
95 else if ( strcmp(argv[arg_index], "-solver") == 0 )
96 {
97 arg_index++;
98 solver_id = atoi(argv[arg_index++]);
99 }
100 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
101 {
102 arg_index++;
103 print_solution = 1;
104 }
105 else if ( strcmp(argv[arg_index], "-help") == 0 )
106 {
107 print_usage = 1;
108 break;
109 }
110 else
111 {
112 arg_index++;
113 }
114 }
115
116 if ((print_usage) && (myid == 0))
117 {
118 printf("\n");
119 printf("Usage: %s [<options>]\n", argv[0]);
120 printf("\n");
121 printf(" -n <n> : problem size in each direction (default: 33)\n");
122 printf(" -solver <ID> : solver ID\n");
123 printf(" 0 - AMG (default) \n");
124 printf(" 1 - AMG-PCG\n");
125 printf(" 8 - ParaSails-PCG\n");
126 printf(" 50 - PCG\n");
127 printf(" -print_solution : print the solution vector\n");
128 printf("\n");
129 }
130
131 if (print_usage)
132 {
133 MPI_Finalize();
134 return (0);
135 }
136 }
137
138 /* Preliminaries: want at least one processor per row */
139 if (n*n < num_procs) n = sqrt(num_procs) + 1;
140 N = n*n; /* global number of rows */
141 h = 1.0/(n+1); /* mesh size*/
142 h2 = h*h;
143
144 /* Each processor knows only of its own rows - the range is denoted by ilower
145 and upper. Here we partition the rows. We account for the fact that
146 N may not divide evenly by the number of processors. */
147 local_size = N/num_procs;
148 extra = N - local_size*num_procs;
149
150 ilower = local_size*myid;
151 ilower += hypre_min(myid, extra);
152
153 iupper = local_size*(myid+1);
154 iupper += hypre_min(myid+1, extra);
155 iupper = iupper - 1;
156
157 /* How many rows do I have? */
158 local_size = iupper - ilower + 1;
159
160 /* Create the matrix.
161 Note that this is a square matrix, so we indicate the row partition
162 size twice (since number of rows = number of cols) */
163 parcsr_A = bHYPRE_IJParCSRMatrix_Create( mpi_comm,
164 ilower, iupper, ilower, iupper, &_ex );
165
166 op_A = bHYPRE_Operator__cast( parcsr_A, &_ex ); /* needed later as a function argument */
167
168 /* Choose a parallel csr format storage (see the User's Manual) */
169 /* Note: Here the HYPRE interface requires a SetObjectType call.
170 I am using the bHYPRE interface in a way which does not because
171 the object type is already specified through the class name. */
172
173 /* Initialize before setting coefficients */
174 bHYPRE_IJParCSRMatrix_Initialize( parcsr_A, &_ex );
175
176 /* Now go through my local rows and set the matrix entries.
177 Each row has at most 5 entries. For example, if n=3:
178
179 A = [M -I 0; -I M -I; 0 -I M]
180 M = [4 -1 0; -1 4 -1; 0 -1 4]
181
182 Note that here we are setting one row at a time, though
183 one could set all the rows together (see the User's Manual).
184 */
185 {
186 int nnz;
187 double values[5];
188 int cols[5];
189
190 for (i = ilower; i <= iupper; i++)
191 {
192 nnz = 0;
193
194 /* The left identity block:position i-n */
195 if ((i-n)>=0)
196 {
197 cols[nnz] = i-n;
198 values[nnz] = -1.0;
199 nnz++;
200 }
201
202 /* The left -1: position i-1 */
203 if (i%n)
204 {
205 cols[nnz] = i-1;
206 values[nnz] = -1.0;
207 nnz++;
208 }
209
210 /* Set the diagonal: position i */
211 cols[nnz] = i;
212 values[nnz] = 4.0;
213 nnz++;
214
215 /* The right -1: position i+1 */
216 if ((i+1)%n)
217 {
218 cols[nnz] = i+1;
219 values[nnz] = -1.0;
220 nnz++;
221 }
222
223 /* The right identity block:position i+n */
224 if ((i+n)< N)
225 {
226 cols[nnz] = i+n;
227 values[nnz] = -1.0;
228 nnz++;
229 }
230
231 /* Set the values for row i */
232 bHYPRE_IJParCSRMatrix_SetValues( parcsr_A, 1, &nnz, &i, cols, values, 5, &_ex );
233 }
234 }
235
236 /* Assemble after setting the coefficients */
237 bHYPRE_IJParCSRMatrix_Assemble( parcsr_A, &_ex );
238
239 /* Create the rhs and solution */
240 par_b = bHYPRE_IJParCSRVector_Create( mpi_comm, ilower, iupper, &_ex );
241
242 vec_b = bHYPRE_Vector__cast( par_b, &_ex ); /* needed later for function arguments */
243
244 bHYPRE_IJParCSRVector_Initialize( par_b, &_ex );
245
246 par_x = bHYPRE_IJParCSRVector_Create( mpi_comm, ilower, iupper, &_ex );
247
248 vec_x = bHYPRE_Vector__cast( par_x, &_ex ); /* needed later for function arguments */
249
250 bHYPRE_IJParCSRVector_Initialize( par_x, &_ex );
251
252 /* Set the rhs values to h^2 and the solution to zero */
253 {
254 double *rhs_values, *x_values;
255 int *rows;
256
257 rhs_values = calloc(local_size, sizeof(double));
258 x_values = calloc(local_size, sizeof(double));
259 rows = calloc(local_size, sizeof(int));
260
261 for (i=0; i<local_size; i++)
262 {
263 rhs_values[i] = h2;
264 x_values[i] = 0.0;
265 rows[i] = ilower + i;
266 }
267
268 bHYPRE_IJParCSRVector_SetValues( par_b, local_size, rows, rhs_values, &_ex );
269 bHYPRE_IJParCSRVector_SetValues( par_x, local_size, rows, x_values, &_ex );
270
271 free(x_values);
272 free(rhs_values);
273 free(rows);
274 }
275
276 bHYPRE_IJParCSRVector_Assemble( par_b, &_ex );
277 bHYPRE_IJParCSRVector_Assemble( par_x, &_ex );
278
279 /* Choose a solver and solve the system */
280
281 /* AMG */
282 if (solver_id == 0)
283 {
284 int num_iterations;
285 double final_res_norm;
286
287 /* Create solver */
288 amg_solver = bHYPRE_BoomerAMG_Create( mpi_comm, parcsr_A, &_ex );
289
290 /* Set some parameters (See Reference Manual for more parameters) */
291 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "PrintLevel", 3, &_ex ); /* print solve info + parameters */
292 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "CoarsenType", 6, &_ex ); /* Falgout coarsening */
293 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "RelaxType", 3, &_ex ); /* G-S/Jacobi hybrid relaxation */
294 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "NumSweeps", 1, &_ex ); /* Sweeeps on each level */
295 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "MaxLevels", 20, &_ex ); /* maximum number of levels */
296 bHYPRE_BoomerAMG_SetDoubleParameter( amg_solver, "Tolerance", 1e-7, &_ex ); /* conv. tolerance */
297
298 /* Now setup and solve! */
299 bHYPRE_BoomerAMG_Setup( amg_solver, vec_b, vec_x, &_ex );
300 bHYPRE_BoomerAMG_Apply( amg_solver, vec_b, &vec_x, &_ex );
301
302 /* Run info - needed logging turned on */
303 ierr += bHYPRE_BoomerAMG_GetIntValue( amg_solver, "NumIterations", &num_iterations, &_ex );
304 bHYPRE_BoomerAMG_GetDoubleValue( amg_solver, "RelResidualNorm", &final_res_norm, &_ex );
305
306 if (myid == 0)
307 {
308 printf("\n");
309 printf("Iterations = %d\n", num_iterations);
310 printf("Final Relative Residual Norm = %e\n", final_res_norm);
311 printf("\n");
312 }
313
314 /* Destroy solver */
315 bHYPRE_BoomerAMG_deleteRef( amg_solver, &_ex );
316 }
317
318 /* PCG */
319 else if (solver_id == 50)
320 {
321 int num_iterations;
322 double final_res_norm;
323
324 /* Create solver */
325 pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A, &_ex );
326
327 /* Set some parameters (See Reference Manual for more parameters) */
328 bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000, &_ex ); /* max iterations */
329 bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7, &_ex ); /* conv. tolerance */
330 bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1, &_ex ); /* use the two norm as the stopping criteria */
331 bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2, &_ex ); /* prints out the iteration info */
332 bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1, &_ex ); /* needed to get run info later */
333
334 identity = bHYPRE_IdentitySolver_Create( mpi_comm, &_ex );
335 precond = bHYPRE_Solver__cast( identity, &_ex );
336 bHYPRE_PCG_SetPreconditioner( pcg_solver, precond, &_ex );
337
338 /* Now setup and solve! */
339 bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x, &_ex );
340 bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x, &_ex );
341
342 /* Run info - needed logging turned on */
343 bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations, &_ex );
344 bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm, &_ex );
345 if (myid == 0)
346 {
347 printf("\n");
348 printf("Iterations = %d\n", num_iterations);
349 printf("Final Relative Residual Norm = %e\n", final_res_norm);
350 printf("\n");
351 }
352
353 /* Destroy solvers */
354 bHYPRE_PCG_deleteRef( pcg_solver, &_ex );
355 bHYPRE_Solver_deleteRef( precond, &_ex );
356 }
357 /* PCG with AMG preconditioner */
358 else if (solver_id == 1)
359 {
360 int num_iterations;
361 double final_res_norm;
362
363 /* Create solver */
364 pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A, &_ex );
365
366 /* Set some parameters (See Reference Manual for more parameters) */
367 bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000, &_ex ); /* max iterations */
368 bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7, &_ex ); /* conv. tolerance */
369 bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1, &_ex ); /* use the two norm as the stopping criteria */
370 bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2, &_ex ); /* prints out the iteration info */
371 bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1, &_ex ); /* needed to get run info later */
372
373 /* Now set up the AMG preconditioner and specify any parameters */
374 amg_solver = bHYPRE_BoomerAMG_Create( mpi_comm, parcsr_A, &_ex );
375 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "PrintLevel", 1, &_ex ); /* print amg solution info*/
376 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "CoarsenType", 6, &_ex ); /* Falgout coarsening */
377 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "RelaxType", 6, &_ex ); /* Sym G-S/Jacobi hybrid relaxation */
378 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "NumSweeps", 1, &_ex ); /* Sweeeps on each level */
379 bHYPRE_BoomerAMG_SetDoubleParameter( amg_solver, "Tolerance", 0, &_ex); /* conv. tolerance */
380 bHYPRE_BoomerAMG_SetIntParameter( amg_solver, "MaxIter", 1, &_ex ); /* do only one iteration! */
381
382 /* Set the PCG preconditioner */
383 precond = bHYPRE_Solver__cast( amg_solver, &_ex );
384 bHYPRE_PCG_SetPreconditioner( pcg_solver, precond, &_ex );
385
386 /* Now setup and solve! */
387 bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x, &_ex );
388 bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x, &_ex );
389
390 /* Run info - needed logging turned on */
391 bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations, &_ex );
392 bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm, &_ex );
393 if (myid == 0)
394 {
395 printf("\n");
396 printf("Iterations = %d\n", num_iterations);
397 printf("Final Relative Residual Norm = %e\n", final_res_norm);
398 printf("\n");
399 }
400
401 /* Destroy solver and preconditioner */
402 bHYPRE_PCG_deleteRef( pcg_solver, &_ex );
403 bHYPRE_Solver_deleteRef( precond, &_ex );
404 }
405
406 /* PCG with Parasails Preconditioner */
407 else if (solver_id == 8)
408 {
409 int num_iterations;
410 double final_res_norm;
411
412 int sai_max_levels = 1;
413 double sai_threshold = 0.1;
414 double sai_filter = 0.05;
415 int sai_sym = 1;
416
417 /* Create solver */
418 pcg_solver = bHYPRE_PCG_Create( mpi_comm, op_A, &_ex );
419
420 /* Set some parameters (See Reference Manual for more parameters) */
421 bHYPRE_PCG_SetIntParameter( pcg_solver, "MaxIter", 1000, &_ex ); /* max iterations */
422 bHYPRE_PCG_SetDoubleParameter( pcg_solver, "Tolerance", 1e-7, &_ex ); /* conv. tolerance */
423 bHYPRE_PCG_SetIntParameter( pcg_solver, "TwoNorm", 1, &_ex ); /* use the two norm as the stopping criteria */
424 bHYPRE_PCG_SetIntParameter( pcg_solver, "PrintLevel", 2, &_ex ); /* prints out the iteration info */
425 bHYPRE_PCG_SetIntParameter( pcg_solver, "Logging", 1, &_ex ); /* needed to get run info later */
426
427 /* Now set up the ParaSails preconditioner and specify any parameters */
428 ps_solver = bHYPRE_ParaSails_Create( mpi_comm, parcsr_A, &_ex );
429
430 /* Set some parameters (See Reference Manual for more parameters) */
431 bHYPRE_ParaSails_SetDoubleParameter( ps_solver, "Thresh", sai_threshold, &_ex );
432 bHYPRE_ParaSails_SetIntParameter( ps_solver, "Nlevels", sai_max_levels, &_ex );
433 bHYPRE_ParaSails_SetDoubleParameter( ps_solver, "Filter", sai_filter, &_ex );
434 bHYPRE_ParaSails_SetIntParameter( ps_solver, "Sym", sai_sym, &_ex );
435 bHYPRE_ParaSails_SetIntParameter( ps_solver, "Logging", 3, &_ex );
436
437 /* Set the PCG preconditioner */
438 precond = bHYPRE_Solver__cast( ps_solver, &_ex );
439 bHYPRE_PCG_SetPreconditioner( pcg_solver, precond, &_ex );
440
441 /* Now setup and solve! */
442 bHYPRE_PCG_Setup( pcg_solver, vec_b, vec_x, &_ex );
443 bHYPRE_PCG_Apply( pcg_solver, vec_b, &vec_x, &_ex );
444
445
446 /* Run info - needed logging turned on */
447 bHYPRE_PCG_GetIntValue( pcg_solver, "NumIterations", &num_iterations, &_ex );
448 bHYPRE_PCG_GetDoubleValue( pcg_solver, "RelResidualNorm", &final_res_norm, &_ex );
449 if (myid == 0)
450 {
451 printf("\n");
452 printf("Iterations = %d\n", num_iterations);
453 printf("Final Relative Residual Norm = %e\n", final_res_norm);
454 printf("\n");
455 }
456
457 /* Destory solver and preconditioner */
458 bHYPRE_PCG_deleteRef( pcg_solver, &_ex );
459 bHYPRE_Solver_deleteRef( precond, &_ex );
460 }
461 else
462 {
463 if (myid ==0) printf("Invalid solver id specified.\n");
464 }
465
466 /* Print the solution */
467 if (print_solution)
468 bHYPRE_IJParCSRVector_Print( par_x, "ij.out.x", &_ex );
469
470 /* Clean up */
471 bHYPRE_Operator_deleteRef( op_A, &_ex );
472 bHYPRE_Vector_deleteRef( vec_b, &_ex );
473 bHYPRE_Vector_deleteRef( vec_x, &_ex );
474 bHYPRE_IJParCSRMatrix_deleteRef( parcsr_A, &_ex );
475 bHYPRE_IJParCSRVector_deleteRef( par_b, &_ex );
476 bHYPRE_IJParCSRVector_deleteRef( par_x, &_ex );
477 bHYPRE_MPICommunicator_deleteRef( mpi_comm, &_ex );
478
479 hypre_assert( ierr == 0 );
480
481 /* Finalize MPI*/
482 MPI_Finalize();
483
484 return(0);
485 }
syntax highlighted by Code2HTML, v. 0.9.1