[petsc-users] PCFieldSplit gives different results for direct and iterative solver

Smith, Barry F. bsmith at mcs.anl.gov
Sat Mar 16 19:05:42 CDT 2019



> On Mar 16, 2019, at 6:50 PM, Y. Shidi via petsc-users <petsc-users at mcs.anl.gov> wrote:
> 
> Hello,
> 
> I am trying to solve the incompressible n-s equations by
> PCFieldSplit.
> 
> The large matrix and vectors are formed by MatCreateNest()
> and VecCreateNest().
> The system is solved directly by the following command:
>    -ksp_type fgmres \
>    -pc_type fieldsplit \
>    -pc_fieldsplit_type schur \
>    -pc_fieldsplit_schur_fact_type full \
>    -ksp_converged_reason \
>    -ksp_monitor_true_residual \
>    -fieldsplit_0_ksp_type preonly \
>    -fieldsplit_0_pc_type cholesky \
>    -fieldsplit_0_pc_factor_mat_solver_package mumps \
>    -mat_mumps_icntl_28 2 \
>    -mat_mumps_icntl_29 2 \
>    -fieldsplit_1_ksp_type preonly \
>    -fieldsplit_1_pc_type jacobi \
> Output:
>  0 KSP unpreconditioned resid norm 1.214252932161e+04 true resid norm 1.214252932161e+04 ||r(i)||/||b|| 1.000000000000e+00
>  1 KSP unpreconditioned resid norm 1.642782495109e-02 true resid norm 1.642782495109e-02 ||r(i)||/||b|| 1.352916226594e-06
> Linear solve converged due to CONVERGED_RTOL iterations 1
> 
> The system is solved iteratively by the following command:
>    -ksp_type fgmres \
>    -pc_type fieldsplit \
>    -pc_fieldsplit_type schur \
>    -pc_fieldsplit_schur_factorization_type diag \
>    -ksp_converged_reason \
>    -ksp_monitor_true_residual \
>    -fieldsplit_0_ksp_type preonly \
>    -fieldsplit_0_pc_type gamg \
>    -fieldsplit_1_ksp_type minres \
>    -fieldsplit_1_pc_type none \
> Output:
>  0 KSP unpreconditioned resid norm 1.214252932161e+04 true resid norm 1.214252932161e+04 ||r(i)||/||b|| 1.000000000000e+00
>  1 KSP unpreconditioned resid norm 2.184037364915e+02 true resid norm 2.184037364915e+02 ||r(i)||/||b|| 1.798667565109e-02
>  2 KSP unpreconditioned resid norm 2.120097409539e+02 true resid norm 2.120097409635e+02 ||r(i)||/||b|| 1.746009709742e-02
>  3 KSP unpreconditioned resid norm 4.364091658268e+01 true resid norm 4.364091658575e+01 ||r(i)||/||b|| 3.594054865332e-03
>  4 KSP unpreconditioned resid norm 2.632671796885e+00 true resid norm 2.632671797020e+00 ||r(i)||/||b|| 2.168141189773e-04
>  5 KSP unpreconditioned resid norm 2.209213998004e+00 true resid norm 2.209213980361e+00 ||r(i)||/||b|| 1.819401808180e-04
>  6 KSP unpreconditioned resid norm 4.683775185840e-01 true resid norm 4.683775085753e-01 ||r(i)||/||b|| 3.857330677735e-05
>  7 KSP unpreconditioned resid norm 3.042503284736e-02 true resid norm 3.042503349258e-02 ||r(i)||/||b|| 2.505658638883e-06
> 
> 
> Both methods give answers, but they are different

   What do you mean the answers are different? Do you mean the solution x from KSPSolve() is different? How are you calculating their difference and how different are they?

    Since the solutions are only approximate; true residual norm is around 1.642782495109e-02 and 3.042503349258e-02  for the two different solvers there will only be a certain number of identical digits in the two solutions (which depends on the condition number of the original matrix). You can run both solvers with -ksp_rtol 1.e-12 and then (assuming everything is working correctly) the two solutions will be much closer to each other.

   Barry

> so I am wondering
> if it is possible that you can help me figure out which part I am
> doing wrong.
> 
> Thank you for your time.
> 
> Kind Regards,
> Shidi



More information about the petsc-users mailing list