xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 82d7ce88f5d4a1951e85e8b7f7dc3fe3cb2854f2)
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 
18226b5f05dSStefano Zampini /* We need to transfer X0 which will be copied into sol_prev */
18326b5f05dSStefano Zampini static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg)
18426b5f05dSStefano Zampini {
18526b5f05dSStefano Zampini   TS_Theta  *th     = (TS_Theta *)ts->data;
18626b5f05dSStefano Zampini   const char name[] = "ts:theta:X0";
18726b5f05dSStefano Zampini 
18826b5f05dSStefano Zampini   PetscFunctionBegin;
18926b5f05dSStefano Zampini   if (reg && th->vec_sol_prev) {
19026b5f05dSStefano Zampini     PetscCall(TSResizeRegisterVec(ts, name, th->X0));
19126b5f05dSStefano Zampini   } else if (!reg) {
19226b5f05dSStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, name, &th->X0));
19326b5f05dSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->X0));
19426b5f05dSStefano Zampini   }
19526b5f05dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
19626b5f05dSStefano Zampini }
19726b5f05dSStefano 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));
222eb284becSJed Brown     }
2239566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
224*82d7ce88SStefano Zampini     PetscCall(TSTheta_SNESSolve(ts, th->endpoint ? th->affine : NULL, th->X));
2259566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
2269566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
227fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
228051f2191SLisandro Dalcin 
229051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2301566a47fSLisandro Dalcin     if (th->endpoint) {
2319566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X, ts->vec_sol));
2321566a47fSLisandro Dalcin     } else {
2339566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
2349566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
2351566a47fSLisandro Dalcin     }
2369566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2371566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2381566a47fSLisandro Dalcin     if (!accept) {
2399566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
240be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
241be5899b3SLisandro Dalcin       goto reject_step;
242051f2191SLisandro Dalcin     }
2433b1890cdSShri Abhyankar 
244715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2451a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2461a352fa8SHong Zhang       th->time_step0 = ts->time_step;
24717145e75SHong Zhang     }
2482b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
249cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
250051f2191SLisandro Dalcin     break;
251051f2191SLisandro Dalcin 
252051f2191SLisandro Dalcin   reject_step:
2539371c9d4SSatish Balay     ts->reject++;
2549371c9d4SSatish Balay     accept = PETSC_FALSE;
2551566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
256051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
25763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
2583b1890cdSShri Abhyankar     }
2593b1890cdSShri Abhyankar   }
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
261316643e7SJed Brown }
262316643e7SJed Brown 
263d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
264d71ae5a4SJacob Faibussowitsch {
265c9681e5dSHong Zhang   TS_Theta      *th           = (TS_Theta *)ts->data;
266cd4cee2dSHong Zhang   TS             quadts       = ts->quadraturets;
267c9681e5dSHong Zhang   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
268c9681e5dSHong Zhang   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
269c9681e5dSHong Zhang   PetscInt       nadj;
2707207e0fdSHong Zhang   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
271c9681e5dSHong Zhang   KSP            ksp;
272c9681e5dSHong Zhang   PetscScalar   *xarr;
273077c4feaSHong Zhang   TSEquationType eqtype;
274077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2751a352fa8SHong Zhang   PetscReal      adjoint_time_step;
276c9681e5dSHong Zhang 
277c9681e5dSHong Zhang   PetscFunctionBegin;
2789566063dSJacob Faibussowitsch   PetscCall(TSGetEquationType(ts, &eqtype));
279077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
280077c4feaSHong Zhang     isexplicitode = PETSC_TRUE;
281077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
282077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
283077c4feaSHong Zhang   }
284c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
2859566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
2869566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
2877207e0fdSHong Zhang   if (quadts) {
2889566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
2899566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
2907207e0fdSHong Zhang   }
291c9681e5dSHong Zhang 
2921a352fa8SHong Zhang   th->stage_time    = ts->ptime;
2931a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
294c9681e5dSHong Zhang 
29533f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2961baa6e33SBarry Smith   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
297cd4cee2dSHong Zhang 
298c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
2999566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
3009566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
301cd4cee2dSHong Zhang     if (quadJ) {
3029566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
3039566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
3049566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
3059566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
3069566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
307c9681e5dSHong Zhang     }
308c9681e5dSHong Zhang   }
309c9681e5dSHong Zhang 
310c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
311c87ba875SIvan Yashchuk   th->shift = 1. / adjoint_time_step;
3129566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3139566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
314c9681e5dSHong Zhang 
315c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
316c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
317c9681e5dSHong Zhang     KSPConvergedReason kspreason;
3189566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
3199566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
320c9681e5dSHong Zhang     if (kspreason < 0) {
321c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32263a3b9bcSJacob 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));
323c9681e5dSHong Zhang     }
324c9681e5dSHong Zhang   }
325c9681e5dSHong Zhang 
326c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
327c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3289566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
3299566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
330c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3319566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
332c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3339566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
334c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
3359566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
3369566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
3379566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
33848a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
339c9681e5dSHong Zhang     }
340c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
341c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
34205c9950eSHong Zhang       KSPConvergedReason kspreason;
3439566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
3449566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
34505c9950eSHong Zhang       if (kspreason < 0) {
34605c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
34763a3b9bcSJacob 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));
34805c9950eSHong Zhang       }
349c9681e5dSHong Zhang     }
350c9681e5dSHong Zhang   }
351c9681e5dSHong Zhang 
352c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
353077c4feaSHong Zhang   if (!isexplicitode) {
3543e05ceb1SHong Zhang     th->shift = 0.0;
3559566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3569566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
357c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
35833f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3599566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
3609566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
3619566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
3629566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
363c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3649566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
3659566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
3669566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
367c9681e5dSHong Zhang       }
368c9681e5dSHong Zhang     }
369077c4feaSHong Zhang   }
370c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3719566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
3729566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
3731baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
374c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
375c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3769566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
37733f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3789566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
379c9681e5dSHong Zhang     }
380cd4cee2dSHong Zhang 
381c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
3829566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
3839566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
384cd4cee2dSHong Zhang       if (quadJp) {
3859566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
3869566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
3879566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
3889566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3899566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
39033f72c5dSHong Zhang       }
391c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
3929566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
3939566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
39448a46eb9SPierre Jolivet         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
39548a46eb9SPierre Jolivet         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
396c9681e5dSHong Zhang       }
397c9681e5dSHong Zhang     }
398c9681e5dSHong Zhang   }
399c9681e5dSHong Zhang 
400c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
4019566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4029566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
403c9681e5dSHong Zhang   }
404c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
406c9681e5dSHong Zhang }
407c9681e5dSHong Zhang 
408d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_Theta(TS ts)
409d71ae5a4SJacob Faibussowitsch {
4102ca6e920SHong Zhang   TS_Theta    *th           = (TS_Theta *)ts->data;
411cd4cee2dSHong Zhang   TS           quadts       = ts->quadraturets;
412b5e0532dSHong Zhang   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
413b5b4017aSHong Zhang   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
4142ca6e920SHong Zhang   PetscInt     nadj;
4157207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
4162ca6e920SHong Zhang   KSP          ksp;
417b5b4017aSHong Zhang   PetscScalar *xarr;
4181a352fa8SHong Zhang   PetscReal    adjoint_time_step;
419da81f932SPierre 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) */
4202ca6e920SHong Zhang 
4212ca6e920SHong Zhang   PetscFunctionBegin;
422c9681e5dSHong Zhang   if (th->Theta == 1.) {
4239566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
4243ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
425c9681e5dSHong Zhang   }
4262ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4279566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
4289566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
4297207e0fdSHong Zhang   if (quadts) {
4309566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
4319566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
4327207e0fdSHong Zhang   }
433b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4341a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
4351a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4361a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
4372ca6e920SHong Zhang 
43882ad9101SHong Zhang   if (!th->endpoint) {
4395072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4409566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
44182ad9101SHong Zhang   }
44282ad9101SHong Zhang 
443b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
44433f72c5dSHong Zhang   /* Cost function has an integral term */
4457207e0fdSHong Zhang   if (quadts) {
44605755b9cSHong Zhang     if (th->endpoint) {
4479566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
44805755b9cSHong Zhang     } else {
4499566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
45005755b9cSHong Zhang     }
4517207e0fdSHong Zhang   }
45233f72c5dSHong Zhang 
453abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4549566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
4559566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
456cd4cee2dSHong Zhang     if (quadJ) {
4579566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
4589566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
4599566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
4609566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4619566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
46236eaed60SHong Zhang     }
4632ca6e920SHong Zhang   }
4643c54f92cSHong Zhang 
465b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4661a352fa8SHong Zhang   th->shift = 1. / (th->Theta * adjoint_time_step);
4673c54f92cSHong Zhang   if (th->endpoint) {
4689566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
4693c54f92cSHong Zhang   } else {
4709566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
4713c54f92cSHong Zhang   }
4729566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
4732ca6e920SHong Zhang 
474b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
475abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4761d14d03bSHong Zhang     KSPConvergedReason kspreason;
4779566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
4789566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
4791d14d03bSHong Zhang     if (kspreason < 0) {
4801d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
48163a3b9bcSJacob 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));
4821d14d03bSHong Zhang     }
4832ca6e920SHong Zhang   }
4843c54f92cSHong Zhang 
4851d14d03bSHong Zhang   /* Second-order adjoint */
486b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4873c633725SBarry Smith     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
488b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4899566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
4909566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
491b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
4929566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
4939566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4949566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
495b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
4969566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
497b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
4989566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
4999566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
5009566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
50148a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
502b5b4017aSHong Zhang     }
503b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
504b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
5051d14d03bSHong Zhang       KSPConvergedReason kspreason;
5069566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
5079566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
5081d14d03bSHong Zhang       if (kspreason < 0) {
5091d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
51063a3b9bcSJacob 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));
5111d14d03bSHong Zhang       }
512b5b4017aSHong Zhang     }
513b5b4017aSHong Zhang   }
514b5b4017aSHong Zhang 
51536eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5161d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5171a352fa8SHong Zhang     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
5181a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
5199566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
5209566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
52133f72c5dSHong Zhang     /* R_U at t_n */
5221baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
523abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
5249566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
525cd4cee2dSHong Zhang       if (quadJ) {
5269566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
5279566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
5289566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
5299566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5309566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
531b5b4017aSHong Zhang       }
5329566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
533b5b4017aSHong Zhang     }
5341d14d03bSHong Zhang 
5351d14d03bSHong Zhang     /* Second-order adjoint */
5361d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
537b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5389566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
5399566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
540b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5419566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
5429566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5439566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
544b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5459566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
546b5b4017aSHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
54733f72c5dSHong 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  */
5489566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
5499566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
55048a46eb9SPierre Jolivet         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
5519566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
552b5b4017aSHong Zhang       }
553b5e0532dSHong Zhang     }
5543fd52205SHong Zhang 
555c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
556c586f2e8SHong Zhang 
5573fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
558b5b4017aSHong Zhang       /* U_{n+1} */
55982ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
5609566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
5619566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5621baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
563abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5649566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
5659566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
5663b711c3fSHong Zhang         if (quadJp) {
5679566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
5689566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
5699566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
5709566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5719566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
5723b711c3fSHong Zhang         }
5733fd52205SHong Zhang       }
57433f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
57533f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5769566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
5779566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
578b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5799566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
5809566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5819566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
58233f72c5dSHong Zhang 
583b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5849566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
585b5b4017aSHong Zhang         for (nadj = 0; nadj < ts->numcost; nadj++) {
58633f72c5dSHong 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)  */
5879566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
5889566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
58948a46eb9SPierre Jolivet           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
59048a46eb9SPierre Jolivet           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
591b5b4017aSHong Zhang         }
592b5b4017aSHong Zhang       }
593b5b4017aSHong Zhang 
594b5b4017aSHong Zhang       /* U_s */
5959566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
5969566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5971baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
598abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5999566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6009566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
6013b711c3fSHong Zhang         if (quadJp) {
6029566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6039566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6049566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
6059566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6069566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
6073b711c3fSHong Zhang         }
60833f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
60933f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
6109566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
6119566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
612b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6139566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
6149566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
6159566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
616b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6179566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
618b5b4017aSHong Zhang           for (nadj = 0; nadj < ts->numcost; nadj++) {
61933f72c5dSHong 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) */
6209566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
6219566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
62248a46eb9SPierre Jolivet             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
62348a46eb9SPierre Jolivet             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
62436eaed60SHong Zhang           }
62536eaed60SHong Zhang         }
6263fd52205SHong Zhang       }
627b5e0532dSHong Zhang     }
6283fd52205SHong Zhang   } else { /* one-stage case */
6293e05ceb1SHong Zhang     th->shift = 0.0;
6309566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
6319566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
6321baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
633abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
6349566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
6359566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
636cd4cee2dSHong Zhang       if (quadJ) {
6379566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
6389566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
6399566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
6409566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6419566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
64236eaed60SHong Zhang       }
6432ca6e920SHong Zhang     }
6443fd52205SHong Zhang     if (ts->vecs_sensip) {
64582ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
6469566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
6479566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
6481baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
649abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6509566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6519566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
652cd4cee2dSHong Zhang         if (quadJp) {
6539566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6549566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6559566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
6569566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6579566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
65836eaed60SHong Zhang         }
65936eaed60SHong Zhang       }
6603fd52205SHong Zhang     }
6612ca6e920SHong Zhang   }
6622ca6e920SHong Zhang 
6632ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6652ca6e920SHong Zhang }
6662ca6e920SHong Zhang 
667d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X)
668d71ae5a4SJacob Faibussowitsch {
669cd652676SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
670be5899b3SLisandro Dalcin   PetscReal dt = t - ts->ptime;
671cd652676SJed Brown 
672cd652676SJed Brown   PetscFunctionBegin;
6739566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X));
674be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6759566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
6763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
677cd652676SJed Brown }
678cd652676SJed Brown 
679d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
680d71ae5a4SJacob Faibussowitsch {
6811566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
6821566a47fSLisandro Dalcin   Vec       X  = ts->vec_sol;      /* X = solution */
6831566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
6847453f775SEmil Constantinescu   PetscReal wltea, wlter;
6851566a47fSLisandro Dalcin 
6861566a47fSLisandro Dalcin   PetscFunctionBegin;
6879371c9d4SSatish Balay   if (!th->vec_sol_prev) {
6889371c9d4SSatish Balay     *wlte = -1;
6893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6909371c9d4SSatish Balay   }
6911566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
6929371c9d4SSatish Balay   if (ts->steprestart) {
6939371c9d4SSatish Balay     *wlte = -1;
6943ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
6959371c9d4SSatish Balay   }
6961566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
697fecfb714SLisandro Dalcin   {
698be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
699be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
7009371c9d4SSatish Balay     PetscScalar scal[3];
7019371c9d4SSatish Balay     Vec         vecs[3];
7029371c9d4SSatish Balay     scal[0] = +1 / a;
7039371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
7049371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
7059371c9d4SSatish Balay     vecs[0] = X;
7069371c9d4SSatish Balay     vecs[1] = th->X0;
7079371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
7089566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
7099566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
7109566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
7111566a47fSLisandro Dalcin   }
7121566a47fSLisandro Dalcin   if (order) *order = 2;
7133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7141566a47fSLisandro Dalcin }
7151566a47fSLisandro Dalcin 
716d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Theta(TS ts)
717d71ae5a4SJacob Faibussowitsch {
7181566a47fSLisandro Dalcin   TS_Theta *th     = (TS_Theta *)ts->data;
7197207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
7201566a47fSLisandro Dalcin 
7211566a47fSLisandro Dalcin   PetscFunctionBegin;
7229566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
72348a46eb9SPierre Jolivet   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
724715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
7251baa6e33SBarry Smith   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
72648a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
728715f1b00SHong Zhang }
729715f1b00SHong Zhang 
730d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_Theta(TS ts)
731d71ae5a4SJacob Faibussowitsch {
732715f1b00SHong Zhang   TS_Theta    *th                   = (TS_Theta *)ts->data;
733cd4cee2dSHong Zhang   TS           quadts               = ts->quadraturets;
73413b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
735b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7367207e0fdSHong Zhang   PetscInt     ntlm;
737715f1b00SHong Zhang   KSP          ksp;
7387207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
73913b1b0ffSHong Zhang   PetscScalar *barr, *xarr;
7401a352fa8SHong Zhang   PetscReal    previous_shift;
741715f1b00SHong Zhang 
742715f1b00SHong Zhang   PetscFunctionBegin;
7431a352fa8SHong Zhang   previous_shift = th->shift;
7449566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));
74513b1b0ffSHong Zhang 
74648a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
7479566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
7489566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
7497207e0fdSHong Zhang   if (quadts) {
7509566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
7519566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
7527207e0fdSHong Zhang   }
753715f1b00SHong Zhang 
754715f1b00SHong Zhang   /* Build RHS */
755715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7561a352fa8SHong Zhang     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
7579566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7589566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7599566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));
760715f1b00SHong Zhang 
761022c081aSHong Zhang     /* Add the f_p forcing terms */
7620e11c21fSHong Zhang     if (ts->Jacp) {
7639566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7649566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7659566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
76682ad9101SHong Zhang       th->shift = previous_shift;
7679566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
7689566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7699566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
7700e11c21fSHong Zhang     }
771715f1b00SHong Zhang   } else { /* 1-stage method */
7723e05ceb1SHong Zhang     th->shift = 0.0;
7739566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7749566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7759566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, -1.));
77613b1b0ffSHong Zhang 
777022c081aSHong Zhang     /* Add the f_p forcing terms */
7780e11c21fSHong Zhang     if (ts->Jacp) {
77982ad9101SHong Zhang       th->shift = previous_shift;
7809566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
7819566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7829566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
783715f1b00SHong Zhang     }
7840e11c21fSHong Zhang   }
785715f1b00SHong Zhang 
786715f1b00SHong Zhang   /* Build LHS */
7871a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
788715f1b00SHong Zhang   if (th->endpoint) {
7899566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
790715f1b00SHong Zhang   } else {
7919566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
792715f1b00SHong Zhang   }
7939566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
794715f1b00SHong Zhang 
795715f1b00SHong Zhang   /*
796715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
797c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
798715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
799715f1b00SHong Zhang   */
800715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8017207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8029566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
8039566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
8049566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8059566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8069566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
807715f1b00SHong Zhang     }
808715f1b00SHong Zhang   }
809715f1b00SHong Zhang 
810715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
811022c081aSHong Zhang   for (ntlm = 0; ntlm < th->num_tlm; ntlm++) {
812b5b4017aSHong Zhang     KSPConvergedReason kspreason;
8139566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
8149566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
815715f1b00SHong Zhang     if (th->endpoint) {
8169566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
8179566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
8189566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
8199566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
8209566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
821715f1b00SHong Zhang     } else {
8229566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
823715f1b00SHong Zhang     }
8249566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
825b5b4017aSHong Zhang     if (kspreason < 0) {
826b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
82763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
828b5b4017aSHong Zhang     }
8299566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
8309566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
831715f1b00SHong Zhang   }
832715f1b00SHong Zhang 
83313b1b0ffSHong Zhang   /*
83413b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
835c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
83613b1b0ffSHong Zhang   */
8377207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
83813b1b0ffSHong Zhang     if (!th->endpoint) {
8399566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
8409566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
8419566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
8429566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8439566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8449566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
8459566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
84613b1b0ffSHong Zhang     } else {
8479566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
8489566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
8499566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8509566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8519566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
85213b1b0ffSHong Zhang     }
85313b1b0ffSHong Zhang   } else {
85448a46eb9SPierre Jolivet     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
855715f1b00SHong Zhang   }
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8571566a47fSLisandro Dalcin }
8581566a47fSLisandro Dalcin 
859d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[])
860d71ae5a4SJacob Faibussowitsch {
8616affb6f8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
8626affb6f8SHong Zhang 
8636affb6f8SHong Zhang   PetscFunctionBegin;
8647409eceaSHong Zhang   if (ns) {
8657409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8667409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
8677409eceaSHong Zhang   }
8685072d23cSHong Zhang   if (stagesensip) {
869cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8707409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
871cc4f23bcSHong Zhang     } else {
872cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
873cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
874cc4f23bcSHong Zhang     }
875cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8765072d23cSHong Zhang   }
8773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8786affb6f8SHong Zhang }
8796affb6f8SHong Zhang 
880316643e7SJed Brown /*------------------------------------------------------------*/
881d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Theta(TS ts)
882d71ae5a4SJacob Faibussowitsch {
883316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
884316643e7SJed Brown 
885316643e7SJed Brown   PetscFunctionBegin;
8869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
8879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
8889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
8899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
8901566a47fSLisandro Dalcin 
8919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
8929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
8931566a47fSLisandro Dalcin 
8949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
896277b19d0SLisandro Dalcin }
897277b19d0SLisandro Dalcin 
898d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_Theta(TS ts)
899d71ae5a4SJacob Faibussowitsch {
900ece46509SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
901ece46509SHong Zhang 
902ece46509SHong Zhang   PetscFunctionBegin;
9039566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
9049566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
9059566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
9069566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
9079566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
9089566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
9093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
910ece46509SHong Zhang }
911ece46509SHong Zhang 
912d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Theta(TS ts)
913d71ae5a4SJacob Faibussowitsch {
914277b19d0SLisandro Dalcin   PetscFunctionBegin;
9159566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
916b3a6b972SJed Brown   if (ts->dm) {
9179566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
9189566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
919b3a6b972SJed Brown   }
9209566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
9229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
9239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
9249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
926316643e7SJed Brown }
927316643e7SJed Brown 
928316643e7SJed Brown /*
929316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9302b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9310fd10508SMatthew Knepley 
9320fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9330fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9340fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
935316643e7SJed Brown */
936d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts)
937d71ae5a4SJacob Faibussowitsch {
938316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9397445fe48SJed Brown   Vec       X0, Xdot;
9407445fe48SJed Brown   DM        dm, dmsave;
9413e05ceb1SHong Zhang   PetscReal shift = th->shift;
942316643e7SJed Brown 
943316643e7SJed Brown   PetscFunctionBegin;
9449566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
9455a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9469566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
947494ce9fcSHong Zhang   if (x != X0) {
9489566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
949494ce9fcSHong Zhang   } else {
9509566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
951494ce9fcSHong Zhang   }
9527445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9537445fe48SJed Brown   dmsave = ts->dm;
9547445fe48SJed Brown   ts->dm = dm;
9559566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
9567445fe48SJed Brown   ts->dm = dmsave;
9579566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
9583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
959316643e7SJed Brown }
960316643e7SJed Brown 
961d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts)
962d71ae5a4SJacob Faibussowitsch {
963316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9647445fe48SJed Brown   Vec       Xdot;
9657445fe48SJed Brown   DM        dm, dmsave;
9663e05ceb1SHong Zhang   PetscReal shift = th->shift;
967316643e7SJed Brown 
968316643e7SJed Brown   PetscFunctionBegin;
9699566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
970be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9719566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));
9727445fe48SJed Brown 
9737445fe48SJed Brown   dmsave = ts->dm;
9747445fe48SJed Brown   ts->dm = dm;
9759566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
9767445fe48SJed Brown   ts->dm = dmsave;
9779566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
9783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
979316643e7SJed Brown }
980316643e7SJed Brown 
981d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_Theta(TS ts)
982d71ae5a4SJacob Faibussowitsch {
983715f1b00SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
9847207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
985715f1b00SHong Zhang 
986715f1b00SHong Zhang   PetscFunctionBegin;
987022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
98813b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
9899566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
9907207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9919566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
9929566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
993715f1b00SHong Zhang   }
9946f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
9959566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));
99613b1b0ffSHong Zhang 
9979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
999715f1b00SHong Zhang }
1000715f1b00SHong Zhang 
1001d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_Theta(TS ts)
1002d71ae5a4SJacob Faibussowitsch {
10037adebddeSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
10047207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
10057adebddeSHong Zhang 
10067adebddeSHong Zhang   PetscFunctionBegin;
10077207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10089566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
10099566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
10107adebddeSHong Zhang   }
10119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
10129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
10139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
10143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10157adebddeSHong Zhang }
10167adebddeSHong Zhang 
1017d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Theta(TS ts)
1018d71ae5a4SJacob Faibussowitsch {
1019316643e7SJed Brown   TS_Theta *th     = (TS_Theta *)ts->data;
1020cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
10212ffb9264SLisandro Dalcin   PetscBool match;
1022316643e7SJed Brown 
1023316643e7SJed Brown   PetscFunctionBegin;
1024cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
10259566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
1026b1cb36f3SHong Zhang   }
102748a46eb9SPierre Jolivet   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
102848a46eb9SPierre Jolivet   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
102948a46eb9SPierre Jolivet   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
103048a46eb9SPierre Jolivet   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
10313b1890cdSShri Abhyankar 
10321566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
103320d49056SMatthew G. Knepley   th->shift = 1 / (th->Theta * ts->time_step);
10341566a47fSLisandro Dalcin 
10359566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &ts->dm));
10369566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
10379566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
10381566a47fSLisandro Dalcin 
10399566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
10409566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
10422ffb9264SLisandro Dalcin   if (!match) {
10439566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
10449566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
10453b1890cdSShri Abhyankar   }
10469566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
1047cc4f23bcSHong Zhang 
1048cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
10493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1050b4dd3bf9SBarry Smith }
10510c7376e5SHong Zhang 
1052b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1053b4dd3bf9SBarry Smith 
1054d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1055d71ae5a4SJacob Faibussowitsch {
1056b4dd3bf9SBarry Smith   TS_Theta *th = (TS_Theta *)ts->data;
1057b4dd3bf9SBarry Smith 
1058b4dd3bf9SBarry Smith   PetscFunctionBegin;
10599566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
10609566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
106148a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1062b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10639566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
10649566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
106567633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
106667633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
106767633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1068b5b4017aSHong Zhang   }
1069c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
10709566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
107167633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
107267633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
107367633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1074b5b4017aSHong Zhang   }
10753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1076316643e7SJed Brown }
1077316643e7SJed Brown 
1078d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject)
1079d71ae5a4SJacob Faibussowitsch {
1080316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1081316643e7SJed Brown 
1082316643e7SJed Brown   PetscFunctionBegin;
1083d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1084316643e7SJed Brown   {
10859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
10869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
10879566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1088316643e7SJed Brown   }
1089d0609cedSBarry Smith   PetscOptionsHeadEnd();
10903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1091316643e7SJed Brown }
1092316643e7SJed Brown 
1093d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer)
1094d71ae5a4SJacob Faibussowitsch {
1095316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1096ace3abfcSBarry Smith   PetscBool iascii;
1097316643e7SJed Brown 
1098316643e7SJed Brown   PetscFunctionBegin;
10999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1100316643e7SJed Brown   if (iascii) {
11019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
11029566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1103316643e7SJed Brown   }
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1105316643e7SJed Brown }
1106316643e7SJed Brown 
1107d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta)
1108d71ae5a4SJacob Faibussowitsch {
11090de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11100de4c49aSJed Brown 
11110de4c49aSJed Brown   PetscFunctionBegin;
11120de4c49aSJed Brown   *theta = th->Theta;
11133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11140de4c49aSJed Brown }
11150de4c49aSJed Brown 
1116d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta)
1117d71ae5a4SJacob Faibussowitsch {
11180de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11190de4c49aSJed Brown 
11200de4c49aSJed Brown   PetscFunctionBegin;
1121cad9d221SBarry Smith   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
11220de4c49aSJed Brown   th->Theta = theta;
11231566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11250de4c49aSJed Brown }
1126eb284becSJed Brown 
1127d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint)
1128d71ae5a4SJacob Faibussowitsch {
112926f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
113026f2ff8fSLisandro Dalcin 
113126f2ff8fSLisandro Dalcin   PetscFunctionBegin;
113226f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113426f2ff8fSLisandro Dalcin }
113526f2ff8fSLisandro Dalcin 
1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg)
1137d71ae5a4SJacob Faibussowitsch {
1138eb284becSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1139eb284becSJed Brown 
1140eb284becSJed Brown   PetscFunctionBegin;
1141eb284becSJed Brown   th->endpoint = flg;
11423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1143eb284becSJed Brown }
11440de4c49aSJed Brown 
1145f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1146d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
1147d71ae5a4SJacob Faibussowitsch {
1148f9c1d6abSBarry Smith   PetscComplex    z   = xr + xi * PETSC_i, f;
1149f9c1d6abSBarry Smith   TS_Theta       *th  = (TS_Theta *)ts->data;
11503fd8ae06SJed Brown   const PetscReal one = 1.0;
1151f9c1d6abSBarry Smith 
1152f9c1d6abSBarry Smith   PetscFunctionBegin;
11533fd8ae06SJed Brown   f   = (one + (one - th->Theta) * z) / (one - th->Theta * z);
1154f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1155f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
11563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1157f9c1d6abSBarry Smith }
1158f9c1d6abSBarry Smith #endif
1159f9c1d6abSBarry Smith 
1160d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[])
1161d71ae5a4SJacob Faibussowitsch {
116242682096SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
116342682096SHong Zhang 
116442682096SHong Zhang   PetscFunctionBegin;
11657409eceaSHong Zhang   if (ns) {
11667409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11677409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
11687409eceaSHong Zhang   }
11695072d23cSHong Zhang   if (Y) {
1170cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
11717409eceaSHong Zhang       th->Stages[0] = th->X;
1172cc4f23bcSHong Zhang     } else {
1173cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1174cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1175cc4f23bcSHong Zhang     }
1176cc4f23bcSHong Zhang     *Y = th->Stages;
11775072d23cSHong Zhang   }
11783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117942682096SHong Zhang }
1180f9c1d6abSBarry Smith 
1181316643e7SJed Brown /* ------------------------------------------------------------ */
1182316643e7SJed Brown /*MC
118396f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1184316643e7SJed Brown 
1185316643e7SJed Brown    Level: beginner
1186316643e7SJed Brown 
1187bcf0153eSBarry Smith    Options Database Keys:
1188c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1189c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
119003842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11914eb428fdSBarry Smith 
1192eb284becSJed Brown    Notes:
119337fdd005SBarry Smith .vb
119437fdd005SBarry Smith   -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
119537fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
119637fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
119737fdd005SBarry Smith .ve
11984eb428fdSBarry Smith 
11997409eceaSHong 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.
1200eb284becSJed Brown 
12017409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1202eb284becSJed Brown 
1203eb284becSJed Brown .vb
1204eb284becSJed Brown   Theta | Theta
1205eb284becSJed Brown   -------------
1206eb284becSJed Brown         |  1
1207eb284becSJed Brown .ve
1208eb284becSJed Brown 
1209eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1210eb284becSJed Brown 
1211eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1212eb284becSJed Brown 
1213eb284becSJed Brown .vb
1214eb284becSJed Brown   0 | 0         0
1215eb284becSJed Brown   1 | 1-Theta   Theta
1216eb284becSJed Brown   -------------------
1217eb284becSJed Brown     | 1-Theta   Theta
1218eb284becSJed Brown .ve
1219eb284becSJed Brown 
1220eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1221eb284becSJed Brown 
1222eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1223eb284becSJed Brown 
1224eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1225eb284becSJed Brown 
12264eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1227eb284becSJed Brown 
12281cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1229316643e7SJed Brown M*/
1230d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1231d71ae5a4SJacob Faibussowitsch {
1232316643e7SJed Brown   TS_Theta *th;
1233316643e7SJed Brown 
1234316643e7SJed Brown   PetscFunctionBegin;
1235277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
1236ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1237316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
1238316643e7SJed Brown   ts->ops->view           = TSView_Theta;
1239316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
124042f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1241ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1242316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
1243cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
12441566a47fSLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
124524655328SShri   ts->ops->rollback       = TSRollBack_Theta;
124626b5f05dSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Theta;
1247316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
12480f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
12490f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1250f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1251f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1252f9c1d6abSBarry Smith #endif
125342682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
125442f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1255b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1256b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12572ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12586affb6f8SHong Zhang 
1259715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12607adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1261715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12626affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1263316643e7SJed Brown 
1264efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1265efd4aadfSBarry Smith 
12664dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
1267316643e7SJed Brown   ts->data = (void *)th;
1268316643e7SJed Brown 
1269715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
1270715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
1271715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
1272b5b4017aSHong Zhang   th->VecsSensi2Temp = NULL;
1273715f1b00SHong Zhang 
12746f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1275316643e7SJed Brown   th->Theta       = 0.5;
12761566a47fSLisandro Dalcin   th->order       = 2;
12779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
12789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
12799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
12809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
12813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1282316643e7SJed Brown }
12830de4c49aSJed Brown 
12840de4c49aSJed Brown /*@
1285bcf0153eSBarry Smith   TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA`
12860de4c49aSJed Brown 
12870de4c49aSJed Brown   Not Collective
12880de4c49aSJed Brown 
12890de4c49aSJed Brown   Input Parameter:
12900de4c49aSJed Brown . ts - timestepping context
12910de4c49aSJed Brown 
12920de4c49aSJed Brown   Output Parameter:
12930de4c49aSJed Brown . theta - stage abscissa
12940de4c49aSJed Brown 
1295bcf0153eSBarry Smith   Level: advanced
1296bcf0153eSBarry Smith 
12970de4c49aSJed Brown   Note:
1298bcf0153eSBarry Smith   Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme.
12990de4c49aSJed Brown 
13001cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA`
13010de4c49aSJed Brown @*/
1302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta)
1303d71ae5a4SJacob Faibussowitsch {
13040de4c49aSJed Brown   PetscFunctionBegin;
1305afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1306dadcf809SJacob Faibussowitsch   PetscValidRealPointer(theta, 2);
1307cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
13083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13090de4c49aSJed Brown }
13100de4c49aSJed Brown 
13110de4c49aSJed Brown /*@
1312bcf0153eSBarry Smith   TSThetaSetTheta - Set the abscissa of the stage in (0,1]  for `TSTHETA`
13130de4c49aSJed Brown 
13140de4c49aSJed Brown   Not Collective
13150de4c49aSJed Brown 
1316d8d19677SJose E. Roman   Input Parameters:
13170de4c49aSJed Brown + ts    - timestepping context
13180de4c49aSJed Brown - theta - stage abscissa
13190de4c49aSJed Brown 
1320bcf0153eSBarry Smith   Options Database Key:
132167b8a455SSatish Balay . -ts_theta_theta <theta> - set theta
13220de4c49aSJed Brown 
1323bcf0153eSBarry Smith   Level: intermediate
13240de4c49aSJed Brown 
13251cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN`
13260de4c49aSJed Brown @*/
1327d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta)
1328d71ae5a4SJacob Faibussowitsch {
13290de4c49aSJed Brown   PetscFunctionBegin;
1330afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1331cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
13323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13330de4c49aSJed Brown }
1334f33bbcb6SJed Brown 
133526f2ff8fSLisandro Dalcin /*@
1336bcf0153eSBarry Smith   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
133726f2ff8fSLisandro Dalcin 
133826f2ff8fSLisandro Dalcin   Not Collective
133926f2ff8fSLisandro Dalcin 
134026f2ff8fSLisandro Dalcin   Input Parameter:
134126f2ff8fSLisandro Dalcin . ts - timestepping context
134226f2ff8fSLisandro Dalcin 
134326f2ff8fSLisandro Dalcin   Output Parameter:
1344bcf0153eSBarry Smith . endpoint - `PETSC_TRUE` when using the endpoint variant
134526f2ff8fSLisandro Dalcin 
1346bcf0153eSBarry Smith   Level: advanced
134726f2ff8fSLisandro Dalcin 
13481cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
134926f2ff8fSLisandro Dalcin @*/
1350d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint)
1351d71ae5a4SJacob Faibussowitsch {
135226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
135326f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1354dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(endpoint, 2);
1355cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
13563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135726f2ff8fSLisandro Dalcin }
135826f2ff8fSLisandro Dalcin 
1359eb284becSJed Brown /*@
1360bcf0153eSBarry Smith   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
1361eb284becSJed Brown 
1362eb284becSJed Brown   Not Collective
1363eb284becSJed Brown 
1364d8d19677SJose E. Roman   Input Parameters:
1365eb284becSJed Brown + ts  - timestepping context
1366bcf0153eSBarry Smith - flg - `PETSC_TRUE` to use the endpoint variant
1367eb284becSJed Brown 
1368bcf0153eSBarry Smith   Options Database Key:
136967b8a455SSatish Balay . -ts_theta_endpoint <flg> - use the endpoint variant
1370eb284becSJed Brown 
1371bcf0153eSBarry Smith   Level: intermediate
1372eb284becSJed Brown 
13731cc06b55SBarry Smith .seealso: [](ch_ts), `TSTHETA`, `TSCN`
1374eb284becSJed Brown @*/
1375d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg)
1376d71ae5a4SJacob Faibussowitsch {
1377eb284becSJed Brown   PetscFunctionBegin;
1378eb284becSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1379cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
13803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1381eb284becSJed Brown }
1382eb284becSJed Brown 
1383f33bbcb6SJed Brown /*
1384f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1385f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1386f33bbcb6SJed Brown  */
1387f33bbcb6SJed Brown 
1388d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BEuler(TS ts)
1389d71ae5a4SJacob Faibussowitsch {
13901566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
13911566a47fSLisandro Dalcin 
13921566a47fSLisandro Dalcin   PetscFunctionBegin;
139308401ef6SPierre 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");
139428b400f6SJacob 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");
13959566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
13963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13971566a47fSLisandro Dalcin }
13981566a47fSLisandro Dalcin 
1399d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer)
1400d71ae5a4SJacob Faibussowitsch {
1401f33bbcb6SJed Brown   PetscFunctionBegin;
14023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1403f33bbcb6SJed Brown }
1404f33bbcb6SJed Brown 
1405f33bbcb6SJed Brown /*MC
1406f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1407f33bbcb6SJed Brown 
1408f33bbcb6SJed Brown   Level: beginner
1409f33bbcb6SJed Brown 
1410bcf0153eSBarry Smith   Note:
141137fdd005SBarry Smith   `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0`
14124eb428fdSBarry Smith 
14131cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1414f33bbcb6SJed Brown M*/
1415d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1416d71ae5a4SJacob Faibussowitsch {
1417f33bbcb6SJed Brown   PetscFunctionBegin;
14189566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14199566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 1.0));
14209566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
14210c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1422f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
14233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1424f33bbcb6SJed Brown }
1425f33bbcb6SJed Brown 
1426d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_CN(TS ts)
1427d71ae5a4SJacob Faibussowitsch {
14281566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
14291566a47fSLisandro Dalcin 
14301566a47fSLisandro Dalcin   PetscFunctionBegin;
143108401ef6SPierre 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");
14323c633725SBarry 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");
14339566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14351566a47fSLisandro Dalcin }
14361566a47fSLisandro Dalcin 
1437d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer)
1438d71ae5a4SJacob Faibussowitsch {
1439f33bbcb6SJed Brown   PetscFunctionBegin;
14403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1441f33bbcb6SJed Brown }
1442f33bbcb6SJed Brown 
1443f33bbcb6SJed Brown /*MC
1444f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1445f33bbcb6SJed Brown 
1446f33bbcb6SJed Brown   Level: beginner
1447f33bbcb6SJed Brown 
1448f33bbcb6SJed Brown   Notes:
1449bcf0153eSBarry Smith   `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e.
1450bcf0153eSBarry Smith .vb
1451bcf0153eSBarry Smith   -ts_type theta
1452bcf0153eSBarry Smith   -ts_theta_theta 0.5
1453bcf0153eSBarry Smith   -ts_theta_endpoint
1454bcf0153eSBarry Smith .ve
14557cf5af47SJed Brown 
14561cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`,
1457f33bbcb6SJed Brown M*/
1458d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1459d71ae5a4SJacob Faibussowitsch {
1460f33bbcb6SJed Brown   PetscFunctionBegin;
14619566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14629566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 0.5));
14639566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
14640c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1465f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
14663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1467f33bbcb6SJed Brown }
1468