TSAdaptHistoryGetStep#

Gets time and time step for a given step number in the history

Synopsis#

#include "petscts.h" 
PetscErrorCode TSAdaptHistoryGetStep(TSAdapt adapt, PetscInt step, PetscReal *t, PetscReal *dt)

Logically Collective

Input Parameters#

  • adapt - the TSAdapt context

  • step - the step number

Output Parameters#

  • t - the time corresponding to the requested step (can be NULL)

  • dt - the time step to be taken at the requested step (can be NULL)

Note#

The time history is internally copied, and the user can free the hist array. The user still needs to specify the termination of the solve via TSSetMaxSteps().

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptHistorySetTrajectory(), TSADAPTHISTORY

Level#

advanced

Location#

src/ts/adapt/impls/history/adapthist.c

Examples#

src/ts/tutorials/ex40.c
src/ts/tutorials/ex41.c


Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages