[Trilinos-Users] Belos PCG with ML preconditioner problem

Berndt, Markus berndt at lanl.gov
Tue Jul 8 16:54:10 MDT 2014


We are having a problem with Belos PCG with ML preconditioning and I am hoping somebody on this list has an idea what might be wrong.

The linear problem we are solving is symmetric positive definite (to be specific, it is the pressure projection equation). We are using the same matrix (an Epetra_FECrsMatrix) for Belos PCG as well as to initialize the ML preconditioner.

During the CG iteration the code stops with the following exception:

terminate called after throwing an instance of 'Belos::PCPGIterateFailure'
  what(): <...>/include/BelosPCPGIter.hpp:794:

Throw number = 1

Throw test that evaluated to true: alpha(0,0) <= zero

Belos::CGIter::iterate(): non-positive value for alpha encountered!

This is for pretty much any setup we choose for ML, except when we set number of levels = 1 and choose symmetric Gauss-Seidel as a coarse grid solver (or other symmetric smoothers). This should result in a symmetric preconditioner and it does indeed. As soon as we set the max number of levels > 1, we see the error above.

In the rest of this email I am attaching some debug output that includes the all the parameters used for ML setup. We think we should have set things up such that the ML preconditioner should be symmetric, but perhaps there is a problem. We don't see it.

We appreciate any help!

Thanks

- Markus

------------------------------------------------------------------------------
Using `cg' for eigen-computations
***
*** ML_Epetra::MultiLevelPreconditioner
***
Matrix has 242 rows and 3.8440e+03 nonzeros, distributed over 1 process(es)
The linear system matrix is an Epetra_CrsMatrix
** Leaving column map of Main linear system matrix unchanged
Default values for `SA'
Maximum number of levels = 10
Using increasing levels. Finest level  = 0, coarsest level = 9
Number of applications of the ML cycle = 1
Number of PDE equations = 1
Aggregation threshold = 0.0000e+00
Max coarse size = 16
R and P smoothing : P = (I-\omega A) P_t, R = P^T
R and P smoothing : \omega = 1.3330e+00/lambda_max
Null space type      = default (constants)
Null space dimension = 1
-----------------------------------------------------------------------
ML_Gen_MultiLevelHierarchy (level 0) : Gen Restriction and Prolongator
-----------------------------------------------------------------------
ML_Aggregate_Coarsen (level 0) begins
ML_Aggregate_CoarsenMIS : current level = 0
ML_Aggregate_CoarsenMIS : current eps = 0.000000e+00
Aggregation(MIS) : Total nonzeros = 3844 (Nrows=242)
Aggregation(MIS) : Phase 1 - nodes aggregated = 242
Aggregation(MIS) : Phase 1 - total aggregates = 16
Aggregation(MIS) : Phase 2a- additional aggregates = 0
Aggregation(MIS) : Phase 2 - total aggregates = 16
Aggregation(MIS) : Phase 2 - boundary nodes   = 0
Aggregation(MIS) : Phase 3 - leftovers = 0 and singletons = 0
Communicating phase 2/3 info
Calling ML_Operator_UnAmalgamateAndDropWeak
Aggregation(MIS) : QR factorization - too small aggregates = 0
Gen_Prolongator (level 0) : Max eigenvalue = 2.1838e+00

Prolongator/Restriction smoother (level 0) : damping factor #1 = 6.1041e-01
Prolongator/Restriction smoother (level 0) : ( = 1.3330e+00 / 2.1838e+00)

Smoothed Aggregation : operator complexity = 1.026015e+00.
Time to build the hierarchy = 2.1720e-03 (s)
Number of actual levels : 2

