[petsc-users] Confusing Schur preconditioner behaviour

Justin Chang jychang48 at gmail.com
Mon Mar 18 16:44:22 CDT 2019


Hi Colin,

Perhaps I may be talking out of my butt at this point, but it looks like
this term in your aP:

dt**2*inner(g*H*div(u)/ai,div(v))

is secretly a problem. When I see guys like this in my schur complement
approximation, I found that scaling this term by some penalty/weighting
value helps especially when dt is really small. Not exactly sure what the
mathematical explanation is, but I found that this heuristical approach
improves my schur complement approximations.

Justin

On Mon, Mar 18, 2019 at 3:33 PM Cotter, Colin J <colin.cotter at imperial.ac.uk>
wrote:

> Hi Justin,
>
>   1) Here is some UFL for my mixed system (u in BDM2, h in DG1)
>
> a = (
>
>
>     inner(ai*u,v) - dt*inner(f*perp(u),v) +
>
>
>     dt*inner(g*h,div(v))
>
>
>     +inner(ai*h,phi) - dt*inner(H*div(u),phi)
>
>
> )*dx
>
>
> Which then has the exact Schur complement (when using an H(div)-L2 element
> pair
>
>
> aP = (inner(ai*u,v) + dt**2*inner(g*H*div(u)/ai,div(v))
>
>
>       - dt*inner(f*perp(u),v))*dx
>
>
> 2) I put FGMRES on the outer iteration and appended the result to the
> bottom of this message.
>
>
> all the best
>
> --Colin
>
>
>     Residual norms for firedrake_0_ solve.
>
>     0 KSP unpreconditioned resid norm 1.086016610848e-03 true resid norm
> 1.086016610848e-03 ||r(i)||/||b|| 1.000000000000e+00
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 8.120721494511e+01 true resid norm
> 1.655197512879e+02 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 9.439702024300e+00 true resid norm
> 3.075997034466e+01 ||r(i)||/||b|| 1.858386694356e-01
>
>       2 KSP preconditioned resid norm 1.137279699060e+00 true resid norm
> 7.503013051026e+00 ||r(i)||/||b|| 4.533001646417e-02
>
>       3 KSP preconditioned resid norm 1.496062340941e-01 true resid norm
> 1.484527899318e+00 ||r(i)||/||b|| 8.968886720568e-03
>
>       4 KSP preconditioned resid norm 2.056482137526e-02 true resid norm
> 3.031663733277e-01 ||r(i)||/||b|| 1.831602397713e-03
>
>       5 KSP preconditioned resid norm 1.745270896489e-03 true resid norm
> 3.767168288302e-02 ||r(i)||/||b|| 2.275962994742e-04
>
>       6 KSP preconditioned resid norm 1.564628458491e-04 true resid norm
> 3.546145962849e-03 ||r(i)||/||b|| 2.142430698002e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     1 KSP unpreconditioned resid norm 1.085836886832e-03 true resid norm
> 1.085836886832e-03 ||r(i)||/||b|| 9.998345108044e-01
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 3.877459187203e-01 true resid norm
> 1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 4.065088730543e-02 true resid norm
> 1.756656283999e-01 ||r(i)||/||b|| 1.756656283999e-01
>
>       2 KSP preconditioned resid norm 5.708009085053e-03 true resid norm
> 4.422956647739e-02 ||r(i)||/||b|| 4.422956647739e-02
>
>       3 KSP preconditioned resid norm 7.309033373261e-04 true resid norm
> 8.666099848662e-03 ||r(i)||/||b|| 8.666099848662e-03
>
>       4 KSP preconditioned resid norm 1.133694903547e-04 true resid norm
> 1.609496315342e-03 ||r(i)||/||b|| 1.609496315342e-03
>
>       5 KSP preconditioned resid norm 1.227169910863e-05 true resid norm
> 2.710100172099e-04 ||r(i)||/||b|| 2.710100172099e-04
>
>       6 KSP preconditioned resid norm 1.047216213052e-06 true resid norm
> 2.489323919135e-05 ||r(i)||/||b|| 2.489323919135e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     2 KSP unpreconditioned resid norm 1.075476649437e-03 true resid norm
> 1.075476649437e-03 ||r(i)||/||b|| 9.902948432781e-01
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 4.856780223694e-01 true resid norm
> 1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 4.745469162937e-02 true resid norm
> 3.302334019978e-01 ||r(i)||/||b|| 3.302334019978e-01
>
>       2 KSP preconditioned resid norm 7.434637319370e-03 true resid norm
> 7.716456937097e-02 ||r(i)||/||b|| 7.716456937097e-02
>
>       3 KSP preconditioned resid norm 9.477403463953e-04 true resid norm
> 1.245956505675e-02 ||r(i)||/||b|| 1.245956505675e-02
>
>       4 KSP preconditioned resid norm 1.197149803247e-04 true resid norm
> 2.445415954040e-03 ||r(i)||/||b|| 2.445415954041e-03
>
>       5 KSP preconditioned resid norm 8.629298003468e-06 true resid norm
> 2.172099238673e-04 ||r(i)||/||b|| 2.172099238673e-04
>
>       6 KSP preconditioned resid norm 8.102365661753e-07 true resid norm
> 2.101591443240e-05 ||r(i)||/||b|| 2.101591443241e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     3 KSP unpreconditioned resid norm 8.228196318802e-04 true resid norm
> 8.228196318802e-04 ||r(i)||/||b|| 7.576492142582e-01
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 5.432442358474e-01 true resid norm
> 1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 5.878289354883e-02 true resid norm
> 3.137642427673e-01 ||r(i)||/||b|| 3.137642427673e-01
>
>       2 KSP preconditioned resid norm 5.788229173110e-03 true resid norm
> 8.897854046325e-02 ||r(i)||/||b|| 8.897854046324e-02
>
>       3 KSP preconditioned resid norm 8.498991030513e-04 true resid norm
> 1.495053243269e-02 ||r(i)||/||b|| 1.495053243269e-02
>
>       4 KSP preconditioned resid norm 1.033515356385e-04 true resid norm
> 1.845694418294e-03 ||r(i)||/||b|| 1.845694418294e-03
>
>       5 KSP preconditioned resid norm 7.847588222485e-06 true resid norm
> 2.466173373852e-04 ||r(i)||/||b|| 2.466173373852e-04
>
>       6 KSP preconditioned resid norm 7.146079509981e-07 true resid norm
> 1.719476126742e-05 ||r(i)||/||b|| 1.719476126742e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     4 KSP unpreconditioned resid norm 1.657400277959e-04 true resid norm
> 1.657400277959e-04 ||r(i)||/||b|| 1.526127926041e-01
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 5.120302137337e-01 true resid norm
> 9.999999999998e-01 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 6.454992757385e-02 true resid norm
> 2.912694339903e-01 ||r(i)||/||b|| 2.912694339903e-01
>
>       2 KSP preconditioned resid norm 6.671724869182e-03 true resid norm
> 6.254491231334e-02 ||r(i)||/||b|| 6.254491231335e-02
>
>       3 KSP preconditioned resid norm 8.212101901860e-04 true resid norm
> 1.773238703254e-02 ||r(i)||/||b|| 1.773238703254e-02
>
>       4 KSP preconditioned resid norm 1.093679881645e-04 true resid norm
> 2.931666247026e-03 ||r(i)||/||b|| 2.931666247027e-03
>
>       5 KSP preconditioned resid norm 1.186578053622e-05 true resid norm
> 3.303748399822e-04 ||r(i)||/||b|| 3.303748399822e-04
>
>       6 KSP preconditioned resid norm 6.503381106855e-07 true resid norm
> 1.912788140490e-05 ||r(i)||/||b|| 1.912788140490e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     5 KSP unpreconditioned resid norm 2.448137330609e-05 true resid norm
> 2.448137330609e-05 ||r(i)||/||b|| 2.254235622324e-02
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 5.179583235178e-01 true resid norm
> 1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 5.254142246569e-02 true resid norm
> 2.516327547106e-01 ||r(i)||/||b|| 2.516327547105e-01
>
>       2 KSP preconditioned resid norm 7.363444586475e-03 true resid norm
> 5.234766639265e-02 ||r(i)||/||b|| 5.234766639263e-02
>
>       3 KSP preconditioned resid norm 9.704643808280e-04 true resid norm
> 1.267956066012e-02 ||r(i)||/||b|| 1.267956066011e-02
>
>       4 KSP preconditioned resid norm 1.086598853731e-04 true resid norm
> 3.391601116038e-03 ||r(i)||/||b|| 3.391601116037e-03
>
>       5 KSP preconditioned resid norm 7.629924161284e-06 true resid norm
> 3.006619653634e-04 ||r(i)||/||b|| 3.006619653632e-04
>
>       6 KSP preconditioned resid norm 7.239392986410e-07 true resid norm
> 2.150547391218e-05 ||r(i)||/||b|| 2.150547391217e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     6 KSP unpreconditioned resid norm 2.402006637155e-06 true resid norm
> 2.402006637153e-06 ||r(i)||/||b|| 2.211758653745e-03
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 3.416765078712e-01 true resid norm
> 9.999999999993e-01 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 3.841617294439e-02 true resid norm
> 1.983060322744e-01 ||r(i)||/||b|| 1.983060322745e-01
>
>       2 KSP preconditioned resid norm 5.192499214944e-03 true resid norm
> 3.903498722784e-02 ||r(i)||/||b|| 3.903498722787e-02
>
>       3 KSP preconditioned resid norm 7.152864973143e-04 true resid norm
> 7.896855714048e-03 ||r(i)||/||b|| 7.896855714054e-03
>
>       4 KSP preconditioned resid norm 9.150050794784e-05 true resid norm
> 1.709682316519e-03 ||r(i)||/||b|| 1.709682316520e-03
>
>       5 KSP preconditioned resid norm 1.369470900870e-05 true resid norm
> 4.755514972773e-04 ||r(i)||/||b|| 4.755514972776e-04
>
>       6 KSP preconditioned resid norm 8.553609674167e-07 true resid norm
> 4.128479317005e-05 ||r(i)||/||b|| 4.128479317008e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     7 KSP unpreconditioned resid norm 1.986730253405e-07 true resid norm
> 1.986730253416e-07 ||r(i)||/||b|| 1.829373725568e-04
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 3.737152121707e-01 true resid norm
> 1.000000000001e+00 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 2.327226787051e-02 true resid norm
> 1.728532258160e-01 ||r(i)||/||b|| 1.728532258158e-01
>
>       2 KSP preconditioned resid norm 2.663283601003e-03 true resid norm
> 2.117224456449e-02 ||r(i)||/||b|| 2.117224456446e-02
>
>       3 KSP preconditioned resid norm 4.000029435799e-04 true resid norm
> 4.527340194475e-03 ||r(i)||/||b|| 4.527340194470e-03
>
>       4 KSP preconditioned resid norm 4.886053047259e-05 true resid norm
> 7.959488416139e-04 ||r(i)||/||b|| 7.959488416128e-04
>
>       5 KSP preconditioned resid norm 7.135914581342e-06 true resid norm
> 1.973225005950e-04 ||r(i)||/||b|| 1.973225005947e-04
>
>       6 KSP preconditioned resid norm 1.059888416230e-06 true resid norm
> 4.026894864745e-05 ||r(i)||/||b|| 4.026894864740e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>     8 KSP unpreconditioned resid norm 1.847329119073e-08 true resid norm
> 1.847329119335e-08 ||r(i)||/||b|| 1.701013686976e-05
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 3.105703514818e-01 true resid norm
> 9.999999999993e-01 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 2.274848341406e-02 true resid norm
> 1.655926612050e-01 ||r(i)||/||b|| 1.655926612051e-01
>
>       2 KSP preconditioned resid norm 1.896323766042e-03 true resid norm
> 2.108620526506e-02 ||r(i)||/||b|| 2.108620526507e-02
>
>       3 KSP preconditioned resid norm 2.098171871374e-04 true resid norm
> 2.562423770591e-03 ||r(i)||/||b|| 2.562423770592e-03
>
>       4 KSP preconditioned resid norm 2.802038862500e-05 true resid norm
> 4.753351438216e-04 ||r(i)||/||b|| 4.753351438219e-04
>
>       5 KSP preconditioned resid norm 2.875826461982e-06 true resid norm
> 8.618109518277e-05 ||r(i)||/||b|| 8.618109518283e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 5
>
>     9 KSP unpreconditioned resid norm 1.093052879829e-09 true resid norm
> 1.093052879480e-09 ||r(i)||/||b|| 1.006478969623e-06
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>
>       0 KSP preconditioned resid norm 2.400304155508e-01 true resid norm
> 9.999999999373e-01 ||r(i)||/||b|| 1.000000000000e+00
>
>       1 KSP preconditioned resid norm 1.375935214529e-02 true resid norm
> 1.268731042794e-01 ||r(i)||/||b|| 1.268731042874e-01
>
>       2 KSP preconditioned resid norm 1.359927799420e-03 true resid norm
> 1.743482323300e-02 ||r(i)||/||b|| 1.743482323409e-02
>
>       3 KSP preconditioned resid norm 1.329143848901e-04 true resid norm
> 2.112118748978e-03 ||r(i)||/||b|| 2.112118749111e-03
>
>       4 KSP preconditioned resid norm 2.299763657223e-05 true resid norm
> 3.661832963830e-04 ||r(i)||/||b|| 3.661832964060e-04
>
>       5 KSP preconditioned resid norm 3.042002639898e-06 true resid norm
> 8.825775436656e-05 ||r(i)||/||b|| 8.825775437210e-05
>
>       6 KSP preconditioned resid norm 2.755403722945e-07 true resid norm
> 1.154450754603e-05 ||r(i)||/||b|| 1.154450754675e-05
>
>     Linear firedrake_0_fieldsplit_1_ solve converged due to
> CONVERGED_RTOL iterations 6
>
>    10 KSP unpreconditioned resid norm 6.865249627765e-11 true resid norm
> 6.865249509492e-11 ||r(i)||/||b|| 6.321495860116e-08
>
> KSP Object: (firedrake_0_) 1 MPI processes
>
>   type: fgmres
>
>     restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>
>     happy breakdown tolerance 1e-30
>
>   maximum iterations=10000, initial guess is zero
>
>   tolerances:  relative=1e-07, absolute=1e-50, divergence=10000.
>
>   right preconditioning
>
>   using UNPRECONDITIONED norm type for convergence test
>
> PC Object: (firedrake_0_) 1 MPI processes
>
>   type: fieldsplit
>
>     FieldSplit with Schur preconditioner, factorization FULL
>
>     Preconditioner for the Schur complement formed from A11
>
>     Split info:
>
>     Split number 0 Defined by IS
>
>     Split number 1 Defined by IS
>
>     KSP solver for A00 block
>
>       KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>
>             Factored matrix follows:
>
>               Mat Object: 1 MPI processes
>
>                 type: seqaij
>
>                 rows=6144, cols=6144
>
>                 package used to perform factorization: petsc
>
>                 total: nonzeros=18432, allocated nonzeros=18432
>
>                 total number of mallocs used during MatSetValues calls =0
>
>                   using I-node routines: found 2048 nodes, limit used is 5
>
>         linear system matrix = precond matrix:
>
>         Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>
>           type: seqaij
>
>           rows=6144, cols=6144
>
>           total: nonzeros=18432, allocated nonzeros=18432
>
>           total number of mallocs used during MatSetValues calls =0
>
>             using I-node routines: found 2048 nodes, limit used is 5
>
>     KSP solver for S = A11 - A10 inv(A00) A01
>
>       KSP Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>
>         type: gmres
>
>           restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>
>           happy breakdown tolerance 1e-30
>
>         maximum iterations=10000, initial guess is zero
>
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>
>         left preconditioning
>
>         using PRECONDITIONED norm type for convergence test
>
>       PC Object: (firedrake_0_fieldsplit_1_) 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 5.09173
>
>             Factored matrix follows:
>
>               Mat Object: 1 MPI processes
>
>                 type: seqaij
>
>                 rows=15360, cols=15360
>
>                 package used to perform factorization: petsc
>
>                 total: nonzeros=1360836, allocated nonzeros=1360836
>
>                 total number of mallocs used during MatSetValues calls =0
>
>                   using I-node routines: found 5120 nodes, limit used is 5
>
>         linear system matrix followed by preconditioner matrix:
>
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>
>           type: schurcomplement
>
>           rows=15360, cols=15360
>
>             Schur complement A11 - A10 inv(A00) A01
>
>             A11
>
>               Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>
>                 type: seqaij
>
>                 rows=15360, cols=15360
>
>                 total: nonzeros=267264, allocated nonzeros=267264
>
>                 total number of mallocs used during MatSetValues calls =0
>
>                   using I-node routines: found 5120 nodes, limit used is 5
>
>             A10
>
>               Mat Object: 1 MPI processes
>
>                 type: seqaij
>
>                 rows=15360, cols=6144
>
>                 total: nonzeros=73728, allocated nonzeros=73728
>
>                 total number of mallocs used during MatSetValues calls =0
>
>                   using I-node routines: found 5120 nodes, limit used is 5
>
>             KSP of A00
>
>               KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>
>                     Factored matrix follows:
>
>                       Mat Object: 1 MPI processes
>
>                         type: seqaij
>
>                         rows=6144, cols=6144
>
>                         package used to perform factorization: petsc
>
>                         total: nonzeros=18432, allocated nonzeros=18432
>
>                         total number of mallocs used during MatSetValues
> calls =0
>
>                           using I-node routines: found 2048 nodes, limit
> used is 5
>
>                 linear system matrix = precond matrix:
>
>                 Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>
>                   type: seqaij
>
>                   rows=6144, cols=6144
>
>                   total: nonzeros=18432, allocated nonzeros=18432
>
>                   total number of mallocs used during MatSetValues calls
> =0
>
>                     using I-node routines: found 2048 nodes, limit used
> is 5
>
>             A01
>
>               Mat Object: 1 MPI processes
>
>                 type: seqaij
>
>                 rows=6144, cols=15360
>
>                 total: nonzeros=73728, allocated nonzeros=73728
>
>                 total number of mallocs used during MatSetValues calls =0
>
>                   using I-node routines: found 2048 nodes, limit used is 5
>
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>
>           type: seqaij
>
>           rows=15360, cols=15360
>
>           total: nonzeros=267264, allocated nonzeros=267264
>
>           total number of mallocs used during MatSetValues calls =0
>
>             using I-node routines: found 5120 nodes, limit used is 5
>
>   linear system matrix followed by preconditioner matrix:
>
>   Mat Object: (firedrake_0_) 1 MPI processes
>
>     type: seqaij
>
>     rows=21504, cols=21504
>
>     total: nonzeros=433152, allocated nonzeros=433152
>
>     total number of mallocs used during MatSetValues calls =0
>
>       using I-node routines: found 7168 nodes, limit used is 5
>
>   Mat Object: (firedrake_0_) 1 MPI processes
>
>     type: seqaij
>
>     rows=21504, cols=21504
>
>     total: nonzeros=433152, allocated nonzeros=433152
>
>     total number of mallocs used during MatSetValues calls =0
>
>       using I-node routines: found 7168 nodes, limit used is 5
> ------------------------------
> *From:* Justin Chang <jychang48 at gmail.com>
> *Sent:* 18 March 2019 21:01:17
> *To:* Cotter, Colin J
> *Cc:* petsc-users at mcs.anl.gov
> *Subject:* Re: [petsc-users] Confusing Schur preconditioner behaviour
>
> Colin,
>
> 1) What equations are you solving?
>
> 2) In your second case, you set the outer ksp to preonly, thus we are
> unable to see the ksp_monitor for the (firedrake_0_) solver. Set it to
> gmres and see if you have a similar output to your first case:
>
> 0 KSP preconditioned resid norm 4.985448866758e+00 true resid norm
> 1.086016610848e-03 ||r(i)||/||b|| 1.000000000000e+00
> 1 KSP preconditioned resid norm 1.245615753306e-13 true resid norm
> 2.082000915439e-14 ||r(i)||/||b|| 1.917098591903e-11
>
> Because according to the first ksp_view output, after one lu sweep for the
> (firedrake_0_fieldsplit_1_) solver. That is, going from:
>
> 0 KSP preconditioned resid norm 8.819238435108e-02 true resid norm
> 1.797571993221e-01 ||r(i)||/||b|| 1.000000000000e+00
>
> to
>
> 1 KSP preconditioned resid norm 1.025167319984e-02 true resid norm
> 3.340583874349e-02 ||r(i)||/||b|| 1.858386694356e-01
>
> appeared to give you an exact schur complement.
>
> Justin
>
> On Mon, Mar 18, 2019 at 2:18 PM Cotter, Colin J via petsc-users <
> petsc-users at mcs.anl.gov> wrote:
>
> Sorry, just to clarify, in the second case I see several *inner*
> iterations, even though I'm using LU on a supposedly exact Schur complement
> as the preconditioner for the Schur system.
> ------------------------------
> *From:* petsc-users <petsc-users-bounces at mcs.anl.gov> on behalf of
> Cotter, Colin J via petsc-users <petsc-users at mcs.anl.gov>
> *Sent:* 18 March 2019 20:14:48
> *To:* petsc-users at mcs.anl.gov
> *Subject:* [petsc-users] Confusing Schur preconditioner behaviour
>
>
> Dear petsc-users,
>
>   I'm solving a 2x2 block system, for which I can construct the Schur
> complement analytically (through compatible FEM stuff),
>
> which I can pass as the preconditioning matrix.
>
>
> When using gmres on the outer iteration, and preonly+lu on the inner
> iterations with a Schur complement preconditioner,
>
> I see convergence in 1 iteration as expected. However, when I set gmres+lu
> on the inner iteration for S, I see several iterations.
>
>
> This seems strange to me, as the first result seems to confirm that I have
> an exact Schur complement, but the second result
>
> implies not.
>
> What could be going on here?
>
>
> I've appended output to the bottom of this message, first the preonly+lu
> and then for gmres+lu.
>
>
> all the best
>
> --Colin
>
>
>    Linear firedrake_0_fieldsplit_1_ solve converged due to CONVERGED_ITS
> iterations 1
>     Residual norms for firedrake_0_ solve.
>     0 KSP preconditioned resid norm 4.985448866758e+00 true resid norm
> 1.086016610848e-03 ||r(i)||/||b|| 1.000000000000e+00
>     Linear firedrake_0_fieldsplit_1_ solve converged due to CONVERGED_ITS
> iterations 1
>     1 KSP preconditioned resid norm 1.245615753306e-13 true resid norm
> 2.082000915439e-14 ||r(i)||/||b|| 1.917098591903e-11
> KSP Object: (firedrake_0_) 1 MPI processes
>   type: gmres
>     restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>     happy breakdown tolerance 1e-30
>   maximum iterations=10000, initial guess is zero
>   tolerances:  relative=1e-07, absolute=1e-50, divergence=10000.
>   left preconditioning
>   using PRECONDITIONED norm type for convergence test
> PC Object: (firedrake_0_) 1 MPI processes
>   type: fieldsplit
>     FieldSplit with Schur preconditioner, factorization FULL
>     Preconditioner for the Schur complement formed from A11
>     Split info:
>     Split number 0 Defined by IS
>     Split number 1 Defined by IS
>     KSP solver for A00 block
>       KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>             Factored matrix follows:
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=6144, cols=6144
>                 package used to perform factorization: petsc
>                 total: nonzeros=18432, allocated nonzeros=18432
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 2048 nodes, limit used is 5
>         linear system matrix = precond matrix:
>         Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>           type: seqaij
>           rows=6144, cols=6144
>           total: nonzeros=18432, allocated nonzeros=18432
>           total number of mallocs used during MatSetValues calls =0
>             using I-node routines: found 2048 nodes, limit used is 5
>     KSP solver for S = A11 - A10 inv(A00) A01
>       KSP Object: (firedrake_0_fieldsplit_1_) 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: (firedrake_0_fieldsplit_1_) 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 5.09173
>             Factored matrix follows:
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=15360
>                 package used to perform factorization: petsc
>                 total: nonzeros=1360836, allocated nonzeros=1360836
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>         linear system matrix followed by preconditioner matrix:
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>           type: schurcomplement
>           rows=15360, cols=15360
>             Schur complement A11 - A10 inv(A00) A01
>             A11
>               Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=15360
>                 total: nonzeros=267264, allocated nonzeros=267264
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>             A10
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=6144
>                 total: nonzeros=73728, allocated nonzeros=73728
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>             KSP of A00
>               KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>                     Factored matrix follows:
>                       Mat Object: 1 MPI processes
>                         type: seqaij
>                         rows=6144, cols=6144
>                         package used to perform factorization: petsc
>                         total: nonzeros=18432, allocated nonzeros=18432
>                         total number of mallocs used during MatSetValues
> calls =0
>                           using I-node routines: found 2048 nodes, limit
> used is 5
>                 linear system matrix = precond matrix:
>                 Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>                   type: seqaij
>                   rows=6144, cols=6144
>                   total: nonzeros=18432, allocated nonzeros=18432
>                   total number of mallocs used during MatSetValues calls =0
>                     using I-node routines: found 2048 nodes, limit used is
> 5
>             A01
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=6144, cols=15360
>                 total: nonzeros=73728, allocated nonzeros=73728
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 2048 nodes, limit used is 5
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>           type: seqaij
>           rows=15360, cols=15360
>           total: nonzeros=267264, allocated nonzeros=267264
>           total number of mallocs used during MatSetValues calls =0
>             using I-node routines: found 5120 nodes, limit used is 5
>   linear system matrix followed by preconditioner matrix:
>   Mat Object: (firedrake_0_) 1 MPI processes
>     type: nest
>     rows=21504, cols=21504
>       Matrix object:
>         type=nest, rows=2, cols=2
>         MatNest structure:
>         (0,0) : type=seqaij, rows=15360, cols=15360
>         (0,1) : type=seqaij, rows=15360, cols=6144
>         (1,0) : type=seqaij, rows=6144, cols=15360
>         (1,1) : type=seqaij, rows=6144, cols=6144
>   Mat Object: (firedrake_0_) 1 MPI processes
>     type: nest
>     rows=21504, cols=21504
>       Matrix object:
>         type=nest, rows=2, cols=2
>         MatNest structure:
>         (0,0) : prefix="firedrake_0_fieldsplit_1_", type=seqaij,
> rows=15360, cols=15360
>         (0,1) : type=seqaij, rows=15360, cols=6144
>         (1,0) : type=seqaij, rows=6144, cols=15360
>         (1,1) : prefix="firedrake_0_fieldsplit_0_", type=seqaij,
> rows=6144, cols=6144
>
> =====
>
>
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>       0 KSP preconditioned resid norm 8.819238435108e-02 true resid norm
> 1.797571993221e-01 ||r(i)||/||b|| 1.000000000000e+00
>       1 KSP preconditioned resid norm 1.025167319984e-02 true resid norm
> 3.340583874349e-02 ||r(i)||/||b|| 1.858386694356e-01
>       2 KSP preconditioned resid norm 1.235104644359e-03 true resid norm
> 8.148396804822e-03 ||r(i)||/||b|| 4.533001646417e-02
>       3 KSP preconditioned resid norm 1.624748553125e-04 true resid norm
> 1.612221957927e-03 ||r(i)||/||b|| 8.968886720573e-03
>       4 KSP preconditioned resid norm 2.233373761266e-05 true resid norm
> 3.292437172839e-04 ||r(i)||/||b|| 1.831602397710e-03
>       5 KSP preconditioned resid norm 1.895393184017e-06 true resid norm
> 4.091207337005e-05 ||r(i)||/||b|| 2.275962994770e-04
>       6 KSP preconditioned resid norm 1.699212495729e-07 true resid norm
> 3.851173419652e-06 ||r(i)||/||b|| 2.142430697728e-05
>     Linear firedrake_0_fieldsplit_1_ solve converged due to CONVERGED_RTOL
> iterations 6
> KSP Object: (firedrake_0_) 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: (firedrake_0_) 1 MPI processes
>   type: fieldsplit
>     FieldSplit with Schur preconditioner, factorization FULL
>     Preconditioner for the Schur complement formed from A11
>     Split info:
>     Split number 0 Defined by IS
>     Split number 1 Defined by IS
>     KSP solver for A00 block
>       KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>             Factored matrix follows:
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=6144, cols=6144
>                 package used to perform factorization: petsc
>                 total: nonzeros=18432, allocated nonzeros=18432
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 2048 nodes, limit used is 5
>         linear system matrix = precond matrix:
>         Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>           type: seqaij
>           rows=6144, cols=6144
>           total: nonzeros=18432, allocated nonzeros=18432
>           total number of mallocs used during MatSetValues calls =0
>             using I-node routines: found 2048 nodes, limit used is 5
>     KSP solver for S = A11 - A10 inv(A00) A01
>       KSP Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>         type: gmres
>           restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>           happy breakdown tolerance 1e-30
>         maximum iterations=10000, initial guess is zero
>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>         left preconditioning
>         using PRECONDITIONED norm type for convergence test
>       PC Object: (firedrake_0_fieldsplit_1_) 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 5.09173
>             Factored matrix follows:
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=15360
>                 package used to perform factorization: petsc
>                 total: nonzeros=1360836, allocated nonzeros=1360836
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>         linear system matrix followed by preconditioner matrix:
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>           type: schurcomplement
>           rows=15360, cols=15360
>             Schur complement A11 - A10 inv(A00) A01
>             A11
>               Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=15360
>                 total: nonzeros=267264, allocated nonzeros=267264
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>             A10
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=15360, cols=6144
>                 total: nonzeros=73728, allocated nonzeros=73728
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 5120 nodes, limit used is 5
>             KSP of A00
>               KSP Object: (firedrake_0_fieldsplit_0_) 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: (firedrake_0_fieldsplit_0_) 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 1.
>                     Factored matrix follows:
>                       Mat Object: 1 MPI processes
>                         type: seqaij
>                         rows=6144, cols=6144
>                         package used to perform factorization: petsc
>                         total: nonzeros=18432, allocated nonzeros=18432
>                         total number of mallocs used during MatSetValues
> calls =0
>                           using I-node routines: found 2048 nodes, limit
> used is 5
>                 linear system matrix = precond matrix:
>                 Mat Object: (firedrake_0_fieldsplit_0_) 1 MPI processes
>                   type: seqaij
>                   rows=6144, cols=6144
>                   total: nonzeros=18432, allocated nonzeros=18432
>                   total number of mallocs used during MatSetValues calls =0
>                     using I-node routines: found 2048 nodes, limit used is
> 5
>             A01
>               Mat Object: 1 MPI processes
>                 type: seqaij
>                 rows=6144, cols=15360
>                 total: nonzeros=73728, allocated nonzeros=73728
>                 total number of mallocs used during MatSetValues calls =0
>                   using I-node routines: found 2048 nodes, limit used is 5
>         Mat Object: (firedrake_0_fieldsplit_1_) 1 MPI processes
>           type: seqaij
>           rows=15360, cols=15360
>           total: nonzeros=267264, allocated nonzeros=267264
>           total number of mallocs used during MatSetValues calls =0
>             using I-node routines: found 5120 nodes, limit used is 5
>   linear system matrix followed by preconditioner matrix:
>   Mat Object: (firedrake_0_) 1 MPI processes
>     type: nest
>     rows=21504, cols=21504
>       Matrix object:
>         type=nest, rows=2, cols=2
>         MatNest structure:
>         (0,0) : type=seqaij, rows=15360, cols=15360
>         (0,1) : type=seqaij, rows=15360, cols=6144
>         (1,0) : type=seqaij, rows=6144, cols=15360
>         (1,1) : type=seqaij, rows=6144, cols=6144
>   Mat Object: (firedrake_0_) 1 MPI processes
>     type: nest
>     rows=21504, cols=21504
>       Matrix object:
>         type=nest, rows=2, cols=2
>         MatNest structure:
>         (0,0) : prefix="firedrake_0_fieldsplit_1_", type=seqaij,
> rows=15360, cols=15360
>         (0,1) : type=seqaij, rows=15360, cols=6144
>         (1,0) : type=seqaij, rows=6144, cols=15360
>         (1,1) : prefix="firedrake_0_fieldsplit_0_", type=seqaij,
> rows=6144, cols=6144
>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190318/5cdd1286/attachment-0001.html>


More information about the petsc-users mailing list