download the original source code.
1 /*
2 Example 5big
3
4 Interface: Linear-Algebraic (IJ)
5
6 Compile with: make ex5big
7
8 Sample run: mpirun -np 4 ex5big
9
10 Description: This example is a slight modification of Example 5 that
11 illustrates the 64-bit integer support in hypre needed to run
12 problems with more than 2B unknowns.
13
14 Specifically, the changes compared to Example 5 are as follows:
15
16 1) All integer arguments to HYPRE functions should be declared
17 of type HYPRE_Int.
18
19 2) Variables of type HYPRE_Int are 64-bit integers, so they
20 should be printed in the %lld format (not %d).
21
22 To enable the 64-bit integer support, you need to build hypre
23 with the --enable-bigint option of the configure script. We
24 recommend comparing this example with Example 5.
25 */
26
27 #include <math.h>
28 #include "_hypre_utilities.h"
29 #include "HYPRE_krylov.h"
30 #include "HYPRE.h"
31 #include "HYPRE_parcsr_ls.h"
32
33
34 int hypre_FlexGMRESModifyPCAMGExample(void *precond_data, int iterations,
35 double rel_residual_norm);
36
37
38 int main (int argc, char *argv[])
39 {
40 HYPRE_Int i;
41 int myid, num_procs;
42 int N, n;
43
44 HYPRE_Int ilower, iupper;
45 HYPRE_Int local_size, extra;
46
47 int solver_id;
48 int print_solution, print_system;
49
50 double h, h2;
51
52 HYPRE_IJMatrix A;
53 HYPRE_ParCSRMatrix parcsr_A;
54 HYPRE_IJVector b;
55 HYPRE_ParVector par_b;
56 HYPRE_IJVector x;
57 HYPRE_ParVector par_x;
58
59 HYPRE_Solver solver, precond;
60
61 /* Initialize MPI */
62 MPI_Init(&argc, &argv);
63 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
64 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
65
66 /* Default problem parameters */
67 n = 33;
68 solver_id = 0;
69 print_solution = 0;
70 print_system = 0;
71
72
73 /* Parse command line */
74 {
75 int arg_index = 0;
76 int print_usage = 0;
77
78 while (arg_index < argc)
79 {
80 if ( strcmp(argv[arg_index], "-n") == 0 )
81 {
82 arg_index++;
83 n = atoi(argv[arg_index++]);
84 }
85 else if ( strcmp(argv[arg_index], "-solver") == 0 )
86 {
87 arg_index++;
88 solver_id = atoi(argv[arg_index++]);
89 }
90 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
91 {
92 arg_index++;
93 print_solution = 1;
94 }
95 else if ( strcmp(argv[arg_index], "-print_system") == 0 )
96 {
97 arg_index++;
98 print_system = 1;
99 }
100
101
102 else if ( strcmp(argv[arg_index], "-help") == 0 )
103 {
104 print_usage = 1;
105 break;
106 }
107 else
108 {
109 arg_index++;
110 }
111 }
112
113 if ((print_usage) && (myid == 0))
114 {
115 printf("\n");
116 printf("Usage: %s [<options>]\n", argv[0]);
117 printf("\n");
118 printf(" -n <n> : problem size in each direction (default: 33)\n");
119 printf(" -solver <ID> : solver ID\n");
120 printf(" 0 - AMG (default) \n");
121 printf(" 1 - AMG-PCG\n");
122 printf(" 8 - ParaSails-PCG\n");
123 printf(" 50 - PCG\n");
124 printf(" 61 - AMG-FlexGMRES\n");
125 printf(" -print_solution : print the solution vector\n");
126 printf(" -print_system : print the matrix and rhs\n");
127 printf("\n");
128 }
129
130 if (print_usage)
131 {
132 MPI_Finalize();
133 return (0);
134 }
135 }
136
137 /* Preliminaries: want at least one processor per row */
138 if (n*n < num_procs) n = sqrt(num_procs) + 1;
139 N = n*n; /* global number of rows */
140 h = 1.0/(n+1); /* mesh size*/
141 h2 = h*h;
142
143 /* Each processor knows only of its own rows - the range is denoted by ilower
144 and upper. Here we partition the rows. We account for the fact that
145 N may not divide evenly by the number of processors. */
146 local_size = N/num_procs;
147 extra = N - local_size*num_procs;
148
149 ilower = local_size*myid;
150 ilower += hypre_min(myid, extra);
151
152 iupper = local_size*(myid+1);
153 iupper += hypre_min(myid+1, extra);
154 iupper = iupper - 1;
155
156 /* How many rows do I have? */
157 local_size = iupper - ilower + 1;
158
159 /* Create the matrix.
160 Note that this is a square matrix, so we indicate the row partition
161 size twice (since number of rows = number of cols) */
162 HYPRE_IJMatrixCreate(MPI_COMM_WORLD, ilower, iupper, ilower, iupper, &A);
163
164 /* Choose a parallel csr format storage (see the User's Manual) */
165 HYPRE_IJMatrixSetObjectType(A, HYPRE_PARCSR);
166
167 /* Initialize before setting coefficients */
168 HYPRE_IJMatrixInitialize(A);
169
170 /* Now go through my local rows and set the matrix entries.
171 Each row has at most 5 entries. For example, if n=3:
172
173 A = [M -I 0; -I M -I; 0 -I M]
174 M = [4 -1 0; -1 4 -1; 0 -1 4]
175
176 Note that here we are setting one row at a time, though
177 one could set all the rows together (see the User's Manual).
178 */
179 {
180 HYPRE_Int nnz;
181 double values[5];
182 HYPRE_Int cols[5];
183
184 for (i = ilower; i <= iupper; i++)
185 {
186 nnz = 0;
187
188 /* The left identity block:position i-n */
189 if ((i-n)>=0)
190 {
191 cols[nnz] = i-n;
192 values[nnz] = -1.0;
193 nnz++;
194 }
195
196 /* The left -1: position i-1 */
197 if (i%n)
198 {
199 cols[nnz] = i-1;
200 values[nnz] = -1.0;
201 nnz++;
202 }
203
204 /* Set the diagonal: position i */
205 cols[nnz] = i;
206 values[nnz] = 4.0;
207 nnz++;
208
209 /* The right -1: position i+1 */
210 if ((i+1)%n)
211 {
212 cols[nnz] = i+1;
213 values[nnz] = -1.0;
214 nnz++;
215 }
216
217 /* The right identity block:position i+n */
218 if ((i+n)< N)
219 {
220 cols[nnz] = i+n;
221 values[nnz] = -1.0;
222 nnz++;
223 }
224
225 /* Set the values for row i */
226 HYPRE_IJMatrixSetValues(A, 1, &nnz, &i, cols, values);
227 }
228 }
229
230 /* Assemble after setting the coefficients */
231 HYPRE_IJMatrixAssemble(A);
232
233 /* Note: for the testing of small problems, one may wish to read
234 in a matrix in IJ format (for the format, see the output files
235 from the -print_system option).
236 In this case, one would use the following routine:
237 HYPRE_IJMatrixRead( <filename>, MPI_COMM_WORLD,
238 HYPRE_PARCSR, &A );
239 <filename> = IJ.A.out to read in what has been printed out
240 by -print_system (processor numbers are omitted).
241 A call to HYPRE_IJMatrixRead is an *alternative* to the
242 following sequence of HYPRE_IJMatrix calls:
243 Create, SetObjectType, Initialize, SetValues, and Assemble
244 */
245
246
247 /* Get the parcsr matrix object to use */
248 HYPRE_IJMatrixGetObject(A, (void**) &parcsr_A);
249
250
251 /* Create the rhs and solution */
252 HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&b);
253 HYPRE_IJVectorSetObjectType(b, HYPRE_PARCSR);
254 HYPRE_IJVectorInitialize(b);
255
256 HYPRE_IJVectorCreate(MPI_COMM_WORLD, ilower, iupper,&x);
257 HYPRE_IJVectorSetObjectType(x, HYPRE_PARCSR);
258 HYPRE_IJVectorInitialize(x);
259
260 /* Set the rhs values to h^2 and the solution to zero */
261 {
262 double *rhs_values, *x_values;
263 HYPRE_Int *rows;
264
265 rhs_values = calloc(local_size, sizeof(double));
266 x_values = calloc(local_size, sizeof(double));
267 rows = calloc(local_size, sizeof(HYPRE_Int));
268
269 for (i=0; i<local_size; i++)
270 {
271 rhs_values[i] = h2;
272 x_values[i] = 0.0;
273 rows[i] = ilower + i;
274 }
275
276 HYPRE_IJVectorSetValues(b, local_size, rows, rhs_values);
277 HYPRE_IJVectorSetValues(x, local_size, rows, x_values);
278
279 free(x_values);
280 free(rhs_values);
281 free(rows);
282 }
283
284
285 HYPRE_IJVectorAssemble(b);
286 /* As with the matrix, for testing purposes, one may wish to read in a rhs:
287 HYPRE_IJVectorRead( <filename>, MPI_COMM_WORLD,
288 HYPRE_PARCSR, &b );
289 as an alternative to the
290 following sequence of HYPRE_IJVectors calls:
291 Create, SetObjectType, Initialize, SetValues, and Assemble
292 */
293 HYPRE_IJVectorGetObject(b, (void **) &par_b);
294
295 HYPRE_IJVectorAssemble(x);
296 HYPRE_IJVectorGetObject(x, (void **) &par_x);
297
298
299 /* Print out the system - files names will be IJ.out.A.XXXXX
300 and IJ.out.b.XXXXX, where XXXXX = processor id */
301 if (print_system)
302 {
303 HYPRE_IJMatrixPrint(A, "IJ.out.A");
304 HYPRE_IJVectorPrint(b, "IJ.out.b");
305 }
306
307
308 /* Choose a solver and solve the system */
309
310 /* AMG */
311 if (solver_id == 0)
312 {
313 HYPRE_Int num_iterations;
314 double final_res_norm;
315
316 /* Create solver */
317 HYPRE_BoomerAMGCreate(&solver);
318
319 /* Set some parameters (See Reference Manual for more parameters) */
320 HYPRE_BoomerAMGSetPrintLevel(solver, 3); /* print solve info + parameters */
321 HYPRE_BoomerAMGSetCoarsenType(solver, 6); /* Falgout coarsening */
322 HYPRE_BoomerAMGSetRelaxType(solver, 3); /* G-S/Jacobi hybrid relaxation */
323 HYPRE_BoomerAMGSetNumSweeps(solver, 1); /* Sweeeps on each level */
324 HYPRE_BoomerAMGSetMaxLevels(solver, 20); /* maximum number of levels */
325 HYPRE_BoomerAMGSetTol(solver, 1e-7); /* conv. tolerance */
326
327 /* Now setup and solve! */
328 HYPRE_BoomerAMGSetup(solver, parcsr_A, par_b, par_x);
329 HYPRE_BoomerAMGSolve(solver, parcsr_A, par_b, par_x);
330
331 /* Run info - needed logging turned on */
332 HYPRE_BoomerAMGGetNumIterations(solver, &num_iterations);
333 HYPRE_BoomerAMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
334 if (myid == 0)
335 {
336 printf("\n");
337 printf("Iterations = %lld\n", num_iterations);
338 printf("Final Relative Residual Norm = %e\n", final_res_norm);
339 printf("\n");
340 }
341
342 /* Destroy solver */
343 HYPRE_BoomerAMGDestroy(solver);
344 }
345 /* PCG */
346 else if (solver_id == 50)
347 {
348 HYPRE_Int num_iterations;
349 double final_res_norm;
350
351 /* Create solver */
352 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
353
354 /* Set some parameters (See Reference Manual for more parameters) */
355 HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
356 HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
357 HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
358 HYPRE_PCGSetPrintLevel(solver, 2); /* prints out the iteration info */
359 HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
360
361 /* Now setup and solve! */
362 HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
363 HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
364
365 /* Run info - needed logging turned on */
366 HYPRE_PCGGetNumIterations(solver, &num_iterations);
367 HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
368 if (myid == 0)
369 {
370 printf("\n");
371 printf("Iterations = %lld\n", num_iterations);
372 printf("Final Relative Residual Norm = %e\n", final_res_norm);
373 printf("\n");
374 }
375
376 /* Destroy solver */
377 HYPRE_ParCSRPCGDestroy(solver);
378 }
379 /* PCG with AMG preconditioner */
380 else if (solver_id == 1)
381 {
382 HYPRE_Int num_iterations;
383 double final_res_norm;
384
385 /* Create solver */
386 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
387
388 /* Set some parameters (See Reference Manual for more parameters) */
389 HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
390 HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
391 HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
392 HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
393 HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
394
395 /* Now set up the AMG preconditioner and specify any parameters */
396 HYPRE_BoomerAMGCreate(&precond);
397 HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
398 HYPRE_BoomerAMGSetCoarsenType(precond, 6);
399 HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
400 HYPRE_BoomerAMGSetNumSweeps(precond, 1);
401 HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
402 HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
403
404 /* Set the PCG preconditioner */
405 HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
406 (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup, precond);
407
408 /* Now setup and solve! */
409 HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
410 HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
411
412 /* Run info - needed logging turned on */
413 HYPRE_PCGGetNumIterations(solver, &num_iterations);
414 HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
415 if (myid == 0)
416 {
417 printf("\n");
418 printf("Iterations = %lld\n", num_iterations);
419 printf("Final Relative Residual Norm = %e\n", final_res_norm);
420 printf("\n");
421 }
422
423 /* Destroy solver and preconditioner */
424 HYPRE_ParCSRPCGDestroy(solver);
425 HYPRE_BoomerAMGDestroy(precond);
426 }
427 /* PCG with Parasails Preconditioner */
428 else if (solver_id == 8)
429 {
430 HYPRE_Int num_iterations;
431 double final_res_norm;
432
433 int sai_max_levels = 1;
434 double sai_threshold = 0.1;
435 double sai_filter = 0.05;
436 int sai_sym = 1;
437
438 /* Create solver */
439 HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, &solver);
440
441 /* Set some parameters (See Reference Manual for more parameters) */
442 HYPRE_PCGSetMaxIter(solver, 1000); /* max iterations */
443 HYPRE_PCGSetTol(solver, 1e-7); /* conv. tolerance */
444 HYPRE_PCGSetTwoNorm(solver, 1); /* use the two norm as the stopping criteria */
445 HYPRE_PCGSetPrintLevel(solver, 2); /* print solve info */
446 HYPRE_PCGSetLogging(solver, 1); /* needed to get run info later */
447
448 /* Now set up the ParaSails preconditioner and specify any parameters */
449 HYPRE_ParaSailsCreate(MPI_COMM_WORLD, &precond);
450
451 /* Set some parameters (See Reference Manual for more parameters) */
452 HYPRE_ParaSailsSetParams(precond, sai_threshold, sai_max_levels);
453 HYPRE_ParaSailsSetFilter(precond, sai_filter);
454 HYPRE_ParaSailsSetSym(precond, sai_sym);
455 HYPRE_ParaSailsSetLogging(precond, 3);
456
457 /* Set the PCG preconditioner */
458 HYPRE_PCGSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSolve,
459 (HYPRE_PtrToSolverFcn) HYPRE_ParaSailsSetup, precond);
460
461 /* Now setup and solve! */
462 HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b, par_x);
463 HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b, par_x);
464
465
466 /* Run info - needed logging turned on */
467 HYPRE_PCGGetNumIterations(solver, &num_iterations);
468 HYPRE_PCGGetFinalRelativeResidualNorm(solver, &final_res_norm);
469 if (myid == 0)
470 {
471 printf("\n");
472 printf("Iterations = %lld\n", num_iterations);
473 printf("Final Relative Residual Norm = %e\n", final_res_norm);
474 printf("\n");
475 }
476
477 /* Destory solver and preconditioner */
478 HYPRE_ParCSRPCGDestroy(solver);
479 HYPRE_ParaSailsDestroy(precond);
480 }
481 /* Flexible GMRES with AMG Preconditioner */
482 else if (solver_id == 61)
483 {
484 HYPRE_Int num_iterations;
485 double final_res_norm;
486 int restart = 30;
487 int modify = 1;
488
489
490 /* Create solver */
491 HYPRE_ParCSRFlexGMRESCreate(MPI_COMM_WORLD, &solver);
492
493 /* Set some parameters (See Reference Manual for more parameters) */
494 HYPRE_FlexGMRESSetKDim(solver, restart);
495 HYPRE_FlexGMRESSetMaxIter(solver, 1000); /* max iterations */
496 HYPRE_FlexGMRESSetTol(solver, 1e-7); /* conv. tolerance */
497 HYPRE_FlexGMRESSetPrintLevel(solver, 2); /* print solve info */
498 HYPRE_FlexGMRESSetLogging(solver, 1); /* needed to get run info later */
499
500
501 /* Now set up the AMG preconditioner and specify any parameters */
502 HYPRE_BoomerAMGCreate(&precond);
503 HYPRE_BoomerAMGSetPrintLevel(precond, 1); /* print amg solution info */
504 HYPRE_BoomerAMGSetCoarsenType(precond, 6);
505 HYPRE_BoomerAMGSetRelaxType(precond, 6); /* Sym G.S./Jacobi hybrid */
506 HYPRE_BoomerAMGSetNumSweeps(precond, 1);
507 HYPRE_BoomerAMGSetTol(precond, 0.0); /* conv. tolerance zero */
508 HYPRE_BoomerAMGSetMaxIter(precond, 1); /* do only one iteration! */
509
510 /* Set the FlexGMRES preconditioner */
511 HYPRE_FlexGMRESSetPrecond(solver, (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSolve,
512 (HYPRE_PtrToSolverFcn) HYPRE_BoomerAMGSetup, precond);
513
514
515 if (modify)
516 /* this is an optional call - if you don't call it, hypre_FlexGMRESModifyPCDefault
517 is used - which does nothing. Otherwise, you can define your own, similar to
518 the one used here */
519 HYPRE_FlexGMRESSetModifyPC( solver,
520 (HYPRE_PtrToModifyPCFcn) hypre_FlexGMRESModifyPCAMGExample);
521
522
523 /* Now setup and solve! */
524 HYPRE_ParCSRFlexGMRESSetup(solver, parcsr_A, par_b, par_x);
525 HYPRE_ParCSRFlexGMRESSolve(solver, parcsr_A, par_b, par_x);
526
527 /* Run info - needed logging turned on */
528 HYPRE_FlexGMRESGetNumIterations(solver, &num_iterations);
529 HYPRE_FlexGMRESGetFinalRelativeResidualNorm(solver, &final_res_norm);
530 if (myid == 0)
531 {
532 printf("\n");
533 printf("Iterations = %lld\n", num_iterations);
534 printf("Final Relative Residual Norm = %e\n", final_res_norm);
535 printf("\n");
536 }
537
538 /* Destory solver and preconditioner */
539 HYPRE_ParCSRFlexGMRESDestroy(solver);
540 HYPRE_BoomerAMGDestroy(precond);
541
542 }
543 else
544 {
545 if (myid ==0) printf("Invalid solver id specified.\n");
546 }
547
548 /* Print the solution */
549 if (print_solution)
550 HYPRE_IJVectorPrint(x, "ij.out.x");
551
552 /* Clean up */
553 HYPRE_IJMatrixDestroy(A);
554 HYPRE_IJVectorDestroy(b);
555 HYPRE_IJVectorDestroy(x);
556
557 /* Finalize MPI*/
558 MPI_Finalize();
559
560 return(0);
561 }
562
563 /*--------------------------------------------------------------------------
564 hypre_FlexGMRESModifyPCAMGExample -
565
566 This is an example (not recommended)
567 of how we can modify things about AMG that
568 affect the solve phase based on how FlexGMRES is doing...For
569 another preconditioner it may make sense to modify the tolerance..
570
571 *--------------------------------------------------------------------------*/
572
573 int hypre_FlexGMRESModifyPCAMGExample(void *precond_data, int iterations,
574 double rel_residual_norm)
575 {
576
577
578 if (rel_residual_norm > .1)
579 {
580 HYPRE_BoomerAMGSetNumSweeps(precond_data, 10);
581 }
582 else
583 {
584 HYPRE_BoomerAMGSetNumSweeps(precond_data, 1);
585 }
586
587
588 return 0;
589 }
syntax highlighted by Code2HTML, v. 0.9.1