[petsc-users] KSP not converge when turn on ksp_monitor_true_residual

Matthew Knepley knepley at gmail.com
Fri Jun 5 10:13:32 CDT 2020


On Fri, Jun 5, 2020 at 10:17 AM Y. Shidi <ys453 at cam.ac.uk> wrote:

> Thank you for your reply.
> I set the relative tolerance of 1e-10, below is the output with
> ksp_monitor_true_residual,
>
> step 1 (converged):
>
>    0 KSP unpreconditioned resid norm 5.853072011788e-01 true resid norm
> 5.853072011788e-01 ||r(i)||/||b|| 1.000000000000e+00
>    1 KSP unpreconditioned resid norm 1.263467051263e-02 true resid norm
> 1.263467051263e-02 ||r(i)||/||b|| 2.158639170539e-02
>    2 KSP unpreconditioned resid norm 1.341051880494e-03 true resid norm
> 1.341051880494e-03 ||r(i)||/||b|| 2.291193202122e-03
>    3 KSP unpreconditioned resid norm 5.281711107599e-04 true resid norm
> 5.281711107599e-04 ||r(i)||/||b|| 9.023827311473e-04
>    4 KSP unpreconditioned resid norm 2.265161267248e-04 true resid norm
> 2.265161267248e-04 ||r(i)||/||b|| 3.870038268256e-04
>    5 KSP unpreconditioned resid norm 6.083108409953e-05 true resid norm
> 6.083108409954e-05 ||r(i)||/||b|| 1.039301822650e-04
>    6 KSP unpreconditioned resid norm 2.049038801905e-05 true resid norm
> 2.049038801906e-05 ||r(i)||/||b|| 3.500792058903e-05
>    7 KSP unpreconditioned resid norm 4.908575378578e-06 true resid norm
> 4.908575378575e-06 ||r(i)||/||b|| 8.386323231099e-06
>    8 KSP unpreconditioned resid norm 1.160371295964e-06 true resid norm
> 1.160371295960e-06 ||r(i)||/||b|| 1.982499606399e-06
>    9 KSP unpreconditioned resid norm 2.993349323170e-07 true resid norm
> 2.993349323154e-07 ||r(i)||/||b|| 5.114150854671e-07
>   10 KSP unpreconditioned resid norm 6.423860393220e-08 true resid norm
> 6.423860392241e-08 ||r(i)||/||b|| 1.097519452913e-07
>   11 KSP unpreconditioned resid norm 1.486343352971e-08 true resid norm
> 1.486343352542e-08 ||r(i)||/||b|| 2.539424339131e-08
>   12 KSP unpreconditioned resid norm 3.388518888745e-09 true resid norm
> 3.388518889406e-09 ||r(i)||/||b|| 5.789299845588e-09
>   13 KSP unpreconditioned resid norm 8.414456498118e-10 true resid norm
> 8.414456551004e-10 ||r(i)||/||b|| 1.437613706795e-09
>   14 KSP unpreconditioned resid norm 2.052670281766e-10 true resid norm
> 2.052670332352e-10 ||r(i)||/||b|| 3.506996545092e-10
>   15 KSP unpreconditioned resid norm 4.739342260680e-11 true resid norm
> 4.739342219912e-11 ||r(i)||/||b|| 8.097187614243e-11
> Linear solve converged due to CONVERGED_RTOL iterations 15
> KSP Object: 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, nonzero initial guess
>    tolerances:  relative=1e-10, absolute=1e-50, divergence=10000.
>    right preconditioning
>    using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>    type: fieldsplit
>      FieldSplit with Schur preconditioner, factorization FULL
>      Preconditioner for the Schur complement formed from Sp, an assembled
> approximation to S, which uses A00's diagonal's inverse
>      Split info:
>      Split number 0 Defined by IS
>      Split number 1 Defined by IS
>      KSP solver for A00 block
>        KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix = precond matrix:
>          Mat Object: (fieldsplit_0_) 1 MPI processes
>            type: seqaij
>            rows=28578, cols=28578
>            total: nonzeros=653916, allocated nonzeros=653916
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>      KSP solver for S = A11 - A10 inv(A00) A01
>        KSP Object: (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: (fieldsplit_1_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix followed by preconditioner matrix:
>          Mat Object: (fieldsplit_1_) 1 MPI processes
>            type: schurcomplement
>            rows=3874, cols=3874
>              Schur complement A11 - A10 inv(A00) A01
>              A11
>                Mat Object: (fieldsplit_1_) 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=3874
>                  total: nonzeros=0, allocated nonzeros=0
>                  total number of mallocs used during MatSetValues calls
> =0
>                    using I-node routines: found 775 nodes, limit used is
> 5
>              A10
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=28578
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>              KSP of A00
>                KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>                  type: hypre
>                    HYPRE BoomerAMG preconditioning
>                      Cycle type V
>                      Maximum number of levels 25
>                      Maximum number of iterations PER hypre call 1
>                      Convergence tolerance PER hypre call 0.
>                      Threshold for strong coupling 0.25
>                      Interpolation truncation factor 0.
>                      Interpolation: max elements per row 0
>                      Number of levels of aggressive coarsening 0
>                      Number of paths for aggressive coarsening 1
>                      Maximum row sums 0.9
>                      Sweeps down         1
>                      Sweeps up           1
>                      Sweeps on coarse    1
>                      Relax down          symmetric-SOR/Jacobi
>                      Relax up            symmetric-SOR/Jacobi
>                      Relax on coarse     Gaussian-elimination
>                      Relax weight  (all)      1.
>                      Outer relax weight (all) 1.
>                      Using CF-relaxation
>                      Not using more complex smoothers.
>                      Measure type        local
>                      Coarsen type        Falgout
>                      Interpolation type  classical
>                  linear system matrix = precond matrix:
>                  Mat Object: (fieldsplit_0_) 1 MPI processes
>                    type: seqaij
>                    rows=28578, cols=28578
>                    total: nonzeros=653916, allocated nonzeros=653916
>                    total number of mallocs used during MatSetValues calls
> =0
>                      not using I-node routines
>              A01
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=28578, cols=3874
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>          Mat Object: 1 MPI processes
>            type: seqaij
>            rows=3874, cols=3874
>            total: nonzeros=77704, allocated nonzeros=77704
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>    linear system matrix = precond matrix:
>    Mat Object: 1 MPI processes
>      type: nest
>      rows=32452, cols=32452
>        Matrix object:
>          type=nest, rows=2, cols=2
>          MatNest structure:
>          (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=28578,
> cols=28578
>          (0,1) : type=seqaij, rows=28578, cols=3874
>          (1,0) : type=seqaij, rows=3874, cols=28578
>          (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=3874,
> cols=3874
>
> step 2 (NOT converged):
>
> Linear solve did not converge due to DIVERGED_NANORINF iterations 0
>

It looks like your matrix is invalid.

  Thanks,

     Matt


> KSP Object: 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, nonzero initial guess
>    tolerances:  relative=1e-10, absolute=1e-50, divergence=10000.
>    right preconditioning
>    using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>    type: fieldsplit
>      FieldSplit with Schur preconditioner, factorization FULL
>      Preconditioner for the Schur complement formed from Sp, an assembled
> approximation to S, which uses A00's diagonal's inverse
>      Split info:
>      Split number 0 Defined by IS
>      Split number 1 Defined by IS
>      KSP solver for A00 block
>        KSP Object: (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 DEFAULT norm type for convergence test
>        PC Object: (fieldsplit_0_) 1 MPI processes
>          type: hypre
>          PC has not been set up so information may be incomplete
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix = precond matrix:
>          Mat Object: (fieldsplit_0_) 1 MPI processes
>            type: seqaij
>            rows=28578, cols=28578
>            total: nonzeros=653916, allocated nonzeros=653916
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>      KSP solver for S = A11 - A10 inv(A00) A01
>        KSP Object: (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 DEFAULT norm type for convergence test
>        PC Object: (fieldsplit_1_) 1 MPI processes
>          type: hypre
>          PC has not been set up so information may be incomplete
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix followed by preconditioner matrix:
>          Mat Object: (fieldsplit_1_) 1 MPI processes
>            type: schurcomplement
>            rows=3874, cols=3874
>              Schur complement A11 - A10 inv(A00) A01
>              A11
>                Mat Object: (fieldsplit_1_) 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=3874
>                  total: nonzeros=0, allocated nonzeros=0
>                  total number of mallocs used during MatSetValues calls
> =0
>                    using I-node routines: found 775 nodes, limit used is
> 5
>              A10
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=28578
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>              KSP of A00
>                KSP Object: (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 DEFAULT norm type for convergence test
>                PC Object: (fieldsplit_0_) 1 MPI processes
>                  type: hypre
>                  PC has not been set up so information may be incomplete
>                    HYPRE BoomerAMG preconditioning
>                      Cycle type V
>                      Maximum number of levels 25
>                      Maximum number of iterations PER hypre call 1
>                      Convergence tolerance PER hypre call 0.
>                      Threshold for strong coupling 0.25
>                      Interpolation truncation factor 0.
>                      Interpolation: max elements per row 0
>                      Number of levels of aggressive coarsening 0
>                      Number of paths for aggressive coarsening 1
>                      Maximum row sums 0.9
>                      Sweeps down         1
>                      Sweeps up           1
>                      Sweeps on coarse    1
>                      Relax down          symmetric-SOR/Jacobi
>                      Relax up            symmetric-SOR/Jacobi
>                      Relax on coarse     Gaussian-elimination
>                      Relax weight  (all)      1.
>                      Outer relax weight (all) 1.
>                      Using CF-relaxation
>                      Not using more complex smoothers.
>                      Measure type        local
>                      Coarsen type        Falgout
>                      Interpolation type  classical
>                  linear system matrix = precond matrix:
>                  Mat Object: (fieldsplit_0_) 1 MPI processes
>                    type: seqaij
>                    rows=28578, cols=28578
>                    total: nonzeros=653916, allocated nonzeros=653916
>                    total number of mallocs used during MatSetValues calls
> =0
>                      not using I-node routines
>              A01
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=28578, cols=3874
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>          Mat Object: 1 MPI processes
>            type: seqaij
>            rows=3874, cols=3874
>            total: nonzeros=77704, allocated nonzeros=77704
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>    linear system matrix = precond matrix:
>    Mat Object: 1 MPI processes
>      type: nest
>      rows=32452, cols=32452
>        Matrix object:
>          type=nest, rows=2, cols=2
>          MatNest structure:
>          (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=28578,
> cols=28578
>          (0,1) : type=seqaij, rows=28578, cols=3874
>          (1,0) : type=seqaij, rows=3874, cols=28578
>          (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=3874,
> cols=3874
>
> HOWEVER, if I do not use -ksp_monitor_true_residual, both two
> steps converge.
>
> step 1:
> Linear solve converged due to CONVERGED_RTOL iterations 15
> KSP Object: 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, nonzero initial guess
>    tolerances:  relative=1e-10, absolute=1e-50, divergence=10000.
>    right preconditioning
>    using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>    type: fieldsplit
>      FieldSplit with Schur preconditioner, factorization FULL
>      Preconditioner for the Schur complement formed from Sp, an assembled
> approximation to S, which uses A00's diagonal's inverse
>      Split info:
>      Split number 0 Defined by IS
>      Split number 1 Defined by IS
>      KSP solver for A00 block
>        KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix = precond matrix:
>          Mat Object: (fieldsplit_0_) 1 MPI processes
>            type: seqaij
>            rows=28578, cols=28578
>            total: nonzeros=653916, allocated nonzeros=653916
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>      KSP solver for S = A11 - A10 inv(A00) A01
>        KSP Object: (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: (fieldsplit_1_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix followed by preconditioner matrix:
>          Mat Object: (fieldsplit_1_) 1 MPI processes
>            type: schurcomplement
>            rows=3874, cols=3874
>              Schur complement A11 - A10 inv(A00) A01
>              A11
>                Mat Object: (fieldsplit_1_) 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=3874
>                  total: nonzeros=0, allocated nonzeros=0
>                  total number of mallocs used during MatSetValues calls
> =0
>                    using I-node routines: found 775 nodes, limit used is
> 5
>              A10
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=28578
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>              KSP of A00
>                KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>                  type: hypre
>                    HYPRE BoomerAMG preconditioning
>                      Cycle type V
>                      Maximum number of levels 25
>                      Maximum number of iterations PER hypre call 1
>                      Convergence tolerance PER hypre call 0.
>                      Threshold for strong coupling 0.25
>                      Interpolation truncation factor 0.
>                      Interpolation: max elements per row 0
>                      Number of levels of aggressive coarsening 0
>                      Number of paths for aggressive coarsening 1
>                      Maximum row sums 0.9
>                      Sweeps down         1
>                      Sweeps up           1
>                      Sweeps on coarse    1
>                      Relax down          symmetric-SOR/Jacobi
>                      Relax up            symmetric-SOR/Jacobi
>                      Relax on coarse     Gaussian-elimination
>                      Relax weight  (all)      1.
>                      Outer relax weight (all) 1.
>                      Using CF-relaxation
>                      Not using more complex smoothers.
>                      Measure type        local
>                      Coarsen type        Falgout
>                      Interpolation type  classical
>                  linear system matrix = precond matrix:
>                  Mat Object: (fieldsplit_0_) 1 MPI processes
>                    type: seqaij
>                    rows=28578, cols=28578
>                    total: nonzeros=653916, allocated nonzeros=653916
>                    total number of mallocs used during MatSetValues calls
> =0
>                      not using I-node routines
>              A01
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=28578, cols=3874
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>          Mat Object: 1 MPI processes
>            type: seqaij
>            rows=3874, cols=3874
>            total: nonzeros=77704, allocated nonzeros=77704
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>    linear system matrix = precond matrix:
>    Mat Object: 1 MPI processes
>      type: nest
>      rows=32452, cols=32452
>        Matrix object:
>          type=nest, rows=2, cols=2
>          MatNest structure:
>          (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=28578,
> cols=28578
>          (0,1) : type=seqaij, rows=28578, cols=3874
>          (1,0) : type=seqaij, rows=3874, cols=28578
>          (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=3874,
> cols=3874
>
> step 2:
>
> Linear solve converged due to CONVERGED_RTOL iterations 1
> KSP Object: 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, nonzero initial guess
>    tolerances:  relative=1e-10, absolute=1e-50, divergence=10000.
>    right preconditioning
>    using UNPRECONDITIONED norm type for convergence test
> PC Object: 1 MPI processes
>    type: fieldsplit
>      FieldSplit with Schur preconditioner, factorization FULL
>      Preconditioner for the Schur complement formed from Sp, an assembled
> approximation to S, which uses A00's diagonal's inverse
>      Split info:
>      Split number 0 Defined by IS
>      Split number 1 Defined by IS
>      KSP solver for A00 block
>        KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix = precond matrix:
>          Mat Object: (fieldsplit_0_) 1 MPI processes
>            type: seqaij
>            rows=28578, cols=28578
>            total: nonzeros=653916, allocated nonzeros=653916
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>      KSP solver for S = A11 - A10 inv(A00) A01
>        KSP Object: (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: (fieldsplit_1_) 1 MPI processes
>          type: hypre
>            HYPRE BoomerAMG preconditioning
>              Cycle type V
>              Maximum number of levels 25
>              Maximum number of iterations PER hypre call 1
>              Convergence tolerance PER hypre call 0.
>              Threshold for strong coupling 0.25
>              Interpolation truncation factor 0.
>              Interpolation: max elements per row 0
>              Number of levels of aggressive coarsening 0
>              Number of paths for aggressive coarsening 1
>              Maximum row sums 0.9
>              Sweeps down         1
>              Sweeps up           1
>              Sweeps on coarse    1
>              Relax down          symmetric-SOR/Jacobi
>              Relax up            symmetric-SOR/Jacobi
>              Relax on coarse     Gaussian-elimination
>              Relax weight  (all)      1.
>              Outer relax weight (all) 1.
>              Using CF-relaxation
>              Not using more complex smoothers.
>              Measure type        local
>              Coarsen type        Falgout
>              Interpolation type  classical
>          linear system matrix followed by preconditioner matrix:
>          Mat Object: (fieldsplit_1_) 1 MPI processes
>            type: schurcomplement
>            rows=3874, cols=3874
>              Schur complement A11 - A10 inv(A00) A01
>              A11
>                Mat Object: (fieldsplit_1_) 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=3874
>                  total: nonzeros=0, allocated nonzeros=0
>                  total number of mallocs used during MatSetValues calls
> =0
>                    using I-node routines: found 775 nodes, limit used is
> 5
>              A10
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=3874, cols=28578
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>              KSP of A00
>                KSP Object: (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: (fieldsplit_0_) 1 MPI processes
>                  type: hypre
>                    HYPRE BoomerAMG preconditioning
>                      Cycle type V
>                      Maximum number of levels 25
>                      Maximum number of iterations PER hypre call 1
>                      Convergence tolerance PER hypre call 0.
>                      Threshold for strong coupling 0.25
>                      Interpolation truncation factor 0.
>                      Interpolation: max elements per row 0
>                      Number of levels of aggressive coarsening 0
>                      Number of paths for aggressive coarsening 1
>                      Maximum row sums 0.9
>                      Sweeps down         1
>                      Sweeps up           1
>                      Sweeps on coarse    1
>                      Relax down          symmetric-SOR/Jacobi
>                      Relax up            symmetric-SOR/Jacobi
>                      Relax on coarse     Gaussian-elimination
>                      Relax weight  (all)      1.
>                      Outer relax weight (all) 1.
>                      Using CF-relaxation
>                      Not using more complex smoothers.
>                      Measure type        local
>                      Coarsen type        Falgout
>                      Interpolation type  classical
>                  linear system matrix = precond matrix:
>                  Mat Object: (fieldsplit_0_) 1 MPI processes
>                    type: seqaij
>                    rows=28578, cols=28578
>                    total: nonzeros=653916, allocated nonzeros=653916
>                    total number of mallocs used during MatSetValues calls
> =0
>                      not using I-node routines
>              A01
>                Mat Object: 1 MPI processes
>                  type: seqaij
>                  rows=28578, cols=3874
>                  total: nonzeros=138180, allocated nonzeros=138180
>                  total number of mallocs used during MatSetValues calls
> =0
>                    not using I-node routines
>          Mat Object: 1 MPI processes
>            type: seqaij
>            rows=3874, cols=3874
>            total: nonzeros=77704, allocated nonzeros=77704
>            total number of mallocs used during MatSetValues calls =0
>              not using I-node routines
>    linear system matrix = precond matrix:
>    Mat Object: 1 MPI processes
>      type: nest
>      rows=32452, cols=32452
>        Matrix object:
>          type=nest, rows=2, cols=2
>          MatNest structure:
>          (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=28578,
> cols=28578
>          (0,1) : type=seqaij, rows=28578, cols=3874
>          (1,0) : type=seqaij, rows=3874, cols=28578
>          (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=3874,
> cols=3874
>
> Thank you for your help and time.
>
> Kind regards,
> Shidi
>
> On 2020-06-05 13:40, Matthew Knepley wrote:
> > On Fri, Jun 5, 2020 at 7:42 AM Y. Shidi <ys453 at cam.ac.uk> wrote:
> >
> >>> This algebraic tolerance is strange. Why put it below unit
> >> roundoff?
> >>
> >> We are testing a static bubble problem, if the tolerance is 1e-10 it
> >> won't approach static solution.
> >
> > I do not understand. Every computation in your problem makes errors of
> > order 1e-15,
> > since that is the roundoff. Moreover, you cannot expect a solution to
> > be more accurate
> > than (condition number * roundoff). I would be surprised if the
> > condition number is near 1,
> > since you would not need a solver. So what is going on?
> >
> >   Thanks,
> >
> >      Matt
> >
> >> Thanks,
> >> Shidi
> >>
> >> On 2020-06-05 11:48, Matthew Knepley wrote:
> >>> On Fri, Jun 5, 2020 at 6:12 AM Y. Shidi <ys453 at cam.ac.uk> wrote:
> >>>
> >>>> Dear developers,
> >>>>
> >>>> We are using filed splitting method solve saddle point problem.
> >>>> Below is our ksp configuration:
> >>>>
> >>>> PetscOptionsSetValue(NULL,"-ksp_type","gmres");
> >>>> PetscOptionsSetValue(NULL,"-ksp_initial_guess_nonzero","");
> >>>> PetscOptionsSetValue(NULL,"-ksp_rtol","1e-20");
> >>>
> >>> This algebraic tolerance is strange. Why put it below unit
> >> roundoff?
> >>>
> >>>> PetscOptionsSetValue(NULL,"-ksp_converged_reason", "");
> >>>> PetscOptionsSetValue(NULL,"-pc_type", "fieldsplit");
> >>>> PetscOptionsSetValue(NULL,"-pc_fieldsplit_type", "schur");
> >>>> PetscOptionsSetValue(NULL,"-pc_fieldsplit_schur_fact_type",
> >>>> "lower");
> >>>> PetscOptionsSetValue(NULL,"-pc_fieldsplit_schur_precondition",
> >>>> "selfp");
> >>>> PetscOptionsSetValue(NULL,"-fieldsplit_0_ksp_type", "preonly");
> >>>> PetscOptionsSetValue(NULL,"-fieldsplit_0_pc_type", "hypre");
> >>>> PetscOptionsSetValue(NULL,"-fieldsplit_1_ksp_type", "preonly");
> >>>> PetscOptionsSetValue(NULL,"-fieldsplit_1_pc_type", "hypre");
> >>>>
> >>>> It can be converged. However, when we turn on the ksp monitor:
> >>>>
> >>>> PetscOptionsSetValue(NULL,"-ksp_monitor_true_residual", "");
> >>>>
> >>>> It won't converge. It seems very odd, I am wondering why it
> >> happens.
> >>>
> >>> Can you try it with a relative tolerance of 1e-10, and send the
> >> output
> >>> (with -ksp_view as well).
> >>>
> >>> Thanks,
> >>>
> >>> Matt
> >>>
> >>>> Thank you for your time and help.
> >>>>
> >>>> Kind regards,
> >>>> Shidi
> >>>
> >>> --
> >>>
> >>> What most experimenters take for granted before they begin their
> >>> experiments is infinitely more interesting than any results to
> >> which
> >>> their experiments lead.
> >>> -- Norbert Wiener
> >>>
> >>> https://www.cse.buffalo.edu/~knepley/ [1]
> >>>
> >>>
> >>> Links:
> >>> ------
> >>> [1] http://www.cse.buffalo.edu/~knepley/
> >
> > --
> >
> > What most experimenters take for granted before they begin their
> > experiments is infinitely more interesting than any results to which
> > their experiments lead.
> > -- Norbert Wiener
> >
> > https://www.cse.buffalo.edu/~knepley/ [1]
> >
> >
> > Links:
> > ------
> > [1] http://www.cse.buffalo.edu/~knepley/
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200605/6d1bf85b/attachment-0001.html>


More information about the petsc-users mailing list