[petsc-users] TS question 1: how to stop explicit methods because they do not use SNES(VI)?

Matthew Knepley knepley at gmail.com
Tue Feb 14 06:53:06 CST 2017


On Mon, Feb 13, 2017 at 3:59 PM, Zhang, Hong <hongzhang at anl.gov> wrote:

>
> On Feb 13, 2017, at 3:08 PM, Ed Bueler <elbueler at alaska.edu> wrote:
>
> Barry --
>
> >   Sounds like a bug to me. The methods should be checking if an
> > IFunction is being provided and error out in that case.
>
> I don't think it is that simple.  I.e. having an IFunction and an explicit
> scheme is not, by itself, a bug.  I think.
>
> If the user has provided IFunction
>    F(t,u,u_t) = u_t - f(t,u)
> which is the convenient form for e.g. TSARKIMEX,
> and RHSFunction
>   G(t,u)
> then I guess I assumed that explicit methods like TSRK should,for
> unconstrained cases, get their RHS by callback like this:
>   g(t,u) = - F(t,u,0) + G(t,u)
> so that the ODE is in form
>   u_t = g(t,u) = f(t,u) + g(t,u)
>
>
> PETSc does not transform IFunction to the RHS this way, and PETSc is not
> supposed to do it automatically. Note that IFunction may involve a mass
> matrix, e.g. F(t,u,u_t) = M*u_t - f(t,u) and sometimes M is not
> invertible.
>
>
> I think that is the behavior I am seeing (i.e. on my problem, using PETSc
> master).
>
>
> Explicit methods use only RHSFunction and ignore IFunction, so in your
> case, if you change TS type to rk and ssp at run time, you are actually
> solving u_t = G(t,u). If RHSFunction is not provided, PETSc will assume the
> RHS is zero (u_t=0).
>

I do not agree with this. I thought our aim was to support users comparing
explicit with implicit with semi-implicit. This breaks that model
completely.

   Matt


>
> The "nonsense" I am referring to is only from the non-enforcement of the
> constraint; it would be o.k. for a pure ODE.
>
> I would love to have some projection mechanism to try, for comparing
> explicit methods with projection to the SNESVI way (i.e. the right way),
> but that is asking for a lot of PETSc refactoring, I think.  For now I just
> want to error-out if the method is not going to call the SNES at each time
> step.
>
>
> Check if TSGetSNES() returns NULL in your code?
>
> Ed
>
>
>
> On Mon, Feb 13, 2017 at 11:26 AM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
>>
>> > On Feb 13, 2017, at 1:53 PM, Ed Bueler <elbueler at alaska.edu> wrote:
>> >
>> > Dear Petsc --
>> >
>> > This is the first of two short TS usage questions.
>> >
>> > My problem is both stiff (diffusive PDE) and constrained, so I require
>> >
>> >    -snes_type vinewton{rs|ss}ls
>> >
>> > *and* I split my ODE system into IFunction and RHSFunction
>> >
>> >    F(t,u,u_t) = G(t,u)
>> >
>> > where F(t,u,u_t) = u_t + f(t,u) in my case (i.e. no mass matrix
>> needed), and the stiff part goes in f(t,u).
>> >
>> > With this arrangement TS types beuler, theta, bdf, arkimex all work
>> quite well.  However, the program runs and produces nonsense with type rk
>> and ssp, that is, explicit methods.
>>
>>    Sounds like a bug to me. The methods should be checking if an
>> IFunction is being provided and error out in that case.
>>
>>   Barry
>>
>> >
>> > So my question is, how do I ask the TS (at run time) whether the chosen
>> TS type will or will not call its SNES at each step?  If SNES is not going
>> to be used then I want to SETERRQ and stop.  That is, I want to error-out
>> if the *method* is fully explicit.
>> >
>> > Note the constraints are enforced by the SNESVI, as they should be, not
>> ad hoc projection.  Also, as a technical matter, I cannot require my
>> iterates to be feasible inside my IFunction evaluation because that would
>> break VINEWTONSSLS.
>> >
>> > Neither TSProblemType (mine is NONLINEAR) nor TSEquationType (mine is
>> IMPLICIT I guess) seem to address this?  My problem is indeed nonlinear and
>> has stiff parts, but it is not a DAE and it *is* constrained.
>> >
>> > Thanks!
>> >
>> > Ed
>> >
>> > PS  I'd prefer not to enumerate the existing TS types and error on the
>> bad ones.  It is not nicely-maintainable.
>> >
>> >
>> >
>> > --
>> > Ed Bueler
>> > Dept of Math and Stat and Geophysical Institute
>> > University of Alaska Fairbanks
>> > Fairbanks, AK 99775-6660
>> > 301C Chapman and 410D Elvey
>> > 907 474-7693 and 907 474-7199  (fax 907 474-5394)
>>
>>
>
>
> --
> Ed Bueler
> Dept of Math and Stat and Geophysical Institute
> University of Alaska Fairbanks
> Fairbanks, AK 99775-6660
> 301C Chapman and 410D Elvey
> 907 474-7693 <(907)%20474-7693> and 907 474-7199 <(907)%20474-7199>  (fax 907
> 474-5394 <(907)%20474-5394>)
>
>
>


-- 
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/20170214/65d93443/attachment.html>


More information about the petsc-users mailing list