[petsc-dev] Using PCFieldSplitSetIS

Jed Brown jed at 59a2.org
Wed Mar 16 09:32:17 CDT 2011


On Wed, Mar 16, 2011 at 14:27, Thomas Witkowski <
thomas.witkowski at tu-dresden.de> wrote:

>
>> I'm going to guess that you still have an outer KSP that (in the global
>> norm, rather than the partitioned norm used inside of splits) has a tighter
>> tolerance, therefore it takes a few outer iterations. If you use loose inner
>> tolerances then the preconditioner becomes nonlinear and you'll need to use
>> FGMRES for the outer. When in doubt, run with -ksp_view and show us the
>> results if you don't understand. Additionally, monitoring inner solves
>> separately can be useful, e.g. -fieldsplit_0_ksp_converged_reason
>> -fieldsplit_1_ksp_monitor -ksp_monitor_true_residual.
>>
> Okay, I run my code with the options "-pc_fieldsplit_type schur
> -fieldsplit_interior_ksp_type preonly -fieldsplit_interior_pc_type bjacobi
> -fieldsplit_interior_sub_pc_type lu -fieldsplit_boundary_ksp_monitor
> -ksp_monitor_true_residual -fieldsplit_interior_ksp_converged_reason". The
> splits are named "interior" and "boundary". The ksp output is as follows:
>
>   Linear solve converged due to CONVERGED_ITS iterations 1
>

This is the initial interior solve.


>   Residual norms for fieldsplit_boundary_ solve.
>   0 KSP Residual norm 1.790059331071e-04
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   1 KSP Residual norm 1.237356212928e-04
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   2 KSP Residual norm 7.952220245101e-05
>   Linear solve converged due to CONVERGED_ITS iterations 1
>

[...]

Each iteration on the boundary problem requires an interior solve. We see
alternating boundary residual and interior solves.


>
> [0]  Petsc-Iteration 0: 0.574504
>

You must be printing this somewhere.


>  0 KSP preconditioned resid norm 5.745043818120e-01 true resid norm
> 1.584249437360e-01 ||r(i)||/||b|| 1.000000000000e+00
>

This is the norm of P^{-1} b, upon which all following norms are computed.
Note that this is not actually a "preconditioned residual" because the
initial guess was zero, however. In fact, since P is an exact preconditioner
for your problem, it is the solution of your problem, but GMRES can't check
that until it does one more iteration (the output below). You can use
-ksp_type richardson to cheaply confirm that the norm is small, it should
show a small residual after one iteration. Once you are happy that the inner
solve really is solving your system, you can switch to -ksp_type preonly.


>   Linear solve converged due to CONVERGED_ITS iterations 1
>   Residual norms for fieldsplit_boundary_ solve.
>   0 KSP Residual norm 1.790059331000e-04
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   1 KSP Residual norm 1.237356237744e-04
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   2 KSP Residual norm 7.952220078415e-05
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   3 KSP Residual norm 3.502859285561e-05
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   4 KSP Residual norm 1.601086485810e-05
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   5 KSP Residual norm 8.491169384185e-06
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   6 KSP Residual norm 4.778229157872e-06
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   7 KSP Residual norm 2.571917112191e-06
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   8 KSP Residual norm 1.353977280183e-06
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   9 KSP Residual norm 7.408352317360e-07
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  10 KSP Residual norm 3.810710583336e-07
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  11 KSP Residual norm 1.955898169394e-07
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  12 KSP Residual norm 9.938892788465e-08
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  13 KSP Residual norm 4.893428546840e-08
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  14 KSP Residual norm 2.395120857135e-08
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  15 KSP Residual norm <1.2464367927> <1.2464367927> <1.2464367927>
> 1.246436792784e-08
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  16 KSP Residual norm 6.316409779118e-09
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  17 KSP Residual norm 3.088649355374e-09
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  18 KSP Residual norm 1.593449602302e-09
>   Linear solve converged due to CONVERGED_ITS iterations 1
>

The way GMRES computes residuals explains why we needed this Schur solve,
but it shouldn't do the one below. Are you using GMRES or some other Krylov
solver (perhaps set in an options file or inside your code, -ksp_view would
help)?


>   Linear solve converged due to CONVERGED_ITS iterations 1
>   Residual norms for fieldsplit_boundary_ solve.
>   0 KSP Residual norm 1.593449603155e-09
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   1 KSP Residual norm 1.011857676662e-09
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   2 KSP Residual norm 5.660892669788e-10
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   3 KSP Residual norm 3.046799343903e-10
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   4 KSP Residual norm <1.6004713993> <1.6004713993> <1.6004713993>
> 1.600471399329e-10
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   5 KSP Residual norm 8.083216601016e-11
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   6 KSP Residual norm 4.308671648345e-11
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   7 KSP Residual norm 2.381959817139e-11
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   8 KSP Residual norm 1.264839790114e-11
>   Linear solve converged due to CONVERGED_ITS iterations 1
>   9 KSP Residual norm 6.435098853282e-12
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  10 KSP Residual norm 3.184769579797e-12
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  11 KSP Residual norm <1.5022019430> <1.5022019430> <1.5022019430>
> 1.502201943016e-12
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  12 KSP Residual norm 7.099486208425e-13
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  13 KSP Residual norm 3.443222176388e-13
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  14 KSP Residual norm <1.6478434989> <1.6478434989> <1.6478434989>
> 1.647843498912e-13
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  15 KSP Residual norm 8.208508173926e-14
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  16 KSP Residual norm 4.376604355287e-14
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  17 KSP Residual norm 2.596069152755e-14
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  18 KSP Residual norm 1.335744434252e-14
>   Linear solve converged due to CONVERGED_ITS iterations 1
>  1 KSP preconditioned resid norm 2.608069232839e-13 true resid norm
> 1.449869994655e-13 ||r(i)||/||b|| 9.151778504469e-13
>
> I have no idea how to interpret this output! Could you help me with it?
>
> Thomas
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20110316/947d8479/attachment.html>


More information about the petsc-dev mailing list