[petsc-dev] TSSolve() shouldn't return a time
Jed Brown
jedbrown at mcs.anl.gov
Sun Oct 23 15:50:04 CDT 2011
On Sun, Oct 23, 2011 at 15:31, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
> On Oct 23, 2011, at 3:26 PM, Jed Brown wrote:
>
> > On Sun, Oct 23, 2011 at 15:18, Barry Smith <bsmith at mcs.anl.gov> wrote:
> > I wasn't thinking of anything so complicated. Just having a separate
> routine TSGetCurrentTime() or something to tell you the time after you
> called TSSolve().
> >
> > But it's not the current time in the sense that it doesn't correspond to
> TSGetSolution or any Vec held by TS and it is not used if you call TSStep()
> or TSSolve(). It's a purely diagnostic point in time with no other
> significance. I can change the TSSolve() interface if we think of a good
> name for this thing.
>
> Wow, that is confusing then :-) So almost every thing in the manual
> page is wrong? It may include "partial time-steps"?
In general, the steppers hold internal state that is disrupted if you change
the solution. Some people like to call TSSolve() in a loop (or the "single
step" version, TSStep()), but they still want high order and error control.
Some methods are unstable if you change the step size too rapidly, so it's
generally better to interpolate to the "exact final time" instead of
adjusting the step sizes to hit it exactly at the end of a step.
I'm not completely opposed to getting rid of the interface for hitting an
exact final time, and instead have users call TSInterpolate themselves, but
I feel like it's a couple unnecessary lines and if we do that, then we can't
make it a command line option.
And x is some weird shit that should not be passed back into TSSolve()?
I think the X argument to TSSolve() should become output-only and require
TSSetSolution() if you really want to throw away internal state and start
stepping from a different starting point (recognizing that this X is not
normally the propagated solution).
> What the ? Is ftime the "interpolated to" time?
>
Yes, on return, X is the state at t=ftime. If -ts_exact_final_time, then
ftime will be the exact final time (assuming the method completed), but if
not -ts_exact_final_time,
>
> Barry
>
>
> TSSolve - Steps the requested number of timesteps.
>
> Collective on TS
>
> Input Parameter:
> + ts - the TS context obtained from TSCreate()
> - x - the solution vector
>
> Output Parameter:
> . ftime - time of the state vector x upon completion
>
> Level: beginner
>
> Notes:
> The final time returned by this function may be different from the time
> of the internally
> held state accessible by TSGetSolution() and TSGetTime() because the
> method may have
> stepped over the final time.
>
>
>
>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20111023/5d33f290/attachment.html>
More information about the petsc-dev
mailing list