TSSetFromOptions#

Sets various TS parameters from the options database

Synopsis#

#include "petscts.h"  
PetscErrorCode TSSetFromOptions(TS ts)

Collective

Input Parameter#

Options Database Keys#

  • -ts_type - EULER, BEULER, SUNDIALS, PSEUDO, CN, RK, THETA, ALPHA, GLLE, SSP, GLEE, BSYMP, IRK, see TSType

  • -ts_save_trajectory - checkpoint the solution at each time-step

  • -ts_max_time maximum time to compute to

  • -ts_time_span <t0,…tf> - sets the time span, solutions are computed and stored for each indicated time

  • -ts_max_steps - maximum number of time-steps to take

  • -ts_init_time initial time to start computation

  • -ts_final_time final time to compute to (deprecated: use -ts_max_time)

  • -ts_dt

    - initial time step

  • -ts_exact_final_time <stepover,interpolate,matchstep> - whether to stop at the exact given final time and how to compute the solution at that time

  • -ts_max_snes_failures - Maximum number of nonlinear solve failures allowed

  • -ts_max_reject - Maximum number of step rejections before step fails

  • -ts_error_if_step_fails <true,false> - Error if no step succeeds

  • -ts_rtol - relative tolerance for local truncation error

  • -ts_atol - Absolute tolerance for local truncation error

  • -ts_rhs_jacobian_test_mult - mat_shell_test_mult_view - test the Jacobian at each iteration against finite difference with RHS function

  • -ts_rhs_jacobian_test_mult_transpose - test the Jacobian at each iteration against finite difference with RHS function

  • -ts_adjoint_solve <yes,no> - After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)

  • -ts_fd_color - Use finite differences with coloring to compute IJacobian

  • -ts_monitor - print information at each timestep

  • -ts_monitor_cancel - Cancel all monitors

  • -ts_monitor_lg_solution - Monitor solution graphically

  • -ts_monitor_lg_error - Monitor error graphically

  • -ts_monitor_error - Monitors norm of error

  • -ts_monitor_lg_timestep - Monitor timestep size graphically

  • -ts_monitor_lg_timestep_log - Monitor log timestep size graphically

  • -ts_monitor_lg_snes_iterations - Monitor number nonlinear iterations for each timestep graphically

  • -ts_monitor_lg_ksp_iterations - Monitor number nonlinear iterations for each timestep graphically

  • -ts_monitor_sp_eig - Monitor eigenvalues of linearized operator graphically

  • -ts_monitor_draw_solution - Monitor solution graphically

  • -ts_monitor_draw_solution_phase <xleft,yleft,xright,yright> - Monitor solution graphically with phase diagram, requires problem with exactly 2 degrees of freedom

  • -ts_monitor_draw_error - Monitor error graphically, requires use to have provided TSSetSolutionFunction()

  • -ts_monitor_solution [ascii binary draw][:filename][:viewerformat] - monitors the solution at each timestep

  • -ts_monitor_solution_interval - output once every interval (default=1) time steps

  • -ts_monitor_solution_vtk <filename.vts,filename.vtu> - Save each time step to a binary file, use filename-%%03” PetscInt_FMT “.vts (filename-%%03” PetscInt_FMT “.vtu)

  • -ts_monitor_envelope - determine maximum and minimum value of each component of the solution over the solution time

Notes#

See SNESSetFromOptions() and KSPSetFromOptions() for how to control the nonlinear and linear solves used by the time-stepper.

Certain SNES options get reset for each new nonlinear solver, for example -snes_lag_jacobian its and -snes_lag_preconditioner its, in order to retain them over the multiple nonlinear solves that TS uses you mush also provide -snes_lag_jacobian_persists true and -snes_lag_preconditioner_persists true

Developer Notes#

We should unify all the -ts_monitor options in the way that -xxx_view has been unified

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSGetType()

Level#

beginner

Location#

src/ts/interface/ts.c

Examples#

src/ts/tutorials/ex13.c
src/ts/tutorials/ex16fwd.c
src/ts/tutorials/ex5.c
src/ts/tutorials/ex41.c
src/ts/tutorials/ex36A.c
src/ts/tutorials/ex20td.c
src/ts/tutorials/ex11.c
src/ts/tutorials/ex9.c
src/ts/tutorials/ex30.c
src/ts/tutorials/ex20opt_ic.c

Implementations#

TSSetFromOptions_ARKIMEX() in src/ts/impls/arkimex/arkimex.c
TSSetFromOptions_BDF() in src/ts/impls/bdf/bdf.c
TSSetFromOptions_EIMEX() in src/ts/impls/eimex/eimex.c
TSSetFromOptions_Euler() in src/ts/impls/explicit/euler/euler.c
TSSetFromOptions_RK() in src/ts/impls/explicit/rk/rk.c
TSSetFromOptions_SSP() in src/ts/impls/explicit/ssp/ssp.c
TSSetFromOptions_GLEE() in src/ts/impls/glee/glee.c
TSSetFromOptions_Alpha() in src/ts/impls/implicit/alpha/alpha1.c
TSSetFromOptions_Alpha() in src/ts/impls/implicit/alpha/alpha2.c
TSSetFromOptions_DiscGrad() in src/ts/impls/implicit/discgrad/tsdiscgrad.c
TSSetFromOptions_GLLE() in src/ts/impls/implicit/glle/glle.c
TSSetFromOptions_IRK() in src/ts/impls/implicit/irk/irk.c
TSSetFromOptions_Sundials() in src/ts/impls/implicit/sundials/sundials.c
TSSetFromOptions_Theta() in src/ts/impls/implicit/theta/theta.c
TSSetFromOptions_Mimex() in src/ts/impls/mimex/mimex.c
TSSetFromOptions_MPRK() in src/ts/impls/multirate/mprk.c
TSSetFromOptions_Pseudo() in src/ts/impls/pseudo/posindep.c
TSSetFromOptions_RosW() in src/ts/impls/rosw/rosw.c
TSSetFromOptions_BasicSymplectic() in src/ts/impls/symplectic/basicsymplectic/basicsymplectic.c


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