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