1*7f296bb3SBarry Smith(ch_ts)= 2*7f296bb3SBarry Smith 3*7f296bb3SBarry Smith# TS: Scalable ODE and DAE Solvers 4*7f296bb3SBarry Smith 5*7f296bb3SBarry SmithThe `TS` library provides a framework for the scalable solution of 6*7f296bb3SBarry SmithODEs and DAEs arising from the discretization of time-dependent PDEs. 7*7f296bb3SBarry Smith 8*7f296bb3SBarry Smith**Simple Example:** Consider the PDE 9*7f296bb3SBarry Smith 10*7f296bb3SBarry Smith$$ 11*7f296bb3SBarry Smithu_t = u_{xx} 12*7f296bb3SBarry Smith$$ 13*7f296bb3SBarry Smith 14*7f296bb3SBarry Smithdiscretized with centered finite differences in space yielding the 15*7f296bb3SBarry Smithsemi-discrete equation 16*7f296bb3SBarry Smith 17*7f296bb3SBarry Smith$$ 18*7f296bb3SBarry Smith\begin{aligned} 19*7f296bb3SBarry Smith (u_i)_t & = & \frac{u_{i+1} - 2 u_{i} + u_{i-1}}{h^2}, \\ 20*7f296bb3SBarry Smith u_t & = & \tilde{A} u;\end{aligned} 21*7f296bb3SBarry Smith$$ 22*7f296bb3SBarry Smith 23*7f296bb3SBarry Smithor with piecewise linear finite elements approximation in space 24*7f296bb3SBarry Smith$u(x,t) \doteq \sum_i \xi_i(t) \phi_i(x)$ yielding the 25*7f296bb3SBarry Smithsemi-discrete equation 26*7f296bb3SBarry Smith 27*7f296bb3SBarry Smith$$ 28*7f296bb3SBarry SmithB {\xi}'(t) = A \xi(t) 29*7f296bb3SBarry Smith$$ 30*7f296bb3SBarry Smith 31*7f296bb3SBarry SmithNow applying the backward Euler method results in 32*7f296bb3SBarry Smith 33*7f296bb3SBarry Smith$$ 34*7f296bb3SBarry Smith( B - dt^n A ) u^{n+1} = B u^n, 35*7f296bb3SBarry Smith$$ 36*7f296bb3SBarry Smith 37*7f296bb3SBarry Smithin which 38*7f296bb3SBarry Smith 39*7f296bb3SBarry Smith$$ 40*7f296bb3SBarry Smith{u^n}_i = \xi_i(t_n) \doteq u(x_i,t_n), 41*7f296bb3SBarry Smith$$ 42*7f296bb3SBarry Smith 43*7f296bb3SBarry Smith$$ 44*7f296bb3SBarry Smith{\xi}'(t_{n+1}) \doteq \frac{{u^{n+1}}_i - {u^{n}}_i }{dt^{n}}, 45*7f296bb3SBarry Smith$$ 46*7f296bb3SBarry Smith 47*7f296bb3SBarry Smith$A$ is the stiffness matrix, and $B$ is the identity for 48*7f296bb3SBarry Smithfinite differences or the mass matrix for the finite element method. 49*7f296bb3SBarry Smith 50*7f296bb3SBarry SmithThe PETSc interface for solving time dependent problems assumes the 51*7f296bb3SBarry Smithproblem is written in the form 52*7f296bb3SBarry Smith 53*7f296bb3SBarry Smith$$ 54*7f296bb3SBarry SmithF(t,u,\dot{u}) = G(t,u), \quad u(t_0) = u_0. 55*7f296bb3SBarry Smith$$ 56*7f296bb3SBarry Smith 57*7f296bb3SBarry SmithIn general, this is a differential algebraic equation (DAE) [^id5]. For 58*7f296bb3SBarry SmithODE with nontrivial mass matrices such as arise in FEM, the implicit/DAE 59*7f296bb3SBarry Smithinterface significantly reduces overhead to prepare the system for 60*7f296bb3SBarry Smithalgebraic solvers (`SNES`/`KSP`) by having the user assemble the 61*7f296bb3SBarry Smithcorrectly shifted matrix. Therefore this interface is also useful for 62*7f296bb3SBarry SmithODE systems. 63*7f296bb3SBarry Smith 64*7f296bb3SBarry SmithTo solve an ODE or DAE one uses: 65*7f296bb3SBarry Smith 66*7f296bb3SBarry Smith- Function $F(t,u,\dot{u})$ 67*7f296bb3SBarry Smith 68*7f296bb3SBarry Smith ``` 69*7f296bb3SBarry Smith TSSetIFunction(TS ts,Vec R,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,Vec,void*),void *funP); 70*7f296bb3SBarry Smith ``` 71*7f296bb3SBarry Smith 72*7f296bb3SBarry Smith The vector `R` is an optional location to store the residual. The 73*7f296bb3SBarry Smith arguments to the function `f()` are the timestep context, current 74*7f296bb3SBarry Smith time, input state $u$, input time derivative $\dot{u}$, 75*7f296bb3SBarry Smith and the (optional) user-provided context `funP`. If 76*7f296bb3SBarry Smith $F(t,u,\dot{u}) = \dot{u}$ then one need not call this 77*7f296bb3SBarry Smith function. 78*7f296bb3SBarry Smith 79*7f296bb3SBarry Smith- Function $G(t,u)$, if it is nonzero, is provided with the 80*7f296bb3SBarry Smith function 81*7f296bb3SBarry Smith 82*7f296bb3SBarry Smith ``` 83*7f296bb3SBarry Smith TSSetRHSFunction(TS ts,Vec R,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *funP); 84*7f296bb3SBarry Smith ``` 85*7f296bb3SBarry Smith 86*7f296bb3SBarry Smith- Jacobian 87*7f296bb3SBarry Smith 88*7f296bb3SBarry Smith 89*7f296bb3SBarry Smith $\sigma F_{\dot{u}}(t^n,u^n,\dot{u}^n) + F_u(t^n,u^n,\dot{u}^n)$ 90*7f296bb3SBarry Smith 91*7f296bb3SBarry Smith If using a fully implicit or semi-implicit (IMEX) method one also 92*7f296bb3SBarry Smith can provide an appropriate (approximate) Jacobian matrix of 93*7f296bb3SBarry Smith 94*7f296bb3SBarry Smith 95*7f296bb3SBarry Smith $F()$ 96*7f296bb3SBarry Smith 97*7f296bb3SBarry Smith . 98*7f296bb3SBarry Smith 99*7f296bb3SBarry Smith ``` 100*7f296bb3SBarry Smith TSSetIJacobian(TS ts,Mat A,Mat B,PetscErrorCode (*fjac)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void *jacP); 101*7f296bb3SBarry Smith ``` 102*7f296bb3SBarry Smith 103*7f296bb3SBarry Smith The arguments for the function `fjac()` are the timestep context, 104*7f296bb3SBarry Smith current time, input state $u$, input derivative 105*7f296bb3SBarry Smith $\dot{u}$, input shift $\sigma$, matrix $A$, 106*7f296bb3SBarry Smith preconditioning matrix $B$, and the (optional) user-provided 107*7f296bb3SBarry Smith context `jacP`. 108*7f296bb3SBarry Smith 109*7f296bb3SBarry Smith The Jacobian needed for the nonlinear system is, by the chain rule, 110*7f296bb3SBarry Smith 111*7f296bb3SBarry Smith $$ 112*7f296bb3SBarry Smith \begin{aligned} 113*7f296bb3SBarry Smith \frac{d F}{d u^n} & = & \frac{\partial F}{\partial \dot{u}}|_{u^n} \frac{\partial \dot{u}}{\partial u}|_{u^n} + \frac{\partial F}{\partial u}|_{u^n}.\end{aligned} 114*7f296bb3SBarry Smith $$ 115*7f296bb3SBarry Smith 116*7f296bb3SBarry Smith For any ODE integration method the approximation of $\dot{u}$ 117*7f296bb3SBarry Smith is linear in $u^n$ hence 118*7f296bb3SBarry Smith $\frac{\partial \dot{u}}{\partial u}|_{u^n} = \sigma$, where 119*7f296bb3SBarry Smith the shift $\sigma$ depends on the ODE integrator and time step 120*7f296bb3SBarry Smith but not on the function being integrated. Thus 121*7f296bb3SBarry Smith 122*7f296bb3SBarry Smith $$ 123*7f296bb3SBarry Smith \begin{aligned} 124*7f296bb3SBarry Smith \frac{d F}{d u^n} & = & \sigma F_{\dot{u}}(t^n,u^n,\dot{u}^n) + F_u(t^n,u^n,\dot{u}^n).\end{aligned} 125*7f296bb3SBarry Smith $$ 126*7f296bb3SBarry Smith 127*7f296bb3SBarry Smith This explains why the user provide Jacobian is in the given form for 128*7f296bb3SBarry Smith all integration methods. An equivalent way to derive the formula is 129*7f296bb3SBarry Smith to note that 130*7f296bb3SBarry Smith 131*7f296bb3SBarry Smith $$ 132*7f296bb3SBarry Smith F(t^n,u^n,\dot{u}^n) = F(t^n,u^n,w+\sigma*u^n) 133*7f296bb3SBarry Smith $$ 134*7f296bb3SBarry Smith 135*7f296bb3SBarry Smith where $w$ is some linear combination of previous time solutions 136*7f296bb3SBarry Smith of $u$ so that 137*7f296bb3SBarry Smith 138*7f296bb3SBarry Smith $$ 139*7f296bb3SBarry Smith \frac{d F}{d u^n} = \sigma F_{\dot{u}}(t^n,u^n,\dot{u}^n) + F_u(t^n,u^n,\dot{u}^n) 140*7f296bb3SBarry Smith $$ 141*7f296bb3SBarry Smith 142*7f296bb3SBarry Smith again by the chain rule. 143*7f296bb3SBarry Smith 144*7f296bb3SBarry Smith For example, consider backward Euler’s method applied to the ODE 145*7f296bb3SBarry Smith $F(t, u, \dot{u}) = \dot{u} - f(t, u)$ with 146*7f296bb3SBarry Smith $\dot{u} = (u^n - u^{n-1})/\delta t$ and 147*7f296bb3SBarry Smith $\frac{\partial \dot{u}}{\partial u}|_{u^n} = 1/\delta t$ 148*7f296bb3SBarry Smith resulting in 149*7f296bb3SBarry Smith 150*7f296bb3SBarry Smith $$ 151*7f296bb3SBarry Smith \begin{aligned} 152*7f296bb3SBarry Smith \frac{d F}{d u^n} & = & (1/\delta t)F_{\dot{u}} + F_u(t^n,u^n,\dot{u}^n).\end{aligned} 153*7f296bb3SBarry Smith $$ 154*7f296bb3SBarry Smith 155*7f296bb3SBarry Smith But $F_{\dot{u}} = 1$, in this special case, resulting in the 156*7f296bb3SBarry Smith expected Jacobian $I/\delta t - f_u(t,u^n)$. 157*7f296bb3SBarry Smith 158*7f296bb3SBarry Smith- Jacobian 159*7f296bb3SBarry Smith 160*7f296bb3SBarry Smith $G_u$ 161*7f296bb3SBarry Smith 162*7f296bb3SBarry Smith If using a fully implicit method and the function 163*7f296bb3SBarry Smith 164*7f296bb3SBarry Smith $G()$ 165*7f296bb3SBarry Smith 166*7f296bb3SBarry Smith is 167*7f296bb3SBarry Smith provided, one also can provide an appropriate (approximate) 168*7f296bb3SBarry Smith Jacobian matrix of 169*7f296bb3SBarry Smith 170*7f296bb3SBarry Smith $G()$ 171*7f296bb3SBarry Smith 172*7f296bb3SBarry Smith . 173*7f296bb3SBarry Smith 174*7f296bb3SBarry Smith ``` 175*7f296bb3SBarry Smith TSSetRHSJacobian(TS ts,Mat A,Mat B, 176*7f296bb3SBarry Smith PetscErrorCode (*fjac)(TS,PetscReal,Vec,Mat,Mat,void*),void *jacP); 177*7f296bb3SBarry Smith ``` 178*7f296bb3SBarry Smith 179*7f296bb3SBarry Smith The arguments for the function `fjac()` are the timestep context, 180*7f296bb3SBarry Smith current time, input state $u$, matrix $A$, 181*7f296bb3SBarry Smith preconditioning matrix $B$, and the (optional) user-provided 182*7f296bb3SBarry Smith context `jacP`. 183*7f296bb3SBarry Smith 184*7f296bb3SBarry SmithProviding appropriate $F()$ and $G()$ for your problem 185*7f296bb3SBarry Smithallows for the easy runtime switching between explicit, semi-implicit 186*7f296bb3SBarry Smith(IMEX), and fully implicit methods. 187*7f296bb3SBarry Smith 188*7f296bb3SBarry Smith(sec_ts_basic)= 189*7f296bb3SBarry Smith 190*7f296bb3SBarry Smith## Basic TS Options 191*7f296bb3SBarry Smith 192*7f296bb3SBarry SmithThe user first creates a `TS` object with the command 193*7f296bb3SBarry Smith 194*7f296bb3SBarry Smith``` 195*7f296bb3SBarry Smithint TSCreate(MPI_Comm comm,TS *ts); 196*7f296bb3SBarry Smith``` 197*7f296bb3SBarry Smith 198*7f296bb3SBarry Smith``` 199*7f296bb3SBarry Smithint TSSetProblemType(TS ts,TSProblemType problemtype); 200*7f296bb3SBarry Smith``` 201*7f296bb3SBarry Smith 202*7f296bb3SBarry SmithThe `TSProblemType` is one of `TS_LINEAR` or `TS_NONLINEAR`. 203*7f296bb3SBarry Smith 204*7f296bb3SBarry SmithTo set up `TS` for solving an ODE, one must set the “initial 205*7f296bb3SBarry Smithconditions” for the ODE with 206*7f296bb3SBarry Smith 207*7f296bb3SBarry Smith``` 208*7f296bb3SBarry SmithTSSetSolution(TS ts, Vec initialsolution); 209*7f296bb3SBarry Smith``` 210*7f296bb3SBarry Smith 211*7f296bb3SBarry SmithOne can set the solution method with the routine 212*7f296bb3SBarry Smith 213*7f296bb3SBarry Smith``` 214*7f296bb3SBarry SmithTSSetType(TS ts,TSType type); 215*7f296bb3SBarry Smith``` 216*7f296bb3SBarry Smith 217*7f296bb3SBarry SmithSome of the currently supported types are `TSEULER`, `TSRK` (Runge-Kutta), `TSBEULER`, `TSCN` (Crank-Nicolson), `TSTHETA`, `TSGLLE` (generalized linear), and `TSPSEUDO`. 218*7f296bb3SBarry SmithThey can also be set with the options database option `-ts_type euler, rk, beuler, cn, theta, gl, pseudo, sundials, eimex, arkimex, rosw`. 219*7f296bb3SBarry SmithA list of available methods is given in {any}`integrator_table`. 220*7f296bb3SBarry Smith 221*7f296bb3SBarry SmithSet the initial time with the command 222*7f296bb3SBarry Smith 223*7f296bb3SBarry Smith``` 224*7f296bb3SBarry SmithTSSetTime(TS ts,PetscReal time); 225*7f296bb3SBarry Smith``` 226*7f296bb3SBarry Smith 227*7f296bb3SBarry SmithOne can change the timestep with the command 228*7f296bb3SBarry Smith 229*7f296bb3SBarry Smith``` 230*7f296bb3SBarry SmithTSSetTimeStep(TS ts,PetscReal dt); 231*7f296bb3SBarry Smith``` 232*7f296bb3SBarry Smith 233*7f296bb3SBarry Smithcan determine the current timestep with the routine 234*7f296bb3SBarry Smith 235*7f296bb3SBarry Smith``` 236*7f296bb3SBarry SmithTSGetTimeStep(TS ts,PetscReal* dt); 237*7f296bb3SBarry Smith``` 238*7f296bb3SBarry Smith 239*7f296bb3SBarry SmithHere, “current” refers to the timestep being used to attempt to promote 240*7f296bb3SBarry Smiththe solution form $u^n$ to $u^{n+1}.$ 241*7f296bb3SBarry Smith 242*7f296bb3SBarry SmithOne sets the total number of timesteps to run or the total time to run 243*7f296bb3SBarry Smith(whatever is first) with the commands 244*7f296bb3SBarry Smith 245*7f296bb3SBarry Smith``` 246*7f296bb3SBarry SmithTSSetMaxSteps(TS ts,PetscInt maxsteps); 247*7f296bb3SBarry SmithTSSetMaxTime(TS ts,PetscReal maxtime); 248*7f296bb3SBarry Smith``` 249*7f296bb3SBarry Smith 250*7f296bb3SBarry Smithand determines the behavior near the final time with 251*7f296bb3SBarry Smith 252*7f296bb3SBarry Smith``` 253*7f296bb3SBarry SmithTSSetExactFinalTime(TS ts,TSExactFinalTimeOption eftopt); 254*7f296bb3SBarry Smith``` 255*7f296bb3SBarry Smith 256*7f296bb3SBarry Smithwhere `eftopt` is one of 257*7f296bb3SBarry Smith`TS_EXACTFINALTIME_STEPOVER`,`TS_EXACTFINALTIME_INTERPOLATE`, or 258*7f296bb3SBarry Smith`TS_EXACTFINALTIME_MATCHSTEP`. One performs the requested number of 259*7f296bb3SBarry Smithtime steps with 260*7f296bb3SBarry Smith 261*7f296bb3SBarry Smith``` 262*7f296bb3SBarry SmithTSSolve(TS ts,Vec U); 263*7f296bb3SBarry Smith``` 264*7f296bb3SBarry Smith 265*7f296bb3SBarry SmithThe solve call implicitly sets up the timestep context; this can be done 266*7f296bb3SBarry Smithexplicitly with 267*7f296bb3SBarry Smith 268*7f296bb3SBarry Smith``` 269*7f296bb3SBarry SmithTSSetUp(TS ts); 270*7f296bb3SBarry Smith``` 271*7f296bb3SBarry Smith 272*7f296bb3SBarry SmithOne destroys the context with 273*7f296bb3SBarry Smith 274*7f296bb3SBarry Smith``` 275*7f296bb3SBarry SmithTSDestroy(TS *ts); 276*7f296bb3SBarry Smith``` 277*7f296bb3SBarry Smith 278*7f296bb3SBarry Smithand views it with 279*7f296bb3SBarry Smith 280*7f296bb3SBarry Smith``` 281*7f296bb3SBarry SmithTSView(TS ts,PetscViewer viewer); 282*7f296bb3SBarry Smith``` 283*7f296bb3SBarry Smith 284*7f296bb3SBarry SmithIn place of `TSSolve()`, a single step can be taken using 285*7f296bb3SBarry Smith 286*7f296bb3SBarry Smith``` 287*7f296bb3SBarry SmithTSStep(TS ts); 288*7f296bb3SBarry Smith``` 289*7f296bb3SBarry Smith 290*7f296bb3SBarry Smith(sec_imex)= 291*7f296bb3SBarry Smith 292*7f296bb3SBarry Smith## DAE Formulations 293*7f296bb3SBarry Smith 294*7f296bb3SBarry SmithYou can find a discussion of DAEs in {cite}`ascherpetzold1998` or [Scholarpedia](http://www.scholarpedia.org/article/Differential-algebraic_equations). In PETSc, TS deals with the semi-discrete form of the equations, so that space has already been discretized. If the DAE depends explicitly on the coordinate $x$, then this will just appear as any other data for the equation, not as an explicit argument. Thus we have 295*7f296bb3SBarry Smith 296*7f296bb3SBarry Smith$$ 297*7f296bb3SBarry SmithF(t, u, \dot{u}) = 0 298*7f296bb3SBarry Smith$$ 299*7f296bb3SBarry Smith 300*7f296bb3SBarry SmithIn this form, only fully implicit solvers are appropriate. However, specialized solvers for restricted forms of DAE are supported by PETSc. Below we consider an ODE which is augmented with algebraic constraints on the variables. 301*7f296bb3SBarry Smith 302*7f296bb3SBarry Smith### Hessenberg Index-1 DAE 303*7f296bb3SBarry Smith 304*7f296bb3SBarry Smith> This is a Semi-Explicit Index-1 DAE which has the form 305*7f296bb3SBarry Smith 306*7f296bb3SBarry Smith$$ 307*7f296bb3SBarry Smith\begin{aligned} 308*7f296bb3SBarry Smith \dot{u} &= f(t, u, z) \\ 309*7f296bb3SBarry Smith 0 &= h(t, u, z) 310*7f296bb3SBarry Smith\end{aligned} 311*7f296bb3SBarry Smith$$ 312*7f296bb3SBarry Smith 313*7f296bb3SBarry Smithwhere $z$ is a new constraint variable, and the Jacobian $\frac{dh}{dz}$ is non-singular everywhere. We have suppressed the $x$ dependence since it plays no role here. Using the non-singularity of the Jacobian and the Implicit Function Theorem, we can solve for $z$ in terms of $u$. This means we could, in principle, plug $z(u)$ into the first equation to obtain a simple ODE, even if this is not the numerical process we use. Below we show that this type of DAE can be used with IMEX schemes. 314*7f296bb3SBarry Smith 315*7f296bb3SBarry Smith### Hessenberg Index-2 DAE 316*7f296bb3SBarry Smith 317*7f296bb3SBarry Smith> This DAE has the form 318*7f296bb3SBarry Smith 319*7f296bb3SBarry Smith$$ 320*7f296bb3SBarry Smith\begin{aligned} 321*7f296bb3SBarry Smith \dot{u} &= f(t, u, z) \\ 322*7f296bb3SBarry Smith 0 &= h(t, u) 323*7f296bb3SBarry Smith\end{aligned} 324*7f296bb3SBarry Smith$$ 325*7f296bb3SBarry Smith 326*7f296bb3SBarry SmithNotice that the constraint equation $h$ is not a function of the constraint variable $z$. This means that we cannot naively invert as we did in the index-1 case. Our strategy will be to convert this into an index-1 DAE using a time derivative, which loosely corresponds to the idea of an index being the number of derivatives necessary to get back to an ODE. If we differentiate the constraint equation with respect to time, we can use the ODE to simplify it, 327*7f296bb3SBarry Smith 328*7f296bb3SBarry Smith$$ 329*7f296bb3SBarry Smith\begin{aligned} 330*7f296bb3SBarry Smith 0 &= \dot{h}(t, u) \\ 331*7f296bb3SBarry Smith &= \frac{dh}{du} \dot{u} + \frac{\partial h}{\partial t} \\ 332*7f296bb3SBarry Smith &= \frac{dh}{du} f(t, u, z) + \frac{\partial h}{\partial t} 333*7f296bb3SBarry Smith\end{aligned} 334*7f296bb3SBarry Smith$$ 335*7f296bb3SBarry Smith 336*7f296bb3SBarry SmithIf the Jacobian $\frac{dh}{du} \frac{df}{dz}$ is non-singular, then we have precisely a semi-explicit index-1 DAE, and we can once again use the PETSc IMEX tools to solve it. A common example of an index-2 DAE is the incompressible Navier-Stokes equations, since the continuity equation $\nabla\cdot u = 0$ does not involve the pressure. Using PETSc IMEX with the above conversion then corresponds to the Segregated Runge-Kutta method applied to this equation {cite}`colomesbadia2016`. 337*7f296bb3SBarry Smith 338*7f296bb3SBarry Smith## Using Implicit-Explicit (IMEX) Methods 339*7f296bb3SBarry Smith 340*7f296bb3SBarry SmithFor “stiff” problems or those with multiple time scales $F()$ will 341*7f296bb3SBarry Smithbe treated implicitly using a method suitable for stiff problems and 342*7f296bb3SBarry Smith$G()$ will be treated explicitly when using an IMEX method like 343*7f296bb3SBarry SmithTSARKIMEX. $F()$ is typically linear or weakly nonlinear while 344*7f296bb3SBarry Smith$G()$ may have very strong nonlinearities such as arise in 345*7f296bb3SBarry Smithnon-oscillatory methods for hyperbolic PDE. The user provides three 346*7f296bb3SBarry Smithpieces of information, the APIs for which have been described above. 347*7f296bb3SBarry Smith 348*7f296bb3SBarry Smith- “Slow” part $G(t,u)$ using `TSSetRHSFunction()`. 349*7f296bb3SBarry Smith- “Stiff” part $F(t,u,\dot u)$ using `TSSetIFunction()`. 350*7f296bb3SBarry Smith- Jacobian $F_u + \sigma F_{\dot u}$ using `TSSetIJacobian()`. 351*7f296bb3SBarry Smith 352*7f296bb3SBarry SmithThe user needs to set `TSSetEquationType()` to `TS_EQ_IMPLICIT` or 353*7f296bb3SBarry Smithhigher if the problem is implicit; e.g., 354*7f296bb3SBarry Smith$F(t,u,\dot u) = M \dot u - f(t,u)$, where $M$ is not the 355*7f296bb3SBarry Smithidentity matrix: 356*7f296bb3SBarry Smith 357*7f296bb3SBarry Smith- the problem is an implicit ODE (defined implicitly through 358*7f296bb3SBarry Smith `TSSetIFunction()`) or 359*7f296bb3SBarry Smith- a DAE is being solved. 360*7f296bb3SBarry Smith 361*7f296bb3SBarry SmithAn IMEX problem representation can be made implicit by setting `TSARKIMEXSetFullyImplicit()`. 362*7f296bb3SBarry SmithNote that multilevel preconditioners (e.g. `PCMG`), won't work in the fully implicit case; the 363*7f296bb3SBarry Smithsame holds true for any other `TS` type requiring a fully implicit formulation in case both 364*7f296bb3SBarry SmithJacobians are specified. 365*7f296bb3SBarry Smith 366*7f296bb3SBarry SmithIn PETSc, DAEs and ODEs are formulated as $F(t,u,\dot{u})=G(t,u)$, where $F()$ is meant to be integrated implicitly and $G()$ explicitly. An IMEX formulation such as $M\dot{u}=f(t,u)+g(t,u)$ requires the user to provide $M^{-1} g(t,u)$ or solve $g(t,u) - M x=0$ in place of $G(t,u)$. General cases such as $F(t,u,\dot{u})=G(t,u)$ are not amenable to IMEX Runge-Kutta, but can be solved by using fully implicit methods. Some use-case examples for `TSARKIMEX` are listed in {numref}`tab_DE_forms` and a list of methods with a summary of their properties is given in {any}`tab_IMEX_RK_PETSc`. 367*7f296bb3SBarry Smith 368*7f296bb3SBarry Smith```{eval-rst} 369*7f296bb3SBarry Smith.. list-table:: Use case examples for ``TSARKIMEX`` 370*7f296bb3SBarry Smith :name: tab_DE_forms 371*7f296bb3SBarry Smith :widths: 40 40 80 372*7f296bb3SBarry Smith 373*7f296bb3SBarry Smith * - :math:`\dot{u} = g(t,u)` 374*7f296bb3SBarry Smith - nonstiff ODE 375*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= \dot{u} \\ G(t,u) &= g(t,u)\end{aligned}` 376*7f296bb3SBarry Smith * - :math:`M \dot{u} = g(t,u)` 377*7f296bb3SBarry Smith - nonstiff ODE with mass matrix 378*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= \dot{u} \\ G(t,u) &= M^{-1} g(t,u)\end{aligned}` 379*7f296bb3SBarry Smith * - :math:`\dot{u} = f(t,u)` 380*7f296bb3SBarry Smith - stiff ODE 381*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= \dot{u} - f(t,u) \\ G(t,u) &= 0\end{aligned}` 382*7f296bb3SBarry Smith * - :math:`M \dot{u} = f(t,u)` 383*7f296bb3SBarry Smith - stiff ODE with mass matrix 384*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= M \dot{u} - f(t,u) \\ G(t,u) &= 0\end{aligned}` 385*7f296bb3SBarry Smith * - :math:`\dot{u} = f(t,u) + g(t,u)` 386*7f296bb3SBarry Smith - stiff-nonstiff ODE 387*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= \dot{u} - f(t,u) \\ G(t,u) &= g(t,u)\end{aligned}` 388*7f296bb3SBarry Smith * - :math:`M \dot{u} = f(t,u) + g(t,u)` 389*7f296bb3SBarry Smith - stiff-nonstiff ODE with mass matrix 390*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= M\dot{u} - f(t,u) \\ G(t,u) &= M^{-1} g(t,u)\end{aligned}` 391*7f296bb3SBarry Smith * - :math:`\begin{aligned}\dot{u} &= f(t,u,z) + g(t,u,z)\\0 &= h(t,y,z)\end{aligned}` 392*7f296bb3SBarry Smith - semi-explicit index-1 DAE 393*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= \begin{pmatrix}\dot{u} - f(t,u,z)\\h(t, u, z)\end{pmatrix}\\G(t,u) &= g(t,u)\end{aligned}` 394*7f296bb3SBarry Smith * - :math:`f(t,u,\dot{u})=0` 395*7f296bb3SBarry Smith - fully implicit ODE/DAE 396*7f296bb3SBarry Smith - :math:`\begin{aligned}F(t,u,\dot{u}) &= f(t,u,\dot{u})\\G(t,u) &= 0\end{aligned}`; the user needs to set ``TSSetEquationType()`` to ``TS_EQ_IMPLICIT`` or higher 397*7f296bb3SBarry Smith``` 398*7f296bb3SBarry Smith 399*7f296bb3SBarry Smith{numref}`tab_IMEX_RK_PETSc` lists of the currently available IMEX Runge-Kutta schemes. For each method, it gives the `-ts_arkimex_type` name, the reference, the total number of stages/implicit stages, the order/stage-order, the implicit stability properties (IM), stiff accuracy (SA), the existence of an embedded scheme, and dense output (DO). 400*7f296bb3SBarry Smith 401*7f296bb3SBarry Smith```{eval-rst} 402*7f296bb3SBarry Smith.. list-table:: IMEX Runge-Kutta schemes 403*7f296bb3SBarry Smith :name: tab_IMEX_RK_PETSc 404*7f296bb3SBarry Smith :header-rows: 1 405*7f296bb3SBarry Smith 406*7f296bb3SBarry Smith * - Name 407*7f296bb3SBarry Smith - Reference 408*7f296bb3SBarry Smith - Stages (IM) 409*7f296bb3SBarry Smith - Order (Stage) 410*7f296bb3SBarry Smith - IM 411*7f296bb3SBarry Smith - SA 412*7f296bb3SBarry Smith - Embed 413*7f296bb3SBarry Smith - DO 414*7f296bb3SBarry Smith - Remarks 415*7f296bb3SBarry Smith * - a2 416*7f296bb3SBarry Smith - based on CN 417*7f296bb3SBarry Smith - 2 (1) 418*7f296bb3SBarry Smith - 2 (2) 419*7f296bb3SBarry Smith - A-Stable 420*7f296bb3SBarry Smith - yes 421*7f296bb3SBarry Smith - yes (1) 422*7f296bb3SBarry Smith - yes (2) 423*7f296bb3SBarry Smith - 424*7f296bb3SBarry Smith * - l2 425*7f296bb3SBarry Smith - SSP2(2,2,2) :cite:`pareschi_2005` 426*7f296bb3SBarry Smith - 2 (2) 427*7f296bb3SBarry Smith - 2 (1) 428*7f296bb3SBarry Smith - L-Stable 429*7f296bb3SBarry Smith - yes 430*7f296bb3SBarry Smith - yes (1) 431*7f296bb3SBarry Smith - yes (2) 432*7f296bb3SBarry Smith - SSP SDIRK 433*7f296bb3SBarry Smith * - ars122 434*7f296bb3SBarry Smith - ARS122 :cite:`ascher_1997` 435*7f296bb3SBarry Smith - 2 (1) 436*7f296bb3SBarry Smith - 3 (1) 437*7f296bb3SBarry Smith - A-Stable 438*7f296bb3SBarry Smith - yes 439*7f296bb3SBarry Smith - yes (1) 440*7f296bb3SBarry Smith - yes (2) 441*7f296bb3SBarry Smith - 442*7f296bb3SBarry Smith * - 2c 443*7f296bb3SBarry Smith - :cite:`giraldo_2013` 444*7f296bb3SBarry Smith - 3 (2) 445*7f296bb3SBarry Smith - 2 (2) 446*7f296bb3SBarry Smith - L-Stable 447*7f296bb3SBarry Smith - yes 448*7f296bb3SBarry Smith - yes (1) 449*7f296bb3SBarry Smith - yes (2) 450*7f296bb3SBarry Smith - SDIRK 451*7f296bb3SBarry Smith * - 2d 452*7f296bb3SBarry Smith - :cite:`giraldo_2013` 453*7f296bb3SBarry Smith - 3 (2) 454*7f296bb3SBarry Smith - 2 (2) 455*7f296bb3SBarry Smith - L-Stable 456*7f296bb3SBarry Smith - yes 457*7f296bb3SBarry Smith - yes (1) 458*7f296bb3SBarry Smith - yes (2) 459*7f296bb3SBarry Smith - SDIRK 460*7f296bb3SBarry Smith * - 2e 461*7f296bb3SBarry Smith - :cite:`giraldo_2013` 462*7f296bb3SBarry Smith - 3 (2) 463*7f296bb3SBarry Smith - 2 (2) 464*7f296bb3SBarry Smith - L-Stable 465*7f296bb3SBarry Smith - yes 466*7f296bb3SBarry Smith - yes (1) 467*7f296bb3SBarry Smith - yes (2) 468*7f296bb3SBarry Smith - SDIRK 469*7f296bb3SBarry Smith * - prssp2 470*7f296bb3SBarry Smith - PRS(3,3,2) :cite:`pareschi_2005` 471*7f296bb3SBarry Smith - 3 (3) 472*7f296bb3SBarry Smith - 3 (1) 473*7f296bb3SBarry Smith - L-Stable 474*7f296bb3SBarry Smith - yes 475*7f296bb3SBarry Smith - no 476*7f296bb3SBarry Smith - no 477*7f296bb3SBarry Smith - SSP 478*7f296bb3SBarry Smith * - 3 479*7f296bb3SBarry Smith - :cite:`kennedy_2003` 480*7f296bb3SBarry Smith - 4 (3) 481*7f296bb3SBarry Smith - 3 (2) 482*7f296bb3SBarry Smith - L-Stable 483*7f296bb3SBarry Smith - yes 484*7f296bb3SBarry Smith - yes (2) 485*7f296bb3SBarry Smith - yes (2) 486*7f296bb3SBarry Smith - SDIRK 487*7f296bb3SBarry Smith * - bpr3 488*7f296bb3SBarry Smith - :cite:`boscarino_tr2011` 489*7f296bb3SBarry Smith - 5 (4) 490*7f296bb3SBarry Smith - 3 (2) 491*7f296bb3SBarry Smith - L-Stable 492*7f296bb3SBarry Smith - yes 493*7f296bb3SBarry Smith - no 494*7f296bb3SBarry Smith - no 495*7f296bb3SBarry Smith - SDIRK 496*7f296bb3SBarry Smith * - ars443 497*7f296bb3SBarry Smith - :cite:`ascher_1997` 498*7f296bb3SBarry Smith - 5 (4) 499*7f296bb3SBarry Smith - 3 (1) 500*7f296bb3SBarry Smith - L-Stable 501*7f296bb3SBarry Smith - yes 502*7f296bb3SBarry Smith - no 503*7f296bb3SBarry Smith - no 504*7f296bb3SBarry Smith - SDIRK 505*7f296bb3SBarry Smith * - 4 506*7f296bb3SBarry Smith - :cite:`kennedy_2003` 507*7f296bb3SBarry Smith - 6 (5) 508*7f296bb3SBarry Smith - 4 (2) 509*7f296bb3SBarry Smith - L-Stable 510*7f296bb3SBarry Smith - yes 511*7f296bb3SBarry Smith - yes (3) 512*7f296bb3SBarry Smith - yes 513*7f296bb3SBarry Smith - SDIRK 514*7f296bb3SBarry Smith * - 5 515*7f296bb3SBarry Smith - :cite:`kennedy_2003` 516*7f296bb3SBarry Smith - 8 (7) 517*7f296bb3SBarry Smith - 5 (2) 518*7f296bb3SBarry Smith - L-Stable 519*7f296bb3SBarry Smith - yes 520*7f296bb3SBarry Smith - yes (4) 521*7f296bb3SBarry Smith - yes (3) 522*7f296bb3SBarry Smith - SDIRK 523*7f296bb3SBarry Smith``` 524*7f296bb3SBarry Smith 525*7f296bb3SBarry SmithROSW are linearized implicit Runge-Kutta methods known as Rosenbrock 526*7f296bb3SBarry SmithW-methods. They can accommodate inexact Jacobian matrices in their 527*7f296bb3SBarry Smithformulation. A series of methods are available in PETSc are listed in 528*7f296bb3SBarry Smith{numref}`tab_IMEX_RosW_PETSc` below. For each method, it gives the reference, the total number of stages and implicit stages, the scheme order and stage order, the implicit stability properties (IM), stiff accuracy (SA), the existence of an embedded scheme, dense output (DO), the capacity to use inexact Jacobian matrices (-W), and high order integration of differential algebraic equations (PDAE). 529*7f296bb3SBarry Smith 530*7f296bb3SBarry Smith```{eval-rst} 531*7f296bb3SBarry Smith.. list-table:: Rosenbrock W-schemes 532*7f296bb3SBarry Smith :name: tab_IMEX_RosW_PETSc 533*7f296bb3SBarry Smith :header-rows: 1 534*7f296bb3SBarry Smith 535*7f296bb3SBarry Smith * - TS 536*7f296bb3SBarry Smith - Reference 537*7f296bb3SBarry Smith - Stages (IM) 538*7f296bb3SBarry Smith - Order (Stage) 539*7f296bb3SBarry Smith - IM 540*7f296bb3SBarry Smith - SA 541*7f296bb3SBarry Smith - Embed 542*7f296bb3SBarry Smith - DO 543*7f296bb3SBarry Smith - -W 544*7f296bb3SBarry Smith - PDAE 545*7f296bb3SBarry Smith - Remarks 546*7f296bb3SBarry Smith * - theta1 547*7f296bb3SBarry Smith - classical 548*7f296bb3SBarry Smith - 1(1) 549*7f296bb3SBarry Smith - 1(1) 550*7f296bb3SBarry Smith - L-Stable 551*7f296bb3SBarry Smith - - 552*7f296bb3SBarry Smith - - 553*7f296bb3SBarry Smith - - 554*7f296bb3SBarry Smith - - 555*7f296bb3SBarry Smith - - 556*7f296bb3SBarry Smith - - 557*7f296bb3SBarry Smith * - theta2 558*7f296bb3SBarry Smith - classical 559*7f296bb3SBarry Smith - 1(1) 560*7f296bb3SBarry Smith - 2(2) 561*7f296bb3SBarry Smith - A-Stable 562*7f296bb3SBarry Smith - - 563*7f296bb3SBarry Smith - - 564*7f296bb3SBarry Smith - - 565*7f296bb3SBarry Smith - - 566*7f296bb3SBarry Smith - - 567*7f296bb3SBarry Smith - - 568*7f296bb3SBarry Smith * - 2m 569*7f296bb3SBarry Smith - Zoltan 570*7f296bb3SBarry Smith - 2(2) 571*7f296bb3SBarry Smith - 2(1) 572*7f296bb3SBarry Smith - L-Stable 573*7f296bb3SBarry Smith - No 574*7f296bb3SBarry Smith - Yes(1) 575*7f296bb3SBarry Smith - Yes(2) 576*7f296bb3SBarry Smith - Yes 577*7f296bb3SBarry Smith - No 578*7f296bb3SBarry Smith - SSP 579*7f296bb3SBarry Smith * - 2p 580*7f296bb3SBarry Smith - Zoltan 581*7f296bb3SBarry Smith - 2(2) 582*7f296bb3SBarry Smith - 2(1) 583*7f296bb3SBarry Smith - L-Stable 584*7f296bb3SBarry Smith - No 585*7f296bb3SBarry Smith - Yes(1) 586*7f296bb3SBarry Smith - Yes(2) 587*7f296bb3SBarry Smith - Yes 588*7f296bb3SBarry Smith - No 589*7f296bb3SBarry Smith - SSP 590*7f296bb3SBarry Smith * - ra3pw 591*7f296bb3SBarry Smith - :cite:`rang_2005` 592*7f296bb3SBarry Smith - 3(3) 593*7f296bb3SBarry Smith - 3(1) 594*7f296bb3SBarry Smith - A-Stable 595*7f296bb3SBarry Smith - No 596*7f296bb3SBarry Smith - Yes 597*7f296bb3SBarry Smith - Yes(2) 598*7f296bb3SBarry Smith - No 599*7f296bb3SBarry Smith - Yes(3) 600*7f296bb3SBarry Smith - - 601*7f296bb3SBarry Smith * - ra34pw2 602*7f296bb3SBarry Smith - :cite:`rang_2005` 603*7f296bb3SBarry Smith - 4(4) 604*7f296bb3SBarry Smith - 3(1) 605*7f296bb3SBarry Smith - L-Stable 606*7f296bb3SBarry Smith - Yes 607*7f296bb3SBarry Smith - Yes 608*7f296bb3SBarry Smith - Yes(3) 609*7f296bb3SBarry Smith - Yes 610*7f296bb3SBarry Smith - Yes(3) 611*7f296bb3SBarry Smith - - 612*7f296bb3SBarry Smith * - rodas3 613*7f296bb3SBarry Smith - :cite:`sandu_1997` 614*7f296bb3SBarry Smith - 4(4) 615*7f296bb3SBarry Smith - 3(1) 616*7f296bb3SBarry Smith - L-Stable 617*7f296bb3SBarry Smith - Yes 618*7f296bb3SBarry Smith - Yes 619*7f296bb3SBarry Smith - No 620*7f296bb3SBarry Smith - No 621*7f296bb3SBarry Smith - Yes 622*7f296bb3SBarry Smith - - 623*7f296bb3SBarry Smith * - sandu3 624*7f296bb3SBarry Smith - :cite:`sandu_1997` 625*7f296bb3SBarry Smith - 3(3) 626*7f296bb3SBarry Smith - 3(1) 627*7f296bb3SBarry Smith - L-Stable 628*7f296bb3SBarry Smith - Yes 629*7f296bb3SBarry Smith - Yes 630*7f296bb3SBarry Smith - Yes(2) 631*7f296bb3SBarry Smith - No 632*7f296bb3SBarry Smith - No 633*7f296bb3SBarry Smith - - 634*7f296bb3SBarry Smith * - assp3p3s1c 635*7f296bb3SBarry Smith - unpub. 636*7f296bb3SBarry Smith - 3(2) 637*7f296bb3SBarry Smith - 3(1) 638*7f296bb3SBarry Smith - A-Stable 639*7f296bb3SBarry Smith - No 640*7f296bb3SBarry Smith - Yes 641*7f296bb3SBarry Smith - Yes(2) 642*7f296bb3SBarry Smith - Yes 643*7f296bb3SBarry Smith - No 644*7f296bb3SBarry Smith - SSP 645*7f296bb3SBarry Smith * - lassp3p4s2c 646*7f296bb3SBarry Smith - unpub. 647*7f296bb3SBarry Smith - 4(3) 648*7f296bb3SBarry Smith - 3(1) 649*7f296bb3SBarry Smith - L-Stable 650*7f296bb3SBarry Smith - No 651*7f296bb3SBarry Smith - Yes 652*7f296bb3SBarry Smith - Yes(3) 653*7f296bb3SBarry Smith - Yes 654*7f296bb3SBarry Smith - No 655*7f296bb3SBarry Smith - SSP 656*7f296bb3SBarry Smith * - lassp3p4s2c 657*7f296bb3SBarry Smith - unpub. 658*7f296bb3SBarry Smith - 4(3) 659*7f296bb3SBarry Smith - 3(1) 660*7f296bb3SBarry Smith - L-Stable 661*7f296bb3SBarry Smith - No 662*7f296bb3SBarry Smith - Yes 663*7f296bb3SBarry Smith - Yes(3) 664*7f296bb3SBarry Smith - Yes 665*7f296bb3SBarry Smith - No 666*7f296bb3SBarry Smith - SSP 667*7f296bb3SBarry Smith * - ark3 668*7f296bb3SBarry Smith - unpub. 669*7f296bb3SBarry Smith - 4(3) 670*7f296bb3SBarry Smith - 3(1) 671*7f296bb3SBarry Smith - L-Stable 672*7f296bb3SBarry Smith - No 673*7f296bb3SBarry Smith - Yes 674*7f296bb3SBarry Smith - Yes(3) 675*7f296bb3SBarry Smith - Yes 676*7f296bb3SBarry Smith - No 677*7f296bb3SBarry Smith - IMEX-RK 678*7f296bb3SBarry Smith``` 679*7f296bb3SBarry Smith 680*7f296bb3SBarry Smith## IMEX Methods for fast-slow systems 681*7f296bb3SBarry Smith 682*7f296bb3SBarry SmithConsider a fast-slow ODE system 683*7f296bb3SBarry Smith 684*7f296bb3SBarry Smith$$ 685*7f296bb3SBarry Smith\begin{aligned} 686*7f296bb3SBarry Smith\dot{u}^{slow} & = f^{slow}(t, u^{slow},u^{fast}) \\ 687*7f296bb3SBarry SmithM \dot{u}^{fast} & = g^{fast}(t, u^{slow},u^{fast}) + f^{fast}(t, u^{slow},u^{fast}) 688*7f296bb3SBarry Smith\end{aligned} 689*7f296bb3SBarry Smith$$ 690*7f296bb3SBarry Smith 691*7f296bb3SBarry Smithwhere $u^{slow}$ is the slow component and $u^{fast}$ is the 692*7f296bb3SBarry Smithfast component. The fast component can be partitioned additively as 693*7f296bb3SBarry Smithdescribed above. Thus we want to treat $f^{slow}()$ and 694*7f296bb3SBarry Smith$f^{fast}()$ explicitly and the other terms implicitly when using 695*7f296bb3SBarry SmithTSARKIMEX. This is achieved by using the following APIs: 696*7f296bb3SBarry Smith 697*7f296bb3SBarry Smith- `TSARKIMEXSetFastSlowSplit()` informs PETSc to use ARKIMEX to solve a fast-slow system. 698*7f296bb3SBarry Smith- `TSRHSSplitSetIS()` specifies the index set for the slow/fast components. 699*7f296bb3SBarry Smith- `TSRHSSplitSetRHSFunction()` specifies the parts to be handled explicitly $f^{slow}()$ and $f^{fast}()$. 700*7f296bb3SBarry Smith- `TSRHSSplitSetIFunction()` and `TSRHSSplitSetIJacobian()` specify the implicit part and its Jacobian. 701*7f296bb3SBarry Smith 702*7f296bb3SBarry SmithNote that this ODE system can also be solved by padding zeros in the implicit part and using the standard IMEX methods. However, one needs to provide the full-dimensional Jacobian whereas only a partial Jacobian is needed for the fast-slow split which is more efficient in storage and speed. 703*7f296bb3SBarry Smith 704*7f296bb3SBarry Smith## GLEE methods 705*7f296bb3SBarry Smith 706*7f296bb3SBarry SmithIn this section, we describe explicit and implicit time stepping methods 707*7f296bb3SBarry Smithwith global error estimation that are introduced in 708*7f296bb3SBarry Smith{cite}`constantinescu_tr2016b`. The solution vector for a 709*7f296bb3SBarry SmithGLEE method is either \[$y$, $\tilde{y}$\] or 710*7f296bb3SBarry Smith\[$y$,$\varepsilon$\], where $y$ is the solution, 711*7f296bb3SBarry Smith$\tilde{y}$ is the “auxiliary solution,” and $\varepsilon$ 712*7f296bb3SBarry Smithis the error. The working vector that `TSGLEE` uses is $Y$ = 713*7f296bb3SBarry Smith\[$y$,$\tilde{y}$\], or \[$y$,$\varepsilon$\]. A 714*7f296bb3SBarry SmithGLEE method is defined by 715*7f296bb3SBarry Smith 716*7f296bb3SBarry Smith- $(p,r,s)$: (order, steps, and stages), 717*7f296bb3SBarry Smith- $\gamma$: factor representing the global error ratio, 718*7f296bb3SBarry Smith- $A, U, B, V$: method coefficients, 719*7f296bb3SBarry Smith- $S$: starting method to compute the working vector from the 720*7f296bb3SBarry Smith solution (say at the beginning of time integration) so that 721*7f296bb3SBarry Smith $Y = Sy$, 722*7f296bb3SBarry Smith- $F$: finalizing method to compute the solution from the working 723*7f296bb3SBarry Smith vector,$y = FY$. 724*7f296bb3SBarry Smith- $F_\text{embed}$: coefficients for computing the auxiliary 725*7f296bb3SBarry Smith solution $\tilde{y}$ from the working vector 726*7f296bb3SBarry Smith ($\tilde{y} = F_\text{embed} Y$), 727*7f296bb3SBarry Smith- $F_\text{error}$: coefficients to compute the estimated error 728*7f296bb3SBarry Smith vector from the working vector 729*7f296bb3SBarry Smith ($\varepsilon = F_\text{error} Y$). 730*7f296bb3SBarry Smith- $S_\text{error}$: coefficients to initialize the auxiliary 731*7f296bb3SBarry Smith solution ($\tilde{y}$ or $\varepsilon$) from a specified 732*7f296bb3SBarry Smith error vector ($\varepsilon$). It is currently implemented only 733*7f296bb3SBarry Smith for $r = 2$. We have $y_\text{aux} = 734*7f296bb3SBarry Smith S_{error}[0]*\varepsilon + S_\text{error}[1]*y$, where 735*7f296bb3SBarry Smith $y_\text{aux}$ is the 2nd component of the working vector 736*7f296bb3SBarry Smith $Y$. 737*7f296bb3SBarry Smith 738*7f296bb3SBarry SmithThe methods can be described in two mathematically equivalent forms: 739*7f296bb3SBarry Smithpropagate two components (“$y\tilde{y}$ form”) and propagating the 740*7f296bb3SBarry Smithsolution and its estimated error (“$y\varepsilon$ form”). The two 741*7f296bb3SBarry Smithforms are not explicitly specified in `TSGLEE`; rather, the specific 742*7f296bb3SBarry Smithvalues of $B, U, S, F, F_{embed}$, and $F_{error}$ 743*7f296bb3SBarry Smithcharacterize whether the method is in $y\tilde{y}$ or 744*7f296bb3SBarry Smith$y\varepsilon$ form. 745*7f296bb3SBarry Smith 746*7f296bb3SBarry SmithThe API used by this `TS` method includes: 747*7f296bb3SBarry Smith 748*7f296bb3SBarry Smith- `TSGetSolutionComponents`: Get all the solution components of the 749*7f296bb3SBarry Smith working vector 750*7f296bb3SBarry Smith 751*7f296bb3SBarry Smith ``` 752*7f296bb3SBarry Smith ierr = TSGetSolutionComponents(TS,int*,Vec*) 753*7f296bb3SBarry Smith ``` 754*7f296bb3SBarry Smith 755*7f296bb3SBarry Smith Call with `NULL` as the last argument to get the total number of 756*7f296bb3SBarry Smith components in the working vector $Y$ (this is $r$ (not 757*7f296bb3SBarry Smith $r-1$)), then call to get the $i$-th solution component. 758*7f296bb3SBarry Smith 759*7f296bb3SBarry Smith- `TSGetAuxSolution`: Returns the auxiliary solution 760*7f296bb3SBarry Smith $\tilde{y}$ (computed as $F_\text{embed} Y$) 761*7f296bb3SBarry Smith 762*7f296bb3SBarry Smith ``` 763*7f296bb3SBarry Smith ierr = TSGetAuxSolution(TS,Vec*) 764*7f296bb3SBarry Smith ``` 765*7f296bb3SBarry Smith 766*7f296bb3SBarry Smith- `TSGetTimeError`: Returns the estimated error vector 767*7f296bb3SBarry Smith $\varepsilon$ (computed as $F_\text{error} Y$ if 768*7f296bb3SBarry Smith $n=0$ or restores the error estimate at the end of the previous 769*7f296bb3SBarry Smith step if $n=-1$) 770*7f296bb3SBarry Smith 771*7f296bb3SBarry Smith ``` 772*7f296bb3SBarry Smith ierr = TSGetTimeError(TS,PetscInt n,Vec*) 773*7f296bb3SBarry Smith ``` 774*7f296bb3SBarry Smith 775*7f296bb3SBarry Smith- `TSSetTimeError`: Initializes the auxiliary solution 776*7f296bb3SBarry Smith ($\tilde{y}$ or $\varepsilon$) for a specified initial 777*7f296bb3SBarry Smith error. 778*7f296bb3SBarry Smith 779*7f296bb3SBarry Smith ``` 780*7f296bb3SBarry Smith ierr = TSSetTimeError(TS,Vec) 781*7f296bb3SBarry Smith ``` 782*7f296bb3SBarry Smith 783*7f296bb3SBarry SmithThe local error is estimated as $\varepsilon(n+1)-\varepsilon(n)$. 784*7f296bb3SBarry SmithThis is to be used in the error control. The error in $y\tilde{y}$ 785*7f296bb3SBarry SmithGLEE is 786*7f296bb3SBarry Smith$\varepsilon(n) = \frac{1}{1-\gamma} * (\tilde{y}(n) - y(n))$. 787*7f296bb3SBarry Smith 788*7f296bb3SBarry SmithNote that $y$ and $\tilde{y}$ are reported to `TSAdapt` 789*7f296bb3SBarry Smith`basic` (`TSADAPTBASIC`), and thus it computes the local error as 790*7f296bb3SBarry Smith$\varepsilon_{loc} = (\tilde{y} - 791*7f296bb3SBarry Smithy)$. However, the actual local error is $\varepsilon_{loc} 792*7f296bb3SBarry Smith= \varepsilon_{n+1} - \varepsilon_n = \frac{1}{1-\gamma} * [(\tilde{y} - 793*7f296bb3SBarry Smithy)_{n+1} - (\tilde{y} - y)_n]$. 794*7f296bb3SBarry Smith 795*7f296bb3SBarry Smith{numref}`tab_IMEX_GLEE_PETSc` lists currently available GL schemes with global error estimation {cite}`constantinescu_tr2016b`. 796*7f296bb3SBarry Smith 797*7f296bb3SBarry Smith```{eval-rst} 798*7f296bb3SBarry Smith.. list-table:: GL schemes with global error estimation 799*7f296bb3SBarry Smith :name: tab_IMEX_GLEE_PETSc 800*7f296bb3SBarry Smith :header-rows: 1 801*7f296bb3SBarry Smith 802*7f296bb3SBarry Smith * - TS 803*7f296bb3SBarry Smith - Reference 804*7f296bb3SBarry Smith - IM/EX 805*7f296bb3SBarry Smith - :math:`(p,r,s)` 806*7f296bb3SBarry Smith - :math:`\gamma` 807*7f296bb3SBarry Smith - Form 808*7f296bb3SBarry Smith - Notes 809*7f296bb3SBarry Smith * - ``TSGLEEi1`` 810*7f296bb3SBarry Smith - ``BE1`` 811*7f296bb3SBarry Smith - IM 812*7f296bb3SBarry Smith - :math:`(1,3,2)` 813*7f296bb3SBarry Smith - :math:`0.5` 814*7f296bb3SBarry Smith - :math:`y\varepsilon` 815*7f296bb3SBarry Smith - Based on backward Euler 816*7f296bb3SBarry Smith * - ``TSGLEE23`` 817*7f296bb3SBarry Smith - ``23`` 818*7f296bb3SBarry Smith - EX 819*7f296bb3SBarry Smith - :math:`(2,3,2)` 820*7f296bb3SBarry Smith - :math:`0` 821*7f296bb3SBarry Smith - :math:`y\varepsilon` 822*7f296bb3SBarry Smith - 823*7f296bb3SBarry Smith * - ``TSGLEE24`` 824*7f296bb3SBarry Smith - ``24`` 825*7f296bb3SBarry Smith - EX 826*7f296bb3SBarry Smith - :math:`(2,4,2)` 827*7f296bb3SBarry Smith - :math:`0` 828*7f296bb3SBarry Smith - :math:`y\tilde{y}` 829*7f296bb3SBarry Smith - 830*7f296bb3SBarry Smith * - ``TSGLEE25I`` 831*7f296bb3SBarry Smith - ``25i`` 832*7f296bb3SBarry Smith - EX 833*7f296bb3SBarry Smith - :math:`(2,5,2)` 834*7f296bb3SBarry Smith - :math:`0` 835*7f296bb3SBarry Smith - :math:`y\tilde{y}` 836*7f296bb3SBarry Smith - 837*7f296bb3SBarry Smith * - ``TSGLEE35`` 838*7f296bb3SBarry Smith - ``35`` 839*7f296bb3SBarry Smith - EX 840*7f296bb3SBarry Smith - :math:`(3,5,2)` 841*7f296bb3SBarry Smith - :math:`0` 842*7f296bb3SBarry Smith - :math:`y\tilde{y}` 843*7f296bb3SBarry Smith - 844*7f296bb3SBarry Smith * - ``TSGLEEEXRK2A`` 845*7f296bb3SBarry Smith - ``exrk2a`` 846*7f296bb3SBarry Smith - EX 847*7f296bb3SBarry Smith - :math:`(2,6,2)` 848*7f296bb3SBarry Smith - :math:`0.25` 849*7f296bb3SBarry Smith - :math:`y\varepsilon` 850*7f296bb3SBarry Smith - 851*7f296bb3SBarry Smith * - ``TSGLEERK32G1`` 852*7f296bb3SBarry Smith - ``rk32g1`` 853*7f296bb3SBarry Smith - EX 854*7f296bb3SBarry Smith - :math:`(3,8,2)` 855*7f296bb3SBarry Smith - :math:`0` 856*7f296bb3SBarry Smith - :math:`y\varepsilon` 857*7f296bb3SBarry Smith - 858*7f296bb3SBarry Smith * - ``TSGLEERK285EX`` 859*7f296bb3SBarry Smith - ``rk285ex`` 860*7f296bb3SBarry Smith - EX 861*7f296bb3SBarry Smith - :math:`(2,9,2)` 862*7f296bb3SBarry Smith - :math:`0.25` 863*7f296bb3SBarry Smith - :math:`y\varepsilon` 864*7f296bb3SBarry Smith - 865*7f296bb3SBarry Smith``` 866*7f296bb3SBarry Smith 867*7f296bb3SBarry Smith## Using fully implicit methods 868*7f296bb3SBarry Smith 869*7f296bb3SBarry SmithTo use a fully implicit method like `TSTHETA`, `TSBDF` or `TSDIRK`, either 870*7f296bb3SBarry Smithprovide the Jacobian of $F()$ (and $G()$ if $G()$ is 871*7f296bb3SBarry Smithprovided) or use a `DM` that provides a coloring so the Jacobian can 872*7f296bb3SBarry Smithbe computed efficiently via finite differences. 873*7f296bb3SBarry Smith 874*7f296bb3SBarry Smith## Using the Explicit Runge-Kutta timestepper with variable timesteps 875*7f296bb3SBarry Smith 876*7f296bb3SBarry SmithThe explicit Euler and Runge-Kutta methods require the ODE be in the 877*7f296bb3SBarry Smithform 878*7f296bb3SBarry Smith 879*7f296bb3SBarry Smith$$ 880*7f296bb3SBarry Smith\dot{u} = G(u,t). 881*7f296bb3SBarry Smith$$ 882*7f296bb3SBarry Smith 883*7f296bb3SBarry SmithThe user can either call `TSSetRHSFunction()` and/or they can call 884*7f296bb3SBarry Smith`TSSetIFunction()` (so long as the function provided to 885*7f296bb3SBarry Smith`TSSetIFunction()` is equivalent to $\dot{u} + \tilde{F}(t,u)$) 886*7f296bb3SBarry Smithbut the Jacobians need not be provided. [^id6] 887*7f296bb3SBarry Smith 888*7f296bb3SBarry SmithThe Explicit Runge-Kutta timestepper with variable timesteps is an 889*7f296bb3SBarry Smithimplementation of the standard Runge-Kutta with an embedded method. The 890*7f296bb3SBarry Smitherror in each timestep is calculated using the solutions from the 891*7f296bb3SBarry SmithRunge-Kutta method and its embedded method (the 2-norm of the difference 892*7f296bb3SBarry Smithis used). The default method is the $3$rd-order Bogacki-Shampine 893*7f296bb3SBarry Smithmethod with a $2$nd-order embedded method (`TSRK3BS`). Other 894*7f296bb3SBarry Smithavailable methods are the $5$th-order Fehlberg RK scheme with a 895*7f296bb3SBarry Smith$4$th-order embedded method (`TSRK5F`), the 896*7f296bb3SBarry Smith$5$th-order Dormand-Prince RK scheme with a $4$th-order 897*7f296bb3SBarry Smithembedded method (`TSRK5DP`), the $5$th-order Bogacki-Shampine 898*7f296bb3SBarry SmithRK scheme with a $4$th-order embedded method (`TSRK5BS`, and 899*7f296bb3SBarry Smiththe $6$th-, $7$th, and $8$th-order robust Verner 900*7f296bb3SBarry SmithRK schemes with a $5$th-, $6$th, and $7$th-order 901*7f296bb3SBarry Smithembedded method, respectively (`TSRK6VR`, `TSRK7VR`, `TSRK8VR`). 902*7f296bb3SBarry SmithVariable timesteps cannot be used with RK schemes that do not have an 903*7f296bb3SBarry Smithembedded method (`TSRK1FE` - $1$st-order, $1$-stage 904*7f296bb3SBarry Smithforward Euler, `TSRK2A` - $2$nd-order, $2$-stage RK 905*7f296bb3SBarry Smithscheme, `TSRK3` - $3$rd-order, $3$-stage RK scheme, 906*7f296bb3SBarry Smith`TSRK4` - $4$-th order, $4$-stage RK scheme). 907*7f296bb3SBarry Smith 908*7f296bb3SBarry Smith## Special Cases 909*7f296bb3SBarry Smith 910*7f296bb3SBarry Smith- $\dot{u} = A u.$ First compute the matrix $A$ then call 911*7f296bb3SBarry Smith 912*7f296bb3SBarry Smith ``` 913*7f296bb3SBarry Smith TSSetProblemType(ts,TS_LINEAR); 914*7f296bb3SBarry Smith TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,NULL); 915*7f296bb3SBarry Smith TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,NULL); 916*7f296bb3SBarry Smith ``` 917*7f296bb3SBarry Smith 918*7f296bb3SBarry Smith or 919*7f296bb3SBarry Smith 920*7f296bb3SBarry Smith ``` 921*7f296bb3SBarry Smith TSSetProblemType(ts,TS_LINEAR); 922*7f296bb3SBarry Smith TSSetIFunction(ts,NULL,TSComputeIFunctionLinear,NULL); 923*7f296bb3SBarry Smith TSSetIJacobian(ts,A,A,TSComputeIJacobianConstant,NULL); 924*7f296bb3SBarry Smith ``` 925*7f296bb3SBarry Smith 926*7f296bb3SBarry Smith- $\dot{u} = A(t) u.$ Use 927*7f296bb3SBarry Smith 928*7f296bb3SBarry Smith ``` 929*7f296bb3SBarry Smith TSSetProblemType(ts,TS_LINEAR); 930*7f296bb3SBarry Smith TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,NULL); 931*7f296bb3SBarry Smith TSSetRHSJacobian(ts,A,A,YourComputeRHSJacobian, &appctx); 932*7f296bb3SBarry Smith ``` 933*7f296bb3SBarry Smith 934*7f296bb3SBarry Smith where `YourComputeRHSJacobian()` is a function you provide that 935*7f296bb3SBarry Smith computes $A$ as a function of time. Or use 936*7f296bb3SBarry Smith 937*7f296bb3SBarry Smith ``` 938*7f296bb3SBarry Smith TSSetProblemType(ts,TS_LINEAR); 939*7f296bb3SBarry Smith TSSetIFunction(ts,NULL,TSComputeIFunctionLinear,NULL); 940*7f296bb3SBarry Smith TSSetIJacobian(ts,A,A,YourComputeIJacobian, &appctx); 941*7f296bb3SBarry Smith ``` 942*7f296bb3SBarry Smith 943*7f296bb3SBarry Smith## Monitoring and visualizing solutions 944*7f296bb3SBarry Smith 945*7f296bb3SBarry Smith- `-ts_monitor` - prints the time and timestep at each iteration. 946*7f296bb3SBarry Smith- `-ts_adapt_monitor` - prints information about the timestep 947*7f296bb3SBarry Smith adaption calculation at each iteration. 948*7f296bb3SBarry Smith- `-ts_monitor_lg_timestep` - plots the size of each timestep, 949*7f296bb3SBarry Smith `TSMonitorLGTimeStep()`. 950*7f296bb3SBarry Smith- `-ts_monitor_lg_solution` - for ODEs with only a few components 951*7f296bb3SBarry Smith (not arising from the discretization of a PDE) plots the solution as 952*7f296bb3SBarry Smith a function of time, `TSMonitorLGSolution()`. 953*7f296bb3SBarry Smith- `-ts_monitor_lg_error` - for ODEs with only a few components plots 954*7f296bb3SBarry Smith the error as a function of time, only if `TSSetSolutionFunction()` 955*7f296bb3SBarry Smith is provided, `TSMonitorLGError()`. 956*7f296bb3SBarry Smith- `-ts_monitor_draw_solution` - plots the solution at each iteration, 957*7f296bb3SBarry Smith `TSMonitorDrawSolution()`. 958*7f296bb3SBarry Smith- `-ts_monitor_draw_error` - plots the error at each iteration only 959*7f296bb3SBarry Smith if `TSSetSolutionFunction()` is provided, 960*7f296bb3SBarry Smith `TSMonitorDrawSolution()`. 961*7f296bb3SBarry Smith- `-ts_monitor_solution binary[:filename]` - saves the solution at each 962*7f296bb3SBarry Smith iteration to a binary file, `TSMonitorSolution()`. Solution viewers work 963*7f296bb3SBarry Smith with other time-aware formats, e.g., `-ts_monitor_solution cgns:sol.cgns`, 964*7f296bb3SBarry Smith and can output one solution every 10 time steps by adding 965*7f296bb3SBarry Smith `-ts_monitor_solution_interval 10`. Use `-ts_monitor_solution_interval -1` 966*7f296bb3SBarry Smith to output data only at then end of a time loop. 967*7f296bb3SBarry Smith- `-ts_monitor_solution_vtk <filename-%03D.vts>` - saves the solution 968*7f296bb3SBarry Smith at each iteration to a file in vtk format, 969*7f296bb3SBarry Smith `TSMonitorSolutionVTK()`. 970*7f296bb3SBarry Smith 971*7f296bb3SBarry Smith## Error control via variable time-stepping 972*7f296bb3SBarry Smith 973*7f296bb3SBarry SmithMost of the time stepping methods available in PETSc have an error 974*7f296bb3SBarry Smithestimation and error control mechanism. This mechanism is implemented by 975*7f296bb3SBarry Smithchanging the step size in order to maintain user specified absolute and 976*7f296bb3SBarry Smithrelative tolerances. The PETSc object responsible with error control is 977*7f296bb3SBarry Smith`TSAdapt`. The available `TSAdapt` types are listed in the following table. 978*7f296bb3SBarry Smith 979*7f296bb3SBarry Smith```{eval-rst} 980*7f296bb3SBarry Smith.. list-table:: ``TSAdapt``: available adaptors 981*7f296bb3SBarry Smith :name: tab_adaptors 982*7f296bb3SBarry Smith :header-rows: 1 983*7f296bb3SBarry Smith 984*7f296bb3SBarry Smith * - ID 985*7f296bb3SBarry Smith - Name 986*7f296bb3SBarry Smith - Notes 987*7f296bb3SBarry Smith * - ``TSADAPTNONE`` 988*7f296bb3SBarry Smith - ``none`` 989*7f296bb3SBarry Smith - no adaptivity 990*7f296bb3SBarry Smith * - ``TSADAPTBASIC`` 991*7f296bb3SBarry Smith - ``basic`` 992*7f296bb3SBarry Smith - the default adaptor 993*7f296bb3SBarry Smith * - ``TSADAPTGLEE`` 994*7f296bb3SBarry Smith - ``glee`` 995*7f296bb3SBarry Smith - extension of the basic adaptor to treat :math:`{\rm Tol}_{\rm A}` and :math:`{\rm Tol}_{\rm R}` as separate criteria. It can also control global errors if the integrator (e.g., ``TSGLEE``) provides this information 996*7f296bb3SBarry Smith * - ``TSADAPTDSP`` 997*7f296bb3SBarry Smith - ``dsp`` 998*7f296bb3SBarry Smith - adaptive controller for time-stepping based on digital signal processing 999*7f296bb3SBarry Smith``` 1000*7f296bb3SBarry Smith 1001*7f296bb3SBarry SmithWhen using `TSADAPTBASIC` (the default), the user typically provides a 1002*7f296bb3SBarry Smithdesired absolute ${\rm Tol}_{\rm A}$ or a relative 1003*7f296bb3SBarry Smith${\rm Tol}_{\rm R}$ error tolerance by invoking 1004*7f296bb3SBarry Smith`TSSetTolerances()` or at the command line with options `-ts_atol` 1005*7f296bb3SBarry Smithand `-ts_rtol`. The error estimate is based on the local truncation 1006*7f296bb3SBarry Smitherror, so for every step the algorithm verifies that the estimated local 1007*7f296bb3SBarry Smithtruncation error satisfies the tolerances provided by the user and 1008*7f296bb3SBarry Smithcomputes a new step size to be taken. For multistage methods, the local 1009*7f296bb3SBarry Smithtruncation is obtained by comparing the solution $y$ to a lower 1010*7f296bb3SBarry Smithorder $\widehat{p}=p-1$ approximation, $\widehat{y}$, where 1011*7f296bb3SBarry Smith$p$ is the order of the method and $\widehat{p}$ the order 1012*7f296bb3SBarry Smithof $\widehat{y}$. 1013*7f296bb3SBarry Smith 1014*7f296bb3SBarry SmithThe adaptive controller at step $n$ computes a tolerance level 1015*7f296bb3SBarry Smith 1016*7f296bb3SBarry Smith$$ 1017*7f296bb3SBarry Smith\begin{aligned} 1018*7f296bb3SBarry SmithTol_n(i)&=&{\rm Tol}_{\rm A}(i) + \max(y_n(i),\widehat{y}_n(i)) {\rm Tol}_{\rm R}(i)\,,\end{aligned} 1019*7f296bb3SBarry Smith$$ 1020*7f296bb3SBarry Smith 1021*7f296bb3SBarry Smithand forms the acceptable error level 1022*7f296bb3SBarry Smith 1023*7f296bb3SBarry Smith$$ 1024*7f296bb3SBarry Smith\begin{aligned} 1025*7f296bb3SBarry Smith\rm wlte_n&=& \frac{1}{m} \sum_{i=1}^{m}\sqrt{\frac{\left\|y_n(i) 1026*7f296bb3SBarry Smith -\widehat{y}_n(i)\right\|}{Tol(i)}}\,,\end{aligned} 1027*7f296bb3SBarry Smith$$ 1028*7f296bb3SBarry Smith 1029*7f296bb3SBarry Smithwhere the errors are computed componentwise, $m$ is the dimension 1030*7f296bb3SBarry Smithof $y$ and `-ts_adapt_wnormtype` is `2` (default). If 1031*7f296bb3SBarry Smith`-ts_adapt_wnormtype` is `infinity` (max norm), then 1032*7f296bb3SBarry Smith 1033*7f296bb3SBarry Smith$$ 1034*7f296bb3SBarry Smith\begin{aligned} 1035*7f296bb3SBarry Smith\rm wlte_n&=& \max_{1\dots m}\frac{\left\|y_n(i) 1036*7f296bb3SBarry Smith -\widehat{y}_n(i)\right\|}{Tol(i)}\,.\end{aligned} 1037*7f296bb3SBarry Smith$$ 1038*7f296bb3SBarry Smith 1039*7f296bb3SBarry SmithThe error tolerances are satisfied when $\rm wlte\le 1.0$. 1040*7f296bb3SBarry Smith 1041*7f296bb3SBarry SmithThe next step size is based on this error estimate, and determined by 1042*7f296bb3SBarry Smith 1043*7f296bb3SBarry Smith$$ 1044*7f296bb3SBarry Smith\begin{aligned} 1045*7f296bb3SBarry Smith \Delta t_{\rm new}(t)&=&\Delta t_{\rm{old}} \min(\alpha_{\max}, 1046*7f296bb3SBarry Smith \max(\alpha_{\min}, \beta (1/\rm wlte)^\frac{1}{\widehat{p}+1}))\,,\end{aligned} 1047*7f296bb3SBarry Smith$$ (hnew) 1048*7f296bb3SBarry Smith 1049*7f296bb3SBarry Smithwhere $\alpha_{\min}=$`-ts_adapt_clip`[0] and 1050*7f296bb3SBarry Smith$\alpha_{\max}$=`-ts_adapt_clip`[1] keep the change in 1051*7f296bb3SBarry Smith$\Delta t$ to within a certain factor, and $\beta<1$ is 1052*7f296bb3SBarry Smithchosen through `-ts_adapt_safety` so that there is some margin to 1053*7f296bb3SBarry Smithwhich the tolerances are satisfied and so that the probability of 1054*7f296bb3SBarry Smithrejection is decreased. 1055*7f296bb3SBarry Smith 1056*7f296bb3SBarry SmithThis adaptive controller works in the following way. After completing 1057*7f296bb3SBarry Smithstep $k$, if $\rm wlte_{k+1} \le 1.0$, then the step is 1058*7f296bb3SBarry Smithaccepted and the next step is modified according to 1059*7f296bb3SBarry Smith{eq}`hnew`; otherwise, the step is rejected and retaken 1060*7f296bb3SBarry Smithwith the step length computed in {eq}`hnew`. 1061*7f296bb3SBarry Smith 1062*7f296bb3SBarry Smith`TSADAPTGLEE` is an extension of the basic 1063*7f296bb3SBarry Smithadaptor to treat ${\rm Tol}_{\rm A}$ and ${\rm Tol}_{\rm R}$ 1064*7f296bb3SBarry Smithas separate criteria. it can also control global errors if the 1065*7f296bb3SBarry Smithintegrator (e.g., `TSGLEE`) provides this information. 1066*7f296bb3SBarry Smith 1067*7f296bb3SBarry Smith## Handling of discontinuities 1068*7f296bb3SBarry Smith 1069*7f296bb3SBarry SmithFor problems that involve discontinuous right-hand sides, one can set an 1070*7f296bb3SBarry Smith“event” function $g(t,u)$ for PETSc to detect and locate the times 1071*7f296bb3SBarry Smithof discontinuities (zeros of $g(t,u)$). Events can be defined 1072*7f296bb3SBarry Smiththrough the event monitoring routine 1073*7f296bb3SBarry Smith 1074*7f296bb3SBarry Smith``` 1075*7f296bb3SBarry SmithTSSetEventHandler(TS ts,PetscInt nevents,PetscInt *direction,PetscBool *terminate,PetscErrorCode (*indicator)(TS,PetscReal,Vec,PetscScalar*,void* eventP),PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void* eventP),void *eventP); 1076*7f296bb3SBarry Smith``` 1077*7f296bb3SBarry Smith 1078*7f296bb3SBarry SmithHere, `nevents` denotes the number of events, `direction` sets the 1079*7f296bb3SBarry Smithtype of zero crossing to be detected for an event (+1 for positive 1080*7f296bb3SBarry Smithzero-crossing, -1 for negative zero-crossing, and 0 for both), 1081*7f296bb3SBarry Smith`terminate` conveys whether the time-stepping should continue or halt 1082*7f296bb3SBarry Smithwhen an event is located, `eventmonitor` is a user- defined routine 1083*7f296bb3SBarry Smiththat specifies the event description, `postevent` is an optional 1084*7f296bb3SBarry Smithuser-defined routine to take specific actions following an event. 1085*7f296bb3SBarry Smith 1086*7f296bb3SBarry SmithThe arguments to `indicator()` are the timestep context, current 1087*7f296bb3SBarry Smithtime, input state $u$, array of event function value, and the 1088*7f296bb3SBarry Smith(optional) user-provided context `eventP`. 1089*7f296bb3SBarry Smith 1090*7f296bb3SBarry SmithThe arguments to `postevent()` routine are the timestep context, 1091*7f296bb3SBarry Smithnumber of events occurred, indices of events occurred, current time, input 1092*7f296bb3SBarry Smithstate $u$, a boolean flag indicating forward solve (1) or adjoint 1093*7f296bb3SBarry Smithsolve (0), and the (optional) user-provided context `eventP`. 1094*7f296bb3SBarry Smith 1095*7f296bb3SBarry Smith(sec_tchem)= 1096*7f296bb3SBarry Smith 1097*7f296bb3SBarry Smith## Explicit integrators with finite element mass matrices 1098*7f296bb3SBarry Smith 1099*7f296bb3SBarry SmithDiscretized finite element problems often have the form $M \dot u = G(t, u)$ where $M$ is the mass matrix. 1100*7f296bb3SBarry SmithSuch problems can be solved using `DMTSSetIFunction()` with implicit integrators. 1101*7f296bb3SBarry SmithWhen $M$ is nonsingular (i.e., the problem is an ODE, not a DAE), explicit integrators can be applied to $\dot u = M^{-1} G(t, u)$ or $\dot u = \hat M^{-1} G(t, u)$, where $\hat M$ is the lumped mass matrix. 1102*7f296bb3SBarry SmithWhile the true mass matrix generally has a dense inverse and thus must be solved iteratively, the lumped mass matrix is diagonal (e.g., computed via collocated quadrature or row sums of $M$). 1103*7f296bb3SBarry SmithTo have PETSc create and apply a (lumped) mass matrix automatically, first use `DMTSSetRHSFunction()` to specify $G$ and set a `PetscFE` using `DMAddField()` and `DMCreateDS()`, then call either `DMTSCreateRHSMassMatrix()` or `DMTSCreateRHSMassMatrixLumped()` to automatically create the mass matrix and a `KSP` that will be used to apply $M^{-1}$. 1104*7f296bb3SBarry SmithThis `KSP` can be customized using the `"mass_"` prefix. 1105*7f296bb3SBarry Smith 1106*7f296bb3SBarry Smith(section_sa)= 1107*7f296bb3SBarry Smith 1108*7f296bb3SBarry Smith## Performing sensitivity analysis with the TS ODE Solvers 1109*7f296bb3SBarry Smith 1110*7f296bb3SBarry SmithThe `TS` library provides a framework based on discrete adjoint models 1111*7f296bb3SBarry Smithfor sensitivity analysis for ODEs and DAEs. The ODE/DAE solution process 1112*7f296bb3SBarry Smith(henceforth called the forward run) can be obtained by using either 1113*7f296bb3SBarry Smithexplicit or implicit solvers in `TS`, depending on the problem 1114*7f296bb3SBarry Smithproperties. Currently supported method types are `TSRK` (Runge-Kutta) 1115*7f296bb3SBarry Smithexplicit methods and `TSTHETA` implicit methods, which include 1116*7f296bb3SBarry Smith`TSBEULER` and `TSCN`. 1117*7f296bb3SBarry Smith 1118*7f296bb3SBarry Smith### Using the discrete adjoint methods 1119*7f296bb3SBarry Smith 1120*7f296bb3SBarry SmithConsider the ODE/DAE 1121*7f296bb3SBarry Smith 1122*7f296bb3SBarry Smith$$ 1123*7f296bb3SBarry SmithF(t,y,\dot{y},p) = 0, \quad y(t_0)=y_0(p) \quad t_0 \le t \le t_F 1124*7f296bb3SBarry Smith$$ 1125*7f296bb3SBarry Smith 1126*7f296bb3SBarry Smithand the cost function(s) 1127*7f296bb3SBarry Smith 1128*7f296bb3SBarry Smith$$ 1129*7f296bb3SBarry Smith\Psi_i(y_0,p) = \Phi_i(y_F,p) + \int_{t_0}^{t_F} r_i(y(t),p,t)dt \quad i=1,...,n_\text{cost}. 1130*7f296bb3SBarry Smith$$ 1131*7f296bb3SBarry Smith 1132*7f296bb3SBarry SmithThe `TSAdjoint` routines of PETSc provide 1133*7f296bb3SBarry Smith 1134*7f296bb3SBarry Smith$$ 1135*7f296bb3SBarry Smith\frac{\partial \Psi_i}{\partial y_0} = \lambda_i 1136*7f296bb3SBarry Smith$$ 1137*7f296bb3SBarry Smith 1138*7f296bb3SBarry Smithand 1139*7f296bb3SBarry Smith 1140*7f296bb3SBarry Smith$$ 1141*7f296bb3SBarry Smith\frac{\partial \Psi_i}{\partial p} = \mu_i + \lambda_i (\frac{\partial y_0}{\partial p}). 1142*7f296bb3SBarry Smith$$ 1143*7f296bb3SBarry Smith 1144*7f296bb3SBarry SmithTo perform the discrete adjoint sensitivity analysis one first sets up 1145*7f296bb3SBarry Smiththe `TS` object for a regular forward run but with one extra function 1146*7f296bb3SBarry Smithcall 1147*7f296bb3SBarry Smith 1148*7f296bb3SBarry Smith``` 1149*7f296bb3SBarry SmithTSSetSaveTrajectory(TS ts), 1150*7f296bb3SBarry Smith``` 1151*7f296bb3SBarry Smith 1152*7f296bb3SBarry Smiththen calls `TSSolve()` in the usual manner. 1153*7f296bb3SBarry Smith 1154*7f296bb3SBarry SmithOne must create two arrays of $n_\text{cost}$ vectors 1155*7f296bb3SBarry Smith$\lambda$ and $\mu$ (if there are no parameters $p$ 1156*7f296bb3SBarry Smiththen one can use `NULL` for the $\mu$ array.) The 1157*7f296bb3SBarry Smith$\lambda$ vectors are the same dimension and parallel layout as 1158*7f296bb3SBarry Smiththe solution vector for the ODE, the $\mu$ vectors are of dimension 1159*7f296bb3SBarry Smith$p$; when $p$ is small usually all its elements are on the 1160*7f296bb3SBarry Smithfirst MPI process, while the vectors have no entries on the other 1161*7f296bb3SBarry Smithprocesses. $\lambda_i$ and $\mu_i$ should be initialized with 1162*7f296bb3SBarry Smiththe values $d\Phi_i/dy|_{t=t_F}$ and $d\Phi_i/dp|_{t=t_F}$ 1163*7f296bb3SBarry Smithrespectively. Then one calls 1164*7f296bb3SBarry Smith 1165*7f296bb3SBarry Smith``` 1166*7f296bb3SBarry SmithTSSetCostGradients(TS ts,PetscInt numcost, Vec *lambda,Vec *mu); 1167*7f296bb3SBarry Smith``` 1168*7f296bb3SBarry Smith 1169*7f296bb3SBarry Smithwhere `numcost` denotes $n_\text{cost}$. 1170*7f296bb3SBarry SmithIf $F()$ is a function of $p$ one needs to also provide the 1171*7f296bb3SBarry SmithJacobian $-F_p$ with 1172*7f296bb3SBarry Smith 1173*7f296bb3SBarry Smith``` 1174*7f296bb3SBarry SmithTSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*fp)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1175*7f296bb3SBarry Smith``` 1176*7f296bb3SBarry Smith 1177*7f296bb3SBarry Smithor 1178*7f296bb3SBarry Smith 1179*7f296bb3SBarry Smith``` 1180*7f296bb3SBarry SmithTSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*fp)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx) 1181*7f296bb3SBarry Smith``` 1182*7f296bb3SBarry Smith 1183*7f296bb3SBarry Smithor both, depending on which form is used to define the ODE. 1184*7f296bb3SBarry Smith 1185*7f296bb3SBarry SmithThe arguments for the function `fp()` are the timestep context, 1186*7f296bb3SBarry Smithcurrent time, $y$, and the (optional) user-provided context. 1187*7f296bb3SBarry Smith 1188*7f296bb3SBarry SmithIf there is an integral term in the cost function, i.e. $r$ is 1189*7f296bb3SBarry Smithnonzero, it can be transformed into another ODE that is augmented to the 1190*7f296bb3SBarry Smithoriginal ODE. To evaluate the integral, one needs to create a child 1191*7f296bb3SBarry Smith`TS` objective by calling 1192*7f296bb3SBarry Smith 1193*7f296bb3SBarry Smith``` 1194*7f296bb3SBarry SmithTSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts); 1195*7f296bb3SBarry Smith``` 1196*7f296bb3SBarry Smith 1197*7f296bb3SBarry Smithand provide the ODE RHS function (which evaluates the integrand 1198*7f296bb3SBarry Smith$r$) with 1199*7f296bb3SBarry Smith 1200*7f296bb3SBarry Smith``` 1201*7f296bb3SBarry SmithTSSetRHSFunction(TS quadts,Vec R,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),void *ctx) 1202*7f296bb3SBarry Smith``` 1203*7f296bb3SBarry Smith 1204*7f296bb3SBarry SmithSimilar to the settings for the original ODE, Jacobians of the integrand 1205*7f296bb3SBarry Smithcan be provided with 1206*7f296bb3SBarry Smith 1207*7f296bb3SBarry Smith``` 1208*7f296bb3SBarry SmithTSSetRHSJacobian(TS quadts,Vec DRDU,Vec DRDU,PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*),void *ctx) 1209*7f296bb3SBarry SmithTSSetRHSJacobianP(TS quadts,Vec DRDU,Vec DRDU,PetscErrorCode (*drdyp)(TS,PetscReal,Vec,Vec*,void*),void *ctx) 1210*7f296bb3SBarry Smith``` 1211*7f296bb3SBarry Smith 1212*7f296bb3SBarry Smithwhere $\mathrm{drdyf}= dr /dy$, $\mathrm{drdpf} = dr /dp$. 1213*7f296bb3SBarry SmithSince the integral term is additive to the cost function, its gradient 1214*7f296bb3SBarry Smithinformation will be included in $\lambda$ and $\mu$. 1215*7f296bb3SBarry Smith 1216*7f296bb3SBarry SmithLastly, one starts the backward run by calling 1217*7f296bb3SBarry Smith 1218*7f296bb3SBarry Smith``` 1219*7f296bb3SBarry SmithTSAdjointSolve(TS ts). 1220*7f296bb3SBarry Smith``` 1221*7f296bb3SBarry Smith 1222*7f296bb3SBarry SmithOne can obtain the value of the integral term by calling 1223*7f296bb3SBarry Smith 1224*7f296bb3SBarry Smith``` 1225*7f296bb3SBarry SmithTSGetCostIntegral(TS ts,Vec *q). 1226*7f296bb3SBarry Smith``` 1227*7f296bb3SBarry Smith 1228*7f296bb3SBarry Smithor accessing directly the solution vector used by `quadts`. 1229*7f296bb3SBarry Smith 1230*7f296bb3SBarry SmithThe second argument of `TSCreateQuadratureTS()` allows one to choose 1231*7f296bb3SBarry Smithif the integral term is evaluated in the forward run (inside 1232*7f296bb3SBarry Smith`TSSolve()`) or in the backward run (inside `TSAdjointSolve()`) when 1233*7f296bb3SBarry Smith`TSSetCostGradients()` and `TSSetCostIntegrand()` are called before 1234*7f296bb3SBarry Smith`TSSolve()`. Note that this also allows for evaluating the integral 1235*7f296bb3SBarry Smithwithout having to use the adjoint solvers. 1236*7f296bb3SBarry Smith 1237*7f296bb3SBarry SmithTo provide a better understanding of the use of the adjoint solvers, we 1238*7f296bb3SBarry Smithintroduce a simple example, corresponding to 1239*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/power_grid/ex3sa.c.html">TS Power Grid Tutorial ex3sa</a>. 1240*7f296bb3SBarry SmithThe problem is to study dynamic security of power system when there are 1241*7f296bb3SBarry Smithcredible contingencies such as short-circuits or loss of generators, 1242*7f296bb3SBarry Smithtransmission lines, or loads. The dynamic security constraints are 1243*7f296bb3SBarry Smithincorporated as equality constraints in the form of discretized 1244*7f296bb3SBarry Smithdifferential equations and inequality constraints for bounds on the 1245*7f296bb3SBarry Smithtrajectory. The governing ODE system is 1246*7f296bb3SBarry Smith 1247*7f296bb3SBarry Smith$$ 1248*7f296bb3SBarry Smith\begin{aligned} 1249*7f296bb3SBarry Smith \phi' &= &\omega_B (\omega - \omega_S) \\ 1250*7f296bb3SBarry Smith 2H/\omega_S \, \omega' & =& p_m - p_{max} sin(\phi) -D (\omega - \omega_S), \quad t_0 \leq t \leq t_F,\end{aligned} 1251*7f296bb3SBarry Smith$$ 1252*7f296bb3SBarry Smith 1253*7f296bb3SBarry Smithwhere $\phi$ is the phase angle and $\omega$ is the 1254*7f296bb3SBarry Smithfrequency. 1255*7f296bb3SBarry Smith 1256*7f296bb3SBarry SmithThe initial conditions at time $t_0$ are 1257*7f296bb3SBarry Smith 1258*7f296bb3SBarry Smith$$ 1259*7f296bb3SBarry Smith\begin{aligned} 1260*7f296bb3SBarry Smith\phi(t_0) &=& \arcsin \left( p_m / p_{max} \right), \\ 1261*7f296bb3SBarry Smithw(t_0) & =& 1.\end{aligned} 1262*7f296bb3SBarry Smith$$ 1263*7f296bb3SBarry Smith 1264*7f296bb3SBarry Smith$p_{max}$ is a positive number when the system operates normally. 1265*7f296bb3SBarry SmithAt an event such as fault incidence/removal, $p_{max}$ will change 1266*7f296bb3SBarry Smithto $0$ temporarily and back to the original value after the fault 1267*7f296bb3SBarry Smithis fixed. The objective is to maximize $p_m$ subject to the above 1268*7f296bb3SBarry SmithODE constraints and $\phi<\phi_S$ during all times. To accommodate 1269*7f296bb3SBarry Smiththe inequality constraint, we want to compute the sensitivity of the 1270*7f296bb3SBarry Smithcost function 1271*7f296bb3SBarry Smith 1272*7f296bb3SBarry Smith$$ 1273*7f296bb3SBarry Smith\Psi(p_m,\phi) = -p_m + c \int_{t_0}^{t_F} \left( \max(0, \phi - \phi_S ) \right)^2 dt 1274*7f296bb3SBarry Smith$$ 1275*7f296bb3SBarry Smith 1276*7f296bb3SBarry Smithwith respect to the parameter $p_m$. $numcost$ is $1$ 1277*7f296bb3SBarry Smithsince it is a scalar function. 1278*7f296bb3SBarry Smith 1279*7f296bb3SBarry SmithFor ODE solution, PETSc requires user-provided functions to evaluate the 1280*7f296bb3SBarry Smithsystem $F(t,y,\dot{y},p)$ (set by `TSSetIFunction()` ) and its 1281*7f296bb3SBarry Smithcorresponding Jacobian $F_y + \sigma F_{\dot y}$ (set by 1282*7f296bb3SBarry Smith`TSSetIJacobian()`). Note that the solution state $y$ is 1283*7f296bb3SBarry Smith$[ \phi \; \omega ]^T$ here. For sensitivity analysis, we need to 1284*7f296bb3SBarry Smithprovide a routine to compute $\mathrm{f}_p=[0 \; 1]^T$ using 1285*7f296bb3SBarry Smith`TSASetRHSJacobianP()`, and three routines corresponding to the 1286*7f296bb3SBarry Smithintegrand $r=c \left( \max(0, \phi - \phi_S ) \right)^2$, 1287*7f296bb3SBarry Smith$r_p = [0 \; 0]^T$ and 1288*7f296bb3SBarry Smith$r_y= [ 2 c \left( \max(0, \phi - \phi_S ) \right) \; 0]^T$ using 1289*7f296bb3SBarry Smith`TSSetCostIntegrand()`. 1290*7f296bb3SBarry Smith 1291*7f296bb3SBarry SmithIn the adjoint run, $\lambda$ and $\mu$ are initialized as 1292*7f296bb3SBarry Smith$[ 0 \; 0 ]^T$ and $[-1]$ at the final time $t_F$. 1293*7f296bb3SBarry SmithAfter `TSAdjointSolve()`, the sensitivity of the cost function w.r.t. 1294*7f296bb3SBarry Smithinitial conditions is given by the sensitivity variable $\lambda$ 1295*7f296bb3SBarry Smith(at time $t_0$) directly. And the sensitivity of the cost function 1296*7f296bb3SBarry Smithw.r.t. the parameter $p_m$ can be computed (by users) as 1297*7f296bb3SBarry Smith 1298*7f296bb3SBarry Smith$$ 1299*7f296bb3SBarry Smith\frac{\mathrm{d} \Psi}{\mathrm{d} p_m} = \mu(t_0) + \lambda(t_0) \frac{\mathrm{d} \left[ \phi(t_0) \; \omega(t_0) \right]^T}{\mathrm{d} p_m} . 1300*7f296bb3SBarry Smith$$ 1301*7f296bb3SBarry Smith 1302*7f296bb3SBarry SmithFor explicit methods where one does not need to provide the Jacobian 1303*7f296bb3SBarry Smith$F_u$ for the forward solve one still does need it for the 1304*7f296bb3SBarry Smithbackward solve and thus must call 1305*7f296bb3SBarry Smith 1306*7f296bb3SBarry Smith``` 1307*7f296bb3SBarry SmithTSSetRHSJacobian(TS ts,Mat Amat, Mat Pmat,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat,Mat,void*),void *fP); 1308*7f296bb3SBarry Smith``` 1309*7f296bb3SBarry Smith 1310*7f296bb3SBarry SmithExamples include: 1311*7f296bb3SBarry Smith 1312*7f296bb3SBarry Smith- discrete adjoint sensitivity using explicit and implicit time stepping methods for an ODE problem 1313*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20adj.c.html">TS Tutorial ex20adj</a>, 1314*7f296bb3SBarry Smith- an optimization problem using the discrete adjoint models of the ERK (for nonstiff ODEs) 1315*7f296bb3SBarry Smith and the Theta methods (for stiff DAEs) 1316*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20opt_ic.c.html">TS Tutorial ex20opt_ic</a> 1317*7f296bb3SBarry Smith and 1318*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20opt_p.c.html">TS Tutorial ex20opt_p</a>, 1319*7f296bb3SBarry Smith- an ODE-constrained optimization using the discrete adjoint models of the 1320*7f296bb3SBarry Smith Theta methods for cost function with an integral term 1321*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/power_grid/ex3opt.c.html">TS Power Grid Tutorial ex3opt</a>, 1322*7f296bb3SBarry Smith- discrete adjoint sensitivity using the Crank-Nicolson methods for DAEs with discontinuities 1323*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/power_grid/stability_9bus/ex9busadj.c.html">TS Power Grid Stability Tutorial ex9busadj</a>, 1324*7f296bb3SBarry Smith- a DAE-constrained optimization problem using the discrete adjoint models of the Crank-Nicolson 1325*7f296bb3SBarry Smith methods for cost function with an integral term 1326*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/power_grid/stability_9bus/ex9busopt.c.html">TS Power Grid Tutorial ex9busopt</a>, 1327*7f296bb3SBarry Smith- discrete adjoint sensitivity using the Crank-Nicolson methods for a PDE problem 1328*7f296bb3SBarry Smith <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/advection-diffusion-reaction/ex5adj.c.html">TS Advection-Diffusion-Reaction Tutorial ex5adj</a>. 1329*7f296bb3SBarry Smith 1330*7f296bb3SBarry Smith### Checkpointing 1331*7f296bb3SBarry Smith 1332*7f296bb3SBarry SmithThe discrete adjoint model requires the states (and stage values in the 1333*7f296bb3SBarry Smithcontext of multistage timestepping methods) to evaluate the Jacobian 1334*7f296bb3SBarry Smithmatrices during the adjoint (backward) run. By default, PETSc stores the 1335*7f296bb3SBarry Smithwhole trajectory to disk as binary files, each of which contains the 1336*7f296bb3SBarry Smithinformation for a single time step including state, time, and stage 1337*7f296bb3SBarry Smithvalues (optional). One can also make PETSc store the trajectory to 1338*7f296bb3SBarry Smithmemory with the option `-ts_trajectory_type memory`. However, there 1339*7f296bb3SBarry Smithmight not be sufficient memory capacity especially for large-scale 1340*7f296bb3SBarry Smithproblems and long-time integration. 1341*7f296bb3SBarry Smith 1342*7f296bb3SBarry SmithA so-called checkpointing scheme is needed to solve this problem. The 1343*7f296bb3SBarry Smithscheme stores checkpoints at selective time steps and recomputes the 1344*7f296bb3SBarry Smithmissing information. The `revolve` library is used by PETSc 1345*7f296bb3SBarry Smith`TSTrajectory` to generate an optimal checkpointing schedule that 1346*7f296bb3SBarry Smithminimizes the recomputations given a limited number of available 1347*7f296bb3SBarry Smithcheckpoints. One can specify the number of available checkpoints with 1348*7f296bb3SBarry Smiththe option 1349*7f296bb3SBarry Smith`-ts_trajectory_max_cps_ram [maximum number of checkpoints in RAM]`. 1350*7f296bb3SBarry SmithNote that one checkpoint corresponds to one time step. 1351*7f296bb3SBarry Smith 1352*7f296bb3SBarry SmithThe `revolve` library also provides an optimal multistage 1353*7f296bb3SBarry Smithcheckpointing scheme that uses both RAM and disk for storage. This 1354*7f296bb3SBarry Smithscheme is automatically chosen if one uses both the option 1355*7f296bb3SBarry Smith`-ts_trajectory_max_cps_ram [maximum number of checkpoints in RAM]` 1356*7f296bb3SBarry Smithand the option 1357*7f296bb3SBarry Smith`-ts_trajectory_max_cps_disk [maximum number of checkpoints on disk]`. 1358*7f296bb3SBarry Smith 1359*7f296bb3SBarry SmithSome other useful options are listed below. 1360*7f296bb3SBarry Smith 1361*7f296bb3SBarry Smith- `-ts_trajectory_view` prints the total number of recomputations, 1362*7f296bb3SBarry Smith- `-ts_monitor` and `-ts_adjoint_monitor` allow users to monitor 1363*7f296bb3SBarry Smith the progress of the adjoint work flow, 1364*7f296bb3SBarry Smith- `-ts_trajectory_type visualization` may be used to save the whole 1365*7f296bb3SBarry Smith trajectory for visualization. It stores the solution and the time, 1366*7f296bb3SBarry Smith but no stage values. The binary files generated can be read into 1367*7f296bb3SBarry Smith MATLAB via the script 1368*7f296bb3SBarry Smith `$PETSC_DIR/share/petsc/matlab/PetscReadBinaryTrajectory.m`. 1369*7f296bb3SBarry Smith 1370*7f296bb3SBarry Smith(sec_sundials)= 1371*7f296bb3SBarry Smith 1372*7f296bb3SBarry Smith## Using Sundials from PETSc 1373*7f296bb3SBarry Smith 1374*7f296bb3SBarry SmithSundials is a parallel ODE solver developed by Hindmarsh et al. at LLNL. 1375*7f296bb3SBarry SmithThe `TS` library provides an interface to use the CVODE component of 1376*7f296bb3SBarry SmithSundials directly from PETSc. (To configure PETSc to use Sundials, see 1377*7f296bb3SBarry Smiththe installation guide, `installation/index.htm`.) 1378*7f296bb3SBarry Smith 1379*7f296bb3SBarry SmithTo use the Sundials integrators, call 1380*7f296bb3SBarry Smith 1381*7f296bb3SBarry Smith``` 1382*7f296bb3SBarry SmithTSSetType(TS ts,TSType TSSUNDIALS); 1383*7f296bb3SBarry Smith``` 1384*7f296bb3SBarry Smith 1385*7f296bb3SBarry Smithor use the command line option `-ts_type` `sundials`. 1386*7f296bb3SBarry Smith 1387*7f296bb3SBarry SmithSundials’ CVODE solver comes with two main integrator families, Adams 1388*7f296bb3SBarry Smithand BDF (backward differentiation formula). One can select these with 1389*7f296bb3SBarry Smith 1390*7f296bb3SBarry Smith``` 1391*7f296bb3SBarry SmithTSSundialsSetType(TS ts,TSSundialsLmmType [SUNDIALS_ADAMS,SUNDIALS_BDF]); 1392*7f296bb3SBarry Smith``` 1393*7f296bb3SBarry Smith 1394*7f296bb3SBarry Smithor the command line option `-ts_sundials_type <adams,bdf>`. BDF is the 1395*7f296bb3SBarry Smithdefault. 1396*7f296bb3SBarry Smith 1397*7f296bb3SBarry SmithSundials does not use the `SNES` library within PETSc for its 1398*7f296bb3SBarry Smithnonlinear solvers, so one cannot change the nonlinear solver options via 1399*7f296bb3SBarry Smith`SNES`. Rather, Sundials uses the preconditioners within the `PC` 1400*7f296bb3SBarry Smithpackage of PETSc, which can be accessed via 1401*7f296bb3SBarry Smith 1402*7f296bb3SBarry Smith``` 1403*7f296bb3SBarry SmithTSSundialsGetPC(TS ts,PC *pc); 1404*7f296bb3SBarry Smith``` 1405*7f296bb3SBarry Smith 1406*7f296bb3SBarry SmithThe user can then directly set preconditioner options; alternatively, 1407*7f296bb3SBarry Smiththe usual runtime options can be employed via `-pc_xxx`. 1408*7f296bb3SBarry Smith 1409*7f296bb3SBarry SmithFinally, one can set the Sundials tolerances via 1410*7f296bb3SBarry Smith 1411*7f296bb3SBarry Smith``` 1412*7f296bb3SBarry SmithTSSundialsSetTolerance(TS ts,double abs,double rel); 1413*7f296bb3SBarry Smith``` 1414*7f296bb3SBarry Smith 1415*7f296bb3SBarry Smithwhere `abs` denotes the absolute tolerance and `rel` the relative 1416*7f296bb3SBarry Smithtolerance. 1417*7f296bb3SBarry Smith 1418*7f296bb3SBarry SmithOther PETSc-Sundials options include 1419*7f296bb3SBarry Smith 1420*7f296bb3SBarry Smith``` 1421*7f296bb3SBarry SmithTSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type); 1422*7f296bb3SBarry Smith``` 1423*7f296bb3SBarry Smith 1424*7f296bb3SBarry Smithwhere `type` is either `SUNDIALS_MODIFIED_GS` or 1425*7f296bb3SBarry Smith`SUNDIALS_UNMODIFIED_GS`. This may be set via the options data base 1426*7f296bb3SBarry Smithwith `-ts_sundials_gramschmidt_type <modifed,unmodified>`. 1427*7f296bb3SBarry Smith 1428*7f296bb3SBarry SmithThe routine 1429*7f296bb3SBarry Smith 1430*7f296bb3SBarry Smith``` 1431*7f296bb3SBarry SmithTSSundialsSetMaxl(TS ts,PetscInt restart); 1432*7f296bb3SBarry Smith``` 1433*7f296bb3SBarry Smith 1434*7f296bb3SBarry Smithsets the number of vectors in the Krylov subpspace used by GMRES. This 1435*7f296bb3SBarry Smithmay be set in the options database with `-ts_sundials_maxl` `maxl`. 1436*7f296bb3SBarry Smith 1437*7f296bb3SBarry Smith## Using TChem from PETSc 1438*7f296bb3SBarry Smith 1439*7f296bb3SBarry SmithTChem [^id7] is a package originally developed at Sandia National 1440*7f296bb3SBarry SmithLaboratory that can read in CHEMKIN [^id8] data files and compute the 1441*7f296bb3SBarry Smithright-hand side function and its Jacobian for a reaction ODE system. To 1442*7f296bb3SBarry Smithutilize PETSc’s ODE solvers for these systems, first install PETSc with 1443*7f296bb3SBarry Smiththe additional `configure` option `--download-tchem`. We currently 1444*7f296bb3SBarry Smithprovide two examples of its use; one for single cell reaction and one 1445*7f296bb3SBarry Smithfor an “artificial” one dimensional problem with periodic boundary 1446*7f296bb3SBarry Smithconditions and diffusion of all species. The self-explanatory examples 1447*7f296bb3SBarry Smithare the 1448*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/extchem.c.html">The TS tutorial extchem</a> 1449*7f296bb3SBarry Smithand 1450*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/extchemfield.c.html">The TS tutorial extchemfield</a>. 1451*7f296bb3SBarry Smith 1452*7f296bb3SBarry Smith[^id5]: If the matrix $F_{\dot{u}}(t) = \partial F 1453*7f296bb3SBarry Smith / \partial \dot{u}$ is nonsingular then it is an ODE and can be 1454*7f296bb3SBarry Smith transformed to the standard explicit form, although this 1455*7f296bb3SBarry Smith transformation may not lead to efficient algorithms. 1456*7f296bb3SBarry Smith 1457*7f296bb3SBarry Smith[^id6]: PETSc will automatically translate the function provided to the 1458*7f296bb3SBarry Smith appropriate form. 1459*7f296bb3SBarry Smith 1460*7f296bb3SBarry Smith[^id7]: [bitbucket.org/jedbrown/tchem](https://bitbucket.org/jedbrown/tchem) 1461*7f296bb3SBarry Smith 1462*7f296bb3SBarry Smith[^id8]: [en.wikipedia.org/wiki/CHEMKIN](https://en.wikipedia.org/wiki/CHEMKIN) 1463*7f296bb3SBarry Smith 1464*7f296bb3SBarry Smith```{raw} html 1465*7f296bb3SBarry Smith<hr> 1466*7f296bb3SBarry Smith``` 1467*7f296bb3SBarry Smith 1468*7f296bb3SBarry Smith# Solving Steady-State Problems with Pseudo-Timestepping 1469*7f296bb3SBarry Smith 1470*7f296bb3SBarry Smith**Simple Example:** `TS` provides a general code for performing pseudo 1471*7f296bb3SBarry Smithtimestepping with a variable timestep at each physical node point. For 1472*7f296bb3SBarry Smithexample, instead of directly attacking the steady-state problem 1473*7f296bb3SBarry Smith 1474*7f296bb3SBarry Smith$$ 1475*7f296bb3SBarry SmithG(u) = 0, 1476*7f296bb3SBarry Smith$$ 1477*7f296bb3SBarry Smith 1478*7f296bb3SBarry Smithwe can use pseudo-transient continuation by solving 1479*7f296bb3SBarry Smith 1480*7f296bb3SBarry Smith$$ 1481*7f296bb3SBarry Smithu_t = G(u). 1482*7f296bb3SBarry Smith$$ 1483*7f296bb3SBarry Smith 1484*7f296bb3SBarry SmithUsing time differencing 1485*7f296bb3SBarry Smith 1486*7f296bb3SBarry Smith$$ 1487*7f296bb3SBarry Smithu_t \doteq \frac{{u^{n+1}} - {u^{n}} }{dt^{n}} 1488*7f296bb3SBarry Smith$$ 1489*7f296bb3SBarry Smith 1490*7f296bb3SBarry Smithwith the backward Euler method, we obtain nonlinear equations at a 1491*7f296bb3SBarry Smithseries of pseudo-timesteps 1492*7f296bb3SBarry Smith 1493*7f296bb3SBarry Smith$$ 1494*7f296bb3SBarry Smith\frac{1}{dt^n} B (u^{n+1} - u^{n} ) = G(u^{n+1}). 1495*7f296bb3SBarry Smith$$ 1496*7f296bb3SBarry Smith 1497*7f296bb3SBarry SmithFor this problem the user must provide $G(u)$, the time steps 1498*7f296bb3SBarry Smith$dt^{n}$ and the left-hand-side matrix $B$ (or optionally, 1499*7f296bb3SBarry Smithif the timestep is position independent and $B$ is the identity 1500*7f296bb3SBarry Smithmatrix, a scalar timestep), as well as optionally the Jacobian of 1501*7f296bb3SBarry Smith$G(u)$. 1502*7f296bb3SBarry Smith 1503*7f296bb3SBarry SmithMore generally, this can be applied to implicit ODE and DAE for which 1504*7f296bb3SBarry Smiththe transient form is 1505*7f296bb3SBarry Smith 1506*7f296bb3SBarry Smith$$ 1507*7f296bb3SBarry SmithF(u,\dot{u}) = 0. 1508*7f296bb3SBarry Smith$$ 1509*7f296bb3SBarry Smith 1510*7f296bb3SBarry SmithFor solving steady-state problems with pseudo-timestepping one proceeds 1511*7f296bb3SBarry Smithas follows. 1512*7f296bb3SBarry Smith 1513*7f296bb3SBarry Smith- Provide the function `G(u)` with the routine 1514*7f296bb3SBarry Smith 1515*7f296bb3SBarry Smith ``` 1516*7f296bb3SBarry Smith TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *fP); 1517*7f296bb3SBarry Smith ``` 1518*7f296bb3SBarry Smith 1519*7f296bb3SBarry Smith The arguments to the function `f()` are the timestep context, the 1520*7f296bb3SBarry Smith current time, the input for the function, the output for the function 1521*7f296bb3SBarry Smith and the (optional) user-provided context variable `fP`. 1522*7f296bb3SBarry Smith 1523*7f296bb3SBarry Smith- Provide the (approximate) Jacobian matrix of `G(u)` and a function 1524*7f296bb3SBarry Smith to compute it at each Newton iteration. This is done with the command 1525*7f296bb3SBarry Smith 1526*7f296bb3SBarry Smith ``` 1527*7f296bb3SBarry Smith TSSetRHSJacobian(TS ts,Mat Amat, Mat Pmat,PetscErrorCode (*f)(TS,PetscReal,Vec,Mat,Mat,void*),void *fP); 1528*7f296bb3SBarry Smith ``` 1529*7f296bb3SBarry Smith 1530*7f296bb3SBarry Smith The arguments for the function `f()` are the timestep context, the 1531*7f296bb3SBarry Smith current time, the location where the Jacobian is to be computed, the 1532*7f296bb3SBarry Smith (approximate) Jacobian matrix, an alternative approximate Jacobian 1533*7f296bb3SBarry Smith matrix used to construct the preconditioner, and the optional 1534*7f296bb3SBarry Smith user-provided context, passed in as `fP`. The user must provide the 1535*7f296bb3SBarry Smith Jacobian as a matrix; thus, if using a matrix-free approach, one must 1536*7f296bb3SBarry Smith create a `MATSHELL` matrix. 1537*7f296bb3SBarry Smith 1538*7f296bb3SBarry SmithIn addition, the user must provide a routine that computes the 1539*7f296bb3SBarry Smithpseudo-timestep. This is slightly different depending on if one is using 1540*7f296bb3SBarry Smitha constant timestep over the entire grid, or it varies with location. 1541*7f296bb3SBarry Smith 1542*7f296bb3SBarry Smith- For location-independent pseudo-timestepping, one uses the routine 1543*7f296bb3SBarry Smith 1544*7f296bb3SBarry Smith ``` 1545*7f296bb3SBarry Smith TSPseudoSetTimeStep(TS ts,PetscInt(*dt)(TS,PetscReal*,void*),void* dtctx); 1546*7f296bb3SBarry Smith ``` 1547*7f296bb3SBarry Smith 1548*7f296bb3SBarry Smith The function `dt` is a user-provided function that computes the 1549*7f296bb3SBarry Smith next pseudo-timestep. As a default one can use 1550*7f296bb3SBarry Smith `TSPseudoTimeStepDefault(TS,PetscReal*,void*)` for `dt`. This 1551*7f296bb3SBarry Smith routine updates the pseudo-timestep with one of two strategies: the 1552*7f296bb3SBarry Smith default 1553*7f296bb3SBarry Smith 1554*7f296bb3SBarry Smith $$ 1555*7f296bb3SBarry Smith dt^{n} = dt_{\mathrm{increment}}*dt^{n-1}*\frac{|| F(u^{n-1}) ||}{|| F(u^{n})||} 1556*7f296bb3SBarry Smith $$ 1557*7f296bb3SBarry Smith 1558*7f296bb3SBarry Smith or, the alternative, 1559*7f296bb3SBarry Smith 1560*7f296bb3SBarry Smith $$ 1561*7f296bb3SBarry Smith dt^{n} = dt_{\mathrm{increment}}*dt^{0}*\frac{|| F(u^{0}) ||}{|| F(u^{n})||} 1562*7f296bb3SBarry Smith $$ 1563*7f296bb3SBarry Smith 1564*7f296bb3SBarry Smith which can be set with the call 1565*7f296bb3SBarry Smith 1566*7f296bb3SBarry Smith ``` 1567*7f296bb3SBarry Smith TSPseudoIncrementDtFromInitialDt(TS ts); 1568*7f296bb3SBarry Smith ``` 1569*7f296bb3SBarry Smith 1570*7f296bb3SBarry Smith or the option `-ts_pseudo_increment_dt_from_initial_dt`. The value 1571*7f296bb3SBarry Smith $dt_{\mathrm{increment}}$ is by default $1.1$, but can be 1572*7f296bb3SBarry Smith reset with the call 1573*7f296bb3SBarry Smith 1574*7f296bb3SBarry Smith ``` 1575*7f296bb3SBarry Smith TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc); 1576*7f296bb3SBarry Smith ``` 1577*7f296bb3SBarry Smith 1578*7f296bb3SBarry Smith or the option `-ts_pseudo_increment <inc>`. 1579*7f296bb3SBarry Smith 1580*7f296bb3SBarry Smith- For location-dependent pseudo-timestepping, the interface function 1581*7f296bb3SBarry Smith has not yet been created. 1582*7f296bb3SBarry Smith 1583*7f296bb3SBarry Smith```{eval-rst} 1584*7f296bb3SBarry Smith.. bibliography:: /petsc.bib 1585*7f296bb3SBarry Smith :filter: docname in docnames 1586*7f296bb3SBarry Smith 1587*7f296bb3SBarry Smith``` 1588