xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 26b5f05de6058eee318c007462368b24cc8d9cbb)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5a055c8caSBarry Smith #include <petscsnes.h>
61e25c274SJed Brown #include <petscdm.h>
73c54f92cSHong Zhang #include <petscmat.h>
8316643e7SJed Brown 
9316643e7SJed Brown typedef struct {
10715f1b00SHong Zhang   /* context for time stepping */
111566a47fSLisandro Dalcin   PetscReal    stage_time;
12cc4f23bcSHong Zhang   Vec          Stages[2];   /* Storage for stage solutions */
13cc4f23bcSHong Zhang   Vec          X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */
141566a47fSLisandro Dalcin   Vec          affine;      /* Affine vector needed for residual at beginning of step in endpoint formulation */
151566a47fSLisandro Dalcin   PetscReal    Theta;
161a352fa8SHong Zhang   PetscReal    shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20715f1b00SHong Zhang   TSStepStatus status;
211a352fa8SHong Zhang   Vec          VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */
221a352fa8SHong Zhang   PetscReal    ptime0;           /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
231a352fa8SHong Zhang   PetscReal    time_step0;       /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/
241566a47fSLisandro Dalcin 
25715f1b00SHong Zhang   /* context for sensitivity analysis */
26022c081aSHong Zhang   PetscInt num_tlm;               /* Total number of tangent linear equations */
27b5e0532dSHong Zhang   Vec     *VecsDeltaLam;          /* Increment of the adjoint sensitivity w.r.t IC at stage */
28b5e0532dSHong Zhang   Vec     *VecsDeltaMu;           /* Increment of the adjoint sensitivity w.r.t P at stage */
2913b1b0ffSHong Zhang   Vec     *VecsSensiTemp;         /* Vector to be multiplied with Jacobian transpose */
30cc4f23bcSHong Zhang   Mat      MatFwdStages[2];       /* TLM Stages */
3113b1b0ffSHong Zhang   Mat      MatDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
32b5b4017aSHong Zhang   Vec      VecDeltaFwdSensipCol;  /* Working vector for holding one column of the sensitivity matrix */
3313b1b0ffSHong Zhang   Mat      MatFwdSensip0;         /* backup for roll-backs due to events */
347207e0fdSHong Zhang   Mat      MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */
357207e0fdSHong Zhang   Mat      MatIntegralSensip0;    /* backup for roll-backs due to events */
36b5b4017aSHong Zhang   Vec     *VecsDeltaLam2;         /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
37b5b4017aSHong Zhang   Vec     *VecsDeltaMu2;          /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
38b5b4017aSHong Zhang   Vec     *VecsSensi2Temp;        /* Working vectors that holds the residual for the second-order adjoint */
39b5b4017aSHong Zhang   Vec     *VecsAffine;            /* Working vectors to store residuals */
40715f1b00SHong Zhang   /* context for error estimation */
411566a47fSLisandro Dalcin   Vec vec_sol_prev;
421566a47fSLisandro Dalcin   Vec vec_lte_work;
43316643e7SJed Brown } TS_Theta;
44316643e7SJed Brown 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
46d71ae5a4SJacob Faibussowitsch {
477445fe48SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
487445fe48SJed Brown 
497445fe48SJed Brown   PetscFunctionBegin;
507445fe48SJed Brown   if (X0) {
517445fe48SJed Brown     if (dm && dm != ts->dm) {
529566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0));
537445fe48SJed Brown     } else *X0 = ts->vec_sol;
547445fe48SJed Brown   }
557445fe48SJed Brown   if (Xdot) {
567445fe48SJed Brown     if (dm && dm != ts->dm) {
579566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
587445fe48SJed Brown     } else *Xdot = th->Xdot;
597445fe48SJed Brown   }
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617445fe48SJed Brown }
627445fe48SJed Brown 
63d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
64d71ae5a4SJacob Faibussowitsch {
650d0b770aSPeter Brune   PetscFunctionBegin;
660d0b770aSPeter Brune   if (X0) {
6748a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0));
680d0b770aSPeter Brune   }
690d0b770aSPeter Brune   if (Xdot) {
7048a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
710d0b770aSPeter Brune   }
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
730d0b770aSPeter Brune }
740d0b770aSPeter Brune 
75d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx)
76d71ae5a4SJacob Faibussowitsch {
777445fe48SJed Brown   PetscFunctionBegin;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
797445fe48SJed Brown }
807445fe48SJed Brown 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
82d71ae5a4SJacob Faibussowitsch {
837445fe48SJed Brown   TS  ts = (TS)ctx;
847445fe48SJed Brown   Vec X0, Xdot, X0_c, Xdot_c;
857445fe48SJed Brown 
867445fe48SJed Brown   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot));
889566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
899566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, X0, X0_c));
909566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
919566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
929566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
939566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot));
949566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
967445fe48SJed Brown }
977445fe48SJed Brown 
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx)
99d71ae5a4SJacob Faibussowitsch {
100258e1594SPeter Brune   PetscFunctionBegin;
1013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102258e1594SPeter Brune }
103258e1594SPeter Brune 
104d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
105d71ae5a4SJacob Faibussowitsch {
106258e1594SPeter Brune   TS  ts = (TS)ctx;
107258e1594SPeter Brune   Vec X0, Xdot, X0_sub, Xdot_sub;
108258e1594SPeter Brune 
109258e1594SPeter Brune   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
1119566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
112258e1594SPeter Brune 
1139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
1149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
115258e1594SPeter Brune 
1169566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
1179566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
118258e1594SPeter Brune 
1199566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
1209566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122258e1594SPeter Brune }
123258e1594SPeter Brune 
124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
125d71ae5a4SJacob Faibussowitsch {
126022c081aSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
127cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
128022c081aSHong Zhang 
129022c081aSHong Zhang   PetscFunctionBegin;
130022c081aSHong Zhang   if (th->endpoint) {
131022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
132022c081aSHong Zhang     if (th->Theta != 1.0) {
1339566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand));
1349566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand));
135022c081aSHong Zhang     }
1369566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand));
1379566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand));
138022c081aSHong Zhang   } else {
1399566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand));
1409566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand));
141022c081aSHong Zhang   }
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143022c081aSHong Zhang }
144022c081aSHong Zhang 
145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
146d71ae5a4SJacob Faibussowitsch {
147b1cb36f3SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
148cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
149b1cb36f3SHong Zhang 
150b1cb36f3SHong Zhang   PetscFunctionBegin;
151b1cb36f3SHong Zhang   /* backup cost integral */
1529566063dSJacob Faibussowitsch   PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0));
1539566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155b1cb36f3SHong Zhang }
156b1cb36f3SHong Zhang 
157d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
158d71ae5a4SJacob Faibussowitsch {
1591a352fa8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
160b1cb36f3SHong Zhang 
161b1cb36f3SHong Zhang   PetscFunctionBegin;
1621a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
1631a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
1641a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
1659566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16724655328SShri }
16824655328SShri 
169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x)
170d71ae5a4SJacob Faibussowitsch {
1711566a47fSLisandro Dalcin   PetscInt nits, lits;
1721566a47fSLisandro Dalcin 
1731566a47fSLisandro Dalcin   PetscFunctionBegin;
1749566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
1759566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1769566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1779371c9d4SSatish Balay   ts->snes_its += nits;
1789371c9d4SSatish Balay   ts->ksp_its += lits;
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1801566a47fSLisandro Dalcin }
1811566a47fSLisandro Dalcin 
182*26b5f05dSStefano Zampini /* We need to transfer X0 which will be copied into sol_prev */
183*26b5f05dSStefano Zampini static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg)
184*26b5f05dSStefano Zampini {
185*26b5f05dSStefano Zampini   TS_Theta  *th     = (TS_Theta *)ts->data;
186*26b5f05dSStefano Zampini   const char name[] = "ts:theta:X0";
187*26b5f05dSStefano Zampini 
188*26b5f05dSStefano Zampini   PetscFunctionBegin;
189*26b5f05dSStefano Zampini   if (reg && th->vec_sol_prev) {
190*26b5f05dSStefano Zampini     PetscCall(TSResizeRegisterVec(ts, name, th->X0));
191*26b5f05dSStefano Zampini   } else if (!reg) {
192*26b5f05dSStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, name, &th->X0));
193*26b5f05dSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->X0));
194*26b5f05dSStefano Zampini   }
195*26b5f05dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
196*26b5f05dSStefano Zampini }
197*26b5f05dSStefano Zampini 
198d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Theta(TS ts)
199d71ae5a4SJacob Faibussowitsch {
200316643e7SJed Brown   TS_Theta *th         = (TS_Theta *)ts->data;
2011566a47fSLisandro Dalcin   PetscInt  rejections = 0;
2024957b756SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
2031566a47fSLisandro Dalcin   PetscReal next_time_step = ts->time_step;
204316643e7SJed Brown 
205316643e7SJed Brown   PetscFunctionBegin;
2061566a47fSLisandro Dalcin   if (!ts->steprollback) {
2079566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
2089566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
2091566a47fSLisandro Dalcin   }
2101566a47fSLisandro Dalcin 
2111566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
21299260352SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
2133e05ceb1SHong Zhang     th->shift      = 1 / (th->Theta * ts->time_step);
2141566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step;
2159566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X));
21648a46eb9SPierre Jolivet     if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot));
217eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
2189566063dSJacob Faibussowitsch       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
2199566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
2209566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE));
2219566063dSJacob Faibussowitsch       PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta));
2221566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2239566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->affine));
224eb284becSJed Brown     }
2259566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
2269566063dSJacob Faibussowitsch     PetscCall(TSTheta_SNESSolve(ts, th->affine, th->X));
2279566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
2289566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
229fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
230051f2191SLisandro Dalcin 
231051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2321566a47fSLisandro Dalcin     if (th->endpoint) {
2339566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X, ts->vec_sol));
2341566a47fSLisandro Dalcin     } else {
2359566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
2369566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
2371566a47fSLisandro Dalcin     }
2389566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2391566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2401566a47fSLisandro Dalcin     if (!accept) {
2419566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
242be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
243be5899b3SLisandro Dalcin       goto reject_step;
244051f2191SLisandro Dalcin     }
2453b1890cdSShri Abhyankar 
246715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2471a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2481a352fa8SHong Zhang       th->time_step0 = ts->time_step;
24917145e75SHong Zhang     }
2502b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
251cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
252051f2191SLisandro Dalcin     break;
253051f2191SLisandro Dalcin 
254051f2191SLisandro Dalcin   reject_step:
2559371c9d4SSatish Balay     ts->reject++;
2569371c9d4SSatish Balay     accept = PETSC_FALSE;
2571566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
258051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
25963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
2603b1890cdSShri Abhyankar     }
2613b1890cdSShri Abhyankar   }
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263316643e7SJed Brown }
264316643e7SJed Brown 
265d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
266d71ae5a4SJacob Faibussowitsch {
267c9681e5dSHong Zhang   TS_Theta      *th           = (TS_Theta *)ts->data;
268cd4cee2dSHong Zhang   TS             quadts       = ts->quadraturets;
269c9681e5dSHong Zhang   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
270c9681e5dSHong Zhang   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
271c9681e5dSHong Zhang   PetscInt       nadj;
2727207e0fdSHong Zhang   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
273c9681e5dSHong Zhang   KSP            ksp;
274c9681e5dSHong Zhang   PetscScalar   *xarr;
275077c4feaSHong Zhang   TSEquationType eqtype;
276077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2771a352fa8SHong Zhang   PetscReal      adjoint_time_step;
278c9681e5dSHong Zhang 
279c9681e5dSHong Zhang   PetscFunctionBegin;
2809566063dSJacob Faibussowitsch   PetscCall(TSGetEquationType(ts, &eqtype));
281077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
282077c4feaSHong Zhang     isexplicitode = PETSC_TRUE;
283077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
284077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
285077c4feaSHong Zhang   }
286c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
2879566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
2889566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
2897207e0fdSHong Zhang   if (quadts) {
2909566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
2919566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
2927207e0fdSHong Zhang   }
293c9681e5dSHong Zhang 
2941a352fa8SHong Zhang   th->stage_time    = ts->ptime;
2951a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
296c9681e5dSHong Zhang 
29733f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2981baa6e33SBarry Smith   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
299cd4cee2dSHong Zhang 
300c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
3019566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
3029566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
303cd4cee2dSHong Zhang     if (quadJ) {
3049566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
3059566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
3069566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
3079566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
3089566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
309c9681e5dSHong Zhang     }
310c9681e5dSHong Zhang   }
311c9681e5dSHong Zhang 
312c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
313c87ba875SIvan Yashchuk   th->shift = 1. / adjoint_time_step;
3149566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3159566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
316c9681e5dSHong Zhang 
317c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
318c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
319c9681e5dSHong Zhang     KSPConvergedReason kspreason;
3209566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
3219566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
322c9681e5dSHong Zhang     if (kspreason < 0) {
323c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32463a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
325c9681e5dSHong Zhang     }
326c9681e5dSHong Zhang   }
327c9681e5dSHong Zhang 
328c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
329c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3309566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
3319566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
332c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3339566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
334c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3359566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
336c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
3379566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
3389566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
3399566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
34048a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
341c9681e5dSHong Zhang     }
342c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
343c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
34405c9950eSHong Zhang       KSPConvergedReason kspreason;
3459566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
3469566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
34705c9950eSHong Zhang       if (kspreason < 0) {
34805c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
34963a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
35005c9950eSHong Zhang       }
351c9681e5dSHong Zhang     }
352c9681e5dSHong Zhang   }
353c9681e5dSHong Zhang 
354c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
355077c4feaSHong Zhang   if (!isexplicitode) {
3563e05ceb1SHong Zhang     th->shift = 0.0;
3579566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3589566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
359c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
36033f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3619566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
3629566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
3639566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
3649566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
365c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3669566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
3679566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
3689566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
369c9681e5dSHong Zhang       }
370c9681e5dSHong Zhang     }
371077c4feaSHong Zhang   }
372c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3739566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
3749566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
3751baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
376c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
377c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3789566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
37933f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3809566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
381c9681e5dSHong Zhang     }
382cd4cee2dSHong Zhang 
383c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
3849566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
3859566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
386cd4cee2dSHong Zhang       if (quadJp) {
3879566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
3889566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
3899566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
3909566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3919566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
39233f72c5dSHong Zhang       }
393c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
3949566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
3959566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
39648a46eb9SPierre Jolivet         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
39748a46eb9SPierre Jolivet         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
398c9681e5dSHong Zhang       }
399c9681e5dSHong Zhang     }
400c9681e5dSHong Zhang   }
401c9681e5dSHong Zhang 
402c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
4039566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4049566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
405c9681e5dSHong Zhang   }
406c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408c9681e5dSHong Zhang }
409c9681e5dSHong Zhang 
410d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_Theta(TS ts)
411d71ae5a4SJacob Faibussowitsch {
4122ca6e920SHong Zhang   TS_Theta    *th           = (TS_Theta *)ts->data;
413cd4cee2dSHong Zhang   TS           quadts       = ts->quadraturets;
414b5e0532dSHong Zhang   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
415b5b4017aSHong Zhang   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
4162ca6e920SHong Zhang   PetscInt     nadj;
4177207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
4182ca6e920SHong Zhang   KSP          ksp;
419b5b4017aSHong Zhang   PetscScalar *xarr;
4201a352fa8SHong Zhang   PetscReal    adjoint_time_step;
421da81f932SPierre Jolivet   PetscReal    adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, usually ts->ptime is larger than adjoint_ptime) */
4222ca6e920SHong Zhang 
4232ca6e920SHong Zhang   PetscFunctionBegin;
424c9681e5dSHong Zhang   if (th->Theta == 1.) {
4259566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
4263ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
427c9681e5dSHong Zhang   }
4282ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4299566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
4309566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
4317207e0fdSHong Zhang   if (quadts) {
4329566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
4339566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
4347207e0fdSHong Zhang   }
435b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4361a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
4371a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4381a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
4392ca6e920SHong Zhang 
44082ad9101SHong Zhang   if (!th->endpoint) {
4415072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4429566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
44382ad9101SHong Zhang   }
44482ad9101SHong Zhang 
445b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
44633f72c5dSHong Zhang   /* Cost function has an integral term */
4477207e0fdSHong Zhang   if (quadts) {
44805755b9cSHong Zhang     if (th->endpoint) {
4499566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
45005755b9cSHong Zhang     } else {
4519566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
45205755b9cSHong Zhang     }
4537207e0fdSHong Zhang   }
45433f72c5dSHong Zhang 
455abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4569566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
4579566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
458cd4cee2dSHong Zhang     if (quadJ) {
4599566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
4609566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
4619566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
4629566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4639566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
46436eaed60SHong Zhang     }
4652ca6e920SHong Zhang   }
4663c54f92cSHong Zhang 
467b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4681a352fa8SHong Zhang   th->shift = 1. / (th->Theta * adjoint_time_step);
4693c54f92cSHong Zhang   if (th->endpoint) {
4709566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
4713c54f92cSHong Zhang   } else {
4729566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
4733c54f92cSHong Zhang   }
4749566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
4752ca6e920SHong Zhang 
476b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
477abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4781d14d03bSHong Zhang     KSPConvergedReason kspreason;
4799566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
4809566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
4811d14d03bSHong Zhang     if (kspreason < 0) {
4821d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
48363a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
4841d14d03bSHong Zhang     }
4852ca6e920SHong Zhang   }
4863c54f92cSHong Zhang 
4871d14d03bSHong Zhang   /* Second-order adjoint */
488b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4893c633725SBarry Smith     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
490b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4919566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
4929566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
493b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
4949566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
4959566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4969566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
497b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
4989566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
499b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
5009566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
5019566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
5029566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
50348a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
504b5b4017aSHong Zhang     }
505b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
506b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
5071d14d03bSHong Zhang       KSPConvergedReason kspreason;
5089566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
5099566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
5101d14d03bSHong Zhang       if (kspreason < 0) {
5111d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
51263a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
5131d14d03bSHong Zhang       }
514b5b4017aSHong Zhang     }
515b5b4017aSHong Zhang   }
516b5b4017aSHong Zhang 
51736eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5181d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5191a352fa8SHong Zhang     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
5201a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
5219566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
5229566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
52333f72c5dSHong Zhang     /* R_U at t_n */
5241baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
525abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
5269566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
527cd4cee2dSHong Zhang       if (quadJ) {
5289566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
5299566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
5309566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
5319566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5329566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
533b5b4017aSHong Zhang       }
5349566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
535b5b4017aSHong Zhang     }
5361d14d03bSHong Zhang 
5371d14d03bSHong Zhang     /* Second-order adjoint */
5381d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
539b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5409566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
5419566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
542b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5439566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
5449566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5459566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
546b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5479566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
548b5b4017aSHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
54933f72c5dSHong Zhang         /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2  */
5509566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
5519566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
55248a46eb9SPierre Jolivet         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
5539566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
554b5b4017aSHong Zhang       }
555b5e0532dSHong Zhang     }
5563fd52205SHong Zhang 
557c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
558c586f2e8SHong Zhang 
5593fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
560b5b4017aSHong Zhang       /* U_{n+1} */
56182ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
5629566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
5639566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5641baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
565abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5669566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
5679566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
5683b711c3fSHong Zhang         if (quadJp) {
5699566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
5709566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
5719566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
5729566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5739566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
5743b711c3fSHong Zhang         }
5753fd52205SHong Zhang       }
57633f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
57733f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5789566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
5799566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
580b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5819566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
5829566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5839566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
58433f72c5dSHong Zhang 
585b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5869566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
587b5b4017aSHong Zhang         for (nadj = 0; nadj < ts->numcost; nadj++) {
58833f72c5dSHong Zhang           /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2)  */
5899566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
5909566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
59148a46eb9SPierre Jolivet           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
59248a46eb9SPierre Jolivet           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
593b5b4017aSHong Zhang         }
594b5b4017aSHong Zhang       }
595b5b4017aSHong Zhang 
596b5b4017aSHong Zhang       /* U_s */
5979566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
5989566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5991baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
600abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6019566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6029566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
6033b711c3fSHong Zhang         if (quadJp) {
6049566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6059566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6069566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
6079566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6089566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
6093b711c3fSHong Zhang         }
61033f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
61133f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
6129566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
6139566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
614b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6159566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
6169566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
6179566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
618b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6199566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
620b5b4017aSHong Zhang           for (nadj = 0; nadj < ts->numcost; nadj++) {
62133f72c5dSHong Zhang             /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */
6229566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
6239566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
62448a46eb9SPierre Jolivet             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
62548a46eb9SPierre Jolivet             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
62636eaed60SHong Zhang           }
62736eaed60SHong Zhang         }
6283fd52205SHong Zhang       }
629b5e0532dSHong Zhang     }
6303fd52205SHong Zhang   } else { /* one-stage case */
6313e05ceb1SHong Zhang     th->shift = 0.0;
6329566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
6339566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
6341baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
635abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
6369566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
6379566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
638cd4cee2dSHong Zhang       if (quadJ) {
6399566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
6409566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
6419566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
6429566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6439566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
64436eaed60SHong Zhang       }
6452ca6e920SHong Zhang     }
6463fd52205SHong Zhang     if (ts->vecs_sensip) {
64782ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
6489566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
6499566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
6501baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
651abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6529566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6539566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
654cd4cee2dSHong Zhang         if (quadJp) {
6559566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6569566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6579566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
6589566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6599566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
66036eaed60SHong Zhang         }
66136eaed60SHong Zhang       }
6623fd52205SHong Zhang     }
6632ca6e920SHong Zhang   }
6642ca6e920SHong Zhang 
6652ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6672ca6e920SHong Zhang }
6682ca6e920SHong Zhang 
669d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X)
670d71ae5a4SJacob Faibussowitsch {
671cd652676SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
672be5899b3SLisandro Dalcin   PetscReal dt = t - ts->ptime;
673cd652676SJed Brown 
674cd652676SJed Brown   PetscFunctionBegin;
6759566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X));
676be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6779566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
6783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
679cd652676SJed Brown }
680cd652676SJed Brown 
681d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
682d71ae5a4SJacob Faibussowitsch {
6831566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
6841566a47fSLisandro Dalcin   Vec       X  = ts->vec_sol;      /* X = solution */
6851566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
6867453f775SEmil Constantinescu   PetscReal wltea, wlter;
6871566a47fSLisandro Dalcin 
6881566a47fSLisandro Dalcin   PetscFunctionBegin;
6899371c9d4SSatish Balay   if (!th->vec_sol_prev) {
6909371c9d4SSatish Balay     *wlte = -1;
6913ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6929371c9d4SSatish Balay   }
6931566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
6949371c9d4SSatish Balay   if (ts->steprestart) {
6959371c9d4SSatish Balay     *wlte = -1;
6963ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6979371c9d4SSatish Balay   }
6981566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
699fecfb714SLisandro Dalcin   {
700be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
701be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
7029371c9d4SSatish Balay     PetscScalar scal[3];
7039371c9d4SSatish Balay     Vec         vecs[3];
7049371c9d4SSatish Balay     scal[0] = +1 / a;
7059371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
7069371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
7079371c9d4SSatish Balay     vecs[0] = X;
7089371c9d4SSatish Balay     vecs[1] = th->X0;
7099371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
7109566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
7119566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
7129566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
7131566a47fSLisandro Dalcin   }
7141566a47fSLisandro Dalcin   if (order) *order = 2;
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7161566a47fSLisandro Dalcin }
7171566a47fSLisandro Dalcin 
718d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Theta(TS ts)
719d71ae5a4SJacob Faibussowitsch {
7201566a47fSLisandro Dalcin   TS_Theta *th     = (TS_Theta *)ts->data;
7217207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
7221566a47fSLisandro Dalcin 
7231566a47fSLisandro Dalcin   PetscFunctionBegin;
7249566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
72548a46eb9SPierre Jolivet   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
726715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
7271baa6e33SBarry Smith   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
72848a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
7293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
730715f1b00SHong Zhang }
731715f1b00SHong Zhang 
732d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_Theta(TS ts)
733d71ae5a4SJacob Faibussowitsch {
734715f1b00SHong Zhang   TS_Theta    *th                   = (TS_Theta *)ts->data;
735cd4cee2dSHong Zhang   TS           quadts               = ts->quadraturets;
73613b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
737b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7387207e0fdSHong Zhang   PetscInt     ntlm;
739715f1b00SHong Zhang   KSP          ksp;
7407207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
74113b1b0ffSHong Zhang   PetscScalar *barr, *xarr;
7421a352fa8SHong Zhang   PetscReal    previous_shift;
743715f1b00SHong Zhang 
744715f1b00SHong Zhang   PetscFunctionBegin;
7451a352fa8SHong Zhang   previous_shift = th->shift;
7469566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));
74713b1b0ffSHong Zhang 
74848a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
7499566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
7509566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
7517207e0fdSHong Zhang   if (quadts) {
7529566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
7539566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
7547207e0fdSHong Zhang   }
755715f1b00SHong Zhang 
756715f1b00SHong Zhang   /* Build RHS */
757715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7581a352fa8SHong Zhang     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
7599566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7609566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7619566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));
762715f1b00SHong Zhang 
763022c081aSHong Zhang     /* Add the f_p forcing terms */
7640e11c21fSHong Zhang     if (ts->Jacp) {
7659566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7669566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7679566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
76882ad9101SHong Zhang       th->shift = previous_shift;
7699566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
7709566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7719566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
7720e11c21fSHong Zhang     }
773715f1b00SHong Zhang   } else { /* 1-stage method */
7743e05ceb1SHong Zhang     th->shift = 0.0;
7759566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7769566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7779566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, -1.));
77813b1b0ffSHong Zhang 
779022c081aSHong Zhang     /* Add the f_p forcing terms */
7800e11c21fSHong Zhang     if (ts->Jacp) {
78182ad9101SHong Zhang       th->shift = previous_shift;
7829566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
7839566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7849566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
785715f1b00SHong Zhang     }
7860e11c21fSHong Zhang   }
787715f1b00SHong Zhang 
788715f1b00SHong Zhang   /* Build LHS */
7891a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
790715f1b00SHong Zhang   if (th->endpoint) {
7919566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
792715f1b00SHong Zhang   } else {
7939566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
794715f1b00SHong Zhang   }
7959566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
796715f1b00SHong Zhang 
797715f1b00SHong Zhang   /*
798715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
799c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
800715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
801715f1b00SHong Zhang   */
802715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8037207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8049566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
8059566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
8069566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8079566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8089566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
809715f1b00SHong Zhang     }
810715f1b00SHong Zhang   }
811715f1b00SHong Zhang 
812715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
813022c081aSHong Zhang   for (ntlm = 0; ntlm < th->num_tlm; ntlm++) {
814b5b4017aSHong Zhang     KSPConvergedReason kspreason;
8159566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
8169566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
817715f1b00SHong Zhang     if (th->endpoint) {
8189566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
8199566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
8209566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
8219566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
8229566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
823715f1b00SHong Zhang     } else {
8249566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
825715f1b00SHong Zhang     }
8269566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
827b5b4017aSHong Zhang     if (kspreason < 0) {
828b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
82963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
830b5b4017aSHong Zhang     }
8319566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
8329566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
833715f1b00SHong Zhang   }
834715f1b00SHong Zhang 
83513b1b0ffSHong Zhang   /*
83613b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
837c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
83813b1b0ffSHong Zhang   */
8397207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
84013b1b0ffSHong Zhang     if (!th->endpoint) {
8419566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
8429566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
8439566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
8449566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8459566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8469566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
8479566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
84813b1b0ffSHong Zhang     } else {
8499566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
8509566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
8519566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8529566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8539566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
85413b1b0ffSHong Zhang     }
85513b1b0ffSHong Zhang   } else {
85648a46eb9SPierre Jolivet     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
857715f1b00SHong Zhang   }
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8591566a47fSLisandro Dalcin }
8601566a47fSLisandro Dalcin 
861d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[])
862d71ae5a4SJacob Faibussowitsch {
8636affb6f8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
8646affb6f8SHong Zhang 
8656affb6f8SHong Zhang   PetscFunctionBegin;
8667409eceaSHong Zhang   if (ns) {
8677409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8687409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
8697409eceaSHong Zhang   }
8705072d23cSHong Zhang   if (stagesensip) {
871cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8727409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
873cc4f23bcSHong Zhang     } else {
874cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
875cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
876cc4f23bcSHong Zhang     }
877cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8785072d23cSHong Zhang   }
8793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8806affb6f8SHong Zhang }
8816affb6f8SHong Zhang 
882316643e7SJed Brown /*------------------------------------------------------------*/
883d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Theta(TS ts)
884d71ae5a4SJacob Faibussowitsch {
885316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
886316643e7SJed Brown 
887316643e7SJed Brown   PetscFunctionBegin;
8889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
8899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
8909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
8919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
8921566a47fSLisandro Dalcin 
8939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
8949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
8951566a47fSLisandro Dalcin 
8969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
8973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
898277b19d0SLisandro Dalcin }
899277b19d0SLisandro Dalcin 
900d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_Theta(TS ts)
901d71ae5a4SJacob Faibussowitsch {
902ece46509SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
903ece46509SHong Zhang 
904ece46509SHong Zhang   PetscFunctionBegin;
9059566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
9069566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
9079566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
9089566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
9099566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
9109566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
9113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
912ece46509SHong Zhang }
913ece46509SHong Zhang 
914d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Theta(TS ts)
915d71ae5a4SJacob Faibussowitsch {
916277b19d0SLisandro Dalcin   PetscFunctionBegin;
9179566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
918b3a6b972SJed Brown   if (ts->dm) {
9199566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
9209566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
921b3a6b972SJed Brown   }
9229566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
9249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
9259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
9269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
928316643e7SJed Brown }
929316643e7SJed Brown 
930316643e7SJed Brown /*
931316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9322b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9330fd10508SMatthew Knepley 
9340fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9350fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9360fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
937316643e7SJed Brown */
938d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts)
939d71ae5a4SJacob Faibussowitsch {
940316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9417445fe48SJed Brown   Vec       X0, Xdot;
9427445fe48SJed Brown   DM        dm, dmsave;
9433e05ceb1SHong Zhang   PetscReal shift = th->shift;
944316643e7SJed Brown 
945316643e7SJed Brown   PetscFunctionBegin;
9469566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
9475a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9489566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
949494ce9fcSHong Zhang   if (x != X0) {
9509566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
951494ce9fcSHong Zhang   } else {
9529566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
953494ce9fcSHong Zhang   }
9547445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9557445fe48SJed Brown   dmsave = ts->dm;
9567445fe48SJed Brown   ts->dm = dm;
9579566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
9587445fe48SJed Brown   ts->dm = dmsave;
9599566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
9603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
961316643e7SJed Brown }
962316643e7SJed Brown 
963d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts)
964d71ae5a4SJacob Faibussowitsch {
965316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9667445fe48SJed Brown   Vec       Xdot;
9677445fe48SJed Brown   DM        dm, dmsave;
9683e05ceb1SHong Zhang   PetscReal shift = th->shift;
969316643e7SJed Brown 
970316643e7SJed Brown   PetscFunctionBegin;
9719566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
972be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9739566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));
9747445fe48SJed Brown 
9757445fe48SJed Brown   dmsave = ts->dm;
9767445fe48SJed Brown   ts->dm = dm;
9779566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
9787445fe48SJed Brown   ts->dm = dmsave;
9799566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
981316643e7SJed Brown }
982316643e7SJed Brown 
983d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_Theta(TS ts)
984d71ae5a4SJacob Faibussowitsch {
985715f1b00SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
9867207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
987715f1b00SHong Zhang 
988715f1b00SHong Zhang   PetscFunctionBegin;
989022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
99013b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
9919566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
9927207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9939566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
9949566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
995715f1b00SHong Zhang   }
9966f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
9979566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));
99813b1b0ffSHong Zhang 
9999566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1001715f1b00SHong Zhang }
1002715f1b00SHong Zhang 
1003d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_Theta(TS ts)
1004d71ae5a4SJacob Faibussowitsch {
10057adebddeSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
10067207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
10077adebddeSHong Zhang 
10087adebddeSHong Zhang   PetscFunctionBegin;
10097207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10109566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
10119566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
10127adebddeSHong Zhang   }
10139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
10149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
10159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
10163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10177adebddeSHong Zhang }
10187adebddeSHong Zhang 
1019d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Theta(TS ts)
1020d71ae5a4SJacob Faibussowitsch {
1021316643e7SJed Brown   TS_Theta *th     = (TS_Theta *)ts->data;
1022cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
10232ffb9264SLisandro Dalcin   PetscBool match;
1024316643e7SJed Brown 
1025316643e7SJed Brown   PetscFunctionBegin;
1026cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
10279566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
1028b1cb36f3SHong Zhang   }
102948a46eb9SPierre Jolivet   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
103048a46eb9SPierre Jolivet   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
103148a46eb9SPierre Jolivet   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
103248a46eb9SPierre Jolivet   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
10333b1890cdSShri Abhyankar 
10341566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
103520d49056SMatthew G. Knepley   th->shift = 1 / (th->Theta * ts->time_step);
10361566a47fSLisandro Dalcin 
10379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &ts->dm));
10389566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
10399566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
10401566a47fSLisandro Dalcin 
10419566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
10429566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
10442ffb9264SLisandro Dalcin   if (!match) {
10459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
10469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
10473b1890cdSShri Abhyankar   }
10489566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
1049cc4f23bcSHong Zhang 
1050cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
10513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1052b4dd3bf9SBarry Smith }
10530c7376e5SHong Zhang 
1054b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1055b4dd3bf9SBarry Smith 
1056d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1057d71ae5a4SJacob Faibussowitsch {
1058b4dd3bf9SBarry Smith   TS_Theta *th = (TS_Theta *)ts->data;
1059b4dd3bf9SBarry Smith 
1060b4dd3bf9SBarry Smith   PetscFunctionBegin;
10619566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
10629566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
106348a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1064b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10659566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
10669566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
106767633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
106867633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
106967633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1070b5b4017aSHong Zhang   }
1071c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
10729566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
107367633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
107467633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
107567633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1076b5b4017aSHong Zhang   }
10773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1078316643e7SJed Brown }
1079316643e7SJed Brown 
1080d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject)
1081d71ae5a4SJacob Faibussowitsch {
1082316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1083316643e7SJed Brown 
1084316643e7SJed Brown   PetscFunctionBegin;
1085d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1086316643e7SJed Brown   {
10879566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
10889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
10899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1090316643e7SJed Brown   }
1091d0609cedSBarry Smith   PetscOptionsHeadEnd();
10923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1093316643e7SJed Brown }
1094316643e7SJed Brown 
1095d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer)
1096d71ae5a4SJacob Faibussowitsch {
1097316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1098ace3abfcSBarry Smith   PetscBool iascii;
1099316643e7SJed Brown 
1100316643e7SJed Brown   PetscFunctionBegin;
11019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1102316643e7SJed Brown   if (iascii) {
11039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
11049566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1105316643e7SJed Brown   }
11063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1107316643e7SJed Brown }
1108316643e7SJed Brown 
1109d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta)
1110d71ae5a4SJacob Faibussowitsch {
11110de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11120de4c49aSJed Brown 
11130de4c49aSJed Brown   PetscFunctionBegin;
11140de4c49aSJed Brown   *theta = th->Theta;
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11160de4c49aSJed Brown }
11170de4c49aSJed Brown 
1118d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta)
1119d71ae5a4SJacob Faibussowitsch {
11200de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11210de4c49aSJed Brown 
11220de4c49aSJed Brown   PetscFunctionBegin;
1123cad9d221SBarry Smith   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
11240de4c49aSJed Brown   th->Theta = theta;
11251566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11270de4c49aSJed Brown }
1128eb284becSJed Brown 
1129d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint)
1130d71ae5a4SJacob Faibussowitsch {
113126f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
113226f2ff8fSLisandro Dalcin 
113326f2ff8fSLisandro Dalcin   PetscFunctionBegin;
113426f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
11353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113626f2ff8fSLisandro Dalcin }
113726f2ff8fSLisandro Dalcin 
1138d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg)
1139d71ae5a4SJacob Faibussowitsch {
1140eb284becSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1141eb284becSJed Brown 
1142eb284becSJed Brown   PetscFunctionBegin;
1143eb284becSJed Brown   th->endpoint = flg;
11443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1145eb284becSJed Brown }
11460de4c49aSJed Brown 
1147f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1148d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
1149d71ae5a4SJacob Faibussowitsch {
1150f9c1d6abSBarry Smith   PetscComplex    z   = xr + xi * PETSC_i, f;
1151f9c1d6abSBarry Smith   TS_Theta       *th  = (TS_Theta *)ts->data;
11523fd8ae06SJed Brown   const PetscReal one = 1.0;
1153f9c1d6abSBarry Smith 
1154f9c1d6abSBarry Smith   PetscFunctionBegin;
11553fd8ae06SJed Brown   f   = (one + (one - th->Theta) * z) / (one - th->Theta * z);
1156f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1157f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
11583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1159f9c1d6abSBarry Smith }
1160f9c1d6abSBarry Smith #endif
1161f9c1d6abSBarry Smith 
1162d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[])
1163d71ae5a4SJacob Faibussowitsch {
116442682096SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
116542682096SHong Zhang 
116642682096SHong Zhang   PetscFunctionBegin;
11677409eceaSHong Zhang   if (ns) {
11687409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11697409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
11707409eceaSHong Zhang   }
11715072d23cSHong Zhang   if (Y) {
1172cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
11737409eceaSHong Zhang       th->Stages[0] = th->X;
1174cc4f23bcSHong Zhang     } else {
1175cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1176cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1177cc4f23bcSHong Zhang     }
1178cc4f23bcSHong Zhang     *Y = th->Stages;
11795072d23cSHong Zhang   }
11803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118142682096SHong Zhang }
1182f9c1d6abSBarry Smith 
1183316643e7SJed Brown /* ------------------------------------------------------------ */
1184316643e7SJed Brown /*MC
118596f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1186316643e7SJed Brown 
1187316643e7SJed Brown    Level: beginner
1188316643e7SJed Brown 
1189bcf0153eSBarry Smith    Options Database Keys:
1190c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1191c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
119203842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11934eb428fdSBarry Smith 
1194eb284becSJed Brown    Notes:
119537fdd005SBarry Smith .vb
119637fdd005SBarry Smith   -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
119737fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
119837fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
119937fdd005SBarry Smith .ve
12004eb428fdSBarry Smith 
12017409eceaSHong Zhang    The endpoint variant of the Theta method and backward Euler can be applied to DAE. The midpoint variant is not suitable for DAEs because it is not stiffly accurate.
1202eb284becSJed Brown 
12037409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1204eb284becSJed Brown 
1205eb284becSJed Brown .vb
1206eb284becSJed Brown   Theta | Theta
1207eb284becSJed Brown   -------------
1208eb284becSJed Brown         |  1
1209eb284becSJed Brown .ve
1210eb284becSJed Brown 
1211eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1212eb284becSJed Brown 
1213eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1214eb284becSJed Brown 
1215eb284becSJed Brown .vb
1216eb284becSJed Brown   0 | 0         0
1217eb284becSJed Brown   1 | 1-Theta   Theta
1218eb284becSJed Brown   -------------------
1219eb284becSJed Brown     | 1-Theta   Theta
1220eb284becSJed Brown .ve
1221eb284becSJed Brown 
1222eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1223eb284becSJed Brown 
1224eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1225eb284becSJed Brown 
1226eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1227eb284becSJed Brown 
12284eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1229eb284becSJed Brown 
12301cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1231316643e7SJed Brown M*/
1232d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1233d71ae5a4SJacob Faibussowitsch {
1234316643e7SJed Brown   TS_Theta *th;
1235316643e7SJed Brown 
1236316643e7SJed Brown   PetscFunctionBegin;
1237277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
1238ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1239316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
1240316643e7SJed Brown   ts->ops->view           = TSView_Theta;
1241316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
124242f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1243ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1244316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
1245cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
12461566a47fSLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
124724655328SShri   ts->ops->rollback       = TSRollBack_Theta;
1248*26b5f05dSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Theta;
1249316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
12500f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
12510f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1252f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1253f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1254f9c1d6abSBarry Smith #endif
125542682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
125642f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1257b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1258b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12592ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12606affb6f8SHong Zhang 
1261715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12627adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1263715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12646affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1265316643e7SJed Brown 
1266efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1267efd4aadfSBarry Smith 
12684dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
1269316643e7SJed Brown   ts->data = (void *)th;
1270316643e7SJed Brown 
1271715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
1272715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
1273715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
1274b5b4017aSHong Zhang   th->VecsSensi2Temp = NULL;
1275715f1b00SHong Zhang 
12766f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1277316643e7SJed Brown   th->Theta       = 0.5;
12781566a47fSLisandro Dalcin   th->order       = 2;
12799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
12809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
12819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
12829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
12833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1284316643e7SJed Brown }
12850de4c49aSJed Brown 
12860de4c49aSJed Brown /*@
1287bcf0153eSBarry Smith   TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA`
12880de4c49aSJed Brown 
12890de4c49aSJed Brown   Not Collective
12900de4c49aSJed Brown 
12910de4c49aSJed Brown   Input Parameter:
12920de4c49aSJed Brown . ts - timestepping context
12930de4c49aSJed Brown 
12940de4c49aSJed Brown   Output Parameter:
12950de4c49aSJed Brown . theta - stage abscissa
12960de4c49aSJed Brown 
1297bcf0153eSBarry Smith   Level: advanced
1298bcf0153eSBarry Smith 
12990de4c49aSJed Brown   Note:
1300bcf0153eSBarry Smith   Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme.
13010de4c49aSJed Brown 
13021cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA`
13030de4c49aSJed Brown @*/
1304d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta)
1305d71ae5a4SJacob Faibussowitsch {
13060de4c49aSJed Brown   PetscFunctionBegin;
1307afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1308dadcf809SJacob Faibussowitsch   PetscValidRealPointer(theta, 2);
1309cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
13103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13110de4c49aSJed Brown }
13120de4c49aSJed Brown 
13130de4c49aSJed Brown /*@
1314bcf0153eSBarry Smith   TSThetaSetTheta - Set the abscissa of the stage in (0,1]  for `TSTHETA`
13150de4c49aSJed Brown 
13160de4c49aSJed Brown   Not Collective
13170de4c49aSJed Brown 
1318d8d19677SJose E. Roman   Input Parameters:
13190de4c49aSJed Brown + ts    - timestepping context
13200de4c49aSJed Brown - theta - stage abscissa
13210de4c49aSJed Brown 
1322bcf0153eSBarry Smith   Options Database Key:
132367b8a455SSatish Balay . -ts_theta_theta <theta> - set theta
13240de4c49aSJed Brown 
1325bcf0153eSBarry Smith   Level: intermediate
13260de4c49aSJed Brown 
13271cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN`
13280de4c49aSJed Brown @*/
1329d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta)
1330d71ae5a4SJacob Faibussowitsch {
13310de4c49aSJed Brown   PetscFunctionBegin;
1332afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1333cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
13343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13350de4c49aSJed Brown }
1336f33bbcb6SJed Brown 
133726f2ff8fSLisandro Dalcin /*@
1338bcf0153eSBarry Smith   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
133926f2ff8fSLisandro Dalcin 
134026f2ff8fSLisandro Dalcin   Not Collective
134126f2ff8fSLisandro Dalcin 
134226f2ff8fSLisandro Dalcin   Input Parameter:
134326f2ff8fSLisandro Dalcin . ts - timestepping context
134426f2ff8fSLisandro Dalcin 
134526f2ff8fSLisandro Dalcin   Output Parameter:
1346bcf0153eSBarry Smith . endpoint - `PETSC_TRUE` when using the endpoint variant
134726f2ff8fSLisandro Dalcin 
1348bcf0153eSBarry Smith   Level: advanced
134926f2ff8fSLisandro Dalcin 
13501cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
135126f2ff8fSLisandro Dalcin @*/
1352d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint)
1353d71ae5a4SJacob Faibussowitsch {
135426f2ff8fSLisandro Dalcin   PetscFunctionBegin;
135526f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1356dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(endpoint, 2);
1357cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
13583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135926f2ff8fSLisandro Dalcin }
136026f2ff8fSLisandro Dalcin 
1361eb284becSJed Brown /*@
1362bcf0153eSBarry Smith   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
1363eb284becSJed Brown 
1364eb284becSJed Brown   Not Collective
1365eb284becSJed Brown 
1366d8d19677SJose E. Roman   Input Parameters:
1367eb284becSJed Brown + ts  - timestepping context
1368bcf0153eSBarry Smith - flg - `PETSC_TRUE` to use the endpoint variant
1369eb284becSJed Brown 
1370bcf0153eSBarry Smith   Options Database Key:
137167b8a455SSatish Balay . -ts_theta_endpoint <flg> - use the endpoint variant
1372eb284becSJed Brown 
1373bcf0153eSBarry Smith   Level: intermediate
1374eb284becSJed Brown 
13751cc06b55SBarry Smith .seealso: [](ch_ts), `TSTHETA`, `TSCN`
1376eb284becSJed Brown @*/
1377d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg)
1378d71ae5a4SJacob Faibussowitsch {
1379eb284becSJed Brown   PetscFunctionBegin;
1380eb284becSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1381cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
13823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1383eb284becSJed Brown }
1384eb284becSJed Brown 
1385f33bbcb6SJed Brown /*
1386f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1387f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1388f33bbcb6SJed Brown  */
1389f33bbcb6SJed Brown 
1390d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BEuler(TS ts)
1391d71ae5a4SJacob Faibussowitsch {
13921566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
13931566a47fSLisandro Dalcin 
13941566a47fSLisandro Dalcin   PetscFunctionBegin;
139508401ef6SPierre Jolivet   PetscCheck(th->Theta == 1.0, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (1) of theta when using backward Euler");
139628b400f6SJacob Faibussowitsch   PetscCheck(!th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the endpoint form of the Theta methods when using backward Euler");
13979566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
13983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13991566a47fSLisandro Dalcin }
14001566a47fSLisandro Dalcin 
1401d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer)
1402d71ae5a4SJacob Faibussowitsch {
1403f33bbcb6SJed Brown   PetscFunctionBegin;
14043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1405f33bbcb6SJed Brown }
1406f33bbcb6SJed Brown 
1407f33bbcb6SJed Brown /*MC
1408f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1409f33bbcb6SJed Brown 
1410f33bbcb6SJed Brown   Level: beginner
1411f33bbcb6SJed Brown 
1412bcf0153eSBarry Smith   Note:
141337fdd005SBarry Smith   `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0`
14144eb428fdSBarry Smith 
14151cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1416f33bbcb6SJed Brown M*/
1417d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1418d71ae5a4SJacob Faibussowitsch {
1419f33bbcb6SJed Brown   PetscFunctionBegin;
14209566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14219566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 1.0));
14229566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
14230c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1424f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
14253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1426f33bbcb6SJed Brown }
1427f33bbcb6SJed Brown 
1428d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_CN(TS ts)
1429d71ae5a4SJacob Faibussowitsch {
14301566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
14311566a47fSLisandro Dalcin 
14321566a47fSLisandro Dalcin   PetscFunctionBegin;
143308401ef6SPierre Jolivet   PetscCheck(th->Theta == 0.5, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (0.5) of theta when using Crank-Nicolson");
14343c633725SBarry Smith   PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the midpoint form of the Theta methods when using Crank-Nicolson");
14359566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14371566a47fSLisandro Dalcin }
14381566a47fSLisandro Dalcin 
1439d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer)
1440d71ae5a4SJacob Faibussowitsch {
1441f33bbcb6SJed Brown   PetscFunctionBegin;
14423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1443f33bbcb6SJed Brown }
1444f33bbcb6SJed Brown 
1445f33bbcb6SJed Brown /*MC
1446f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1447f33bbcb6SJed Brown 
1448f33bbcb6SJed Brown   Level: beginner
1449f33bbcb6SJed Brown 
1450f33bbcb6SJed Brown   Notes:
1451bcf0153eSBarry Smith   `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e.
1452bcf0153eSBarry Smith .vb
1453bcf0153eSBarry Smith   -ts_type theta
1454bcf0153eSBarry Smith   -ts_theta_theta 0.5
1455bcf0153eSBarry Smith   -ts_theta_endpoint
1456bcf0153eSBarry Smith .ve
14577cf5af47SJed Brown 
14581cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`,
1459f33bbcb6SJed Brown M*/
1460d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1461d71ae5a4SJacob Faibussowitsch {
1462f33bbcb6SJed Brown   PetscFunctionBegin;
14639566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14649566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 0.5));
14659566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
14660c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1467f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
14683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1469f33bbcb6SJed Brown }
1470