TSMonitorEnvelope#
Monitors the maximum and minimum value of each component of the solution
Synopsis#
#include "petscts.h"
PetscErrorCode TSMonitorEnvelope(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
Collective
Input Parameters#
ts - the
TS
contextstep - current time-step
ptime - current time
u - current solution
dctx - the envelope context
Options Database Key#
-ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time
Notes#
After a solve you can use TSMonitorEnvelopeGetBounds()
to access the envelope
This is not called directly by users, rather one calls TSMonitorSet()
, with this function as an argument, to cause the monitor
to be used during the TS
integration.
See Also#
TS: Scalable ODE and DAE Solvers, TSMonitorSet()
, TSMonitorDefault()
, VecView()
, TSMonitorEnvelopeGetBounds()
, TSMonitorEnvelopeCtxCreate()
Level#
intermediate
Location#
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages