[petsc-users] Strongly nonlinear equation solved within the framework of PETSc

Matthew Knepley knepley at gmail.com
Tue Jan 10 13:26:07 CST 2012


On Tue, Jan 10, 2012 at 2:27 AM, Fatcharm <wolfshow at gmail.com> wrote:

> Dear all,
>
> First, I would like to re-describe my problem. I want to numerically
> solve a strongly nonlinear fourth-order equation, which is used to
> describe the dynamics of a liquid film.Please find the form of the
> equation below,
>
> u_t = -(1/3)*C*(u^3*u_xxx)_x + (A*(u_x/u))_x
>
> "u" is the thickness of the film to be solved, is a function of "x"
> and time "t" , "C" and "A" are constant parameters.
> u_xxx is the 3-th order derivative.
>
> I wrote a PETSc programs for this problem, using the central finite
> difference scheme in the space and CN method in time .
>
> I start my PETSc program from the
> /petsc-3.2-p5/src/ts/examples/tutorials/ex13.c
>
> I wrote in my "RHSFunction" function:
>        u       = uarray[i];
>        ux      = (uarray[i+1] - uarray[i-1]);
>        uxx     = (-2.0*u + uarray[i-1] + uarray[i+1]);
>        uxxx    = (uarray[i+2] - 2.0*uarray[i+1] + 2.0*uarray[i-1] -
> uarray[i-2]);
>        uxxxx   = (uarray[i+2] - 4.0*uarray[i+1] + 6.0*u -
> 4.0*uarray[i-1] + uarray[i-2]);
>        ucx    = -(user->c/3.0)*sx*(0.75*sx*u*u*ux*uxxx + sx*u*u*u*uxxxx);
>        uax    = (user->a)*sx*(-0.25*ux*ux/(u*u+l_res) + uxx/(u+l_res));
>        f[i]   = ucx + uax;
>
> Also I provided the "RHSJacobian" to evaluate the changing Jacobian.
>
> Followed Jed's advice, I run the program with -snes_monitor
> -snes_converged_reason -ksp_converged_reason
>
> I found that "Nonlinear solve did not converge due to
> DIVERGED_LINE_SEARCH", at the same time "Linear solve converged due to
> CONVERGED_RTOL iterations 5".
>

This sounds like an error in your Jacobian. Run a small problem with
-snes_fd. This
is in the FAQ.


> I tested to solve this problem with -ts_type beuler, I found:
>
> timestep 0: time 0, solution norm 0.00899198, max 0.057735, min 0.001
>    Linear solve converged due to CONVERGED_RTOL iterations 30
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
> timestep 1: time 0.001, solution norm 0.00899198, max 0.057735, min 0.001
>  Nonlinear solve converged due to CONVERGED_FNORM_ABS
> timestep 2: time 0.002, solution norm 0.00899198, max 0.057735, min 0.001
>    Linear solve converged due to CONVERGED_RTOL iterations 30
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
> timestep 3: time 0.003, solution norm 0.00899198, max 0.057735, min 0.001
>  Nonlinear solve converged due to CONVERGED_FNORM_ABS
> timestep 4: time 0.004, solution norm 0.00899198, max 0.057735, min 0.001
>    Linear solve converged due to CONVERGED_RTOL iterations 30
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
> timestep 5: time 0.005, solution norm 0.00899198, max 0.057735, min 0.001
>  Nonlinear solve converged due to CONVERGED_FNORM_ABS
> timestep 6: time 0.006, solution norm 0.00899198, max 0.057735, min 0.001
>    Linear solve converged due to CONVERGED_RTOL iterations 30
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
>
> I read
> http://scicomp.stackexchange.com/questions/30/why-is-newtons-method-not-converging
>
> I run with -pc_type lu, it was told that
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: No support for this operation for this object type!
> [0]PETSC ERROR: Matrix format mpiaij does not have a built-in PETSc LU!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
>

Always run o n1 process first, and never never never never set the type
explicitly to
MPIAIJ. Use AIJ, which falls back to SEQAIJ on 1 processes.

   Matt


