Hi all,
I'm working on a CPR-AMG Two-Stage preconditioner implemented as multiplicative PCComposite with outer FGMRES, where the first PC is Hypre AMG (PCGalerkin + KSPRichardson + PCHYPRE) and the second stage is Block Jacobi with LU. The pde's describe two-phase subsurface flow, and I kept the problem small at 8000 x 8000 dofs. The first stage is hard-wired because of the PCGalerkin part and the second stage Block Jacobi is configured via command line (with pflotran prefix flow_): -flow_sub_1_pc_type bjacobi \ -flow_sub_1_sub_pc_type lu \ With this configuration I see occasionally that Hypre struggles to converge fast: Step 16 0 2r: 3.95E-03 2x: 0.00E+00 2u: 0.00E+00 ir: 2.53E-03 iu: 0.00E+00 rsn: 0 Residual norms for flow_ solve. 0 KSP unpreconditioned resid norm 3.945216988332e-03 true resid norm 3.945216988332e-03 ||r(i)||/||b|| 1.000000000000e+00 Residual norms for flow_sub_0_galerkin_ solve. 0 KSP preconditioned resid norm 2.495457360562e+08 true resid norm 9.213492769259e-01 ||r(i)||/||b|| 1.000000000000e+00 1 KSP preconditioned resid norm 3.900401635809e+07 true resid norm 1.211813734614e-01 ||r(i)||/||b|| 1.315259874797e-01 2 KSP preconditioned resid norm 7.264015944695e+06 true resid norm 2.127154159346e-02 ||r(i)||/||b|| 2.308738078618e-02 3 KSP preconditioned resid norm 1.523934370189e+06 true resid norm 4.507204888834e-03 ||r(i)||/||b|| 4.891961172285e-03 4 KSP preconditioned resid norm 3.456355485206e+05 true resid norm 1.017486337883e-03 ||r(i)||/||b|| 1.104343774250e-03 5 KSP preconditioned resid norm 8.215494701640e+04 true resid norm 2.386758602821e-04 ||r(i)||/||b|| 2.590503582729e-04 6 KSP preconditioned resid norm 2.006221595869e+04 true resid norm 5.806707975375e-05 ||r(i)||/||b|| 6.302395975986e-05 7 KSP preconditioned resid norm 4.975749682114e+03 true resid norm 1.457831681999e-05 ||r(i)||/||b|| 1.582279075383e-05 8 KSP preconditioned resid norm 1.245359749620e+03 true resid norm 3.746721600730e-06 ||r(i)||/||b|| 4.066559441204e-06 9 KSP preconditioned resid norm 3.134373137075e+02 true resid norm 9.784665277082e-07 ||r(i)||/||b|| 1.061993048904e-06 10 KSP preconditioned resid norm 7.917076489741e+01 true resid norm 2.582765351245e-07 ||r(i)||/||b|| 2.803242392356e-07 11 KSP preconditioned resid norm 2.004702594193e+01 true resid norm 6.867609287185e-08 ||r(i)||/||b|| 7.453860831257e-08 1 KSP unpreconditioned resid norm 3.022346103074e-11 true resid norm 3.022346103592e-11 ||r(i)||/||b|| 7.660785484121e-09 1 2r: 2.87E-04 2x: 3.70E+09 2u: 3.36E+02 ir: 1.67E-04 iu: 2.19E+01 rsn: stol Nonlinear flow_ solve converged due to CONVERGED_SNORM_RELATIVE iterations 1 Step 17 0 2r: 3.85E-03 2x: 0.00E+00 2u: 0.00E+00 ir: 2.69E-03 iu: 0.00E+00 rsn: 0 Residual norms for flow_ solve. 0 KSP unpreconditioned resid norm 3.846677237838e-03 true resid norm 3.846677237838e-03 ||r(i)||/||b|| 1.000000000000e+00 Residual norms for flow_sub_0_galerkin_ solve. 0 KSP preconditioned resid norm 8.359592959751e+07 true resid norm 8.919381920269e-01 ||r(i)||/||b|| 1.000000000000e+00 1 KSP preconditioned resid norm 2.046474217608e+07 true resid norm 1.356172589724e+00 ||r(i)||/||b|| 1.520478214574e+00 2 KSP preconditioned resid norm 5.534610937223e+06 true resid norm 1.361527715124e+00 ||r(i)||/||b|| 1.526482134406e+00 3 KSP preconditioned resid norm 1.642592089665e+06 true resid norm 1.359990274368e+00 ||r(i)||/||b|| 1.524758426677e+00 4 KSP preconditioned resid norm 6.869446528993e+05 true resid norm 1.357740694885e+00 ||r(i)||/||b|| 1.522236301823e+00 5 KSP preconditioned resid norm 5.245968674991e+05 true resid norm 1.355364470917e+00 ||r(i)||/||b|| 1.519572189007e+00 6 KSP preconditioned resid norm 5.042030663187e+05 true resid norm 1.352962944308e+00 ||r(i)||/||b|| 1.516879708036e+00 7 KSP preconditioned resid norm 5.007302249221e+05 true resid norm 1.350558656878e+00 ||r(i)||/||b|| 1.514184131760e+00 8 KSP preconditioned resid norm 4.994105316949e+05 true resid norm 1.348156961110e+00 ||r(i)||/||b|| 1.511491461137e+00 9 KSP preconditioned resid norm 4.984373051647e+05 true resid norm 1.345759135434e+00 ||r(i)||/||b|| 1.508803129481e+00 10 KSP preconditioned resid norm 4.975323739321e+05 true resid norm 1.343365479502e+00 ||r(i)||/||b|| 1.506119472750e+00 11 KSP preconditioned resid norm 4.966432959339e+05 true resid norm 1.340976058673e+00 ||r(i)||/||b|| 1.503440564224e+00 [...] 193 KSP preconditioned resid norm 3.591931201817e+05 true resid norm 9.698521332569e-01 ||r(i)||/||b|| 1.087353520599e+00 194 KSP preconditioned resid norm 3.585542278288e+05 true resid norm 9.681270691497e-01 ||r(i)||/||b|| 1.085419458213e+00 195 KSP preconditioned resid norm 3.579164717745e+05 true resid norm 9.664050733935e-01 ||r(i)||/||b|| 1.083488835922e+00 196 KSP preconditioned resid norm 3.572798501551e+05 true resid norm 9.646861405301e-01 ||r(i)||/||b|| 1.081561647605e+00 197 KSP preconditioned resid norm 3.566443608646e+05 true resid norm 9.629702651108e-01 ||r(i)||/||b|| 1.079637887153e+00 198 KSP preconditioned resid norm 3.560100018703e+05 true resid norm 9.612574416991e-01 ||r(i)||/||b|| 1.077717548471e+00 199 KSP preconditioned resid norm 3.553767713002e+05 true resid norm 9.595476648643e-01 ||r(i)||/||b|| 1.075800625471e+00 200 KSP preconditioned resid norm 3.547446669197e+05 true resid norm 9.578409291897e-01 ||r(i)||/||b|| 1.073887112080e+00 1 KSP unpreconditioned resid norm 3.816569407795e-11 true resid norm 3.816569407353e-11 ||r(i)||/||b|| 9.921730291825e-09 1 2r: 2.74E-02 2x: 3.70E+09 2u: 1.23E+02 ir: 1.99E-02 iu: 2.71E+01 rsn: stol Nonlinear flow_ solve converged due to CONVERGED_SNORM_RELATIVE iterations 1 Step 18 0 2r: 2.73E-02 2x: 0.00E+00 2u: 0.00E+00 ir: 2.02E-02 iu: 0.00E+00 rsn: 0 Residual norms for flow_ solve. 0 KSP unpreconditioned resid norm 2.734891161446e-02 true resid norm 2.734891161446e-02 ||r(i)||/||b|| 1.000000000000e+00 Residual norms for flow_sub_0_galerkin_ solve. 0 KSP preconditioned resid norm 3.550345478098e+07 true resid norm 1.048585361984e+00 ||r(i)||/||b|| 1.000000000000e+00 1 KSP preconditioned resid norm 6.139218831613e+06 true resid norm 1.797822962324e-02 ||r(i)||/||b|| 1.714522276873e-02 2 KSP preconditioned resid norm 1.301871956838e+06 true resid norm 3.761355992926e-03 ||r(i)||/||b|| 3.587076578878e-03 3 KSP preconditioned resid norm 3.070518418113e+05 true resid norm 9.283056182563e-04 ||r(i)||/||b|| 8.852933217570e-04 4 KSP preconditioned resid norm 7.639640178912e+04 true resid norm 2.348078927331e-04 ||r(i)||/||b|| 2.239282572941e-04 5 KSP preconditioned resid norm 1.953032767966e+04 true resid norm 5.930230662989e-05 ||r(i)||/||b|| 5.655458180124e-05 6 KSP preconditioned resid norm 5.066937883132e+03 true resid norm 1.497534370201e-05 ||r(i)||/||b|| 1.428147315892e-05 7 KSP preconditioned resid norm 1.326441080568e+03 true resid norm 3.793872760594e-06 ||r(i)||/||b|| 3.618086708188e-06 8 KSP preconditioned resid norm 3.494353490063e+02 true resid norm 9.659536247849e-07 ||r(i)||/||b|| 9.211969380896e-07 9 KSP preconditioned resid norm 9.251497983280e+01 true resid norm 2.472922526467e-07 ||r(i)||/||b|| 2.358341644011e-07 10 KSP preconditioned resid norm 2.459917675189e+01 true resid norm 6.364691902290e-08 ||r(i)||/||b|| 6.069789006257e-08 11 KSP preconditioned resid norm 6.566117552226e+00 true resid norm 1.646205416458e-08 ||r(i)||/||b|| 1.569929808426e-08 12 KSP preconditioned resid norm 1.758927386308e+00 true resid norm 4.277033775892e-09 ||r(i)||/||b|| 4.078860845245e-09 1 KSP unpreconditioned resid norm 2.831146511164e-10 true resid norm 2.831146511142e-10 ||r(i)||/||b|| 1.035195312725e-08 1 2r: 1.31E-02 2x: 3.70E+09 2u: 3.66E+02 ir: 9.77E-03 iu: 6.03E+01 rsn: stol Nonlinear flow_ solve converged due to CONVERGED_SNORM_RELATIVE iterations 1 SNES_view: SNES Object: (flow_) 2 MPI processes type: newtonls maximum iterations=8, maximum function evaluations=10000 tolerances: relative=1e-05, absolute=1e-05, solution=1e-05 total number of linear solver iterations=1 total number of function evaluations=2 norm schedule ALWAYS SNESLineSearch Object: (flow_) 2 MPI processes type: basic maxstep=1.000000e+08, minlambda=1.000000e-05 tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08 maximum iterations=40 using user-defined precheck step KSP Object: (flow_) 2 MPI processes type: fgmres GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement GMRES: happy breakdown tolerance 1e-30 maximum iterations=200, initial guess is zero tolerances: relative=1e-07, absolute=1e-50, divergence=10000. right preconditioning using UNPRECONDITIONED norm type for convergence test PC Object: (flow_) 2 MPI processes type: composite Composite PC type - MULTIPLICATIVE PCs on composite preconditioner follow --------------------------------- PC Object: (flow_sub_0_) 2 MPI processes type: galerkin Galerkin PC KSP on Galerkin follow --------------------------------- KSP Object: (flow_sub_0_galerkin_) 2 MPI processes type: richardson Richardson: damping factor=1. maximum iterations=200, initial guess is zero tolerances: relative=1e-07, absolute=1e-50, divergence=10000. left preconditioning using PRECONDITIONED norm type for convergence test PC Object: (flow_sub_0_galerkin_) 2 MPI processes type: hypre HYPRE BoomerAMG preconditioning HYPRE BoomerAMG: Cycle type V HYPRE BoomerAMG: Maximum number of levels 25 HYPRE BoomerAMG: Maximum number of iterations PER hypre call 1 HYPRE BoomerAMG: Convergence tolerance PER hypre call 0. HYPRE BoomerAMG: Threshold for strong coupling 0.25 HYPRE BoomerAMG: Interpolation truncation factor 0. HYPRE BoomerAMG: Interpolation: max elements per row 0 HYPRE BoomerAMG: Number of levels of aggressive coarsening 0 HYPRE BoomerAMG: Number of paths for aggressive coarsening 1 HYPRE BoomerAMG: Maximum row sums 0.9 HYPRE BoomerAMG: Sweeps down 1 HYPRE BoomerAMG: Sweeps up 1 HYPRE BoomerAMG: Sweeps on coarse 1 HYPRE BoomerAMG: Relax down symmetric-SOR/Jacobi HYPRE BoomerAMG: Relax up symmetric-SOR/Jacobi HYPRE BoomerAMG: Relax on coarse Gaussian-elimination HYPRE BoomerAMG: Relax weight (all) 1. HYPRE BoomerAMG: Outer relax weight (all) 1. HYPRE BoomerAMG: Using CF-relaxation HYPRE BoomerAMG: Not using more complex smoothers. HYPRE BoomerAMG: Measure type local HYPRE BoomerAMG: Coarsen type Falgout HYPRE BoomerAMG: Interpolation type classical linear system matrix = precond matrix: Mat Object: 2 MPI processes type: mpiaij rows=8000, cols=8000 total: nonzeros=53600, allocated nonzeros=53600 total number of mallocs used during MatSetValues calls =0 not using I-node (on process 0) routines linear system matrix = precond matrix: Mat Object: (flow_) 2 MPI processes type: mpibaij rows=24000, cols=24000, bs=3 total: nonzeros=482400, allocated nonzeros=482400 total number of mallocs used during MatSetValues calls =0 PC Object: (flow_sub_1_) 2 MPI processes type: bjacobi block Jacobi: number of blocks = 2 Local solve is same for all blocks, in the following KSP and PC objects: KSP Object: (flow_sub_1_sub_) 1 MPI processes type: preonly maximum iterations=10000, initial guess is zero tolerances: relative=1e-05, absolute=1e-50, divergence=10000. left preconditioning using NONE norm type for convergence test PC Object: (flow_sub_1_sub_) 1 MPI processes type: lu out-of-place factorization tolerance for zero pivot 2.22045e-14 matrix ordering: nd factor fill ratio given 5., needed 18.3108 Factored matrix follows: Mat Object: 1 MPI processes type: seqbaij rows=12000, cols=12000, bs=3 package used to perform factorization: petsc total: nonzeros=4350654, allocated nonzeros=4350654 total number of mallocs used during MatSetValues calls =0 block size is 3 linear system matrix = precond matrix: Mat Object: (flow_) 1 MPI processes type: seqbaij rows=12000, cols=12000, bs=3 total: nonzeros=237600, allocated nonzeros=237600 total number of mallocs used during MatSetValues calls =0 block size is 3 linear system matrix = precond matrix: Mat Object: (flow_) 2 MPI processes type: mpibaij rows=24000, cols=24000, bs=3 total: nonzeros=482400, allocated nonzeros=482400 total number of mallocs used during MatSetValues calls =0 --------------------------------- linear system matrix = precond matrix: Mat Object: (flow_) 2 MPI processes type: mpibaij rows=24000, cols=24000, bs=3 total: nonzeros=482400, allocated nonzeros=482400 total number of mallocs used during MatSetValues calls =0 Is there a way to improve on the AMG part? Do I have to adjust the tolerances (make the inner tighter)? Which Hypre AMG parameters are worth tuning? This problem occurs for 1 MPI process as well, and solving the problem in Standard PFLOTRAN (i.e. Block Jacobi + ILU) is without any issue. Grateful for any help! Robert