On 2/4/17 11:00 AM, Matthew Knepley wrote:
On Sat, Feb 4, 2017 at 9:47 AM, Gideon Simpson <[email protected]
<mailto:[email protected]>> wrote:
Would setting it up as a DAE in petsc be algorithmically euivalent
to a projected method (i.e., step of standard RK followed by
nonlinear projection)?
Not quite, if you set it as a DAE, then both the time stepping part and
the algebraic constraint are solved at the same time (you are going to
use an implicit time stepping method.
So for backward Euler and DAE, PETSc solves something like:
y_{n+1} - y{n} - dt* f(y_{n+1}) = 0
g(y_{n+1}) = 0
By projection (assuming index 1 or 2): you are solving
y* - y{n} - dt* f(y*) = 0
then initialize SNES to y* and solve:
g(y_{n+1}) = 0
In some cases g(y_{n+1}) may not have a value (you may not land on the
manifold that is defined by the conservation property). Whereas, if you
define it as a DAE in almost all cases (there are crazy exceptions) you
will.
Emil
I am not sure, as I do not understand those solvers. However, I wrote my
own solver that does exactly that MIMEX.
Matt
-gideon
On Feb 3, 2017, at 11:47 PM, Matthew Knepley <[email protected]
<mailto:[email protected]>> wrote:
That is one answer. Another one is that this particular system is
a DAE and we have methods for that.
Matt
On Fri, Feb 3, 2017 at 8:40 PM, Barry Smith <[email protected]
<mailto:[email protected]>> wrote:
TSSetPostStep(); in your function use TSGetSolution() to get
the current solution.
Please let us know how it works out
Barry
> On Feb 3, 2017, at 7:14 PM, Gideon Simpson
<[email protected] <mailto:[email protected]>>
wrote:
>
> I’m interested in implementing a projection method for an
ODE of the form:
>
> y’ = f(y),
>
> such that g(y) = 0 for all time (i.e., g is conserved).
Note that in a projection method, a standard time step is made
to produce y* from y_{n}, and then this is corrected to obtain
y_{n+1} satisfying g(y) = 0.
>
> There were two ways I was thinking of doing this, and I was
hoping to get some input:
>
> Idea 1: Manually loop through using taking a time step and
then implementing the projection routine. I see that there is
a TSStep command, but this doesn’t seem to be much
documentation on how to use it in this scenario. Does anyone
have any guidance?
>
> Idea 2: Is there some analog to TSMonitor that allows me to
modify the solution after each time step, instead of just
allowing for some computation of a statistic?
>
>
--
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
--
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