snes querie
Barry Smith
bsmith at mcs.anl.gov
Wed Mar 26 05:54:18 CDT 2008
The number one reason Newton's method does not converge is due to
an error
in a supplied Jacobian (how is your Jacobian computed?). You can
*first run with -snes_type test and it will attempt to determine
if your Jacobian is correct
* if that works then run with -snes_mf_operator, what happens
then, does the
Newton converge?
The number two reason Newton's method does not converge is due to
an error in
the function evaluation. This is harder for us to help track down
because it is pure
user code. What I suggest doing is trying to run on the smallest size
problem
possible (use a tiny grid), what happens?
The number three reason Newton's method does not converge is due
to too poor
an initial guess. If doing time-stepping you can try cutting the time-
step back until
you get Newton converged.
Good luck,
Barry
On Mar 25, 2008, at 10:39 PM, Dave Lee wrote:
> Hello PETSc,
>
> I'm trying to solve a non-linear problem using the SNES routines,
> and i've been
> finding that for every even time step the solution converges as I'd
> expect, however
> for every odd time step, the solution diverges. I've run an
> additional test problem
> for which the solution diverges on the second time step only (ie:
> converges for every
> successive time step), so I'm pretty sure its not a odd/even logic
> problem in my
> code. The output (for the initial guess and first two time steps)
> when I run with
> the flags -info -snes_monitor is as follows:
>
>
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 9.34921e-14 is less than relative tolerance 1e-05 times initial
> right hand side norm 48.7135 at iteration 1
> [0] SNESSolve_LS(): iter=0, linear solve iterations=1
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 112.285
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=3.0170114380870206e+01,
> gnorm=1.0142108262591121e-12, ynorm=4.8713530783010192e+01,
> lssucceed=1
> 1 SNES Function norm 1.014210826259e-12
> [0] SNESDefaultConverged(): Converged due to function norm
> 1.01421e-12 < 3.01701e-07 (relative tolerance)
> Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE
> [0] SNESSolve_LS(): iter=0, linear solve iterations=8
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 0.925505
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430941e-02
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430942e-03
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430942e-04
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430944e-05
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430946e-06
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430949e-07
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430952e-08
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430954e-09
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430954e-10
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430954e-11
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430954e-12
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430954e-13
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430956e-14
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430958e-15
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430960e-16
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430961e-17
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430961e-18
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430962e-19
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430962e-20
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430964e-21
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430967e-22
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430969e-23
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430971e-24
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430971e-25
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430974e-26
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430976e-27
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430979e-28
> [0] SNESLineSearchCubic(): Cubic step no good, shrinking lambda,
> lambda=2.6904093002430979e-29
> [0] SNESLineSearchCubic(): Unable to find good step length! 29
> [0] SNESLineSearchCubic(): fnorm=9.6575221913187848e-01,
> gnorm=2.0934668664952425e+00, ynorm=5.1772194669712386e+00,
> lambda=2.6904093002430979e-29, initial slope=-9.3267734652255951e-01
> [0] SNESSolve_LS(): fnorm=9.6575221913187848e-01,
> gnorm=2.0934668664952425e+00, ynorm=5.1772194669712386e+00,
> lssucceed=0
> 1 SNES Function norm 2.093466866495e+00
> [0] SNESLSCheckLocalMin_Private(): || J^T F|| 107.04 near zero
> implies found a local minimum
> Nonlinear solve did not converge due to DIVERGED_LS_FAILURE
> 0 SNES Function norm 2.185496756061e+00
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 2.29435e-06 is less than relative tolerance 1e-05 times initial
> right hand side norm 3.74384 at iteration 8
> [0] SNESSolve_LS(): iter=0, linear solve iterations=8
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 119.216
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=2.1854967560609060e+00,
> gnorm=1.1219829200815006e+00, ynorm=6.0664603876861261e+00,
> lssucceed=1
> 1 SNES Function norm 1.121982920082e+00
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 1.49325e-06 is less than relative tolerance 1e-05 times initial
> right hand side norm 1.54195 at iteration 7
> [0] SNESSolve_LS(): iter=1, linear solve iterations=7
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 1.58793
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=1.1219829200815006e+00,
> gnorm=2.6934556136855723e-02, ynorm=1.5405251272601992e+00,
> lssucceed=1
> 2 SNES Function norm 2.693455613686e-02
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 3.31729e-08 is less than relative tolerance 1e-05 times initial
> right hand side norm 0.0334934 at iteration 7
> [0] SNESSolve_LS(): iter=2, linear solve iterations=7
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 1.36601
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=2.6934556136855723e-02,
> gnorm=6.6982884106362914e-04, ynorm=3.3460040092452174e-02,
> lssucceed=1
> 3 SNES Function norm 6.698288410636e-04
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 7.22802e-10 is less than relative tolerance 1e-05 times initial
> right hand side norm 0.000832349 at iteration 7
> [0] SNESSolve_LS(): iter=3, linear solve iterations=7
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 1.30347
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=6.6982884106362914e-04,
> gnorm=1.7317637296224453e-05, ynorm=8.3135074665742882e-04,
> lssucceed=1
> 4 SNES Function norm 1.731763729622e-05
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 1.77439e-11 is less than relative tolerance 1e-05 times initial
> right hand side norm 2.04787e-05 at iteration 7
> [0] SNESSolve_LS(): iter=4, linear solve iterations=7
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 1.14134
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=1.7317637296224453e-05,
> gnorm=4.5421644145512881e-07, ynorm=2.0457661052256775e-05,
> lssucceed=1
> 5 SNES Function norm 4.542164414551e-07
> [0] KSPDefaultConverged(): Linear solver has converged. Residual
> norm 4.11491e-13 is less than relative tolerance 1e-05 times initial
> right hand side norm 5.14482e-07 at iteration 7
> [0] SNESSolve_LS(): iter=5, linear solve iterations=7
> [0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 1.24632
> near zero implies inconsistent rhs
> [0] SNESLineSearchCubic(): Using full step
> [0] SNESSolve_LS(): fnorm=4.5421644145512881e-07,
> gnorm=1.1954313416752086e-08, ynorm=5.1403324849206651e-07,
> lssucceed=1
> 6 SNES Function norm 1.195431341675e-08
> [0] SNESDefaultConverged(): Converged due to function norm
> 1.19543e-08 < 2.1855e-08 (relative tolerance)
> Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE
>
>
> I'm not sure what to make of the statement:
> SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 112.285 near
> zero implies inconsistent rhs
> which arises for both the converged time steps and the unresolved
> time steps - should I be concerned
> about this?
>
> Also at the point where the SNESComputeFunction call is repeatedly
> made (ls.c:645) during the unresolved
> time steps, the system is not actually being solved, as the
> condition on lambda is not being met, so
> the SNESLineSearchCubic function just keeps calling the residual
> assembly function.
>
> Any ideas what all of this means?
>
> cheers, dave.
>
More information about the petsc-users
mailing list