> I run with -viewJacobian, the Jacobian looks reasonable. Only the
> value in the Jacobian is on the order of 10^7 or 10^8, Jed told me
> that "That may cause ill-conditioning"?
>
> If I run with -snes_type test -snes_test_display, I got
>
> timestep 0: time 0, solution norm 0.00899198, max 0.057735, min 0.001
> Testing hand-coded Jacobian, if the ratio is
> O(1.e-8), the hand-coded Jacobian is probably correct.
> Run with -snes_test_display to show difference
> of hand-coded and finite difference Jacobian.
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Floating point exception!
> [0]PETSC ERROR: Infinite or not-a-number generated in norm!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
>
> If I run with -snes_ls_monitor, I got
>
> timestep 0: time 0, solution norm 0.02004, max 0.0299803, min 0.0100197
>      Line search: gnorm after quadratic fit 6.651216893041e+04
>      Line search: Quadratically determined step,
> lambda=4.8079877662732573e-01
>      Line search: gnorm after quadratic fit 5.508857435695e+07
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989545851e+07 lambda=1.0000000000000002e-02
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855229151e+07 lambda=1.0000000000000002e-03
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989770948e+07 lambda=1.0000000000000003e-04
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855206653e+07 lambda=1.0000000000000004e-05
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989773203e+07 lambda=1.0000000000000004e-06
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855206428e+07 lambda=1.0000000000000005e-07
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989773226e+07 lambda=1.0000000000000005e-08
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855206425e+07 lambda=1.0000000000000005e-09
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989773226e+07 lambda=1.0000000000000006e-10
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855206425e+07 lambda=1.0000000000000006e-11
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508989773226e+07 lambda=1.0000000000000006e-12
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.508855206425e+07 lambda=1.0000000000000007e-13
>      Line search: unable to find good step length! After 13 tries
>      Line search: fnorm=6.6512168930411834e+04,
> gnorm=5.5088552064252406e+07, ynorm=2.0443556209235136e-01,
> minlambda=6.8680778683552649e-13, lambda=1.0000000000000007e-13,
> initial slope=-4.4238725909386473e+09
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
> timestep 1: time 0.001, solution norm 0.0192692, max 0.021293, min
> 0.0161292
>      Line search: gnorm after quadratic fit 5.509183229235e+07
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262678284e+07 lambda=1.0000000000000002e-02
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184917687e+07 lambda=1.0000000000000002e-03
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262521680e+07 lambda=1.0000000000000003e-04
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184934034e+07 lambda=1.0000000000000004e-05
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262520119e+07 lambda=1.0000000000000004e-06
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184934198e+07 lambda=1.0000000000000005e-07
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262520104e+07 lambda=1.0000000000000005e-08
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184934199e+07 lambda=1.0000000000000005e-09
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262520104e+07 lambda=1.0000000000000006e-10
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184934199e+07 lambda=1.0000000000000006e-11
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509262520104e+07 lambda=1.0000000000000006e-12
>    Line search: Cubic step no good, shrinking lambda, current gnorm
> 5.509184934199e+07 lambda=1.0000000000000007e-13
>      Line search: unable to find good step length! After 13 tries
>      Line search: fnorm=9.5088545581815284e+04,
> gnorm=5.5091849341994494e+07, ynorm=1.4128453377933892e-01,
> minlambda=8.9929091375150137e-13, lambda=1.0000000000000007e-13,
> initial slope=-9.0419264362675228e+09
>  Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
> timestep 2: time 0.002, solution norm 0.0192692, max 0.021293, min
> 0.0161292
>
>
>
>
> It is weird that if I run with -snes_converged_reason, I got
>
> Nonlinear solve did not converge due to DIVERGED_LINE_SEARCH
>
> If I run with  -snes_converged_reason and -info I found that
>
> Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE
>
>
> Could anyone please give me some suggestions?
>
> Thank you very much.
>
> Feng-Chao
>
>
>
>
> >
> > Message: 6
> > Date: Thu, 29 Dec 2011 08:21:08 -0600
> > From: Jed Brown <jedbrown at mcs.anl.gov>
> > Subject: Re: [petsc-users] Strongly nonlinear equation solved within
> >       the framework of PETSc
> > To: PETSc users list <petsc-users at mcs.anl.gov>
> > Message-ID:
> >       <CAM9tzSkb1yrduqHqwtuRhN+bztCixxBjh_KcCV=
> JXWuNc1PNJA at mail.gmail.com>
> > Content-Type: text/plain; charset="utf-8"
> >
> > On Thu, Dec 29, 2011 at 08:10, Fatcharm <wolfshow at gmail.com> wrote:
> >
> >> We can see that the SNES Function norm is extremely large. I think it
> >> is because that the initial value for the unknown function H(X,T) is
> >> quite small and there is  some terms like (1/H)(dH/dX) or
> >> (1/H^2)(dH/dX) in my equations.
> >>
> >
> > That may cause ill-conditioning, but you could still scale the equations
> so
> > that the initial norm was of order 1. It shouldn't matter here though,
> > because most methods are unaffected by scaling.
> >
> > Are you computing an analytic Jacobian or using finite differencing?
> >
> >
> >>
> >> For "Linear solve did not converge due to DIVERGED_DTOL iterations
> >> 3270", does this mean I should change the ksp_type?
> >>
> >
> > It's important to solve the linear systems before worrying about
> > convergence rates for Newton methods. Try a direct solve on a small
> problem
> > first, then read this
> >
> >
> http://scicomp.stackexchange.com/questions/513/why-is-my-iterative-linear-solver-not-converging
> >
> > If you fix the linear solve issues, but SNES is still not converging,
> read
> >
> >
> http://scicomp.stackexchange.com/questions/30/why-is-newtons-method-not-converging
> > -------------- next part --------------
> > An HTML attachment was scrubbed...
> > URL:
> > <
> http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20111229/0e223952/attachment.htm
> >
> >
> > ------------------------------
> >
> > _______________________________________________
> > petsc-users mailing list
> > petsc-users at mcs.anl.gov
> > https://lists.mcs.anl.gov/mailman/listinfo/petsc-users
> >
> >
> > End of petsc-users Digest, Vol 36, Issue 84
> > *******************************************
> >
>



-- 
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
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120110/94cd208f/attachment.htm>


More information about the petsc-users mailing list