[petsc-users] SNES with constraint diverges

Josh L ysjosh.lo at gmail.com
Wed Sep 12 16:31:16 CDT 2018


Barry,

0 represents fully damaged and 1 is intact. if I unload it, the crack
shouldn't heal(solution goes back to 1), so I set the constraint as
    0 <= u^{i+1} <= u^{i}
    i is either load step or time step.



Matt,

Yes, that is what I am going to do next. a linear penalty gives a quadratic
distribution instead of the exponential distribution given by quadratic
penalty term.


Thanks,
Josh



2018-09-12 16:23 GMT-05:00 Matthew Knepley <knepley at gmail.com>:

> On Wed, Sep 12, 2018 at 5:16 PM Smith, Barry F. <bsmith at mcs.anl.gov>
> wrote:
>
>>
>>    The function norm given by the -snes_monitor is the function norm for
>> the DVI problem (which depends on where the solution is constrained) while
>> FormFunction() norm you compute is for the original unconstrained problem.
>>
>>     I am still not sure why you solve this as a DVI. If you solve it
>> without the constraints do you get the "wrong" answer?
>>
>
> If this is a variational crack problem, you can solve it unconstrained if
> you use a quadratic penalty term. However,
> that induces spurious long range communication between crack tips. If you
> bound the phase field between 0 and 1,
> you can use a linear penalty which has physical behavior.
>
>    Matt
>
>
>>     Barry
>>
>>
>> > On Sep 12, 2018, at 3:32 PM, Josh L <ysjosh.lo at gmail.com> wrote:
>> >
>> > Hi,
>> >
>> > my solution is mostly very close to 1. only  for a very small area
>> where solution goes from 0 to 1(a smeared crack).
>> >
>> > I set -snes_atol 1e-7 and it is converging.
>> >
>> > I've noticed the following:
>> >
>> > There is a difference between the function norm.
>> >
>> > I calculate the function norm in FormFunction, so every time it is
>> called it gives the function norm
>> > , and the result is different from the function norm given by
>> -snes_monitor if i set
>> >
>> >      SNESSetType(snes,SNESVINEWTONRSLS,ierr)
>> >      SNESVISetVariableBounds(snes,xl,xu,ierr)
>> >
>> > The function norm calculated in FromFunction is NOT reducing, however,
>> the function norm given by -snes_monitor is reducing
>> > They are the same if I just use regular SNES without setting variable
>> bounds.
>> >
>> >
>> > Thanks,
>> > Josh
>> >
>> > 2018-09-12 12:02 GMT-05:00 Smith, Barry F. <bsmith at mcs.anl.gov>:
>> >
>> >      You have too tight a convergence tolerance for your problem.  You
>> can't expect to get more than 1.e-12 as the minimum residual norm or even
>> less.
>> >
>> >       How close is your solution to 1 and -1?
>> >
>> >       If you really need much higher convergence you can try
>> ./configure --with-precision=__float128
>> >
>> >
>> >       Barry
>> >
>> > > On Sep 11, 2018, at 11:53 PM, Josh L <ysjosh.lo at gmail.com> wrote:
>> > >
>> > > Yes, I initialize all u_i to 1.0
>> > >
>> > >
>> > >
>> > > 2018-09-11 23:37 GMT-05:00 Smith, Barry F. <bsmith at mcs.anl.gov>:
>> > >
>> > >    Do you start with initial conditions of  0 <= u_i <= 1 ?
>> > >
>> > >     Run with -snes_monitor -snes_converged_reason
>> -ksp_monitor_true_residual -info -snes_linesearch_monitor and send all the
>> output
>> > >
>> > >   Barry
>> > >
>> > >
>> > > > On Sep 11, 2018, at 11:33 PM, Josh L <ysjosh.lo at gmail.com> wrote:
>> > > >
>> > > > Hi,
>> > > >
>> > > > I am using SNES to solve an nonlinear equation f(u), and I know all
>> the u_i should be 0 and 1.
>> > > >
>> > > > First, I use SNES without constraint, and it converges.
>> > > >
>> > > > But, If I set
>> > > >      SNESSetType(snes,SNESVINEWTONRSLS,ierr)
>> > > >      SNESVISetVariableBounds(snes,xl,xu,ierr)
>> > > >
>> > > > where xl and xu is vector, and xl_i=0 and xu_i=1
>> > > >
>> > > > then SNES fails to converge, because linesearch fails(snes reason =
>> -6), and the norm of residual is not reducing(the norm of incremental
>> solution is reducing)
>> > > >
>> > > > The reason to add constraint is that I want to implement some
>> irreversibility.
>> > > >
>> > > >
>> > > > Thanks,
>> > > > Josh
>> > > >
>> > >
>> > >
>> >
>> >
>>
>>
>
> --
> 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/20180912/b6bfcaf1/attachment.html>


More information about the petsc-users mailing list