xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision b44f4de476a10dd52a698c0f9a8f47d015bd1e19)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5a055c8caSBarry Smith #include <petscsnes.h>
61e25c274SJed Brown #include <petscdm.h>
73c54f92cSHong Zhang #include <petscmat.h>
8316643e7SJed Brown 
9316643e7SJed Brown typedef struct {
10715f1b00SHong Zhang   /* context for time stepping */
111566a47fSLisandro Dalcin   PetscReal    stage_time;
12cc4f23bcSHong Zhang   Vec          Stages[2];   /* Storage for stage solutions */
13cc4f23bcSHong Zhang   Vec          X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */
141566a47fSLisandro Dalcin   Vec          affine;      /* Affine vector needed for residual at beginning of step in endpoint formulation */
151566a47fSLisandro Dalcin   PetscReal    Theta;
161a352fa8SHong Zhang   PetscReal    shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20715f1b00SHong Zhang   TSStepStatus status;
211a352fa8SHong Zhang   Vec          VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */
221a352fa8SHong Zhang   PetscReal    ptime0;           /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
231a352fa8SHong Zhang   PetscReal    time_step0;       /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/
241566a47fSLisandro Dalcin 
25715f1b00SHong Zhang   /* context for sensitivity analysis */
26022c081aSHong Zhang   PetscInt num_tlm;               /* Total number of tangent linear equations */
27b5e0532dSHong Zhang   Vec     *VecsDeltaLam;          /* Increment of the adjoint sensitivity w.r.t IC at stage */
28b5e0532dSHong Zhang   Vec     *VecsDeltaMu;           /* Increment of the adjoint sensitivity w.r.t P at stage */
2913b1b0ffSHong Zhang   Vec     *VecsSensiTemp;         /* Vector to be multiplied with Jacobian transpose */
30cc4f23bcSHong Zhang   Mat      MatFwdStages[2];       /* TLM Stages */
3113b1b0ffSHong Zhang   Mat      MatDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
32b5b4017aSHong Zhang   Vec      VecDeltaFwdSensipCol;  /* Working vector for holding one column of the sensitivity matrix */
3313b1b0ffSHong Zhang   Mat      MatFwdSensip0;         /* backup for roll-backs due to events */
347207e0fdSHong Zhang   Mat      MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */
357207e0fdSHong Zhang   Mat      MatIntegralSensip0;    /* backup for roll-backs due to events */
36b5b4017aSHong Zhang   Vec     *VecsDeltaLam2;         /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
37b5b4017aSHong Zhang   Vec     *VecsDeltaMu2;          /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
38b5b4017aSHong Zhang   Vec     *VecsSensi2Temp;        /* Working vectors that holds the residual for the second-order adjoint */
39b5b4017aSHong Zhang   Vec     *VecsAffine;            /* Working vectors to store residuals */
40715f1b00SHong Zhang   /* context for error estimation */
411566a47fSLisandro Dalcin   Vec vec_sol_prev;
421566a47fSLisandro Dalcin   Vec vec_lte_work;
43316643e7SJed Brown } TS_Theta;
44316643e7SJed Brown 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
46d71ae5a4SJacob Faibussowitsch {
477445fe48SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
487445fe48SJed Brown 
497445fe48SJed Brown   PetscFunctionBegin;
507445fe48SJed Brown   if (X0) {
517445fe48SJed Brown     if (dm && dm != ts->dm) {
529566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0));
537445fe48SJed Brown     } else *X0 = ts->vec_sol;
547445fe48SJed Brown   }
557445fe48SJed Brown   if (Xdot) {
567445fe48SJed Brown     if (dm && dm != ts->dm) {
579566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
587445fe48SJed Brown     } else *Xdot = th->Xdot;
597445fe48SJed Brown   }
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
617445fe48SJed Brown }
627445fe48SJed Brown 
63d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
64d71ae5a4SJacob Faibussowitsch {
650d0b770aSPeter Brune   PetscFunctionBegin;
660d0b770aSPeter Brune   if (X0) {
6748a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0));
680d0b770aSPeter Brune   }
690d0b770aSPeter Brune   if (Xdot) {
7048a46eb9SPierre Jolivet     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
710d0b770aSPeter Brune   }
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
730d0b770aSPeter Brune }
740d0b770aSPeter Brune 
75d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx)
76d71ae5a4SJacob Faibussowitsch {
777445fe48SJed Brown   PetscFunctionBegin;
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
797445fe48SJed Brown }
807445fe48SJed Brown 
81d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
82d71ae5a4SJacob Faibussowitsch {
837445fe48SJed Brown   TS  ts = (TS)ctx;
847445fe48SJed Brown   Vec X0, Xdot, X0_c, Xdot_c;
857445fe48SJed Brown 
867445fe48SJed Brown   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot));
889566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
899566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, X0, X0_c));
909566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
919566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
929566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
939566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot));
949566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
967445fe48SJed Brown }
977445fe48SJed Brown 
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx)
99d71ae5a4SJacob Faibussowitsch {
100258e1594SPeter Brune   PetscFunctionBegin;
1013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102258e1594SPeter Brune }
103258e1594SPeter Brune 
104d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
105d71ae5a4SJacob Faibussowitsch {
106258e1594SPeter Brune   TS  ts = (TS)ctx;
107258e1594SPeter Brune   Vec X0, Xdot, X0_sub, Xdot_sub;
108258e1594SPeter Brune 
109258e1594SPeter Brune   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
1119566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
112258e1594SPeter Brune 
1139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
1149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
115258e1594SPeter Brune 
1169566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
1179566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
118258e1594SPeter Brune 
1199566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
1209566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122258e1594SPeter Brune }
123258e1594SPeter Brune 
124d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
125d71ae5a4SJacob Faibussowitsch {
126022c081aSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
127cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
128022c081aSHong Zhang 
129022c081aSHong Zhang   PetscFunctionBegin;
130022c081aSHong Zhang   if (th->endpoint) {
131022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
132022c081aSHong Zhang     if (th->Theta != 1.0) {
1339566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand));
1349566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand));
135022c081aSHong Zhang     }
1369566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand));
1379566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand));
138022c081aSHong Zhang   } else {
1399566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand));
1409566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand));
141022c081aSHong Zhang   }
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
143022c081aSHong Zhang }
144022c081aSHong Zhang 
145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
146d71ae5a4SJacob Faibussowitsch {
147b1cb36f3SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
148cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
149b1cb36f3SHong Zhang 
150b1cb36f3SHong Zhang   PetscFunctionBegin;
151b1cb36f3SHong Zhang   /* backup cost integral */
1529566063dSJacob Faibussowitsch   PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0));
1539566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155b1cb36f3SHong Zhang }
156b1cb36f3SHong Zhang 
157d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
158d71ae5a4SJacob Faibussowitsch {
1591a352fa8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
160b1cb36f3SHong Zhang 
161b1cb36f3SHong Zhang   PetscFunctionBegin;
1621a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
1631a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
1641a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
1659566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16724655328SShri }
16824655328SShri 
169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x)
170d71ae5a4SJacob Faibussowitsch {
1711566a47fSLisandro Dalcin   PetscInt nits, lits;
1721566a47fSLisandro Dalcin 
1731566a47fSLisandro Dalcin   PetscFunctionBegin;
1749566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
1759566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1769566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1779371c9d4SSatish Balay   ts->snes_its += nits;
1789371c9d4SSatish Balay   ts->ksp_its += lits;
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1801566a47fSLisandro Dalcin }
1811566a47fSLisandro Dalcin 
18226b5f05dSStefano Zampini static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg)
18326b5f05dSStefano Zampini {
18426b5f05dSStefano Zampini   TS_Theta *th = (TS_Theta *)ts->data;
18526b5f05dSStefano Zampini 
18626b5f05dSStefano Zampini   PetscFunctionBegin;
187ecc87898SStefano Zampini   if (reg) {
188ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev));
189ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0));
190ecc87898SStefano Zampini   } else {
191ecc87898SStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev));
192ecc87898SStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev));
193ecc87898SStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:X0", &th->X0));
19426b5f05dSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->X0));
19526b5f05dSStefano Zampini   }
19626b5f05dSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
19726b5f05dSStefano Zampini }
19826b5f05dSStefano Zampini 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Theta(TS ts)
200d71ae5a4SJacob Faibussowitsch {
201316643e7SJed Brown   TS_Theta *th         = (TS_Theta *)ts->data;
2021566a47fSLisandro Dalcin   PetscInt  rejections = 0;
2034957b756SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
2041566a47fSLisandro Dalcin   PetscReal next_time_step = ts->time_step;
205316643e7SJed Brown 
206316643e7SJed Brown   PetscFunctionBegin;
2071566a47fSLisandro Dalcin   if (!ts->steprollback) {
2089566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
2099566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
2101566a47fSLisandro Dalcin   }
2111566a47fSLisandro Dalcin 
2121566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
21399260352SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
2143e05ceb1SHong Zhang     th->shift      = 1 / (th->Theta * ts->time_step);
2151566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step;
2169566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X));
21748a46eb9SPierre Jolivet     if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot));
218eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
2199566063dSJacob Faibussowitsch       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
2209566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
2219566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE));
2229566063dSJacob Faibussowitsch       PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta));
223eb284becSJed Brown     }
2249566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
22582d7ce88SStefano Zampini     PetscCall(TSTheta_SNESSolve(ts, th->endpoint ? th->affine : NULL, th->X));
2269566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
2279566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
228fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
229051f2191SLisandro Dalcin 
230b7795fd7SStefano Zampini     if (th->endpoint) {
2319566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X, ts->vec_sol));
2321566a47fSLisandro Dalcin     } else {
233b7795fd7SStefano Zampini       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); /* th->Xdot is needed by TSInterpolate_Theta */
234b7795fd7SStefano Zampini       if (th->Theta == 1.0) PetscCall(VecCopy(th->X, ts->vec_sol));              /* BEULER, stage already checked */
235b7795fd7SStefano Zampini       else {
2369566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
237b7795fd7SStefano Zampini         PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, ts->vec_sol, &stageok));
238b7795fd7SStefano Zampini         if (!stageok) {
239b7795fd7SStefano Zampini           PetscCall(VecCopy(th->X0, ts->vec_sol));
240b7795fd7SStefano Zampini           goto reject_step;
2411566a47fSLisandro Dalcin         }
242b7795fd7SStefano Zampini       }
243b7795fd7SStefano Zampini     }
244b7795fd7SStefano Zampini 
245b7795fd7SStefano Zampini     th->status = TS_STEP_PENDING;
2469566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2471566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2481566a47fSLisandro Dalcin     if (!accept) {
2499566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
250be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
251be5899b3SLisandro Dalcin       goto reject_step;
252051f2191SLisandro Dalcin     }
2533b1890cdSShri Abhyankar 
254715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2551a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2561a352fa8SHong Zhang       th->time_step0 = ts->time_step;
25717145e75SHong Zhang     }
2582b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
259cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
260051f2191SLisandro Dalcin     break;
261051f2191SLisandro Dalcin 
262051f2191SLisandro Dalcin   reject_step:
2639371c9d4SSatish Balay     ts->reject++;
2649371c9d4SSatish Balay     accept = PETSC_FALSE;
2651566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
266051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
26763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
2683b1890cdSShri Abhyankar     }
2693b1890cdSShri Abhyankar   }
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
271316643e7SJed Brown }
272316643e7SJed Brown 
273d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
274d71ae5a4SJacob Faibussowitsch {
275c9681e5dSHong Zhang   TS_Theta      *th           = (TS_Theta *)ts->data;
276cd4cee2dSHong Zhang   TS             quadts       = ts->quadraturets;
277c9681e5dSHong Zhang   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
278c9681e5dSHong Zhang   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
279c9681e5dSHong Zhang   PetscInt       nadj;
2807207e0fdSHong Zhang   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
281c9681e5dSHong Zhang   KSP            ksp;
282c9681e5dSHong Zhang   PetscScalar   *xarr;
283077c4feaSHong Zhang   TSEquationType eqtype;
284077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2851a352fa8SHong Zhang   PetscReal      adjoint_time_step;
286c9681e5dSHong Zhang 
287c9681e5dSHong Zhang   PetscFunctionBegin;
2889566063dSJacob Faibussowitsch   PetscCall(TSGetEquationType(ts, &eqtype));
289077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
290077c4feaSHong Zhang     isexplicitode = PETSC_TRUE;
291077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
292077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
293077c4feaSHong Zhang   }
294c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
2959566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
2969566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
2977207e0fdSHong Zhang   if (quadts) {
2989566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
2999566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
3007207e0fdSHong Zhang   }
301c9681e5dSHong Zhang 
3021a352fa8SHong Zhang   th->stage_time    = ts->ptime;
3031a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
304c9681e5dSHong Zhang 
30533f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
3061baa6e33SBarry Smith   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
307cd4cee2dSHong Zhang 
308c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
3099566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
3109566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
311cd4cee2dSHong Zhang     if (quadJ) {
3129566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
3139566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
3149566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
3159566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
3169566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
317c9681e5dSHong Zhang     }
318c9681e5dSHong Zhang   }
319c9681e5dSHong Zhang 
320c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
321c87ba875SIvan Yashchuk   th->shift = 1. / adjoint_time_step;
3229566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3239566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
324c9681e5dSHong Zhang 
325c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
326c9681e5dSHong Zhang   for (nadj = 0; nadj < ts->numcost; nadj++) {
327c9681e5dSHong Zhang     KSPConvergedReason kspreason;
3289566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
3299566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
330c9681e5dSHong Zhang     if (kspreason < 0) {
331c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
33263a3b9bcSJacob 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));
333c9681e5dSHong Zhang     }
334c9681e5dSHong Zhang   }
335c9681e5dSHong Zhang 
336c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
337c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3389566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
3399566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
340c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3419566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
342c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3439566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
344c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
3459566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
3469566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
3479566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
34848a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
349c9681e5dSHong Zhang     }
350c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
351c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
35205c9950eSHong Zhang       KSPConvergedReason kspreason;
3539566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
3549566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
35505c9950eSHong Zhang       if (kspreason < 0) {
35605c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
35763a3b9bcSJacob 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));
35805c9950eSHong Zhang       }
359c9681e5dSHong Zhang     }
360c9681e5dSHong Zhang   }
361c9681e5dSHong Zhang 
362c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
363077c4feaSHong Zhang   if (!isexplicitode) {
3643e05ceb1SHong Zhang     th->shift = 0.0;
3659566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
3669566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
367c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
36833f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3699566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
3709566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
3719566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
3729566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
373c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3749566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
3759566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
3769566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
377c9681e5dSHong Zhang       }
378c9681e5dSHong Zhang     }
379077c4feaSHong Zhang   }
380c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3819566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
3829566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
3831baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
384c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
385c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3869566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
38733f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3889566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
389c9681e5dSHong Zhang     }
390cd4cee2dSHong Zhang 
391c9681e5dSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
3929566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
3939566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
394cd4cee2dSHong Zhang       if (quadJp) {
3959566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
3969566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
3979566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
3989566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3999566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
40033f72c5dSHong Zhang       }
401c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
4029566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
4039566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
40448a46eb9SPierre Jolivet         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
40548a46eb9SPierre Jolivet         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
406c9681e5dSHong Zhang       }
407c9681e5dSHong Zhang     }
408c9681e5dSHong Zhang   }
409c9681e5dSHong Zhang 
410c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
4119566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4129566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
413c9681e5dSHong Zhang   }
414c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
4153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
416c9681e5dSHong Zhang }
417c9681e5dSHong Zhang 
418d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointStep_Theta(TS ts)
419d71ae5a4SJacob Faibussowitsch {
4202ca6e920SHong Zhang   TS_Theta    *th           = (TS_Theta *)ts->data;
421cd4cee2dSHong Zhang   TS           quadts       = ts->quadraturets;
422b5e0532dSHong Zhang   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
423b5b4017aSHong Zhang   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
4242ca6e920SHong Zhang   PetscInt     nadj;
4257207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
4262ca6e920SHong Zhang   KSP          ksp;
427b5b4017aSHong Zhang   PetscScalar *xarr;
4281a352fa8SHong Zhang   PetscReal    adjoint_time_step;
429da81f932SPierre Jolivet   PetscReal    adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, usually ts->ptime is larger than adjoint_ptime) */
4302ca6e920SHong Zhang 
4312ca6e920SHong Zhang   PetscFunctionBegin;
432c9681e5dSHong Zhang   if (th->Theta == 1.) {
4339566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
4343ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
435c9681e5dSHong Zhang   }
4362ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4379566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
4389566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
4397207e0fdSHong Zhang   if (quadts) {
4409566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
4419566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
4427207e0fdSHong Zhang   }
443b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4441a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
4451a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4461a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
4472ca6e920SHong Zhang 
44882ad9101SHong Zhang   if (!th->endpoint) {
4495072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4509566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
45182ad9101SHong Zhang   }
45282ad9101SHong Zhang 
453b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
45433f72c5dSHong Zhang   /* Cost function has an integral term */
4557207e0fdSHong Zhang   if (quadts) {
45605755b9cSHong Zhang     if (th->endpoint) {
4579566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
45805755b9cSHong Zhang     } else {
4599566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
46005755b9cSHong Zhang     }
4617207e0fdSHong Zhang   }
46233f72c5dSHong Zhang 
463abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4649566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
4659566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
466cd4cee2dSHong Zhang     if (quadJ) {
4679566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
4689566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
4699566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
4709566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4719566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
47236eaed60SHong Zhang     }
4732ca6e920SHong Zhang   }
4743c54f92cSHong Zhang 
475b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4761a352fa8SHong Zhang   th->shift = 1. / (th->Theta * adjoint_time_step);
4773c54f92cSHong Zhang   if (th->endpoint) {
4789566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
4793c54f92cSHong Zhang   } else {
4809566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
4813c54f92cSHong Zhang   }
4829566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
4832ca6e920SHong Zhang 
484b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
485abc2977eSBarry Smith   for (nadj = 0; nadj < ts->numcost; nadj++) {
4861d14d03bSHong Zhang     KSPConvergedReason kspreason;
4879566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
4889566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
4891d14d03bSHong Zhang     if (kspreason < 0) {
4901d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
49163a3b9bcSJacob 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));
4921d14d03bSHong Zhang     }
4932ca6e920SHong Zhang   }
4943c54f92cSHong Zhang 
4951d14d03bSHong Zhang   /* Second-order adjoint */
496b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4973c633725SBarry Smith     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
498b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4999566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
5009566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
501b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
5029566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
5039566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
5049566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
505b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
5069566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
507b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
5089566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
5099566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
5109566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
51148a46eb9SPierre Jolivet       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
512b5b4017aSHong Zhang     }
513b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
514b5b4017aSHong Zhang     for (nadj = 0; nadj < ts->numcost; nadj++) {
5151d14d03bSHong Zhang       KSPConvergedReason kspreason;
5169566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
5179566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
5181d14d03bSHong Zhang       if (kspreason < 0) {
5191d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
52063a3b9bcSJacob 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));
5211d14d03bSHong Zhang       }
522b5b4017aSHong Zhang     }
523b5b4017aSHong Zhang   }
524b5b4017aSHong Zhang 
52536eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5261d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5271a352fa8SHong Zhang     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
5281a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
5299566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
5309566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
53133f72c5dSHong Zhang     /* R_U at t_n */
5321baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
533abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
5349566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
535cd4cee2dSHong Zhang       if (quadJ) {
5369566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
5379566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
5389566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
5399566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5409566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
541b5b4017aSHong Zhang       }
5429566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
543b5b4017aSHong Zhang     }
5441d14d03bSHong Zhang 
5451d14d03bSHong Zhang     /* Second-order adjoint */
5461d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
547b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5489566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
5499566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
550b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5519566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
5529566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5539566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
554b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5559566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
556b5b4017aSHong Zhang       for (nadj = 0; nadj < ts->numcost; nadj++) {
55733f72c5dSHong 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  */
5589566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
5599566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
56048a46eb9SPierre Jolivet         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
5619566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
562b5b4017aSHong Zhang       }
563b5e0532dSHong Zhang     }
5643fd52205SHong Zhang 
565c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
566c586f2e8SHong Zhang 
5673fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
568b5b4017aSHong Zhang       /* U_{n+1} */
56982ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
5709566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
5719566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
5721baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
573abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
5749566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
5759566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
5763b711c3fSHong Zhang         if (quadJp) {
5779566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
5789566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
5799566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
5809566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5819566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
5823b711c3fSHong Zhang         }
5833fd52205SHong Zhang       }
58433f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
58533f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5869566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
5879566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
588b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5899566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
5909566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5919566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
59233f72c5dSHong Zhang 
593b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5949566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
595b5b4017aSHong Zhang         for (nadj = 0; nadj < ts->numcost; nadj++) {
59633f72c5dSHong 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)  */
5979566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
5989566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
59948a46eb9SPierre Jolivet           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
60048a46eb9SPierre Jolivet           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
601b5b4017aSHong Zhang         }
602b5b4017aSHong Zhang       }
603b5b4017aSHong Zhang 
604b5b4017aSHong Zhang       /* U_s */
6059566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
6069566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
6071baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
608abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6099566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6109566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
6113b711c3fSHong Zhang         if (quadJp) {
6129566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6139566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6149566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
6159566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6169566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
6173b711c3fSHong Zhang         }
61833f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
61933f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
6209566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
6219566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
622b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6239566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
6249566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
6259566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
626b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6279566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
628b5b4017aSHong Zhang           for (nadj = 0; nadj < ts->numcost; nadj++) {
62933f72c5dSHong 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) */
6309566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
6319566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
63248a46eb9SPierre Jolivet             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
63348a46eb9SPierre Jolivet             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
63436eaed60SHong Zhang           }
63536eaed60SHong Zhang         }
6363fd52205SHong Zhang       }
637b5e0532dSHong Zhang     }
6383fd52205SHong Zhang   } else { /* one-stage case */
6393e05ceb1SHong Zhang     th->shift = 0.0;
6409566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
6419566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, J, Jpre));
6421baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
643abc2977eSBarry Smith     for (nadj = 0; nadj < ts->numcost; nadj++) {
6449566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
6459566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
646cd4cee2dSHong Zhang       if (quadJ) {
6479566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
6489566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
6499566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
6509566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6519566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
65236eaed60SHong Zhang       }
6532ca6e920SHong Zhang     }
6543fd52205SHong Zhang     if (ts->vecs_sensip) {
65582ad9101SHong Zhang       th->shift = 1.0 / (adjoint_time_step * th->Theta);
6569566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
6579566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
6581baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
659abc2977eSBarry Smith       for (nadj = 0; nadj < ts->numcost; nadj++) {
6609566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
6619566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
662cd4cee2dSHong Zhang         if (quadJp) {
6639566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
6649566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
6659566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
6669566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6679566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
66836eaed60SHong Zhang         }
66936eaed60SHong Zhang       }
6703fd52205SHong Zhang     }
6712ca6e920SHong Zhang   }
6722ca6e920SHong Zhang 
6732ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6752ca6e920SHong Zhang }
6762ca6e920SHong Zhang 
677d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X)
678d71ae5a4SJacob Faibussowitsch {
679cd652676SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
680be5899b3SLisandro Dalcin   PetscReal dt = t - ts->ptime;
681cd652676SJed Brown 
682cd652676SJed Brown   PetscFunctionBegin;
6839566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X));
684be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6859566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
6863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
687cd652676SJed Brown }
688cd652676SJed Brown 
689d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
690d71ae5a4SJacob Faibussowitsch {
6911566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
6921566a47fSLisandro Dalcin   Vec       X  = ts->vec_sol;      /* X = solution */
6931566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
6947453f775SEmil Constantinescu   PetscReal wltea, wlter;
6951566a47fSLisandro Dalcin 
6961566a47fSLisandro Dalcin   PetscFunctionBegin;
6979371c9d4SSatish Balay   if (!th->vec_sol_prev) {
6989371c9d4SSatish Balay     *wlte = -1;
6993ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
7009371c9d4SSatish Balay   }
7011566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
7029371c9d4SSatish Balay   if (ts->steprestart) {
7039371c9d4SSatish Balay     *wlte = -1;
7043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
7059371c9d4SSatish Balay   }
7061566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
707fecfb714SLisandro Dalcin   {
708be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
709be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
7109371c9d4SSatish Balay     PetscScalar scal[3];
7119371c9d4SSatish Balay     Vec         vecs[3];
712225f7d0bSStefano Zampini     scal[0] = -1 / a;
713225f7d0bSStefano Zampini     scal[1] = +1 / (a - 1);
714225f7d0bSStefano Zampini     scal[2] = -1 / (a * (a - 1));
7159371c9d4SSatish Balay     vecs[0] = X;
7169371c9d4SSatish Balay     vecs[1] = th->X0;
7179371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
7189566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
7199566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
7209566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
7211566a47fSLisandro Dalcin   }
7221566a47fSLisandro Dalcin   if (order) *order = 2;
7233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7241566a47fSLisandro Dalcin }
7251566a47fSLisandro Dalcin 
726d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Theta(TS ts)
727d71ae5a4SJacob Faibussowitsch {
7281566a47fSLisandro Dalcin   TS_Theta *th     = (TS_Theta *)ts->data;
7297207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
7301566a47fSLisandro Dalcin 
7311566a47fSLisandro Dalcin   PetscFunctionBegin;
73248a46eb9SPierre Jolivet   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
733715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
7341baa6e33SBarry Smith   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
73548a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
7363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
737715f1b00SHong Zhang }
738715f1b00SHong Zhang 
739d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardStep_Theta(TS ts)
740d71ae5a4SJacob Faibussowitsch {
741715f1b00SHong Zhang   TS_Theta    *th                   = (TS_Theta *)ts->data;
742cd4cee2dSHong Zhang   TS           quadts               = ts->quadraturets;
74313b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
744b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7457207e0fdSHong Zhang   PetscInt     ntlm;
746715f1b00SHong Zhang   KSP          ksp;
7477207e0fdSHong Zhang   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
74813b1b0ffSHong Zhang   PetscScalar *barr, *xarr;
7491a352fa8SHong Zhang   PetscReal    previous_shift;
750715f1b00SHong Zhang 
751715f1b00SHong Zhang   PetscFunctionBegin;
7521a352fa8SHong Zhang   previous_shift = th->shift;
7539566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));
75413b1b0ffSHong Zhang 
75548a46eb9SPierre Jolivet   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
7569566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes, &ksp));
7579566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
7587207e0fdSHong Zhang   if (quadts) {
7599566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
7609566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
7617207e0fdSHong Zhang   }
762715f1b00SHong Zhang 
763715f1b00SHong Zhang   /* Build RHS */
764715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7651a352fa8SHong Zhang     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
7669566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
767fb842aefSJose E. Roman     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip));
7689566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));
769715f1b00SHong Zhang 
770022c081aSHong Zhang     /* Add the f_p forcing terms */
7710e11c21fSHong Zhang     if (ts->Jacp) {
7729566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7739566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7749566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
77582ad9101SHong Zhang       th->shift = previous_shift;
7769566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
7779566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7789566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
7790e11c21fSHong Zhang     }
780715f1b00SHong Zhang   } else { /* 1-stage method */
7813e05ceb1SHong Zhang     th->shift = 0.0;
7829566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
783fb842aefSJose E. Roman     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip));
7849566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip, -1.));
78513b1b0ffSHong Zhang 
786022c081aSHong Zhang     /* Add the f_p forcing terms */
7870e11c21fSHong Zhang     if (ts->Jacp) {
78882ad9101SHong Zhang       th->shift = previous_shift;
7899566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
7909566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
7919566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
792715f1b00SHong Zhang     }
7930e11c21fSHong Zhang   }
794715f1b00SHong Zhang 
795715f1b00SHong Zhang   /* Build LHS */
7961a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
797715f1b00SHong Zhang   if (th->endpoint) {
7989566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
799715f1b00SHong Zhang   } else {
8009566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
801715f1b00SHong Zhang   }
8029566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, J, Jpre));
803715f1b00SHong Zhang 
804715f1b00SHong Zhang   /*
805715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
806c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
807715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
808715f1b00SHong Zhang   */
809715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8107207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8119566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
8129566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
813fb842aefSJose E. Roman       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
8149566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8159566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
816715f1b00SHong Zhang     }
817715f1b00SHong Zhang   }
818715f1b00SHong Zhang 
819715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
820022c081aSHong Zhang   for (ntlm = 0; ntlm < th->num_tlm; ntlm++) {
821b5b4017aSHong Zhang     KSPConvergedReason kspreason;
8229566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
8239566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
824715f1b00SHong Zhang     if (th->endpoint) {
8259566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
8269566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
8279566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
8289566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
8299566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
830715f1b00SHong Zhang     } else {
8319566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
832715f1b00SHong Zhang     }
8339566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
834b5b4017aSHong Zhang     if (kspreason < 0) {
835b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
83663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
837b5b4017aSHong Zhang     }
8389566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
8399566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
840715f1b00SHong Zhang   }
841715f1b00SHong Zhang 
84213b1b0ffSHong Zhang   /*
84313b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
844c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
84513b1b0ffSHong Zhang   */
8467207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
84713b1b0ffSHong Zhang     if (!th->endpoint) {
8489566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
8499566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
8509566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
851fb842aefSJose E. Roman       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
8529566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8539566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
8549566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
85513b1b0ffSHong Zhang     } else {
8569566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
8579566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
858fb842aefSJose E. Roman       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
8599566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
8609566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
86113b1b0ffSHong Zhang     }
86213b1b0ffSHong Zhang   } else {
86348a46eb9SPierre Jolivet     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
864715f1b00SHong Zhang   }
8653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8661566a47fSLisandro Dalcin }
8671566a47fSLisandro Dalcin 
868d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[])
869d71ae5a4SJacob Faibussowitsch {
8706affb6f8SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
8716affb6f8SHong Zhang 
8726affb6f8SHong Zhang   PetscFunctionBegin;
8737409eceaSHong Zhang   if (ns) {
8747409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8757409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
8767409eceaSHong Zhang   }
8775072d23cSHong Zhang   if (stagesensip) {
878cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8797409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
880cc4f23bcSHong Zhang     } else {
881cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
882cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
883cc4f23bcSHong Zhang     }
884cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8855072d23cSHong Zhang   }
8863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8876affb6f8SHong Zhang }
8886affb6f8SHong Zhang 
889316643e7SJed Brown /*------------------------------------------------------------*/
890d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Theta(TS ts)
891d71ae5a4SJacob Faibussowitsch {
892316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
893316643e7SJed Brown 
894316643e7SJed Brown   PetscFunctionBegin;
8959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
8969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
8979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
8989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
8991566a47fSLisandro Dalcin 
9009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
9019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
9021566a47fSLisandro Dalcin 
9039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
9043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
905277b19d0SLisandro Dalcin }
906277b19d0SLisandro Dalcin 
907d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointReset_Theta(TS ts)
908d71ae5a4SJacob Faibussowitsch {
909ece46509SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
910ece46509SHong Zhang 
911ece46509SHong Zhang   PetscFunctionBegin;
9129566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
9139566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
9149566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
9159566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
9169566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
9179566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
9183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
919ece46509SHong Zhang }
920ece46509SHong Zhang 
921d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Theta(TS ts)
922d71ae5a4SJacob Faibussowitsch {
923277b19d0SLisandro Dalcin   PetscFunctionBegin;
9249566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
925b3a6b972SJed Brown   if (ts->dm) {
9269566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
9279566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
928b3a6b972SJed Brown   }
9299566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
9319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
9329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
9339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935316643e7SJed Brown }
936316643e7SJed Brown 
937316643e7SJed Brown /*
938316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9392b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9400fd10508SMatthew Knepley 
9410fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9420fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9430fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
944316643e7SJed Brown */
945d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts)
946d71ae5a4SJacob Faibussowitsch {
947316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9487445fe48SJed Brown   Vec       X0, Xdot;
9497445fe48SJed Brown   DM        dm, dmsave;
9503e05ceb1SHong Zhang   PetscReal shift = th->shift;
951316643e7SJed Brown 
952316643e7SJed Brown   PetscFunctionBegin;
9539566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
9545a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9559566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
956494ce9fcSHong Zhang   if (x != X0) {
9579566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
958494ce9fcSHong Zhang   } else {
9599566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
960494ce9fcSHong Zhang   }
9617445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9627445fe48SJed Brown   dmsave = ts->dm;
9637445fe48SJed Brown   ts->dm = dm;
9649566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
9657445fe48SJed Brown   ts->dm = dmsave;
9669566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
9673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
968316643e7SJed Brown }
969316643e7SJed Brown 
970d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts)
971d71ae5a4SJacob Faibussowitsch {
972316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
9737445fe48SJed Brown   Vec       Xdot;
9747445fe48SJed Brown   DM        dm, dmsave;
9753e05ceb1SHong Zhang   PetscReal shift = th->shift;
976316643e7SJed Brown 
977316643e7SJed Brown   PetscFunctionBegin;
9789566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
979be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9809566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));
9817445fe48SJed Brown 
9827445fe48SJed Brown   dmsave = ts->dm;
9837445fe48SJed Brown   ts->dm = dm;
9849566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
9857445fe48SJed Brown   ts->dm = dmsave;
9869566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
9873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
988316643e7SJed Brown }
989316643e7SJed Brown 
990d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardSetUp_Theta(TS ts)
991d71ae5a4SJacob Faibussowitsch {
992715f1b00SHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
9937207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
994715f1b00SHong Zhang 
995715f1b00SHong Zhang   PetscFunctionBegin;
996022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
99713b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
9989566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
9997207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10009566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
10019566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
1002715f1b00SHong Zhang   }
10036f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
10049566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));
100513b1b0ffSHong Zhang 
10069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
10073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1008715f1b00SHong Zhang }
1009715f1b00SHong Zhang 
1010d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSForwardReset_Theta(TS ts)
1011d71ae5a4SJacob Faibussowitsch {
10127adebddeSHong Zhang   TS_Theta *th     = (TS_Theta *)ts->data;
10137207e0fdSHong Zhang   TS        quadts = ts->quadraturets;
10147adebddeSHong Zhang 
10157adebddeSHong Zhang   PetscFunctionBegin;
10167207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10179566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
10189566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
10197adebddeSHong Zhang   }
10209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
10219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
10229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
10233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10247adebddeSHong Zhang }
10257adebddeSHong Zhang 
1026d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Theta(TS ts)
1027d71ae5a4SJacob Faibussowitsch {
1028316643e7SJed Brown   TS_Theta *th     = (TS_Theta *)ts->data;
1029cd4cee2dSHong Zhang   TS        quadts = ts->quadraturets;
10302ffb9264SLisandro Dalcin   PetscBool match;
1031316643e7SJed Brown 
1032316643e7SJed Brown   PetscFunctionBegin;
1033cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
10349566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
1035b1cb36f3SHong Zhang   }
103648a46eb9SPierre Jolivet   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
103748a46eb9SPierre Jolivet   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
103848a46eb9SPierre Jolivet   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
103948a46eb9SPierre Jolivet   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
10403b1890cdSShri Abhyankar 
10411566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
104220d49056SMatthew G. Knepley   th->shift = 1 / (th->Theta * ts->time_step);
10431566a47fSLisandro Dalcin 
10449566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &ts->dm));
10459566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
10469566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
10471566a47fSLisandro Dalcin 
10489566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
10499566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
10512ffb9264SLisandro Dalcin   if (!match) {
1052ecc87898SStefano Zampini     if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
1053ecc87898SStefano Zampini     if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
10543b1890cdSShri Abhyankar   }
10559566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
1056cc4f23bcSHong Zhang 
1057cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
10583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1059b4dd3bf9SBarry Smith }
10600c7376e5SHong Zhang 
1061b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1062b4dd3bf9SBarry Smith 
1063d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1064d71ae5a4SJacob Faibussowitsch {
1065b4dd3bf9SBarry Smith   TS_Theta *th = (TS_Theta *)ts->data;
1066b4dd3bf9SBarry Smith 
1067b4dd3bf9SBarry Smith   PetscFunctionBegin;
10689566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
10699566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
107048a46eb9SPierre Jolivet   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1071b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10729566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
10739566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
107467633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
107567633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
107667633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1077b5b4017aSHong Zhang   }
1078c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
10799566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
108067633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
108167633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
108267633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1083b5b4017aSHong Zhang   }
10843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1085316643e7SJed Brown }
1086316643e7SJed Brown 
1087ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems PetscOptionsObject)
1088d71ae5a4SJacob Faibussowitsch {
1089316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1090316643e7SJed Brown 
1091316643e7SJed Brown   PetscFunctionBegin;
1092d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1093316643e7SJed Brown   {
10949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
10959566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
10969566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1097316643e7SJed Brown   }
1098d0609cedSBarry Smith   PetscOptionsHeadEnd();
10993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1100316643e7SJed Brown }
1101316643e7SJed Brown 
1102d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer)
1103d71ae5a4SJacob Faibussowitsch {
1104316643e7SJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1105ace3abfcSBarry Smith   PetscBool iascii;
1106316643e7SJed Brown 
1107316643e7SJed Brown   PetscFunctionBegin;
11089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1109316643e7SJed Brown   if (iascii) {
11109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
11119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1112316643e7SJed Brown   }
11133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1114316643e7SJed Brown }
1115316643e7SJed Brown 
1116d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta)
1117d71ae5a4SJacob Faibussowitsch {
11180de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11190de4c49aSJed Brown 
11200de4c49aSJed Brown   PetscFunctionBegin;
11210de4c49aSJed Brown   *theta = th->Theta;
11223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11230de4c49aSJed Brown }
11240de4c49aSJed Brown 
1125d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta)
1126d71ae5a4SJacob Faibussowitsch {
11270de4c49aSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
11280de4c49aSJed Brown 
11290de4c49aSJed Brown   PetscFunctionBegin;
1130cad9d221SBarry Smith   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
11310de4c49aSJed Brown   th->Theta = theta;
11321566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11340de4c49aSJed Brown }
1135eb284becSJed Brown 
1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint)
1137d71ae5a4SJacob Faibussowitsch {
113826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
113926f2ff8fSLisandro Dalcin 
114026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
114126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
11423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114326f2ff8fSLisandro Dalcin }
114426f2ff8fSLisandro Dalcin 
1145d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg)
1146d71ae5a4SJacob Faibussowitsch {
1147eb284becSJed Brown   TS_Theta *th = (TS_Theta *)ts->data;
1148eb284becSJed Brown 
1149eb284becSJed Brown   PetscFunctionBegin;
1150eb284becSJed Brown   th->endpoint = flg;
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152eb284becSJed Brown }
11530de4c49aSJed Brown 
1154f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
1156d71ae5a4SJacob Faibussowitsch {
1157f9c1d6abSBarry Smith   PetscComplex z  = xr + xi * PETSC_i, f;
1158f9c1d6abSBarry Smith   TS_Theta    *th = (TS_Theta *)ts->data;
1159f9c1d6abSBarry Smith 
1160f9c1d6abSBarry Smith   PetscFunctionBegin;
1161520136c7SJose E. Roman   f   = (1.0 + (1.0 - th->Theta) * z) / (1.0 - th->Theta * z);
1162f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1163f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
11643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1165f9c1d6abSBarry Smith }
1166f9c1d6abSBarry Smith #endif
1167f9c1d6abSBarry Smith 
1168d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[])
1169d71ae5a4SJacob Faibussowitsch {
117042682096SHong Zhang   TS_Theta *th = (TS_Theta *)ts->data;
117142682096SHong Zhang 
117242682096SHong Zhang   PetscFunctionBegin;
11737409eceaSHong Zhang   if (ns) {
11747409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11757409eceaSHong Zhang     else *ns = 2;                                   /* endpoint form */
11767409eceaSHong Zhang   }
11775072d23cSHong Zhang   if (Y) {
1178cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
11797409eceaSHong Zhang       th->Stages[0] = th->X;
1180cc4f23bcSHong Zhang     } else {
1181cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1182cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1183cc4f23bcSHong Zhang     }
1184cc4f23bcSHong Zhang     *Y = th->Stages;
11855072d23cSHong Zhang   }
11863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118742682096SHong Zhang }
1188f9c1d6abSBarry Smith 
1189316643e7SJed Brown /* ------------------------------------------------------------ */
1190316643e7SJed Brown /*MC
119196f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1192316643e7SJed Brown 
1193316643e7SJed Brown    Level: beginner
1194316643e7SJed Brown 
1195bcf0153eSBarry Smith    Options Database Keys:
1196c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1197c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
119803842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11994eb428fdSBarry Smith 
1200eb284becSJed Brown    Notes:
120137fdd005SBarry Smith .vb
120237fdd005SBarry Smith   -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
120337fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
120437fdd005SBarry Smith   -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
120537fdd005SBarry Smith .ve
12064eb428fdSBarry Smith 
12077409eceaSHong 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.
1208eb284becSJed Brown 
12097409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1210eb284becSJed Brown 
1211eb284becSJed Brown .vb
1212eb284becSJed Brown   Theta | Theta
1213eb284becSJed Brown   -------------
1214eb284becSJed Brown         |  1
1215eb284becSJed Brown .ve
1216eb284becSJed Brown 
1217eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1218eb284becSJed Brown 
1219eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1220eb284becSJed Brown 
1221eb284becSJed Brown .vb
1222eb284becSJed Brown   0 | 0         0
1223eb284becSJed Brown   1 | 1-Theta   Theta
1224eb284becSJed Brown   -------------------
1225eb284becSJed Brown     | 1-Theta   Theta
1226eb284becSJed Brown .ve
1227eb284becSJed Brown 
1228eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1229eb284becSJed Brown 
1230eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1231*b44f4de4SBarry Smith .vb
1232*b44f4de4SBarry Smith   Y_i = X + h sum_j a_ij Y'_j
1233*b44f4de4SBarry Smith .ve
12344eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1235eb284becSJed Brown 
12361cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1237316643e7SJed Brown M*/
1238d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1239d71ae5a4SJacob Faibussowitsch {
1240316643e7SJed Brown   TS_Theta *th;
1241316643e7SJed Brown 
1242316643e7SJed Brown   PetscFunctionBegin;
1243277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
1244ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1245316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
1246316643e7SJed Brown   ts->ops->view           = TSView_Theta;
1247316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
124842f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1249ece46509SHong Zhang   ts->ops->adjointreset   = TSAdjointReset_Theta;
1250316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
1251cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
12521566a47fSLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
125324655328SShri   ts->ops->rollback       = TSRollBack_Theta;
125426b5f05dSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Theta;
1255316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
12560f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
12570f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1258f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1259f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1260f9c1d6abSBarry Smith #endif
126142682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
126242f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1263b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1264b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12652ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12666affb6f8SHong Zhang 
1267715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12687adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1269715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12706affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1271316643e7SJed Brown 
1272efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1273efd4aadfSBarry Smith 
12744dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
1275316643e7SJed Brown   ts->data = (void *)th;
1276316643e7SJed Brown 
1277715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
1278715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
1279715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
1280b5b4017aSHong Zhang   th->VecsSensi2Temp = NULL;
1281715f1b00SHong Zhang 
12826f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1283316643e7SJed Brown   th->Theta       = 0.5;
12841566a47fSLisandro Dalcin   th->order       = 2;
12859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
12869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
12879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
12889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
12893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1290316643e7SJed Brown }
12910de4c49aSJed Brown 
12920de4c49aSJed Brown /*@
1293bcf0153eSBarry Smith   TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA`
12940de4c49aSJed Brown 
12950de4c49aSJed Brown   Not Collective
12960de4c49aSJed Brown 
12970de4c49aSJed Brown   Input Parameter:
12980de4c49aSJed Brown . ts - timestepping context
12990de4c49aSJed Brown 
13000de4c49aSJed Brown   Output Parameter:
13010de4c49aSJed Brown . theta - stage abscissa
13020de4c49aSJed Brown 
1303bcf0153eSBarry Smith   Level: advanced
1304bcf0153eSBarry Smith 
13050de4c49aSJed Brown   Note:
1306bcf0153eSBarry Smith   Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme.
13070de4c49aSJed Brown 
13081cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA`
13090de4c49aSJed Brown @*/
1310d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta)
1311d71ae5a4SJacob Faibussowitsch {
13120de4c49aSJed Brown   PetscFunctionBegin;
1313afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
13144f572ea9SToby Isaac   PetscAssertPointer(theta, 2);
1315cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
13163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13170de4c49aSJed Brown }
13180de4c49aSJed Brown 
13190de4c49aSJed Brown /*@
1320bcf0153eSBarry Smith   TSThetaSetTheta - Set the abscissa of the stage in (0,1]  for `TSTHETA`
13210de4c49aSJed Brown 
13220de4c49aSJed Brown   Not Collective
13230de4c49aSJed Brown 
1324d8d19677SJose E. Roman   Input Parameters:
13250de4c49aSJed Brown + ts    - timestepping context
13260de4c49aSJed Brown - theta - stage abscissa
13270de4c49aSJed Brown 
1328bcf0153eSBarry Smith   Options Database Key:
132967b8a455SSatish Balay . -ts_theta_theta <theta> - set theta
13300de4c49aSJed Brown 
1331bcf0153eSBarry Smith   Level: intermediate
13320de4c49aSJed Brown 
13331cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN`
13340de4c49aSJed Brown @*/
1335d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta)
1336d71ae5a4SJacob Faibussowitsch {
13370de4c49aSJed Brown   PetscFunctionBegin;
1338afb20b64SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1339cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
13403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13410de4c49aSJed Brown }
1342f33bbcb6SJed Brown 
134326f2ff8fSLisandro Dalcin /*@
1344bcf0153eSBarry Smith   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
134526f2ff8fSLisandro Dalcin 
134626f2ff8fSLisandro Dalcin   Not Collective
134726f2ff8fSLisandro Dalcin 
134826f2ff8fSLisandro Dalcin   Input Parameter:
134926f2ff8fSLisandro Dalcin . ts - timestepping context
135026f2ff8fSLisandro Dalcin 
135126f2ff8fSLisandro Dalcin   Output Parameter:
1352bcf0153eSBarry Smith . endpoint - `PETSC_TRUE` when using the endpoint variant
135326f2ff8fSLisandro Dalcin 
1354bcf0153eSBarry Smith   Level: advanced
135526f2ff8fSLisandro Dalcin 
13561cc06b55SBarry Smith .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
135726f2ff8fSLisandro Dalcin @*/
1358d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint)
1359d71ae5a4SJacob Faibussowitsch {
136026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
136126f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
13624f572ea9SToby Isaac   PetscAssertPointer(endpoint, 2);
1363cac4c232SBarry Smith   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
13643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136526f2ff8fSLisandro Dalcin }
136626f2ff8fSLisandro Dalcin 
1367eb284becSJed Brown /*@
1368bcf0153eSBarry Smith   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`
1369eb284becSJed Brown 
1370eb284becSJed Brown   Not Collective
1371eb284becSJed Brown 
1372d8d19677SJose E. Roman   Input Parameters:
1373eb284becSJed Brown + ts  - timestepping context
1374bcf0153eSBarry Smith - flg - `PETSC_TRUE` to use the endpoint variant
1375eb284becSJed Brown 
1376bcf0153eSBarry Smith   Options Database Key:
137767b8a455SSatish Balay . -ts_theta_endpoint <flg> - use the endpoint variant
1378eb284becSJed Brown 
1379bcf0153eSBarry Smith   Level: intermediate
1380eb284becSJed Brown 
13811cc06b55SBarry Smith .seealso: [](ch_ts), `TSTHETA`, `TSCN`
1382eb284becSJed Brown @*/
1383d71ae5a4SJacob Faibussowitsch PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg)
1384d71ae5a4SJacob Faibussowitsch {
1385eb284becSJed Brown   PetscFunctionBegin;
1386eb284becSJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1387cac4c232SBarry Smith   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
13883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1389eb284becSJed Brown }
1390eb284becSJed Brown 
1391f33bbcb6SJed Brown /*
1392f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1393f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1394f33bbcb6SJed Brown  */
1395f33bbcb6SJed Brown 
1396d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BEuler(TS ts)
1397d71ae5a4SJacob Faibussowitsch {
13981566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
13991566a47fSLisandro Dalcin 
14001566a47fSLisandro Dalcin   PetscFunctionBegin;
140108401ef6SPierre 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");
140228b400f6SJacob 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");
14039566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14051566a47fSLisandro Dalcin }
14061566a47fSLisandro Dalcin 
1407d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer)
1408d71ae5a4SJacob Faibussowitsch {
1409f33bbcb6SJed Brown   PetscFunctionBegin;
14103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1411f33bbcb6SJed Brown }
1412f33bbcb6SJed Brown 
1413f33bbcb6SJed Brown /*MC
1414f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1415f33bbcb6SJed Brown 
1416f33bbcb6SJed Brown   Level: beginner
1417f33bbcb6SJed Brown 
1418bcf0153eSBarry Smith   Note:
141937fdd005SBarry Smith   `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0`
14204eb428fdSBarry Smith 
14211cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1422f33bbcb6SJed Brown M*/
1423d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1424d71ae5a4SJacob Faibussowitsch {
1425f33bbcb6SJed Brown   PetscFunctionBegin;
14269566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14279566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 1.0));
14289566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
14290c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1430f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
14313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1432f33bbcb6SJed Brown }
1433f33bbcb6SJed Brown 
1434d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_CN(TS ts)
1435d71ae5a4SJacob Faibussowitsch {
14361566a47fSLisandro Dalcin   TS_Theta *th = (TS_Theta *)ts->data;
14371566a47fSLisandro Dalcin 
14381566a47fSLisandro Dalcin   PetscFunctionBegin;
143908401ef6SPierre 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");
14403c633725SBarry 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");
14419566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14431566a47fSLisandro Dalcin }
14441566a47fSLisandro Dalcin 
1445d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer)
1446d71ae5a4SJacob Faibussowitsch {
1447f33bbcb6SJed Brown   PetscFunctionBegin;
14483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1449f33bbcb6SJed Brown }
1450f33bbcb6SJed Brown 
1451f33bbcb6SJed Brown /*MC
1452f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1453f33bbcb6SJed Brown 
1454f33bbcb6SJed Brown   Level: beginner
1455f33bbcb6SJed Brown 
1456f33bbcb6SJed Brown   Notes:
1457bcf0153eSBarry Smith   `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e.
1458bcf0153eSBarry Smith .vb
1459bcf0153eSBarry Smith   -ts_type theta
1460bcf0153eSBarry Smith   -ts_theta_theta 0.5
1461bcf0153eSBarry Smith   -ts_theta_endpoint
1462bcf0153eSBarry Smith .ve
14637cf5af47SJed Brown 
14641cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`,
1465f33bbcb6SJed Brown M*/
1466d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1467d71ae5a4SJacob Faibussowitsch {
1468f33bbcb6SJed Brown   PetscFunctionBegin;
14699566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14709566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts, 0.5));
14719566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
14720c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1473f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
14743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1475f33bbcb6SJed Brown }
1476