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