TSTrajectoryGetVecs#

Reconstructs the vector of state and its time derivative using information from the TSTrajectory and, possibly, from the TS

Synopsis#

#include "petscts.h"  
PetscErrorCode TSTrajectoryGetVecs(TSTrajectory tj, TS ts, PetscInt stepnum, PetscReal *time, Vec U, Vec Udot)

Collective

Input Parameters#

  • tj - the trajectory object

  • ts - the time stepper object (optional)

  • stepnum - the requested step number

Output Parameters#

  • time - On input time for the step if step number is PETSC_DECIDE, on output the time associated with the step number

  • U - state vector (can be NULL)

  • Udot - time derivative of state vector (can be NULL)

Notes#

If the step number is PETSC_DECIDE, the time argument is used to inquire the trajectory. If the requested time does not match any in the trajectory, Lagrangian interpolations are returned.

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSTrajectory, TSTrajectorySetUp(), TSTrajectoryDestroy(), TSTrajectorySetType(), TSTrajectorySetVariableNames(), TSGetTrajectory(), TSTrajectorySet(), TSTrajectoryGet()

Level#

developer

Location#

src/ts/trajectory/interface/traj.c


Edit on GitLab

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