xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
459371c9d4SSatish Balay static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
467445fe48SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
477445fe48SJed Brown 
487445fe48SJed Brown   PetscFunctionBegin;
497445fe48SJed Brown   if (X0) {
507445fe48SJed Brown     if (dm && dm != ts->dm) {
519566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0));
527445fe48SJed Brown     } else *X0 = ts->vec_sol;
537445fe48SJed Brown   }
547445fe48SJed Brown   if (Xdot) {
557445fe48SJed Brown     if (dm && dm != ts->dm) {
569566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
577445fe48SJed Brown     } else *Xdot = th->Xdot;
587445fe48SJed Brown   }
597445fe48SJed Brown   PetscFunctionReturn(0);
607445fe48SJed Brown }
617445fe48SJed Brown 
629371c9d4SSatish Balay static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
630d0b770aSPeter Brune   PetscFunctionBegin;
640d0b770aSPeter Brune   if (X0) {
65*48a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0));
660d0b770aSPeter Brune   }
670d0b770aSPeter Brune   if (Xdot) {
68*48a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
690d0b770aSPeter Brune   }
700d0b770aSPeter Brune   PetscFunctionReturn(0);
710d0b770aSPeter Brune }
720d0b770aSPeter Brune 
739371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx) {
747445fe48SJed Brown   PetscFunctionBegin;
757445fe48SJed Brown   PetscFunctionReturn(0);
767445fe48SJed Brown }
777445fe48SJed Brown 
789371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) {
797445fe48SJed Brown   TS  ts = (TS)ctx;
807445fe48SJed Brown   Vec X0, Xdot, X0_c, Xdot_c;
817445fe48SJed Brown 
827445fe48SJed Brown   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot));
849566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
859566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, X0, X0_c));
869566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
879566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
889566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
899566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot));
909566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
917445fe48SJed Brown   PetscFunctionReturn(0);
927445fe48SJed Brown }
937445fe48SJed Brown 
949371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx) {
95258e1594SPeter Brune   PetscFunctionBegin;
96258e1594SPeter Brune   PetscFunctionReturn(0);
97258e1594SPeter Brune }
98258e1594SPeter Brune 
999371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
100258e1594SPeter Brune   TS  ts = (TS)ctx;
101258e1594SPeter Brune   Vec X0, Xdot, X0_sub, Xdot_sub;
102258e1594SPeter Brune 
103258e1594SPeter Brune   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
1059566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
106258e1594SPeter Brune 
1079566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
1089566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
109258e1594SPeter Brune 
1109566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
1119566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
112258e1594SPeter Brune 
1139566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
1149566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
115258e1594SPeter Brune   PetscFunctionReturn(0);
116258e1594SPeter Brune }
117258e1594SPeter Brune 
1189371c9d4SSatish Balay static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) {
119022c081aSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
120cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
121022c081aSHong Zhang 
122022c081aSHong Zhang   PetscFunctionBegin;
123022c081aSHong Zhang   if (th->endpoint) {
124022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
125022c081aSHong Zhang     if (th->Theta != 1.0) {
1269566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand));
1279566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand));
128022c081aSHong Zhang     }
1299566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand));
1309566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand));
131022c081aSHong Zhang   } else {
1329566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand));
1339566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand));
134022c081aSHong Zhang   }
135022c081aSHong Zhang   PetscFunctionReturn(0);
136022c081aSHong Zhang }
137022c081aSHong Zhang 
1389371c9d4SSatish Balay static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) {
139b1cb36f3SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
140cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
141b1cb36f3SHong Zhang 
142b1cb36f3SHong Zhang   PetscFunctionBegin;
143b1cb36f3SHong Zhang   /* backup cost integral */
1449566063dSJacob Faibussowitsch   PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0));
1459566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
146b1cb36f3SHong Zhang   PetscFunctionReturn(0);
147b1cb36f3SHong Zhang }
148b1cb36f3SHong Zhang 
1499371c9d4SSatish Balay static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) {
1501a352fa8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
151b1cb36f3SHong Zhang 
152b1cb36f3SHong Zhang   PetscFunctionBegin;
1531a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
1541a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
1551a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
1569566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
15724655328SShri   PetscFunctionReturn(0);
15824655328SShri }
15924655328SShri 
1609371c9d4SSatish Balay static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x) {
1611566a47fSLisandro Dalcin   PetscInt nits, lits;
1621566a47fSLisandro Dalcin 
1631566a47fSLisandro Dalcin   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
1659566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1669566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1679371c9d4SSatish Balay   ts->snes_its += nits;
1689371c9d4SSatish Balay   ts->ksp_its += lits;
1691566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1701566a47fSLisandro Dalcin }
1711566a47fSLisandro Dalcin 
1729371c9d4SSatish Balay static PetscErrorCode TSStep_Theta(TS ts) {
173316643e7SJed Brown   TS_Theta *th         = (TS_Theta *)ts->data;
1741566a47fSLisandro Dalcin   PetscInt  rejections = 0;
1754957b756SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
1761566a47fSLisandro Dalcin   PetscReal next_time_step = ts->time_step;
177316643e7SJed Brown 
178316643e7SJed Brown   PetscFunctionBegin;
1791566a47fSLisandro Dalcin   if (!ts->steprollback) {
1809566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1819566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
1821566a47fSLisandro Dalcin   }
1831566a47fSLisandro Dalcin 
1841566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
18599260352SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
1863e05ceb1SHong Zhang     th->shift      = 1 / (th->Theta * ts->time_step);
1871566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step;
1889566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X));
189*48a46eb9SPierre Jolivet     if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot));
190eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
1919566063dSJacob Faibussowitsch       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
1929566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
1939566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE));
1949566063dSJacob Faibussowitsch       PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta));
1951566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
1969566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->affine));
197eb284becSJed Brown     }
1989566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
1999566063dSJacob Faibussowitsch     PetscCall(TSTheta_SNESSolve(ts, th->affine, th->X));
2009566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
2019566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
202fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
203051f2191SLisandro Dalcin 
204051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2051566a47fSLisandro Dalcin     if (th->endpoint) {
2069566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X, ts->vec_sol));
2071566a47fSLisandro Dalcin     } else {
2089566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
2099566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
2101566a47fSLisandro Dalcin     }
2119566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2121566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2131566a47fSLisandro Dalcin     if (!accept) {
2149566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
215be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
216be5899b3SLisandro Dalcin       goto reject_step;
217051f2191SLisandro Dalcin     }
2183b1890cdSShri Abhyankar 
219715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2201a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2211a352fa8SHong Zhang       th->time_step0 = ts->time_step;
22217145e75SHong Zhang     }
2232b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
224cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
225051f2191SLisandro Dalcin     break;
226051f2191SLisandro Dalcin 
227051f2191SLisandro Dalcin   reject_step:
2289371c9d4SSatish Balay     ts->reject++;
2299371c9d4SSatish Balay     accept = PETSC_FALSE;
2301566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
231051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
23263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
2333b1890cdSShri Abhyankar     }
2343b1890cdSShri Abhyankar   }
235316643e7SJed Brown   PetscFunctionReturn(0);
236316643e7SJed Brown }
237316643e7SJed Brown 
2389371c9d4SSatish Balay static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) {
239c9681e5dSHong Zhang   TS_Theta      *th           = (TS_Theta *)ts->data;
240cd4cee2dSHong Zhang   TS             quadts       = ts->quadraturets;
241c9681e5dSHong Zhang   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
242c9681e5dSHong Zhang   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
243c9681e5dSHong Zhang   PetscInt       nadj;
2447207e0fdSHong Zhang   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
245c9681e5dSHong Zhang   KSP            ksp;
246c9681e5dSHong Zhang   PetscScalar   *xarr;
247077c4feaSHong Zhang   TSEquationType eqtype;
248077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2491a352fa8SHong Zhang   PetscReal      adjoint_time_step;
250c9681e5dSHong Zhang 
251c9681e5dSHong Zhang   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(TSGetEquationType(ts, &eqtype));
253077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
254077c4feaSHong Zhang     isexplicitode = PETSC_TRUE;
255077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
256077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
257077c4feaSHong Zhang   }
258c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
2599566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
2609566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
2617207e0fdSHong Zhang   if (quadts) {
2629566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
2639566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
2647207e0fdSHong Zhang   }
265c9681e5dSHong Zhang 
2661a352fa8SHong Zhang   th->stage_time    = ts->ptime;
2671a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
268c9681e5dSHong Zhang 
26933f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2701baa6e33SBarry Smith   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
271cd4cee2dSHong Zhang 
272c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
2739566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
2749566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
275cd4cee2dSHong Zhang     if (quadJ) {
2769566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
2779566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
2789566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
2799566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
2809566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
281c9681e5dSHong Zhang     }
282c9681e5dSHong Zhang   }
283c9681e5dSHong Zhang 
284c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
285c87ba875SIvan Yashchuk   th->shift = 1. / adjoint_time_step;
2869566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
2879566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
288c9681e5dSHong Zhang 
289c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
290c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
291c9681e5dSHong Zhang     KSPConvergedReason kspreason;
2929566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
2939566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
294c9681e5dSHong Zhang     if (kspreason < 0) {
295c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
29663a3b9bcSJacob 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));
297c9681e5dSHong Zhang     }
298c9681e5dSHong Zhang   }
299c9681e5dSHong Zhang 
300c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
301c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3029566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
3039566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
304c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3059566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
306c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3079566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
308c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
3099566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
3109566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
3119566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
312*48a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
313c9681e5dSHong Zhang     }
314c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
315c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
31605c9950eSHong Zhang       KSPConvergedReason kspreason;
3179566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
3189566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
31905c9950eSHong Zhang       if (kspreason < 0) {
32005c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32163a3b9bcSJacob 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));
32205c9950eSHong Zhang       }
323c9681e5dSHong Zhang     }
324c9681e5dSHong Zhang   }
325c9681e5dSHong Zhang 
326c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
327077c4feaSHong Zhang   if (!isexplicitode) {
3283e05ceb1SHong Zhang     th->shift = 0.0;
3299566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3309566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
331c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
33233f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3339566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
3349566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
3359566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
3369566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
337c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3389566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
3399566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
3409566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
341c9681e5dSHong Zhang       }
342c9681e5dSHong Zhang     }
343077c4feaSHong Zhang   }
344c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3459566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
3469566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
3471baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
348c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
349c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3509566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
35133f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3529566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
353c9681e5dSHong Zhang     }
354cd4cee2dSHong Zhang 
355c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
3569566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
3579566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
358cd4cee2dSHong Zhang       if (quadJp) {
3599566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
3609566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
3619566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
3629566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3639566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
36433f72c5dSHong Zhang       }
365c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
3669566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
3679566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
368*48a46eb9SPierre Jolivet         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
369*48a46eb9SPierre Jolivet         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
370c9681e5dSHong Zhang       }
371c9681e5dSHong Zhang     }
372c9681e5dSHong Zhang   }
373c9681e5dSHong Zhang 
374c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
3759566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
3769566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
377c9681e5dSHong Zhang   }
378c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
379c9681e5dSHong Zhang   PetscFunctionReturn(0);
380c9681e5dSHong Zhang }
381c9681e5dSHong Zhang 
3829371c9d4SSatish Balay static PetscErrorCode TSAdjointStep_Theta(TS ts) {
3832ca6e920SHong Zhang   TS_Theta    *th           = (TS_Theta *)ts->data;
384cd4cee2dSHong Zhang   TS           quadts       = ts->quadraturets;
385b5e0532dSHong Zhang   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
386b5b4017aSHong Zhang   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
3872ca6e920SHong Zhang   PetscInt     nadj;
3887207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
3892ca6e920SHong Zhang   KSP          ksp;
390b5b4017aSHong Zhang   PetscScalar *xarr;
3911a352fa8SHong Zhang   PetscReal    adjoint_time_step;
3921a352fa8SHong Zhang   PetscReal    adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, ususally ts->ptime is larger than adjoint_ptime) */
3932ca6e920SHong Zhang 
3942ca6e920SHong Zhang   PetscFunctionBegin;
395c9681e5dSHong Zhang   if (th->Theta == 1.) {
3969566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
397c9681e5dSHong Zhang     PetscFunctionReturn(0);
398c9681e5dSHong Zhang   }
3992ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4009566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
4019566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
4027207e0fdSHong Zhang   if (quadts) {
4039566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
4049566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
4057207e0fdSHong Zhang   }
406b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4071a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
4081a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4091a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
4102ca6e920SHong Zhang 
41182ad9101SHong Zhang   if (!th->endpoint) {
4125072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4139566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
41482ad9101SHong Zhang   }
41582ad9101SHong Zhang 
416b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
41733f72c5dSHong Zhang   /* Cost function has an integral term */
4187207e0fdSHong Zhang   if (quadts) {
41905755b9cSHong Zhang     if (th->endpoint) {
4209566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
42105755b9cSHong Zhang     } else {
4229566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
42305755b9cSHong Zhang     }
4247207e0fdSHong Zhang   }
42533f72c5dSHong Zhang 
426abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4279566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
4289566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
429cd4cee2dSHong Zhang     if (quadJ) {
4309566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
4319566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
4329566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
4339566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4349566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
43536eaed60SHong Zhang     }
4362ca6e920SHong Zhang   }
4373c54f92cSHong Zhang 
438b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4391a352fa8SHong Zhang   th->shift = 1. / (th->Theta * adjoint_time_step);
4403c54f92cSHong Zhang   if (th->endpoint) {
4419566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
4423c54f92cSHong Zhang   } else {
4439566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
4443c54f92cSHong Zhang   }
4459566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
4462ca6e920SHong Zhang 
447b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
448abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4491d14d03bSHong Zhang     KSPConvergedReason kspreason;
4509566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
4519566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
4521d14d03bSHong Zhang     if (kspreason < 0) {
4531d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
45463a3b9bcSJacob 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));
4551d14d03bSHong Zhang     }
4562ca6e920SHong Zhang   }
4573c54f92cSHong Zhang 
4581d14d03bSHong Zhang   /* Second-order adjoint */
459b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4603c633725SBarry Smith     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
461b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4629566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
4639566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
464b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
4659566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
4669566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4679566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
468b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
4699566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
470b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
4719566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
4729566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
4739566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
474*48a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
475b5b4017aSHong Zhang     }
476b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
477b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
4781d14d03bSHong Zhang       KSPConvergedReason kspreason;
4799566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
4809566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
4811d14d03bSHong Zhang       if (kspreason < 0) {
4821d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
48363a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
4841d14d03bSHong Zhang       }
485b5b4017aSHong Zhang     }
486b5b4017aSHong Zhang   }
487b5b4017aSHong Zhang 
48836eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
4891d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
4901a352fa8SHong Zhang     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
4911a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
4929566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
4939566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
49433f72c5dSHong Zhang     /* R_U at t_n */
4951baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
496abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
4979566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
498cd4cee2dSHong Zhang       if (quadJ) {
4999566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
5009566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
5019566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
5029566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5039566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
504b5b4017aSHong Zhang       }
5059566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
506b5b4017aSHong Zhang     }
5071d14d03bSHong Zhang 
5081d14d03bSHong Zhang     /* Second-order adjoint */
5091d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
510b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5119566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
5129566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
513b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5149566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
5159566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5169566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
517b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5189566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
519b5b4017aSHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
52033f72c5dSHong 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  */
5219566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
5229566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
523*48a46eb9SPierre Jolivet         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
5249566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
525b5b4017aSHong Zhang       }
526b5e0532dSHong Zhang     }
5273fd52205SHong Zhang 
528c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
529c586f2e8SHong Zhang 
5303fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
531b5b4017aSHong Zhang       /* U_{n+1} */
53282ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
5339566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
5349566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5351baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
536abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5379566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
5389566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
5393b711c3fSHong Zhang         if (quadJp) {
5409566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
5419566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
5429566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
5439566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5449566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
5453b711c3fSHong Zhang         }
5463fd52205SHong Zhang       }
54733f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
54833f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5499566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
5509566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
551b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5529566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
5539566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5549566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
55533f72c5dSHong Zhang 
556b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5579566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
558b5b4017aSHong Zhang         for (nadj = 0; nadj < ts->numcost; nadj++) {
55933f72c5dSHong 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)  */
5609566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
5619566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
562*48a46eb9SPierre Jolivet           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
563*48a46eb9SPierre Jolivet           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
564b5b4017aSHong Zhang         }
565b5b4017aSHong Zhang       }
566b5b4017aSHong Zhang 
567b5b4017aSHong Zhang       /* U_s */
5689566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
5699566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5701baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
571abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5729566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
5739566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
5743b711c3fSHong Zhang         if (quadJp) {
5759566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
5769566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
5779566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
5789566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5799566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
5803b711c3fSHong Zhang         }
58133f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
58233f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
5839566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
5849566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
585b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
5869566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
5879566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
5889566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
589b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
5909566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
591b5b4017aSHong Zhang           for (nadj = 0; nadj < ts->numcost; nadj++) {
59233f72c5dSHong 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) */
5939566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
5949566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
595*48a46eb9SPierre Jolivet             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
596*48a46eb9SPierre Jolivet             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
59736eaed60SHong Zhang           }
59836eaed60SHong Zhang         }
5993fd52205SHong Zhang       }
600b5e0532dSHong Zhang     }
6013fd52205SHong Zhang   } else { /* one-stage case */
6023e05ceb1SHong Zhang     th->shift = 0.0;
6039566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
6049566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
6051baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
606abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
6079566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
6089566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
609cd4cee2dSHong Zhang       if (quadJ) {
6109566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
6119566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
6129566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
6139566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6149566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
61536eaed60SHong Zhang       }
6162ca6e920SHong Zhang     }
6173fd52205SHong Zhang     if (ts->vecs_sensip) {
61882ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
6199566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
6209566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
6211baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
622abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6239566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6249566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
625cd4cee2dSHong Zhang         if (quadJp) {
6269566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6279566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6289566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
6299566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6309566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
63136eaed60SHong Zhang         }
63236eaed60SHong Zhang       }
6333fd52205SHong Zhang     }
6342ca6e920SHong Zhang   }
6352ca6e920SHong Zhang 
6362ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6372ca6e920SHong Zhang   PetscFunctionReturn(0);
6382ca6e920SHong Zhang }
6392ca6e920SHong Zhang 
6409371c9d4SSatish Balay static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X) {
641cd652676SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
642be5899b3SLisandro Dalcin   PetscReal dt = t - ts->ptime;
643cd652676SJed Brown 
644cd652676SJed Brown   PetscFunctionBegin;
6459566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X));
646be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6479566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
648cd652676SJed Brown   PetscFunctionReturn(0);
649cd652676SJed Brown }
650cd652676SJed Brown 
6519371c9d4SSatish Balay static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) {
6521566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
6531566a47fSLisandro Dalcin   Vec       X  = ts->vec_sol;      /* X = solution */
6541566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
6557453f775SEmil Constantinescu   PetscReal wltea, wlter;
6561566a47fSLisandro Dalcin 
6571566a47fSLisandro Dalcin   PetscFunctionBegin;
6589371c9d4SSatish Balay   if (!th->vec_sol_prev) {
6599371c9d4SSatish Balay     *wlte = -1;
6609371c9d4SSatish Balay     PetscFunctionReturn(0);
6619371c9d4SSatish Balay   }
6621566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
6639371c9d4SSatish Balay   if (ts->steprestart) {
6649371c9d4SSatish Balay     *wlte = -1;
6659371c9d4SSatish Balay     PetscFunctionReturn(0);
6669371c9d4SSatish Balay   }
6671566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
668fecfb714SLisandro Dalcin   {
669be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
670be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
6719371c9d4SSatish Balay     PetscScalar scal[3];
6729371c9d4SSatish Balay     Vec         vecs[3];
6739371c9d4SSatish Balay     scal[0] = +1 / a;
6749371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
6759371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
6769371c9d4SSatish Balay     vecs[0] = X;
6779371c9d4SSatish Balay     vecs[1] = th->X0;
6789371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
6799566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
6809566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
6819566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
6821566a47fSLisandro Dalcin   }
6831566a47fSLisandro Dalcin   if (order) *order = 2;
6841566a47fSLisandro Dalcin   PetscFunctionReturn(0);
6851566a47fSLisandro Dalcin }
6861566a47fSLisandro Dalcin 
6879371c9d4SSatish Balay static PetscErrorCode TSRollBack_Theta(TS ts) {
6881566a47fSLisandro Dalcin   TS_Theta *th     = (TS_Theta *)ts->data;
6897207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
6901566a47fSLisandro Dalcin 
6911566a47fSLisandro Dalcin   PetscFunctionBegin;
6929566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
693*48a46eb9SPierre Jolivet   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
694715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
6951baa6e33SBarry Smith   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
696*48a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
697715f1b00SHong Zhang   PetscFunctionReturn(0);
698715f1b00SHong Zhang }
699715f1b00SHong Zhang 
7009371c9d4SSatish Balay static PetscErrorCode TSForwardStep_Theta(TS ts) {
701715f1b00SHong Zhang   TS_Theta    *th                   = (TS_Theta *)ts->data;
702cd4cee2dSHong Zhang   TS           quadts               = ts->quadraturets;
70313b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
704b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7057207e0fdSHong Zhang   PetscInt     ntlm;
706715f1b00SHong Zhang   KSP          ksp;
7077207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
70813b1b0ffSHong Zhang   PetscScalar *barr, *xarr;
7091a352fa8SHong Zhang   PetscReal    previous_shift;
710715f1b00SHong Zhang 
711715f1b00SHong Zhang   PetscFunctionBegin;
7121a352fa8SHong Zhang   previous_shift = th->shift;
7139566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));
71413b1b0ffSHong Zhang 
715*48a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
7169566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
7179566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
7187207e0fdSHong Zhang   if (quadts) {
7199566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
7209566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
7217207e0fdSHong Zhang   }
722715f1b00SHong Zhang 
723715f1b00SHong Zhang   /* Build RHS */
724715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7251a352fa8SHong Zhang     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
7269566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7279566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7289566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));
729715f1b00SHong Zhang 
730022c081aSHong Zhang     /* Add the f_p forcing terms */
7310e11c21fSHong Zhang     if (ts->Jacp) {
7329566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7339566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7349566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
73582ad9101SHong Zhang       th->shift = previous_shift;
7369566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
7379566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7389566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
7390e11c21fSHong Zhang     }
740715f1b00SHong Zhang   } else { /* 1-stage method */
7413e05ceb1SHong Zhang     th->shift = 0.0;
7429566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
7439566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
7449566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, -1.));
74513b1b0ffSHong Zhang 
746022c081aSHong Zhang     /* Add the f_p forcing terms */
7470e11c21fSHong Zhang     if (ts->Jacp) {
74882ad9101SHong Zhang       th->shift = previous_shift;
7499566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
7509566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7519566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
752715f1b00SHong Zhang     }
7530e11c21fSHong Zhang   }
754715f1b00SHong Zhang 
755715f1b00SHong Zhang   /* Build LHS */
7561a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
757715f1b00SHong Zhang   if (th->endpoint) {
7589566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
759715f1b00SHong Zhang   } else {
7609566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
761715f1b00SHong Zhang   }
7629566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
763715f1b00SHong Zhang 
764715f1b00SHong Zhang   /*
765715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
766c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
767715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
768715f1b00SHong Zhang   */
769715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
7707207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
7719566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
7729566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
7739566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
7749566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
7759566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
776715f1b00SHong Zhang     }
777715f1b00SHong Zhang   }
778715f1b00SHong Zhang 
779715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
780022c081aSHong Zhang   for (ntlm = 0; ntlm < th->num_tlm; ntlm++) {
781b5b4017aSHong Zhang     KSPConvergedReason kspreason;
7829566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
7839566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
784715f1b00SHong Zhang     if (th->endpoint) {
7859566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
7869566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
7879566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
7889566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
7899566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
790715f1b00SHong Zhang     } else {
7919566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
792715f1b00SHong Zhang     }
7939566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
794b5b4017aSHong Zhang     if (kspreason < 0) {
795b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
79663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
797b5b4017aSHong Zhang     }
7989566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
7999566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
800715f1b00SHong Zhang   }
801715f1b00SHong Zhang 
80213b1b0ffSHong Zhang   /*
80313b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
804c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
80513b1b0ffSHong Zhang   */
8067207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
80713b1b0ffSHong Zhang     if (!th->endpoint) {
8089566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
8099566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
8109566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
8119566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8129566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8139566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
8149566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
81513b1b0ffSHong Zhang     } else {
8169566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
8179566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
8189566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
8199566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8209566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
82113b1b0ffSHong Zhang     }
82213b1b0ffSHong Zhang   } else {
823*48a46eb9SPierre Jolivet     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
824715f1b00SHong Zhang   }
8251566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8261566a47fSLisandro Dalcin }
8271566a47fSLisandro Dalcin 
8289371c9d4SSatish Balay static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[]) {
8296affb6f8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
8306affb6f8SHong Zhang 
8316affb6f8SHong Zhang   PetscFunctionBegin;
8327409eceaSHong Zhang   if (ns) {
8337409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8347409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
8357409eceaSHong Zhang   }
8365072d23cSHong Zhang   if (stagesensip) {
837cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8387409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
839cc4f23bcSHong Zhang     } else {
840cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
841cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
842cc4f23bcSHong Zhang     }
843cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8445072d23cSHong Zhang   }
8456affb6f8SHong Zhang   PetscFunctionReturn(0);
8466affb6f8SHong Zhang }
8476affb6f8SHong Zhang 
848316643e7SJed Brown /*------------------------------------------------------------*/
8499371c9d4SSatish Balay static PetscErrorCode TSReset_Theta(TS ts) {
850316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
851316643e7SJed Brown 
852316643e7SJed Brown   PetscFunctionBegin;
8539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
8549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
8559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
8569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
8571566a47fSLisandro Dalcin 
8589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
8599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
8601566a47fSLisandro Dalcin 
8619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
862277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
863277b19d0SLisandro Dalcin }
864277b19d0SLisandro Dalcin 
8659371c9d4SSatish Balay static PetscErrorCode TSAdjointReset_Theta(TS ts) {
866ece46509SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
867ece46509SHong Zhang 
868ece46509SHong Zhang   PetscFunctionBegin;
8699566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
8709566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
8719566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
8729566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
8739566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
8749566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
875ece46509SHong Zhang   PetscFunctionReturn(0);
876ece46509SHong Zhang }
877ece46509SHong Zhang 
8789371c9d4SSatish Balay static PetscErrorCode TSDestroy_Theta(TS ts) {
879277b19d0SLisandro Dalcin   PetscFunctionBegin;
8809566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
881b3a6b972SJed Brown   if (ts->dm) {
8829566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
8839566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
884b3a6b972SJed Brown   }
8859566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
8869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
8879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
8889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
8899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
890316643e7SJed Brown   PetscFunctionReturn(0);
891316643e7SJed Brown }
892316643e7SJed Brown 
893316643e7SJed Brown /*
894316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
8952b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
8960fd10508SMatthew Knepley 
8970fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
8980fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
8990fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
900316643e7SJed Brown */
9019371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts) {
902316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9037445fe48SJed Brown   Vec       X0, Xdot;
9047445fe48SJed Brown   DM        dm, dmsave;
9053e05ceb1SHong Zhang   PetscReal shift = th->shift;
906316643e7SJed Brown 
907316643e7SJed Brown   PetscFunctionBegin;
9089566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
9095a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9109566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
911494ce9fcSHong Zhang   if (x != X0) {
9129566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
913494ce9fcSHong Zhang   } else {
9149566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
915494ce9fcSHong Zhang   }
9167445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9177445fe48SJed Brown   dmsave = ts->dm;
9187445fe48SJed Brown   ts->dm = dm;
9199566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
9207445fe48SJed Brown   ts->dm = dmsave;
9219566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
922316643e7SJed Brown   PetscFunctionReturn(0);
923316643e7SJed Brown }
924316643e7SJed Brown 
9259371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts) {
926316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9277445fe48SJed Brown   Vec       Xdot;
9287445fe48SJed Brown   DM        dm, dmsave;
9293e05ceb1SHong Zhang   PetscReal shift = th->shift;
930316643e7SJed Brown 
931316643e7SJed Brown   PetscFunctionBegin;
9329566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
933be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9349566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));
9357445fe48SJed Brown 
9367445fe48SJed Brown   dmsave = ts->dm;
9377445fe48SJed Brown   ts->dm = dm;
9389566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
9397445fe48SJed Brown   ts->dm = dmsave;
9409566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
941316643e7SJed Brown   PetscFunctionReturn(0);
942316643e7SJed Brown }
943316643e7SJed Brown 
9449371c9d4SSatish Balay static PetscErrorCode TSForwardSetUp_Theta(TS ts) {
945715f1b00SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
9467207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
947715f1b00SHong Zhang 
948715f1b00SHong Zhang   PetscFunctionBegin;
949022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
95013b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
9519566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
9527207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9539566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
9549566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
955715f1b00SHong Zhang   }
9566f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
9579566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));
95813b1b0ffSHong Zhang 
9599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
960715f1b00SHong Zhang   PetscFunctionReturn(0);
961715f1b00SHong Zhang }
962715f1b00SHong Zhang 
9639371c9d4SSatish Balay static PetscErrorCode TSForwardReset_Theta(TS ts) {
9647adebddeSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
9657207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
9667adebddeSHong Zhang 
9677adebddeSHong Zhang   PetscFunctionBegin;
9687207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9699566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
9709566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
9717adebddeSHong Zhang   }
9729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
9739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
9749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
9757adebddeSHong Zhang   PetscFunctionReturn(0);
9767adebddeSHong Zhang }
9777adebddeSHong Zhang 
9789371c9d4SSatish Balay static PetscErrorCode TSSetUp_Theta(TS ts) {
979316643e7SJed Brown   TS_Theta *th     = (TS_Theta *)ts->data;
980cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
9812ffb9264SLisandro Dalcin   PetscBool match;
982316643e7SJed Brown 
983316643e7SJed Brown   PetscFunctionBegin;
984cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
9859566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
986b1cb36f3SHong Zhang   }
987*48a46eb9SPierre Jolivet   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
988*48a46eb9SPierre Jolivet   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
989*48a46eb9SPierre Jolivet   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
990*48a46eb9SPierre Jolivet   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
9913b1890cdSShri Abhyankar 
9921566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
99320d49056SMatthew G. Knepley   th->shift = 1 / (th->Theta * ts->time_step);
9941566a47fSLisandro Dalcin 
9959566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &ts->dm));
9969566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
9979566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
9981566a47fSLisandro Dalcin 
9999566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
10009566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
10022ffb9264SLisandro Dalcin   if (!match) {
10039566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
10049566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
10053b1890cdSShri Abhyankar   }
10069566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
1007cc4f23bcSHong Zhang 
1008cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1009b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1010b4dd3bf9SBarry Smith }
10110c7376e5SHong Zhang 
1012b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1013b4dd3bf9SBarry Smith 
10149371c9d4SSatish Balay static PetscErrorCode TSAdjointSetUp_Theta(TS ts) {
1015b4dd3bf9SBarry Smith   TS_Theta *th = (TS_Theta *)ts->data;
1016b4dd3bf9SBarry Smith 
1017b4dd3bf9SBarry Smith   PetscFunctionBegin;
10189566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
10199566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
1020*48a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1021b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10229566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
10239566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
102467633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
102567633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
102667633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1027b5b4017aSHong Zhang   }
1028c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
10299566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
103067633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
103167633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
103267633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1033b5b4017aSHong Zhang   }
1034316643e7SJed Brown   PetscFunctionReturn(0);
1035316643e7SJed Brown }
1036316643e7SJed Brown 
10379371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject) {
1038316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1039316643e7SJed Brown 
1040316643e7SJed Brown   PetscFunctionBegin;
1041d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1042316643e7SJed Brown   {
10439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
10449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
10459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1046316643e7SJed Brown   }
1047d0609cedSBarry Smith   PetscOptionsHeadEnd();
1048316643e7SJed Brown   PetscFunctionReturn(0);
1049316643e7SJed Brown }
1050316643e7SJed Brown 
10519371c9d4SSatish Balay static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer) {
1052316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1053ace3abfcSBarry Smith   PetscBool iascii;
1054316643e7SJed Brown 
1055316643e7SJed Brown   PetscFunctionBegin;
10569566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1057316643e7SJed Brown   if (iascii) {
10589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
10599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1060316643e7SJed Brown   }
1061316643e7SJed Brown   PetscFunctionReturn(0);
1062316643e7SJed Brown }
1063316643e7SJed Brown 
10649371c9d4SSatish Balay static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta) {
10650de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
10660de4c49aSJed Brown 
10670de4c49aSJed Brown   PetscFunctionBegin;
10680de4c49aSJed Brown   *theta = th->Theta;
10690de4c49aSJed Brown   PetscFunctionReturn(0);
10700de4c49aSJed Brown }
10710de4c49aSJed Brown 
10729371c9d4SSatish Balay static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta) {
10730de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
10740de4c49aSJed Brown 
10750de4c49aSJed Brown   PetscFunctionBegin;
1076cad9d221SBarry Smith   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
10770de4c49aSJed Brown   th->Theta = theta;
10781566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
10790de4c49aSJed Brown   PetscFunctionReturn(0);
10800de4c49aSJed Brown }
1081eb284becSJed Brown 
10829371c9d4SSatish Balay static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint) {
108326f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
108426f2ff8fSLisandro Dalcin 
108526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
108626f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
108726f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
108826f2ff8fSLisandro Dalcin }
108926f2ff8fSLisandro Dalcin 
10909371c9d4SSatish Balay static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg) {
1091eb284becSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1092eb284becSJed Brown 
1093eb284becSJed Brown   PetscFunctionBegin;
1094eb284becSJed Brown   th->endpoint = flg;
1095eb284becSJed Brown   PetscFunctionReturn(0);
1096eb284becSJed Brown }
10970de4c49aSJed Brown 
1098f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
10999371c9d4SSatish Balay static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) {
1100f9c1d6abSBarry Smith   PetscComplex    z   = xr + xi * PETSC_i, f;
1101f9c1d6abSBarry Smith   TS_Theta       *th  = (TS_Theta *)ts->data;
11023fd8ae06SJed Brown   const PetscReal one = 1.0;
1103f9c1d6abSBarry Smith 
1104f9c1d6abSBarry Smith   PetscFunctionBegin;
11053fd8ae06SJed Brown   f   = (one + (one - th->Theta) * z) / (one - th->Theta * z);
1106f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1107f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1108f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1109f9c1d6abSBarry Smith }
1110f9c1d6abSBarry Smith #endif
1111f9c1d6abSBarry Smith 
11129371c9d4SSatish Balay static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[]) {
111342682096SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
111442682096SHong Zhang 
111542682096SHong Zhang   PetscFunctionBegin;
11167409eceaSHong Zhang   if (ns) {
11177409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11187409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
11197409eceaSHong Zhang   }
11205072d23cSHong Zhang   if (Y) {
1121cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
11227409eceaSHong Zhang       th->Stages[0] = th->X;
1123cc4f23bcSHong Zhang     } else {
1124cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1125cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1126cc4f23bcSHong Zhang     }
1127cc4f23bcSHong Zhang     *Y = th->Stages;
11285072d23cSHong Zhang   }
112942682096SHong Zhang   PetscFunctionReturn(0);
113042682096SHong Zhang }
1131f9c1d6abSBarry Smith 
1132316643e7SJed Brown /* ------------------------------------------------------------ */
1133316643e7SJed Brown /*MC
113496f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1135316643e7SJed Brown 
1136316643e7SJed Brown    Level: beginner
1137316643e7SJed Brown 
11384eb428fdSBarry Smith    Options Database:
1139c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1140c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
114103842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11424eb428fdSBarry Smith 
1143eb284becSJed Brown    Notes:
11440c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
11450c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
11464eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
11474eb428fdSBarry Smith 
11487409eceaSHong 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.
1149eb284becSJed Brown 
11507409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1151eb284becSJed Brown 
1152eb284becSJed Brown .vb
1153eb284becSJed Brown   Theta | Theta
1154eb284becSJed Brown   -------------
1155eb284becSJed Brown         |  1
1156eb284becSJed Brown .ve
1157eb284becSJed Brown 
1158eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1159eb284becSJed Brown 
1160eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1161eb284becSJed Brown 
1162eb284becSJed Brown .vb
1163eb284becSJed Brown   0 | 0         0
1164eb284becSJed Brown   1 | 1-Theta   Theta
1165eb284becSJed Brown   -------------------
1166eb284becSJed Brown     | 1-Theta   Theta
1167eb284becSJed Brown .ve
1168eb284becSJed Brown 
1169eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1170eb284becSJed Brown 
1171eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1172eb284becSJed Brown 
1173eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1174eb284becSJed Brown 
11754eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1176eb284becSJed Brown 
1177db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1178316643e7SJed Brown 
1179316643e7SJed Brown M*/
11809371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) {
1181316643e7SJed Brown   TS_Theta *th;
1182316643e7SJed Brown 
1183316643e7SJed Brown   PetscFunctionBegin;
1184277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
1185ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1186316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
1187316643e7SJed Brown   ts->ops->view           = TSView_Theta;
1188316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
118942f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1190ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1191316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
1192cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
11931566a47fSLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
119424655328SShri   ts->ops->rollback       = TSRollBack_Theta;
1195316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
11960f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
11970f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1198f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1199f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1200f9c1d6abSBarry Smith #endif
120142682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
120242f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1203b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1204b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12052ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12066affb6f8SHong Zhang 
1207715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12087adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1209715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12106affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1211316643e7SJed Brown 
1212efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1213efd4aadfSBarry Smith 
12149566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &th));
1215316643e7SJed Brown   ts->data = (void *)th;
1216316643e7SJed Brown 
1217715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
1218715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
1219715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
1220b5b4017aSHong Zhang   th->VecsSensi2Temp = NULL;
1221715f1b00SHong Zhang 
12226f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1223316643e7SJed Brown   th->Theta       = 0.5;
12241566a47fSLisandro Dalcin   th->order       = 2;
12259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
12269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
12279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
12289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
1229316643e7SJed Brown   PetscFunctionReturn(0);
1230316643e7SJed Brown }
12310de4c49aSJed Brown 
12320de4c49aSJed Brown /*@
12330de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12340de4c49aSJed Brown 
12350de4c49aSJed Brown   Not Collective
12360de4c49aSJed Brown 
12370de4c49aSJed Brown   Input Parameter:
12380de4c49aSJed Brown .  ts - timestepping context
12390de4c49aSJed Brown 
12400de4c49aSJed Brown   Output Parameter:
12410de4c49aSJed Brown .  theta - stage abscissa
12420de4c49aSJed Brown 
12430de4c49aSJed Brown   Note:
12440de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
12450de4c49aSJed Brown 
12460de4c49aSJed Brown   Level: Advanced
12470de4c49aSJed Brown 
1248db781477SPatrick Sanan .seealso: `TSThetaSetTheta()`
12490de4c49aSJed Brown @*/
12509371c9d4SSatish Balay PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta) {
12510de4c49aSJed Brown   PetscFunctionBegin;
1252afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1253dadcf809SJacob Faibussowitsch   PetscValidRealPointer(theta, 2);
1254cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
12550de4c49aSJed Brown   PetscFunctionReturn(0);
12560de4c49aSJed Brown }
12570de4c49aSJed Brown 
12580de4c49aSJed Brown /*@
12590de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
12600de4c49aSJed Brown 
12610de4c49aSJed Brown   Not Collective
12620de4c49aSJed Brown 
1263d8d19677SJose E. Roman   Input Parameters:
12640de4c49aSJed Brown +  ts - timestepping context
12650de4c49aSJed Brown -  theta - stage abscissa
12660de4c49aSJed Brown 
12670de4c49aSJed Brown   Options Database:
126867b8a455SSatish Balay .  -ts_theta_theta <theta> - set theta
12690de4c49aSJed Brown 
12700de4c49aSJed Brown   Level: Intermediate
12710de4c49aSJed Brown 
1272db781477SPatrick Sanan .seealso: `TSThetaGetTheta()`
12730de4c49aSJed Brown @*/
12749371c9d4SSatish Balay PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta) {
12750de4c49aSJed Brown   PetscFunctionBegin;
1276afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1277cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
12780de4c49aSJed Brown   PetscFunctionReturn(0);
12790de4c49aSJed Brown }
1280f33bbcb6SJed Brown 
128126f2ff8fSLisandro Dalcin /*@
128226f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
128326f2ff8fSLisandro Dalcin 
128426f2ff8fSLisandro Dalcin   Not Collective
128526f2ff8fSLisandro Dalcin 
128626f2ff8fSLisandro Dalcin   Input Parameter:
128726f2ff8fSLisandro Dalcin .  ts - timestepping context
128826f2ff8fSLisandro Dalcin 
128926f2ff8fSLisandro Dalcin   Output Parameter:
129026f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
129126f2ff8fSLisandro Dalcin 
129226f2ff8fSLisandro Dalcin   Level: Advanced
129326f2ff8fSLisandro Dalcin 
1294db781477SPatrick Sanan .seealso: `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
129526f2ff8fSLisandro Dalcin @*/
12969371c9d4SSatish Balay PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint) {
129726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
129826f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1299dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(endpoint, 2);
1300cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
130126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
130226f2ff8fSLisandro Dalcin }
130326f2ff8fSLisandro Dalcin 
1304eb284becSJed Brown /*@
1305eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1306eb284becSJed Brown 
1307eb284becSJed Brown   Not Collective
1308eb284becSJed Brown 
1309d8d19677SJose E. Roman   Input Parameters:
1310eb284becSJed Brown +  ts - timestepping context
1311eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1312eb284becSJed Brown 
1313eb284becSJed Brown   Options Database:
131467b8a455SSatish Balay .  -ts_theta_endpoint <flg> - use the endpoint variant
1315eb284becSJed Brown 
1316eb284becSJed Brown   Level: Intermediate
1317eb284becSJed Brown 
1318db781477SPatrick Sanan .seealso: `TSTHETA`, `TSCN`
1319eb284becSJed Brown @*/
13209371c9d4SSatish Balay PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg) {
1321eb284becSJed Brown   PetscFunctionBegin;
1322eb284becSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1323cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
1324eb284becSJed Brown   PetscFunctionReturn(0);
1325eb284becSJed Brown }
1326eb284becSJed Brown 
1327f33bbcb6SJed Brown /*
1328f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1329f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1330f33bbcb6SJed Brown  */
1331f33bbcb6SJed Brown 
13329371c9d4SSatish Balay static PetscErrorCode TSSetUp_BEuler(TS ts) {
13331566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
13341566a47fSLisandro Dalcin 
13351566a47fSLisandro Dalcin   PetscFunctionBegin;
133608401ef6SPierre 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");
133728b400f6SJacob 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");
13389566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
13391566a47fSLisandro Dalcin   PetscFunctionReturn(0);
13401566a47fSLisandro Dalcin }
13411566a47fSLisandro Dalcin 
13429371c9d4SSatish Balay static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer) {
1343f33bbcb6SJed Brown   PetscFunctionBegin;
1344f33bbcb6SJed Brown   PetscFunctionReturn(0);
1345f33bbcb6SJed Brown }
1346f33bbcb6SJed Brown 
1347f33bbcb6SJed Brown /*MC
1348f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1349f33bbcb6SJed Brown 
1350f33bbcb6SJed Brown   Level: beginner
1351f33bbcb6SJed Brown 
13524eb428fdSBarry Smith   Notes:
1353c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
13544eb428fdSBarry Smith 
13551566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
13564eb428fdSBarry Smith 
1357db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1358f33bbcb6SJed Brown 
1359f33bbcb6SJed Brown M*/
13609371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) {
1361f33bbcb6SJed Brown   PetscFunctionBegin;
13629566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
13639566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 1.0));
13649566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
13650c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1366f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1367f33bbcb6SJed Brown   PetscFunctionReturn(0);
1368f33bbcb6SJed Brown }
1369f33bbcb6SJed Brown 
13709371c9d4SSatish Balay static PetscErrorCode TSSetUp_CN(TS ts) {
13711566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
13721566a47fSLisandro Dalcin 
13731566a47fSLisandro Dalcin   PetscFunctionBegin;
137408401ef6SPierre 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");
13753c633725SBarry 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");
13769566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
13771566a47fSLisandro Dalcin   PetscFunctionReturn(0);
13781566a47fSLisandro Dalcin }
13791566a47fSLisandro Dalcin 
13809371c9d4SSatish Balay static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer) {
1381f33bbcb6SJed Brown   PetscFunctionBegin;
1382f33bbcb6SJed Brown   PetscFunctionReturn(0);
1383f33bbcb6SJed Brown }
1384f33bbcb6SJed Brown 
1385f33bbcb6SJed Brown /*MC
1386f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1387f33bbcb6SJed Brown 
1388f33bbcb6SJed Brown   Level: beginner
1389f33bbcb6SJed Brown 
1390f33bbcb6SJed Brown   Notes:
13917cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
13927cf5af47SJed Brown 
13937cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1394f33bbcb6SJed Brown 
1395db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`
1396f33bbcb6SJed Brown 
1397f33bbcb6SJed Brown M*/
13989371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) {
1399f33bbcb6SJed Brown   PetscFunctionBegin;
14009566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14019566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 0.5));
14029566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
14030c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1404f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1405f33bbcb6SJed Brown   PetscFunctionReturn(0);
1406f33bbcb6SJed Brown }
1407