Smoother (level 0) : # global rows = 2.e+02, # estim. global nnz = 4.e+03, # nnz
per row = 16.
Smoother (level 0) : IFPACK, type=`IC',
Smoother (level 0) : both,overlap=0
Smoother (level 0) : level-of-fill=0.0000e+00,rel. threshold=1.0000e+00,abs.
threshold=0.0000e+00
Smoother (level 0) : Setup time : 6.2299e-04 (s)

WARNING **********************************************************************
You have requested to keep the finest level smoother from a previous
solve.  In order to reuse part or all of the incomplete factorization
from smoother "IC", you must also specify the smoother options
"reuse numeric factorization" and/or "reuse symbolic factorization".
Not doing so  may result in memory leaks or crashes.

Amesos (level 1) : NumGlobalRows = 16
Amesos (level 1) : NumGlobalNonzeros = 100
Amesos (level 1) : Fill-in = 3.9062e+01 %
Amesos (level 1) : Building KLU
Amesos (level 1) : Time for factorization  = 8.2016e-05 (s)

WARNING: Parameter "ML print final list" -1   [unused] is unused
WARNING: Parameter "coarse: list"    [unused] is unused
WARNING: Parameter "aggregation: list (level 0)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 1)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 2)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 3)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 4)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 5)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 6)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 7)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 8)"    [unused] is unused
WARNING: Parameter "aggregation: list (level 9)"    [unused] is unused
WARNING: Parameter "smoother: list (level 0)"    [unused] is unused
++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+++ Printing ML parameter list "SA default values" on pid 0
++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
default values = SA
max levels = 10
prec type = MGV
increasing or decreasing = increasing
aggregation: type = Uncoupled-MIS
aggregation: damping factor = 1.3330e+00
eigen-analysis: type = cg
eigen-analysis: iterations = 10
smoother: sweeps = 1
smoother: damping factor = 1.0000e+00
smoother: pre or post = both
smoother: type = IC
ML output = 10
print unused = -1
ML print initial list = -1
ML print final list = -1
read XML = 0   [default]
ML debug mode = 0   [default]
output prefix =    [default]
analyze memory = 0   [default]
cycle applications = 1   [default]
zero starting solution = 1   [default]
ML label = not-set   [default]
low memory usage = 0   [default]
PDE equations = 1   [default]
ML validate depth = 5   [default]
ML validate parameter list = 1   [default]
profile: operator iterations = 0   [default]
filter: type = 0   [default]
filter: mask = 0   [default]
use crs matrix storage = 0   [default]
aggregation: aggressive = 5.0000e-01   [default]
aggregation: threshold = 0.0000e+00   [default]
aggregation: next-level aggregates per process = 128   [default]
aggregation: symmetrize = 0   [default]
energy minimization: enable = 0   [default]
aggregation: smoothing sweeps = 1   [default]
R and P smoothing: type = classic   [default]
aggregation: aux: enable = 0   [default]
aggregation: block scaling = 0   [default]
aggregation: use tentative restriction = 0   [default]
null space: type = default vectors   [default]
null space: scaling = 0   [default]
x-coordinates = 0   [default]
y-coordinates = 0   [default]
z-coordinates = 0   [default]
repartition: enable = 0   [default]
smoother: Aztec options = Teuchos::RCP<std::vector<int, std::allocator<int> > >{ptr=0x3935640,node=0x3935660,strong_count=2,weak_count=0}   [default]
smoother: Aztec params = Teuchos::RCP<std::vector<double, std::allocator<double> > >{ptr=0x39356a0,node=0x39356c0,strong_count=2,weak_count=0}   [default]
smoother: Aztec as solver = 0   [default]
smoother: MLS polynomial order = -97   [default]
smoother: polynomial order = -97   [default]
smoother: MLS alpha = -2.0000e+00   [default]
smoother: Chebyshev alpha = -2.0000e+00   [default]
smoother: ParaSails levels = 0   [default]
smoother: ParaSails matrix = 0   [default]
smoother: ParaSails threshold = 1.0000e-02   [default]
smoother: ParaSails filter = 5.0000e-02   [default]
smoother: ParaSails load balancing = 0.0000e+00   [default]
smoother: ParaSails factorized = 0   [default]
smoother: ifpack type = Amesos   [default]
smoother: ifpack overlap = 0   [default]
smoother: ifpack level-of-fill = 0.0000e+00   [default]
smoother: ifpack relative threshold = 1.0000e+00   [default]
smoother: ifpack absolute threshold = 0.0000e+00   [default]
smoother: Block Chebyshev number of blocks = -1   [default]
smoother: Block Chebyshev block list = 0   [default]
smoother: Block Chebyshev block starts = 0   [default]
smoother: chebyshev solve normal equations = 0   [default]
filtering: enable = 0   [default]
dump matrix: enable = 0   [default]
print hierarchy = -2   [default]
reuse: enable = 0   [default]
coarse: list ->
 smoother: type = Amesos-KLU
 smoother: max size = 16
 smoother: pre or post = post
 smoother: sweeps = 1
 smoother: damping factor = 1   [default]
 smoother: add to diag = 1e-12   [default]
aggregation: list (level 0) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 1) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 2) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 3) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 4) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 5) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 6) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 7) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 8) ->
 aggregation: type = Uncoupled-MIS   [default]
 aggregation: smoothing sweeps = 1   [default]
aggregation: list (level 9) ->
 aggregation: smoothing sweeps = 1   [default]
smoother: list (level 0) ->
 smoother: sweeps = 1   [default]
 smoother: damping factor = 1   [default]
 smoother: pre or post = both   [default]
 smoother: type = IC   [default]
 smoother: add to diag = 1e-12   [default]
 smoother: ifpack level-of-fill = 0   [default]
 smoother: ifpack overlap = 0   [default]
 smoother: ifpack relative threshold = 1   [default]
 smoother: ifpack absolute threshold = 0   [default]
 smoother: ifpack list ->
  ILU: sweeps = 1   [unused]
  partitioner: map = 0x39432d8   [unused]
  fact: level-of-fill = 0
  fact: relative threshold = 1   [unused]
  fact: absolute threshold = 0   [unused]
-------------------------------------------------
----------- end of ML parameter list ------------
-------------------------------------------------

Cumulative timing for construction so far:
- for initial setup   = 5.0602e-03 (s)
- for hierarchy setup = 2.1951e-03 (s)
- for smoothers setup = 6.2299e-04 (s)
- for coarse setup    = 2.3699e-04 (s)
- for final setup     = 6.2203e-04 (s)
Total for this setup  = 8.9622e-03 (s)
------------------------------------------------------------------------------
============ ML ============

============ Belos ============
Block Size = 1   [unused]
Maximum Iterations = 250   [unused]
Convergence Tolerance = 1.0000e-12   [unused]
Verbosity = 59   [unused]
Output Frequency = 1   [unused]

Belos::StatusTestGeneralOutput: Failed
  (Num calls,Mod test,State test): (1, 1, Passed Failed Undefined)
   Failed.......OR Combination ->
     OK...........Number of Iterations = 0 < 250
     Unconverged..(2-Norm Imp Res Vec) / (2-Norm Res0)
                  residual [ 0 ] = 1.0000e+00 > 1.0000e-12


Belos::StatusTestGeneralOutput: Failed
  (Num calls,Mod test,State test): (2, 1, Passed Failed Undefined)
   Failed.......OR Combination ->
     OK...........Number of Iterations = 1 < 250
     Unconverged..(2-Norm Imp Res Vec) / (2-Norm Res0)
                  residual [ 0 ] = 9.5492e-01 > 1.0000e-12

terminate called after throwing an instance of 'Belos::PCPGIterateFailure'
  what(): <...>/include/BelosPCPGIter.hpp:794:

Throw number = 1

Throw test that evaluated to true: alpha(0,0) <= zero

Belos::CGIter::iterate(): non-positive value for alpha encountered!
============ Belos ============



--
Markus Berndt -- LANL CCS-2 -- 505-665-4711
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://software.sandia.gov/pipermail/trilinos-users/attachments/20140708/e11df12d/attachment.html>


More information about the Trilinos-Users mailing list