TSSetEvaluationTimes#

sets the evaluation points. The solution will be computed and stored for each time requested

Synopsis#

#include "petscts.h"  
PetscErrorCode TSSetEvaluationTimes(TS ts, PetscInt n, PetscReal *time_points)

Collective

Input Parameters#

  • ts - the time-stepper

  • n - number of the time points

  • time_points - array of the time points

Options Database Key#

  • -ts_eval_times <t0,…tn> - Sets the evaluation times

Notes#

The elements in time_points must be all increasing. They correspond to the intermediate points for time integration.

TS_EXACTFINALTIME_MATCHSTEP must be used to make the last time step in each sub-interval match the intermediate points specified.

The intermediate solutions are saved in a vector array that can be accessed with TSGetEvaluationSolutions(). Thus using evaluation times may pressure the memory system when using a large number of time points.

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSGetEvaluationTimes(), TSGetEvaluationSolutions(), TSSetTimeSpan()

Level#

intermediate

Location#

src/ts/interface/ts.c


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