download the original source code.
1 /*
2 Example 4
3
4 Interface: Structured interface (Struct)
5
6 Compile with: make ex4
7
8 Sample run: mpirun -np 16 ex4 -n 33 -solver 10 -K 3 -B 0 -C 1 -U0 2 -F 4
9
10 To see options: ex4 -help
11
12 Description: This example differs from the previous structured example
13 (Example 3) in that a more sophisticated stencil and
14 boundary conditions are implemented. The method illustrated
15 here to implement the boundary conditions is much more general
16 than that in the previous example. Also symmetric storage is
17 utilized when applicable.
18
19 This code solves the convection-reaction-diffusion problem
20 div (-K grad u + B u) + C u = F in the unit square with
21 boundary condition u = U0. The domain is split into N x N
22 processor grid. Thus, the given number of processors should
23 be a perfect square. Each processor has a n x n grid, with
24 nodes connected by a 5-point stencil. Note that the struct
25 interface assumes a cell-centered grid, and, therefore, the
26 nodes are not shared.
27
28 To incorporate the boundary conditions, we do the following:
29 Let x_i and x_b be the interior and boundary parts of the
30 solution vector x. If we split the matrix A as
31 A = [A_ii A_ib; A_bi A_bb],
32 then we solve
33 [A_ii 0; 0 I] [x_i ; x_b] = [b_i - A_ib u_0; u_0].
34 Note that this differs from the previous example in that we
35 are actually solving for the boundary conditions (so they
36 may not be exact as in ex3, where we only solved for the
37 interior). This approach is useful for more general types
38 of b.c.
39
40 A number of solvers are available. More information can be
41 found in the Solvers and Preconditioners chapter of the
42 User's Manual.
43
44 We recommend viewing examples 1, 2, and 3 before viewing this
45 example.
46 */
47
48 #include <math.h>
49 #include "_hypre_utilities.h"
50 #include "HYPRE_krylov.h"
51 #include "HYPRE_struct_ls.h"
52
53 #ifdef M_PI
54 #define PI M_PI
55 #else
56 #define PI 3.14159265358979
57 #endif
58
59 /* Macro to evaluate a function F in the grid point (i,j) */
60 #define Eval(F,i,j) (F( (ilower[0]+(i))*h, (ilower[1]+(j))*h ))
61 #define bcEval(F,i,j) (F( (bc_ilower[0]+(i))*h, (bc_ilower[1]+(j))*h ))
62
63 int optionK, optionB, optionC, optionU0, optionF;
64
65 /* Diffusion coefficient */
66 double K(double x, double y)
67 {
68 switch (optionK)
69 {
70 case 0:
71 return 1.0;
72 case 1:
73 return x*x+exp(y);
74 case 2:
75 if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25))
76 return 100.0;
77 else
78 return 1.0;
79 case 3:
80 if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)) < 0.0625)
81 return 10.0;
82 else
83 return 1.0;
84 default:
85 return 1.0;
86 }
87 }
88
89 /* Convection vector, first component */
90 double B1(double x, double y)
91 {
92 switch (optionB)
93 {
94 case 0:
95 return 0.0;
96 case 1:
97 return -0.1;
98 case 2:
99 return 0.25;
100 case 3:
101 return 1.0;
102 default:
103 return 0.0;
104 }
105 }
106
107 /* Convection vector, second component */
108 double B2(double x, double y)
109 {
110 switch (optionB)
111 {
112 case 0:
113 return 0.0;
114 case 1:
115 return 0.1;
116 case 2:
117 return -0.25;
118 case 3:
119 return 1.0;
120 default:
121 return 0.0;
122 }
123 }
124
125 /* Reaction coefficient */
126 double C(double x, double y)
127 {
128 switch (optionC)
129 {
130 case 0:
131 return 0.0;
132 case 1:
133 return 10.0;
134 case 2:
135 return 100.0;
136 default:
137 return 0.0;
138 }
139 }
140
141 /* Boundary condition */
142 double U0(double x, double y)
143 {
144 switch (optionU0)
145 {
146 case 0:
147 return 0.0;
148 case 1:
149 return (x+y)/100;
150 case 2:
151 return (sin(5*PI*x)+sin(5*PI*y))/1000;
152 default:
153 return 0.0;
154 }
155 }
156
157 /* Right-hand side */
158 double F(double x, double y)
159 {
160 switch (optionF)
161 {
162 case 0:
163 return 1.0;
164 case 1:
165 return 0.0;
166 case 2:
167 return 2*PI*PI*sin(PI*x)*sin(PI*y);
168 case 3:
169 if ((fabs(x-0.5) < 0.25) && (fabs(y-0.5) < 0.25))
170 return -1.0;
171 else
172 return 1.0;
173 case 4:
174 if (((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5)) < 0.0625)
175 return -1.0;
176 else
177 return 1.0;
178 default:
179 return 1.0;
180 }
181 }
182
183 int main (int argc, char *argv[])
184 {
185 int i, j, k;
186
187 int myid, num_procs;
188
189 int n, N, pi, pj;
190 double h, h2;
191 int ilower[2], iupper[2];
192
193 int solver_id;
194 int n_pre, n_post;
195 int rap, relax, skip, sym;
196 int time_index;
197
198 int num_iterations;
199 double final_res_norm;
200
201 int print_solution;
202
203 HYPRE_StructGrid grid;
204 HYPRE_StructStencil stencil;
205 HYPRE_StructMatrix A;
206 HYPRE_StructVector b;
207 HYPRE_StructVector x;
208 HYPRE_StructSolver solver;
209 HYPRE_StructSolver precond;
210
211 /* Initialize MPI */
212 MPI_Init(&argc, &argv);
213 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
214 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
215
216 /* Set default parameters */
217 n = 33;
218 optionK = 0;
219 optionB = 0;
220 optionC = 0;
221 optionU0 = 0;
222 optionF = 0;
223 solver_id = 10;
224 n_pre = 1;
225 n_post = 1;
226 rap = 0;
227 relax = 1;
228 skip = 0;
229 sym = 0;
230
231 print_solution = 0;
232
233 /* Parse command line */
234 {
235 int arg_index = 0;
236 int print_usage = 0;
237
238 while (arg_index < argc)
239 {
240 if ( strcmp(argv[arg_index], "-n") == 0 )
241 {
242 arg_index++;
243 n = atoi(argv[arg_index++]);
244 }
245 else if ( strcmp(argv[arg_index], "-K") == 0 )
246 {
247 arg_index++;
248 optionK = atoi(argv[arg_index++]);
249 }
250 else if ( strcmp(argv[arg_index], "-B") == 0 )
251 {
252 arg_index++;
253 optionB = atoi(argv[arg_index++]);
254 }
255 else if ( strcmp(argv[arg_index], "-C") == 0 )
256 {
257 arg_index++;
258 optionC = atoi(argv[arg_index++]);
259 }
260 else if ( strcmp(argv[arg_index], "-U0") == 0 )
261 {
262 arg_index++;
263 optionU0 = atoi(argv[arg_index++]);
264 }
265 else if ( strcmp(argv[arg_index], "-F") == 0 )
266 {
267 arg_index++;
268 optionF = atoi(argv[arg_index++]);
269 }
270 else if ( strcmp(argv[arg_index], "-solver") == 0 )
271 {
272 arg_index++;
273 solver_id = atoi(argv[arg_index++]);
274 }
275 else if ( strcmp(argv[arg_index], "-v") == 0 )
276 {
277 arg_index++;
278 n_pre = atoi(argv[arg_index++]);
279 n_post = atoi(argv[arg_index++]);
280 }
281 else if ( strcmp(argv[arg_index], "-rap") == 0 )
282 {
283 arg_index++;
284 rap = atoi(argv[arg_index++]);
285 }
286 else if ( strcmp(argv[arg_index], "-relax") == 0 )
287 {
288 arg_index++;
289 relax = atoi(argv[arg_index++]);
290 }
291 else if ( strcmp(argv[arg_index], "-skip") == 0 )
292 {
293 arg_index++;
294 skip = atoi(argv[arg_index++]);
295 }
296 else if ( strcmp(argv[arg_index], "-sym") == 0 )
297 {
298 arg_index++;
299 sym = atoi(argv[arg_index++]);
300 }
301 else if ( strcmp(argv[arg_index], "-print_solution") == 0 )
302 {
303 arg_index++;
304 print_solution = 1;
305 }
306 else if ( strcmp(argv[arg_index], "-help") == 0 )
307 {
308 print_usage = 1;
309 break;
310 }
311 else
312 {
313 arg_index++;
314 }
315 }
316
317 if ((print_usage) && (myid == 0))
318 {
319 printf("\n");
320 printf("Usage: %s [<options>]\n", argv[0]);
321 printf("\n");
322 printf(" -n <n> : problem size per processor (default: 8)\n");
323 printf(" -K <K> : choice for the diffusion coefficient (default: 1.0)\n");
324 printf(" -B <B> : choice for the convection vector (default: 0.0)\n");
325 printf(" -C <C> : choice for the reaction coefficient (default: 0.0)\n");
326 printf(" -U0 <U0> : choice for the boundary condition (default: 0.0)\n");
327 printf(" -F <F> : choice for the right-hand side (default: 1.0) \n");
328 printf(" -solver <ID> : solver ID\n");
329 printf(" 0 - SMG \n");
330 printf(" 1 - PFMG\n");
331 printf(" 10 - CG with SMG precond (default)\n");
332 printf(" 11 - CG with PFMG precond\n");
333 printf(" 17 - CG with 2-step Jacobi\n");
334 printf(" 18 - CG with diagonal scaling\n");
335 printf(" 19 - CG\n");
336 printf(" 30 - GMRES with SMG precond\n");
337 printf(" 31 - GMRES with PFMG precond\n");
338 printf(" 37 - GMRES with 2-step Jacobi\n");
339 printf(" 38 - GMRES with diagonal scaling\n");
340 printf(" 39 - GMRES\n");
341 printf(" -v <n_pre> <n_post> : number of pre and post relaxations\n");
342 printf(" -rap <r> : coarse grid operator type\n");
343 printf(" 0 - Galerkin (default)\n");
344 printf(" 1 - non-Galerkin ParFlow operators\n");
345 printf(" 2 - Galerkin, general operators\n");
346 printf(" -relax <r> : relaxation type\n");
347 printf(" 0 - Jacobi\n");
348 printf(" 1 - Weighted Jacobi (default)\n");
349 printf(" 2 - R/B Gauss-Seidel\n");
350 printf(" 3 - R/B Gauss-Seidel (nonsymmetric)\n");
351 printf(" -skip <s> : skip levels in PFMG (0 or 1)\n");
352 printf(" -sym <s> : symmetric storage (1) or not (0)\n");
353 printf(" -print_solution : print the solution vector\n");
354 printf("\n");
355 }
356
357 if (print_usage)
358 {
359 MPI_Finalize();
360 return (0);
361 }
362 }
363
364 /* Convection produces non-symmetric matrices */
365 if (optionB && sym)
366 optionB = 0;
367
368 /* Figure out the processor grid (N x N). The local
369 problem size is indicated by n (n x n). pi and pj
370 indicate position in the processor grid. */
371 N = sqrt(num_procs);
372 h = 1.0 / (N*n-1);
373 h2 = h*h;
374 pj = myid / N;
375 pi = myid - pj*N;
376
377 /* Define the nodes owned by the current processor (each processor's
378 piece of the global grid) */
379 ilower[0] = pi*n;
380 ilower[1] = pj*n;
381 iupper[0] = ilower[0] + n-1;
382 iupper[1] = ilower[1] + n-1;
383
384 /* 1. Set up a grid */
385 {
386 /* Create an empty 2D grid object */
387 HYPRE_StructGridCreate(MPI_COMM_WORLD, 2, &grid);
388
389 /* Add a new box to the grid */
390 HYPRE_StructGridSetExtents(grid, ilower, iupper);
391
392 /* This is a collective call finalizing the grid assembly.
393 The grid is now ``ready to be used'' */
394 HYPRE_StructGridAssemble(grid);
395 }
396
397 /* 2. Define the discretization stencil */
398 if (sym == 0)
399 {
400 /* Define the geometry of the stencil */
401 int offsets[5][2] = {{0,0}, {-1,0}, {1,0}, {0,-1}, {0,1}};
402
403 /* Create an empty 2D, 5-pt stencil object */
404 HYPRE_StructStencilCreate(2, 5, &stencil);
405
406 /* Assign stencil entries */
407 for (i = 0; i < 5; i++)
408 HYPRE_StructStencilSetElement(stencil, i, offsets[i]);
409 }
410 else /* Symmetric storage */
411 {
412 /* Define the geometry of the stencil */
413 int offsets[3][2] = {{0,0}, {1,0}, {0,1}};
414
415 /* Create an empty 2D, 3-pt stencil object */
416 HYPRE_StructStencilCreate(2, 3, &stencil);
417
418 /* Assign stencil entries */
419 for (i = 0; i < 3; i++)
420 HYPRE_StructStencilSetElement(stencil, i, offsets[i]);
421 }
422
423 /* 3. Set up Struct Vectors for b and x */
424 {
425 double *values;
426
427 /* Create an empty vector object */
428 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &b);
429 HYPRE_StructVectorCreate(MPI_COMM_WORLD, grid, &x);
430
431 /* Indicate that the vector coefficients are ready to be set */
432 HYPRE_StructVectorInitialize(b);
433 HYPRE_StructVectorInitialize(x);
434
435 values = calloc((n*n), sizeof(double));
436
437 /* Set the values of b in left-to-right, bottom-to-top order */
438 for (k = 0, j = 0; j < n; j++)
439 for (i = 0; i < n; i++, k++)
440 values[k] = h2 * Eval(F,i,j);
441 HYPRE_StructVectorSetBoxValues(b, ilower, iupper, values);
442
443 /* Set x = 0 */
444 for (i = 0; i < (n*n); i ++)
445 values[i] = 0.0;
446 HYPRE_StructVectorSetBoxValues(x, ilower, iupper, values);
447
448 free(values);
449
450 /* Assembling is postponed since the vectors will be further modified */
451 }
452
453 /* 4. Set up a Struct Matrix */
454 {
455 /* Create an empty matrix object */
456 HYPRE_StructMatrixCreate(MPI_COMM_WORLD, grid, stencil, &A);
457
458 /* Use symmetric storage? */
459 HYPRE_StructMatrixSetSymmetric(A, sym);
460
461 /* Indicate that the matrix coefficients are ready to be set */
462 HYPRE_StructMatrixInitialize(A);
463
464 /* Set the stencil values in the interior. Here we set the values
465 at every node. We will modify the boundary nodes later. */
466 if (sym == 0)
467 {
468 int stencil_indices[5] = {0, 1, 2, 3, 4}; /* labels correspond
469 to the offsets */
470 double *values;
471
472 values = calloc(5*(n*n), sizeof(double));
473
474 /* The order is left-to-right, bottom-to-top */
475 for (k = 0, j = 0; j < n; j++)
476 for (i = 0; i < n; i++, k+=5)
477 {
478 values[k+1] = - Eval(K,i-0.5,j) - Eval(B1,i-0.5,j);
479
480 values[k+2] = - Eval(K,i+0.5,j) + Eval(B1,i+0.5,j);
481
482 values[k+3] = - Eval(K,i,j-0.5) - Eval(B2,i,j-0.5);
483
484 values[k+4] = - Eval(K,i,j+0.5) + Eval(B2,i,j+0.5);
485
486 values[k] = h2 * Eval(C,i,j)
487 + Eval(K ,i-0.5,j) + Eval(K ,i+0.5,j)
488 + Eval(K ,i,j-0.5) + Eval(K ,i,j+0.5)
489 - Eval(B1,i-0.5,j) + Eval(B1,i+0.5,j)
490 - Eval(B2,i,j-0.5) + Eval(B2,i,j+0.5);
491 }
492
493 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 5,
494 stencil_indices, values);
495
496 free(values);
497 }
498 else /* Symmetric storage */
499 {
500 int stencil_indices[3] = {0, 1, 2};
501 double *values;
502
503 values = calloc(3*(n*n), sizeof(double));
504
505 /* The order is left-to-right, bottom-to-top */
506 for (k = 0, j = 0; j < n; j++)
507 for (i = 0; i < n; i++, k+=3)
508 {
509 values[k+1] = - Eval(K,i+0.5,j);
510 values[k+2] = - Eval(K,i,j+0.5);
511 values[k] = h2 * Eval(C,i,j)
512 + Eval(K,i+0.5,j) + Eval(K,i,j+0.5)
513 + Eval(K,i-0.5,j) + Eval(K,i,j-0.5);
514 }
515
516 HYPRE_StructMatrixSetBoxValues(A, ilower, iupper, 3,
517 stencil_indices, values);
518
519 free(values);
520 }
521 }
522
523 /* 5. Set the boundary conditions, while eliminating the coefficients
524 reaching ouside of the domain boundary. We must modify the matrix
525 stencil and the corresponding rhs entries. */
526 {
527 int bc_ilower[2];
528 int bc_iupper[2];
529
530 int stencil_indices[5] = {0, 1, 2, 3, 4};
531 double *values, *bvalues;
532
533 int nentries;
534 if (sym == 0)
535 nentries = 5;
536 else
537 nentries = 3;
538
539 values = calloc(nentries*n, sizeof(double));
540 bvalues = calloc(n, sizeof(double));
541
542 /* The stencil at the boundary nodes is 1-0-0-0-0. Because
543 we have I x_b = u_0; */
544 for (i = 0; i < nentries*n; i += nentries)
545 {
546 values[i] = 1.0;
547 for (j = 1; j < nentries; j++)
548 values[i+j] = 0.0;
549 }
550
551 /* Processors at y = 0 */
552 if (pj == 0)
553 {
554 bc_ilower[0] = pi*n;
555 bc_ilower[1] = pj*n;
556
557 bc_iupper[0] = bc_ilower[0] + n-1;
558 bc_iupper[1] = bc_ilower[1];
559
560 /* Modify the matrix */
561 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
562 stencil_indices, values);
563
564 /* Put the boundary conditions in b */
565 for (i = 0; i < n; i++)
566 bvalues[i] = bcEval(U0,i,0);
567
568 HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
569 }
570
571 /* Processors at y = 1 */
572 if (pj == N-1)
573 {
574 bc_ilower[0] = pi*n;
575 bc_ilower[1] = pj*n + n-1;
576
577 bc_iupper[0] = bc_ilower[0] + n-1;
578 bc_iupper[1] = bc_ilower[1];
579
580 /* Modify the matrix */
581 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
582 stencil_indices, values);
583
584 /* Put the boundary conditions in b */
585 for (i = 0; i < n; i++)
586 bvalues[i] = bcEval(U0,i,0);
587
588 HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
589 }
590
591 /* Processors at x = 0 */
592 if (pi == 0)
593 {
594 bc_ilower[0] = pi*n;
595 bc_ilower[1] = pj*n;
596
597 bc_iupper[0] = bc_ilower[0];
598 bc_iupper[1] = bc_ilower[1] + n-1;
599
600 /* Modify the matrix */
601 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
602 stencil_indices, values);
603
604 /* Put the boundary conditions in b */
605 for (j = 0; j < n; j++)
606 bvalues[j] = bcEval(U0,0,j);
607
608 HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
609 }
610
611 /* Processors at x = 1 */
612 if (pi == N-1)
613 {
614 bc_ilower[0] = pi*n + n-1;
615 bc_ilower[1] = pj*n;
616
617 bc_iupper[0] = bc_ilower[0];
618 bc_iupper[1] = bc_ilower[1] + n-1;
619
620 /* Modify the matrix */
621 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, nentries,
622 stencil_indices, values);
623
624 /* Put the boundary conditions in b */
625 for (j = 0; j < n; j++)
626 bvalues[j] = bcEval(U0,0,j);
627
628 HYPRE_StructVectorSetBoxValues(b, bc_ilower, bc_iupper, bvalues);
629 }
630
631 /* Recall that the system we are solving is:
632 [A_ii 0; 0 I] [x_i ; x_b] = [b_i - A_ib u_0; u_0].
633 This requires removing the connections between the interior
634 and boundary nodes that we have set up when we set the
635 5pt stencil at each node. We adjust for removing
636 these connections by appropriately modifying the rhs.
637 For the symm ordering scheme, just do the top and right
638 boundary */
639
640 /* Processors at y = 0, neighbors of boundary nodes */
641 if (pj == 0)
642 {
643 bc_ilower[0] = pi*n;
644 bc_ilower[1] = pj*n + 1;
645
646 bc_iupper[0] = bc_ilower[0] + n-1;
647 bc_iupper[1] = bc_ilower[1];
648
649 stencil_indices[0] = 3;
650
651 /* Modify the matrix */
652 for (i = 0; i < n; i++)
653 bvalues[i] = 0.0;
654
655 if (sym == 0)
656 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
657 stencil_indices, bvalues);
658
659 /* Eliminate the boundary conditions in b */
660 for (i = 0; i < n; i++)
661 bvalues[i] = bcEval(U0,i,-1) * (bcEval(K,i,-0.5)+bcEval(B2,i,-0.5));
662
663 if (pi == 0)
664 bvalues[0] = 0.0;
665
666 if (pi == N-1)
667 bvalues[n-1] = 0.0;
668
669 /* Note the use of AddToBoxValues (because we have already set values
670 at these nodes) */
671 HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
672 }
673
674 /* Processors at x = 0, neighbors of boundary nodes */
675 if (pi == 0)
676 {
677 bc_ilower[0] = pi*n + 1;
678 bc_ilower[1] = pj*n;
679
680 bc_iupper[0] = bc_ilower[0];
681 bc_iupper[1] = bc_ilower[1] + n-1;
682
683 stencil_indices[0] = 1;
684
685 /* Modify the matrix */
686 for (j = 0; j < n; j++)
687 bvalues[j] = 0.0;
688
689 if (sym == 0)
690 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
691 stencil_indices, bvalues);
692
693 /* Eliminate the boundary conditions in b */
694 for (j = 0; j < n; j++)
695 bvalues[j] = bcEval(U0,-1,j) * (bcEval(K,-0.5,j)+bcEval(B1,-0.5,j));
696
697 if (pj == 0)
698 bvalues[0] = 0.0;
699
700 if (pj == N-1)
701 bvalues[n-1] = 0.0;
702
703 HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
704 }
705
706 /* Processors at y = 1, neighbors of boundary nodes */
707 if (pj == N-1)
708 {
709 bc_ilower[0] = pi*n;
710 bc_ilower[1] = pj*n + (n-1) -1;
711
712 bc_iupper[0] = bc_ilower[0] + n-1;
713 bc_iupper[1] = bc_ilower[1];
714
715 if (sym == 0)
716 stencil_indices[0] = 4;
717 else
718 stencil_indices[0] = 2;
719
720 /* Modify the matrix */
721 for (i = 0; i < n; i++)
722 bvalues[i] = 0.0;
723
724 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
725 stencil_indices, bvalues);
726
727 /* Eliminate the boundary conditions in b */
728 for (i = 0; i < n; i++)
729 bvalues[i] = bcEval(U0,i,1) * (bcEval(K,i,0.5)+bcEval(B2,i,0.5));
730
731 if (pi == 0)
732 bvalues[0] = 0.0;
733
734 if (pi == N-1)
735 bvalues[n-1] = 0.0;
736
737 HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
738 }
739
740 /* Processors at x = 1, neighbors of boundary nodes */
741 if (pi == N-1)
742 {
743 bc_ilower[0] = pi*n + (n-1) - 1;
744 bc_ilower[1] = pj*n;
745
746 bc_iupper[0] = bc_ilower[0];
747 bc_iupper[1] = bc_ilower[1] + n-1;
748
749 if (sym == 0)
750 stencil_indices[0] = 2;
751 else
752 stencil_indices[0] = 1;
753
754 /* Modify the matrix */
755 for (j = 0; j < n; j++)
756 bvalues[j] = 0.0;
757
758 HYPRE_StructMatrixSetBoxValues(A, bc_ilower, bc_iupper, 1,
759 stencil_indices, bvalues);
760
761 /* Eliminate the boundary conditions in b */
762 for (j = 0; j < n; j++)
763 bvalues[j] = bcEval(U0,1,j) * (bcEval(K,0.5,j)+bcEval(B1,0.5,j));
764
765 if (pj == 0)
766 bvalues[0] = 0.0;
767
768 if (pj == N-1)
769 bvalues[n-1] = 0.0;
770
771 HYPRE_StructVectorAddToBoxValues(b, bc_ilower, bc_iupper, bvalues);
772 }
773
774 free(values);
775 free(bvalues);
776 }
777
778 /* Finalize the vector and matrix assembly */
779 HYPRE_StructMatrixAssemble(A);
780 HYPRE_StructVectorAssemble(b);
781 HYPRE_StructVectorAssemble(x);
782
783 /* 6. Set up and use a solver */
784 if (solver_id == 0) /* SMG */
785 {
786 /* Start timing */
787 time_index = hypre_InitializeTiming("SMG Setup");
788 hypre_BeginTiming(time_index);
789
790 /* Options and setup */
791 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &solver);
792 HYPRE_StructSMGSetMemoryUse(solver, 0);
793 HYPRE_StructSMGSetMaxIter(solver, 50);
794 HYPRE_StructSMGSetTol(solver, 1.0e-06);
795 HYPRE_StructSMGSetRelChange(solver, 0);
796 HYPRE_StructSMGSetNumPreRelax(solver, n_pre);
797 HYPRE_StructSMGSetNumPostRelax(solver, n_post);
798 HYPRE_StructSMGSetPrintLevel(solver, 1);
799 HYPRE_StructSMGSetLogging(solver, 1);
800 HYPRE_StructSMGSetup(solver, A, b, x);
801
802 /* Finalize current timing */
803 hypre_EndTiming(time_index);
804 hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
805 hypre_FinalizeTiming(time_index);
806 hypre_ClearTiming();
807
808 /* Start timing again */
809 time_index = hypre_InitializeTiming("SMG Solve");
810 hypre_BeginTiming(time_index);
811
812 /* Solve */
813 HYPRE_StructSMGSolve(solver, A, b, x);
814
815 /* Finalize current timing */
816 hypre_EndTiming(time_index);
817 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
818 hypre_FinalizeTiming(time_index);
819 hypre_ClearTiming();
820
821 /* Get info and release memory */
822 HYPRE_StructSMGGetNumIterations(solver, &num_iterations);
823 HYPRE_StructSMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
824 HYPRE_StructSMGDestroy(solver);
825 }
826
827 if (solver_id == 1) /* PFMG */
828 {
829 /* Start timing */
830 time_index = hypre_InitializeTiming("PFMG Setup");
831 hypre_BeginTiming(time_index);
832
833 /* Options and setup */
834 HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &solver);
835 HYPRE_StructPFMGSetMaxIter(solver, 50);
836 HYPRE_StructPFMGSetTol(solver, 1.0e-06);
837 HYPRE_StructPFMGSetRelChange(solver, 0);
838 HYPRE_StructPFMGSetRAPType(solver, rap);
839 HYPRE_StructPFMGSetRelaxType(solver, relax);
840 HYPRE_StructPFMGSetNumPreRelax(solver, n_pre);
841 HYPRE_StructPFMGSetNumPostRelax(solver, n_post);
842 HYPRE_StructPFMGSetSkipRelax(solver, skip);
843 HYPRE_StructPFMGSetPrintLevel(solver, 1);
844 HYPRE_StructPFMGSetLogging(solver, 1);
845 HYPRE_StructPFMGSetup(solver, A, b, x);
846
847 /* Finalize current timing */
848 hypre_EndTiming(time_index);
849 hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
850 hypre_FinalizeTiming(time_index);
851 hypre_ClearTiming();
852
853 /* Start timing again */
854 time_index = hypre_InitializeTiming("PFMG Solve");
855 hypre_BeginTiming(time_index);
856
857 /* Solve */
858 HYPRE_StructPFMGSolve(solver, A, b, x);
859
860 /* Finalize current timing */
861 hypre_EndTiming(time_index);
862 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
863 hypre_FinalizeTiming(time_index);
864 hypre_ClearTiming();
865
866 /* Get info and release memory */
867 HYPRE_StructPFMGGetNumIterations(solver, &num_iterations);
868 HYPRE_StructPFMGGetFinalRelativeResidualNorm(solver, &final_res_norm);
869 HYPRE_StructPFMGDestroy(solver);
870 }
871
872 /* Preconditioned CG */
873 if ((solver_id > 9) && (solver_id < 20))
874 {
875 time_index = hypre_InitializeTiming("PCG Setup");
876 hypre_BeginTiming(time_index);
877
878 HYPRE_StructPCGCreate(MPI_COMM_WORLD, &solver);
879 HYPRE_StructPCGSetMaxIter(solver, 200 );
880 HYPRE_StructPCGSetTol(solver, 1.0e-06 );
881 HYPRE_StructPCGSetTwoNorm(solver, 1 );
882 HYPRE_StructPCGSetRelChange(solver, 0 );
883 HYPRE_StructPCGSetPrintLevel(solver, 2 );
884
885 if (solver_id == 10)
886 {
887 /* use symmetric SMG as preconditioner */
888 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
889 HYPRE_StructSMGSetMemoryUse(precond, 0);
890 HYPRE_StructSMGSetMaxIter(precond, 1);
891 HYPRE_StructSMGSetTol(precond, 0.0);
892 HYPRE_StructSMGSetZeroGuess(precond);
893 HYPRE_StructSMGSetNumPreRelax(precond, n_pre);
894 HYPRE_StructSMGSetNumPostRelax(precond, n_post);
895 HYPRE_StructSMGSetPrintLevel(precond, 0);
896 HYPRE_StructSMGSetLogging(precond, 0);
897 HYPRE_StructPCGSetPrecond(solver,
898 HYPRE_StructSMGSolve,
899 HYPRE_StructSMGSetup,
900 precond);
901 }
902
903 else if (solver_id == 11)
904 {
905 /* use symmetric PFMG as preconditioner */
906 HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
907 HYPRE_StructPFMGSetMaxIter(precond, 1);
908 HYPRE_StructPFMGSetTol(precond, 0.0);
909 HYPRE_StructPFMGSetZeroGuess(precond);
910 HYPRE_StructPFMGSetRAPType(precond, rap);
911 HYPRE_StructPFMGSetRelaxType(precond, relax);
912 HYPRE_StructPFMGSetNumPreRelax(precond, n_pre);
913 HYPRE_StructPFMGSetNumPostRelax(precond, n_post);
914 HYPRE_StructPFMGSetSkipRelax(precond, skip);
915 HYPRE_StructPFMGSetPrintLevel(precond, 0);
916 HYPRE_StructPFMGSetLogging(precond, 0);
917 HYPRE_StructPCGSetPrecond(solver,
918 HYPRE_StructPFMGSolve,
919 HYPRE_StructPFMGSetup,
920 precond);
921 }
922
923 else if (solver_id == 17)
924 {
925 /* use two-step Jacobi as preconditioner */
926 HYPRE_StructJacobiCreate(MPI_COMM_WORLD, &precond);
927 HYPRE_StructJacobiSetMaxIter(precond, 2);
928 HYPRE_StructJacobiSetTol(precond, 0.0);
929 HYPRE_StructJacobiSetZeroGuess(precond);
930 HYPRE_StructPCGSetPrecond( solver,
931 HYPRE_StructJacobiSolve,
932 HYPRE_StructJacobiSetup,
933 precond);
934 }
935
936 else if (solver_id == 18)
937 {
938 /* use diagonal scaling as preconditioner */
939 precond = NULL;
940 HYPRE_StructPCGSetPrecond(solver,
941 HYPRE_StructDiagScale,
942 HYPRE_StructDiagScaleSetup,
943 precond);
944 }
945
946 /* PCG Setup */
947 HYPRE_StructPCGSetup(solver, A, b, x );
948
949 hypre_EndTiming(time_index);
950 hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
951 hypre_FinalizeTiming(time_index);
952 hypre_ClearTiming();
953
954 time_index = hypre_InitializeTiming("PCG Solve");
955 hypre_BeginTiming(time_index);
956
957 /* PCG Solve */
958 HYPRE_StructPCGSolve(solver, A, b, x);
959
960 hypre_EndTiming(time_index);
961 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
962 hypre_FinalizeTiming(time_index);
963 hypre_ClearTiming();
964
965 /* Get info and release memory */
966 HYPRE_StructPCGGetNumIterations( solver, &num_iterations );
967 HYPRE_StructPCGGetFinalRelativeResidualNorm( solver, &final_res_norm );
968 HYPRE_StructPCGDestroy(solver);
969
970 if (solver_id == 10)
971 {
972 HYPRE_StructSMGDestroy(precond);
973 }
974 else if (solver_id == 11 )
975 {
976 HYPRE_StructPFMGDestroy(precond);
977 }
978 else if (solver_id == 17)
979 {
980 HYPRE_StructJacobiDestroy(precond);
981 }
982 }
983
984 /* Preconditioned GMRES */
985 if ((solver_id > 29) && (solver_id < 40))
986 {
987 time_index = hypre_InitializeTiming("GMRES Setup");
988 hypre_BeginTiming(time_index);
989
990 HYPRE_StructGMRESCreate(MPI_COMM_WORLD, &solver);
991
992 /* Note that GMRES can be used with all the interfaces - not
993 just the struct. So here we demonstrate the
994 more generic GMRES interface functions. Since we have chosen
995 a struct solver then we must type cast to the more generic
996 HYPRE_Solver when setting options with these generic functions.
997 Note that one could declare the solver to be
998 type HYPRE_Solver, and then the casting would not be necessary.*/
999
1000 HYPRE_GMRESSetMaxIter((HYPRE_Solver) solver, 500 );
1001 HYPRE_GMRESSetKDim((HYPRE_Solver) solver,30);
1002 HYPRE_GMRESSetTol((HYPRE_Solver) solver, 1.0e-06 );
1003 HYPRE_GMRESSetPrintLevel((HYPRE_Solver) solver, 2 );
1004 HYPRE_GMRESSetLogging((HYPRE_Solver) solver, 1 );
1005
1006 if (solver_id == 30)
1007 {
1008 /* use symmetric SMG as preconditioner */
1009 HYPRE_StructSMGCreate(MPI_COMM_WORLD, &precond);
1010 HYPRE_StructSMGSetMemoryUse(precond, 0);
1011 HYPRE_StructSMGSetMaxIter(precond, 1);
1012 HYPRE_StructSMGSetTol(precond, 0.0);
1013 HYPRE_StructSMGSetZeroGuess(precond);
1014 HYPRE_StructSMGSetNumPreRelax(precond, n_pre);
1015 HYPRE_StructSMGSetNumPostRelax(precond, n_post);
1016 HYPRE_StructSMGSetPrintLevel(precond, 0);
1017 HYPRE_StructSMGSetLogging(precond, 0);
1018 HYPRE_StructGMRESSetPrecond(solver,
1019 HYPRE_StructSMGSolve,
1020 HYPRE_StructSMGSetup,
1021 precond);
1022 }
1023
1024 else if (solver_id == 31)
1025 {
1026 /* use symmetric PFMG as preconditioner */
1027 HYPRE_StructPFMGCreate(MPI_COMM_WORLD, &precond);
1028 HYPRE_StructPFMGSetMaxIter(precond, 1);
1029 HYPRE_StructPFMGSetTol(precond, 0.0);
1030 HYPRE_StructPFMGSetZeroGuess(precond);
1031 HYPRE_StructPFMGSetRAPType(precond, rap);
1032 HYPRE_StructPFMGSetRelaxType(precond, relax);
1033 HYPRE_StructPFMGSetNumPreRelax(precond, n_pre);
1034 HYPRE_StructPFMGSetNumPostRelax(precond, n_post);
1035 HYPRE_StructPFMGSetSkipRelax(precond, skip);
1036 HYPRE_StructPFMGSetPrintLevel(precond, 0);
1037 HYPRE_StructPFMGSetLogging(precond, 0);
1038 HYPRE_StructGMRESSetPrecond( solver,
1039 HYPRE_StructPFMGSolve,
1040 HYPRE_StructPFMGSetup,
1041 precond);
1042 }
1043
1044 else if (solver_id == 37)
1045 {
1046 /* use two-step Jacobi as preconditioner */
1047 HYPRE_StructJacobiCreate(MPI_COMM_WORLD, &precond);
1048 HYPRE_StructJacobiSetMaxIter(precond, 2);
1049 HYPRE_StructJacobiSetTol(precond, 0.0);
1050 HYPRE_StructJacobiSetZeroGuess(precond);
1051 HYPRE_StructGMRESSetPrecond( solver,
1052 HYPRE_StructJacobiSolve,
1053 HYPRE_StructJacobiSetup,
1054 precond);
1055 }
1056
1057 else if (solver_id == 38)
1058 {
1059 /* use diagonal scaling as preconditioner */
1060 precond = NULL;
1061 HYPRE_StructGMRESSetPrecond( solver,
1062 HYPRE_StructDiagScale,
1063 HYPRE_StructDiagScaleSetup,
1064 precond);
1065 }
1066
1067 /* GMRES Setup */
1068 HYPRE_StructGMRESSetup(solver, A, b, x );
1069
1070 hypre_EndTiming(time_index);
1071 hypre_PrintTiming("Setup phase times", MPI_COMM_WORLD);
1072 hypre_FinalizeTiming(time_index);
1073 hypre_ClearTiming();
1074
1075 time_index = hypre_InitializeTiming("GMRES Solve");
1076 hypre_BeginTiming(time_index);
1077
1078 /* GMRES Solve */
1079 HYPRE_StructGMRESSolve(solver, A, b, x);
1080
1081 hypre_EndTiming(time_index);
1082 hypre_PrintTiming("Solve phase times", MPI_COMM_WORLD);
1083 hypre_FinalizeTiming(time_index);
1084 hypre_ClearTiming();
1085
1086 /* Get info and release memory */
1087 HYPRE_StructGMRESGetNumIterations(solver, &num_iterations);
1088 HYPRE_StructGMRESGetFinalRelativeResidualNorm(solver, &final_res_norm);
1089 HYPRE_StructGMRESDestroy(solver);
1090
1091 if (solver_id == 30)
1092 {
1093 HYPRE_StructSMGDestroy(precond);
1094 }
1095 else if (solver_id == 31)
1096 {
1097 HYPRE_StructPFMGDestroy(precond);
1098 }
1099 else if (solver_id == 37)
1100 {
1101 HYPRE_StructJacobiDestroy(precond);
1102 }
1103 }
1104
1105 /* Print the solution and other info */
1106 if (print_solution)
1107 HYPRE_StructVectorPrint("struct.out.x", x, 0);
1108
1109 if (myid == 0)
1110 {
1111 printf("\n");
1112 printf("Iterations = %d\n", num_iterations);
1113 printf("Final Relative Residual Norm = %e\n", final_res_norm);
1114 printf("\n");
1115 }
1116
1117 /* Free memory */
1118 HYPRE_StructGridDestroy(grid);
1119 HYPRE_StructStencilDestroy(stencil);
1120 HYPRE_StructMatrixDestroy(A);
1121 HYPRE_StructVectorDestroy(b);
1122 HYPRE_StructVectorDestroy(x);
1123
1124 /* Finalize MPI */
1125 MPI_Finalize();
1126
1127 return (0);
1128 }
syntax highlighted by Code2HTML, v. 0.9.1