TSMonitorSPSwarmSolution#
Graphically displays phase plots of DMSWARM
particles on a scatter plot
Synopsis#
#include "petscts.h"
PetscErrorCode TSMonitorSPSwarmSolution(TS ts, PetscInt step, PetscReal ptime, Vec u, void *dctx)
Input Parameters#
ts - the
TS
contextstep - current time-step
ptime - current time
u - current solution
dctx - the
TSMonitorSPCtx
object that contains all the options for the monitoring, this is created withTSMonitorSPCtxCreate()
Options Database Keys#
-ts_monitor_sp_swarm
- Monitor the solution every n steps, or -1 for plotting only the final solution-ts_monitor_sp_swarm_retain
- Retain n old points so we can see the history, or -1 for all points-ts_monitor_sp_swarm_multi_species
- Color each species differently-ts_monitor_sp_swarm_phase
- Plot in phase space, as opposed to coordinate space
Notes#
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, TS
, TSMonitoSet()
, DMSWARM
, TSMonitorSPCtxCreate()
Level#
intermediate
Location#
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages