xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 /*
2   Code for timestepping with implicit Theta method
3 */
4 #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5 #include <petscsnes.h>
6 #include <petscdm.h>
7 #include <petscmat.h>
8 
9 typedef struct {
10   /* context for time stepping */
11   PetscReal    stage_time;
12   Vec          Stages[2];   /* Storage for stage solutions */
13   Vec          X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */
14   Vec          affine;      /* Affine vector needed for residual at beginning of step in endpoint formulation */
15   PetscReal    Theta;
16   PetscReal    shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
17   PetscInt     order;
18   PetscBool    endpoint;
19   PetscBool    extrapolate;
20   TSStepStatus status;
21   Vec          VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */
22   PetscReal    ptime0;           /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
23   PetscReal    time_step0;       /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/
24 
25   /* context for sensitivity analysis */
26   PetscInt num_tlm;               /* Total number of tangent linear equations */
27   Vec     *VecsDeltaLam;          /* Increment of the adjoint sensitivity w.r.t IC at stage */
28   Vec     *VecsDeltaMu;           /* Increment of the adjoint sensitivity w.r.t P at stage */
29   Vec     *VecsSensiTemp;         /* Vector to be multiplied with Jacobian transpose */
30   Mat      MatFwdStages[2];       /* TLM Stages */
31   Mat      MatDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
32   Vec      VecDeltaFwdSensipCol;  /* Working vector for holding one column of the sensitivity matrix */
33   Mat      MatFwdSensip0;         /* backup for roll-backs due to events */
34   Mat      MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */
35   Mat      MatIntegralSensip0;    /* backup for roll-backs due to events */
36   Vec     *VecsDeltaLam2;         /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
37   Vec     *VecsDeltaMu2;          /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
38   Vec     *VecsSensi2Temp;        /* Working vectors that holds the residual for the second-order adjoint */
39   Vec     *VecsAffine;            /* Working vectors to store residuals */
40   /* context for error estimation */
41   Vec      vec_sol_prev;
42   Vec      vec_lte_work;
43 } TS_Theta;
44 
45 static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
46   TS_Theta *th = (TS_Theta *)ts->data;
47 
48   PetscFunctionBegin;
49   if (X0) {
50     if (dm && dm != ts->dm) {
51       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0));
52     } else *X0 = ts->vec_sol;
53   }
54   if (Xdot) {
55     if (dm && dm != ts->dm) {
56       PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
57     } else *Xdot = th->Xdot;
58   }
59   PetscFunctionReturn(0);
60 }
61 
62 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
63   PetscFunctionBegin;
64   if (X0) {
65     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0));
66   }
67   if (Xdot) {
68     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
69   }
70   PetscFunctionReturn(0);
71 }
72 
73 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, void *ctx) {
74   PetscFunctionBegin;
75   PetscFunctionReturn(0);
76 }
77 
78 static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) {
79   TS  ts = (TS)ctx;
80   Vec X0, Xdot, X0_c, Xdot_c;
81 
82   PetscFunctionBegin;
83   PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot));
84   PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
85   PetscCall(MatRestrict(restrct, X0, X0_c));
86   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
87   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
88   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
89   PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot));
90   PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
91   PetscFunctionReturn(0);
92 }
93 
94 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, void *ctx) {
95   PetscFunctionBegin;
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
100   TS  ts = (TS)ctx;
101   Vec X0, Xdot, X0_sub, Xdot_sub;
102 
103   PetscFunctionBegin;
104   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
105   PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
106 
107   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
108   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
109 
110   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
111   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
112 
113   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
114   PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
115   PetscFunctionReturn(0);
116 }
117 
118 static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) {
119   TS_Theta *th     = (TS_Theta *)ts->data;
120   TS        quadts = ts->quadraturets;
121 
122   PetscFunctionBegin;
123   if (th->endpoint) {
124     /* Evolve ts->vec_costintegral to compute integrals */
125     if (th->Theta != 1.0) {
126       PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand));
127       PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand));
128     }
129     PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand));
130     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand));
131   } else {
132     PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand));
133     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand));
134   }
135   PetscFunctionReturn(0);
136 }
137 
138 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) {
139   TS_Theta *th     = (TS_Theta *)ts->data;
140   TS        quadts = ts->quadraturets;
141 
142   PetscFunctionBegin;
143   /* backup cost integral */
144   PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0));
145   PetscCall(TSThetaEvaluateCostIntegral(ts));
146   PetscFunctionReturn(0);
147 }
148 
149 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) {
150   TS_Theta *th = (TS_Theta *)ts->data;
151 
152   PetscFunctionBegin;
153   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
154   th->ptime0     = ts->ptime + ts->time_step;
155   th->time_step0 = -ts->time_step;
156   PetscCall(TSThetaEvaluateCostIntegral(ts));
157   PetscFunctionReturn(0);
158 }
159 
160 static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x) {
161   PetscInt nits, lits;
162 
163   PetscFunctionBegin;
164   PetscCall(SNESSolve(ts->snes, b, x));
165   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
166   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
167   ts->snes_its += nits;
168   ts->ksp_its += lits;
169   PetscFunctionReturn(0);
170 }
171 
172 static PetscErrorCode TSStep_Theta(TS ts) {
173   TS_Theta *th         = (TS_Theta *)ts->data;
174   PetscInt  rejections = 0;
175   PetscBool stageok, accept = PETSC_TRUE;
176   PetscReal next_time_step = ts->time_step;
177 
178   PetscFunctionBegin;
179   if (!ts->steprollback) {
180     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
181     PetscCall(VecCopy(ts->vec_sol, th->X0));
182   }
183 
184   th->status = TS_STEP_INCOMPLETE;
185   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
186     th->shift      = 1 / (th->Theta * ts->time_step);
187     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step;
188     PetscCall(VecCopy(th->X0, th->X));
189     if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot));
190     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
191       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
192       PetscCall(VecZeroEntries(th->Xdot));
193       PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE));
194       PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta));
195     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
196       PetscCall(VecZeroEntries(th->affine));
197     }
198     PetscCall(TSPreStage(ts, th->stage_time));
199     PetscCall(TSTheta_SNESSolve(ts, th->affine, th->X));
200     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
201     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
202     if (!stageok) goto reject_step;
203 
204     th->status = TS_STEP_PENDING;
205     if (th->endpoint) {
206       PetscCall(VecCopy(th->X, ts->vec_sol));
207     } else {
208       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
209       PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
210     }
211     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
212     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
213     if (!accept) {
214       PetscCall(VecCopy(th->X0, ts->vec_sol));
215       ts->time_step = next_time_step;
216       goto reject_step;
217     }
218 
219     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
220       th->ptime0     = ts->ptime;
221       th->time_step0 = ts->time_step;
222     }
223     ts->ptime += ts->time_step;
224     ts->time_step = next_time_step;
225     break;
226 
227   reject_step:
228     ts->reject++;
229     accept = PETSC_FALSE;
230     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
231       ts->reason = TS_DIVERGED_STEP_REJECTED;
232       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
233     }
234   }
235   PetscFunctionReturn(0);
236 }
237 
238 static PetscErrorCode TSAdjointStepBEuler_Private(TS ts) {
239   TS_Theta      *th           = (TS_Theta *)ts->data;
240   TS             quadts       = ts->quadraturets;
241   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
242   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
243   PetscInt       nadj;
244   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
245   KSP            ksp;
246   PetscScalar   *xarr;
247   TSEquationType eqtype;
248   PetscBool      isexplicitode = PETSC_FALSE;
249   PetscReal      adjoint_time_step;
250 
251   PetscFunctionBegin;
252   PetscCall(TSGetEquationType(ts, &eqtype));
253   if (eqtype == TS_EQ_ODE_EXPLICIT) {
254     isexplicitode = PETSC_TRUE;
255     VecsDeltaLam  = ts->vecs_sensi;
256     VecsDeltaLam2 = ts->vecs_sensi2;
257   }
258   th->status = TS_STEP_INCOMPLETE;
259   PetscCall(SNESGetKSP(ts->snes, &ksp));
260   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
261   if (quadts) {
262     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
263     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
264   }
265 
266   th->stage_time    = ts->ptime;
267   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
268 
269   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
270   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
271 
272   for (nadj = 0; nadj < ts->numcost; nadj++) {
273     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
274     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
275     if (quadJ) {
276       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
277       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
278       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
279       PetscCall(VecResetArray(ts->vec_drdu_col));
280       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
281     }
282   }
283 
284   /* Build LHS for first-order adjoint */
285   th->shift = 1. / adjoint_time_step;
286   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
287   PetscCall(KSPSetOperators(ksp, J, Jpre));
288 
289   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
290   for (nadj = 0; nadj < ts->numcost; nadj++) {
291     KSPConvergedReason kspreason;
292     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
293     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
294     if (kspreason < 0) {
295       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
296       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
297     }
298   }
299 
300   if (ts->vecs_sensi2) { /* U_{n+1} */
301     /* Get w1 at t_{n+1} from TLM matrix */
302     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
303     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
304     /* lambda_s^T F_UU w_1 */
305     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
306     /* lambda_s^T F_UP w_2 */
307     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
308     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
309       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
310       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
311       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
312       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
313     }
314     /* Solve stage equation LHS X = RHS for second-order adjoint */
315     for (nadj = 0; nadj < ts->numcost; nadj++) {
316       KSPConvergedReason kspreason;
317       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
318       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
319       if (kspreason < 0) {
320         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
321         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
322       }
323     }
324   }
325 
326   /* Update sensitivities, and evaluate integrals if there is any */
327   if (!isexplicitode) {
328     th->shift = 0.0;
329     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
330     PetscCall(KSPSetOperators(ksp, J, Jpre));
331     for (nadj = 0; nadj < ts->numcost; nadj++) {
332       /* Add f_U \lambda_s to the original RHS */
333       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
334       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
335       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
336       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
337       if (ts->vecs_sensi2) {
338         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
339         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
340         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
341       }
342     }
343   }
344   if (ts->vecs_sensip) {
345     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
346     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
347     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
348     if (ts->vecs_sensi2p) {
349       /* lambda_s^T F_PU w_1 */
350       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
351       /* lambda_s^T F_PP w_2 */
352       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
353     }
354 
355     for (nadj = 0; nadj < ts->numcost; nadj++) {
356       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
357       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
358       if (quadJp) {
359         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
360         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
361         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
362         PetscCall(VecResetArray(ts->vec_drdp_col));
363         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
364       }
365       if (ts->vecs_sensi2p) {
366         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
367         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
368         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
369         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
370       }
371     }
372   }
373 
374   if (ts->vecs_sensi2) {
375     PetscCall(VecResetArray(ts->vec_sensip_col));
376     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
377   }
378   th->status = TS_STEP_COMPLETE;
379   PetscFunctionReturn(0);
380 }
381 
382 static PetscErrorCode TSAdjointStep_Theta(TS ts) {
383   TS_Theta    *th           = (TS_Theta *)ts->data;
384   TS           quadts       = ts->quadraturets;
385   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
386   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
387   PetscInt     nadj;
388   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
389   KSP          ksp;
390   PetscScalar *xarr;
391   PetscReal    adjoint_time_step;
392   PetscReal    adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, ususally ts->ptime is larger than adjoint_ptime) */
393 
394   PetscFunctionBegin;
395   if (th->Theta == 1.) {
396     PetscCall(TSAdjointStepBEuler_Private(ts));
397     PetscFunctionReturn(0);
398   }
399   th->status = TS_STEP_INCOMPLETE;
400   PetscCall(SNESGetKSP(ts->snes, &ksp));
401   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
402   if (quadts) {
403     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
404     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
405   }
406   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
407   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
408   adjoint_ptime     = ts->ptime + ts->time_step;
409   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
410 
411   if (!th->endpoint) {
412     /* recover th->X0 using vec_sol and the stage value th->X */
413     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
414   }
415 
416   /* Build RHS for first-order adjoint */
417   /* Cost function has an integral term */
418   if (quadts) {
419     if (th->endpoint) {
420       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
421     } else {
422       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
423     }
424   }
425 
426   for (nadj = 0; nadj < ts->numcost; nadj++) {
427     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
428     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
429     if (quadJ) {
430       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
431       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
432       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
433       PetscCall(VecResetArray(ts->vec_drdu_col));
434       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
435     }
436   }
437 
438   /* Build LHS for first-order adjoint */
439   th->shift = 1. / (th->Theta * adjoint_time_step);
440   if (th->endpoint) {
441     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
442   } else {
443     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
444   }
445   PetscCall(KSPSetOperators(ksp, J, Jpre));
446 
447   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
448   for (nadj = 0; nadj < ts->numcost; nadj++) {
449     KSPConvergedReason kspreason;
450     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
451     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
452     if (kspreason < 0) {
453       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
454       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
455     }
456   }
457 
458   /* Second-order adjoint */
459   if (ts->vecs_sensi2) { /* U_{n+1} */
460     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
461     /* Get w1 at t_{n+1} from TLM matrix */
462     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
463     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
464     /* lambda_s^T F_UU w_1 */
465     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
466     PetscCall(VecResetArray(ts->vec_sensip_col));
467     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
468     /* lambda_s^T F_UP w_2 */
469     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
470     for (nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
471       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
472       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
473       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
474       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
475     }
476     /* Solve stage equation LHS X = RHS for second-order adjoint */
477     for (nadj = 0; nadj < ts->numcost; nadj++) {
478       KSPConvergedReason kspreason;
479       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
480       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
481       if (kspreason < 0) {
482         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
483         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
484       }
485     }
486   }
487 
488   /* Update sensitivities, and evaluate integrals if there is any */
489   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
490     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
491     th->stage_time = adjoint_ptime;
492     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
493     PetscCall(KSPSetOperators(ksp, J, Jpre));
494     /* R_U at t_n */
495     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
496     for (nadj = 0; nadj < ts->numcost; nadj++) {
497       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
498       if (quadJ) {
499         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
500         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
501         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
502         PetscCall(VecResetArray(ts->vec_drdu_col));
503         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
504       }
505       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
506     }
507 
508     /* Second-order adjoint */
509     if (ts->vecs_sensi2) { /* U_n */
510       /* Get w1 at t_n from TLM matrix */
511       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
512       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
513       /* lambda_s^T F_UU w_1 */
514       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
515       PetscCall(VecResetArray(ts->vec_sensip_col));
516       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
517       /* lambda_s^T F_UU w_2 */
518       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
519       for (nadj = 0; nadj < ts->numcost; nadj++) {
520         /* 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  */
521         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
522         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
523         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
524         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
525       }
526     }
527 
528     th->stage_time = ts->ptime; /* recover the old value */
529 
530     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
531       /* U_{n+1} */
532       th->shift = 1.0 / (adjoint_time_step * th->Theta);
533       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
534       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
535       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
536       for (nadj = 0; nadj < ts->numcost; nadj++) {
537         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
538         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
539         if (quadJp) {
540           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
541           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
542           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
543           PetscCall(VecResetArray(ts->vec_drdp_col));
544           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
545         }
546       }
547       if (ts->vecs_sensi2p) { /* second-order */
548         /* Get w1 at t_{n+1} from TLM matrix */
549         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
550         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
551         /* lambda_s^T F_PU w_1 */
552         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
553         PetscCall(VecResetArray(ts->vec_sensip_col));
554         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
555 
556         /* lambda_s^T F_PP w_2 */
557         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
558         for (nadj = 0; nadj < ts->numcost; nadj++) {
559           /* 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)  */
560           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
561           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
562           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
563           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
564         }
565       }
566 
567       /* U_s */
568       PetscCall(VecZeroEntries(th->Xdot));
569       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
570       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
571       for (nadj = 0; nadj < ts->numcost; nadj++) {
572         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
573         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
574         if (quadJp) {
575           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
576           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
577           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
578           PetscCall(VecResetArray(ts->vec_drdp_col));
579           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
580         }
581         if (ts->vecs_sensi2p) { /* second-order */
582           /* Get w1 at t_n from TLM matrix */
583           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
584           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
585           /* lambda_s^T F_PU w_1 */
586           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
587           PetscCall(VecResetArray(ts->vec_sensip_col));
588           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
589           /* lambda_s^T F_PP w_2 */
590           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
591           for (nadj = 0; nadj < ts->numcost; nadj++) {
592             /* 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) */
593             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
594             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
595             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
596             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
597           }
598         }
599       }
600     }
601   } else { /* one-stage case */
602     th->shift = 0.0;
603     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
604     PetscCall(KSPSetOperators(ksp, J, Jpre));
605     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
606     for (nadj = 0; nadj < ts->numcost; nadj++) {
607       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
608       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
609       if (quadJ) {
610         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
611         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
612         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
613         PetscCall(VecResetArray(ts->vec_drdu_col));
614         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
615       }
616     }
617     if (ts->vecs_sensip) {
618       th->shift = 1.0 / (adjoint_time_step * th->Theta);
619       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
620       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
621       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
622       for (nadj = 0; nadj < ts->numcost; nadj++) {
623         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
624         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
625         if (quadJp) {
626           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
627           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
628           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
629           PetscCall(VecResetArray(ts->vec_drdp_col));
630           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
631         }
632       }
633     }
634   }
635 
636   th->status = TS_STEP_COMPLETE;
637   PetscFunctionReturn(0);
638 }
639 
640 static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X) {
641   TS_Theta *th = (TS_Theta *)ts->data;
642   PetscReal dt = t - ts->ptime;
643 
644   PetscFunctionBegin;
645   PetscCall(VecCopy(ts->vec_sol, th->X));
646   if (th->endpoint) dt *= th->Theta;
647   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
648   PetscFunctionReturn(0);
649 }
650 
651 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) {
652   TS_Theta *th = (TS_Theta *)ts->data;
653   Vec       X  = ts->vec_sol;      /* X = solution */
654   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
655   PetscReal wltea, wlter;
656 
657   PetscFunctionBegin;
658   if (!th->vec_sol_prev) {
659     *wlte = -1;
660     PetscFunctionReturn(0);
661   }
662   /* Cannot compute LTE in first step or in restart after event */
663   if (ts->steprestart) {
664     *wlte = -1;
665     PetscFunctionReturn(0);
666   }
667   /* Compute LTE using backward differences with non-constant time step */
668   {
669     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
670     PetscReal   a = 1 + h_prev / h;
671     PetscScalar scal[3];
672     Vec         vecs[3];
673     scal[0] = +1 / a;
674     scal[1] = -1 / (a - 1);
675     scal[2] = +1 / (a * (a - 1));
676     vecs[0] = X;
677     vecs[1] = th->X0;
678     vecs[2] = th->vec_sol_prev;
679     PetscCall(VecCopy(X, Y));
680     PetscCall(VecMAXPY(Y, 3, scal, vecs));
681     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
682   }
683   if (order) *order = 2;
684   PetscFunctionReturn(0);
685 }
686 
687 static PetscErrorCode TSRollBack_Theta(TS ts) {
688   TS_Theta *th     = (TS_Theta *)ts->data;
689   TS        quadts = ts->quadraturets;
690 
691   PetscFunctionBegin;
692   PetscCall(VecCopy(th->X0, ts->vec_sol));
693   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
694   th->status = TS_STEP_INCOMPLETE;
695   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
696   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
697   PetscFunctionReturn(0);
698 }
699 
700 static PetscErrorCode TSForwardStep_Theta(TS ts) {
701   TS_Theta    *th                   = (TS_Theta *)ts->data;
702   TS           quadts               = ts->quadraturets;
703   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
704   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
705   PetscInt     ntlm;
706   KSP          ksp;
707   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
708   PetscScalar *barr, *xarr;
709   PetscReal    previous_shift;
710 
711   PetscFunctionBegin;
712   previous_shift = th->shift;
713   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));
714 
715   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
716   PetscCall(SNESGetKSP(ts->snes, &ksp));
717   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
718   if (quadts) {
719     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
720     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
721   }
722 
723   /* Build RHS */
724   if (th->endpoint) { /* 2-stage method*/
725     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
726     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
727     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
728     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));
729 
730     /* Add the f_p forcing terms */
731     if (ts->Jacp) {
732       PetscCall(VecZeroEntries(th->Xdot));
733       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
734       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
735       th->shift = previous_shift;
736       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
737       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
738       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
739     }
740   } else { /* 1-stage method */
741     th->shift = 0.0;
742     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
743     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DEFAULT, &MatDeltaFwdSensip));
744     PetscCall(MatScale(MatDeltaFwdSensip, -1.));
745 
746     /* Add the f_p forcing terms */
747     if (ts->Jacp) {
748       th->shift = previous_shift;
749       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
750       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
751       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
752     }
753   }
754 
755   /* Build LHS */
756   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
757   if (th->endpoint) {
758     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
759   } else {
760     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
761   }
762   PetscCall(KSPSetOperators(ksp, J, Jpre));
763 
764   /*
765     Evaluate the first stage of integral gradients with the 2-stage method:
766     drdu|t_n*S(t_n) + drdp|t_n
767     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
768   */
769   if (th->endpoint) { /* 2-stage method only */
770     if (quadts && quadts->mat_sensip) {
771       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
772       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
773       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
774       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
775       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
776     }
777   }
778 
779   /* Solve the tangent linear equation for forward sensitivities to parameters */
780   for (ntlm = 0; ntlm < th->num_tlm; ntlm++) {
781     KSPConvergedReason kspreason;
782     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
783     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
784     if (th->endpoint) {
785       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
786       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
787       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
788       PetscCall(VecResetArray(ts->vec_sensip_col));
789       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
790     } else {
791       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
792     }
793     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
794     if (kspreason < 0) {
795       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
796       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
797     }
798     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
799     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
800   }
801 
802   /*
803     Evaluate the second stage of integral gradients with the 2-stage method:
804     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
805   */
806   if (quadts && quadts->mat_sensip) {
807     if (!th->endpoint) {
808       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
809       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
810       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
811       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
812       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
813       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
814       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
815     } else {
816       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
817       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
818       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DEFAULT, &th->MatIntegralSensipTemp));
819       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
820       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
821     }
822   } else {
823     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
824   }
825   PetscFunctionReturn(0);
826 }
827 
828 static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[]) {
829   TS_Theta *th = (TS_Theta *)ts->data;
830 
831   PetscFunctionBegin;
832   if (ns) {
833     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
834     else *ns = 2;                                   /* endpoint form */
835   }
836   if (stagesensip) {
837     if (!th->endpoint && th->Theta != 1.0) {
838       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
839     } else {
840       th->MatFwdStages[0] = th->MatFwdSensip0;
841       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
842     }
843     *stagesensip = th->MatFwdStages;
844   }
845   PetscFunctionReturn(0);
846 }
847 
848 /*------------------------------------------------------------*/
849 static PetscErrorCode TSReset_Theta(TS ts) {
850   TS_Theta *th = (TS_Theta *)ts->data;
851 
852   PetscFunctionBegin;
853   PetscCall(VecDestroy(&th->X));
854   PetscCall(VecDestroy(&th->Xdot));
855   PetscCall(VecDestroy(&th->X0));
856   PetscCall(VecDestroy(&th->affine));
857 
858   PetscCall(VecDestroy(&th->vec_sol_prev));
859   PetscCall(VecDestroy(&th->vec_lte_work));
860 
861   PetscCall(VecDestroy(&th->VecCostIntegral0));
862   PetscFunctionReturn(0);
863 }
864 
865 static PetscErrorCode TSAdjointReset_Theta(TS ts) {
866   TS_Theta *th = (TS_Theta *)ts->data;
867 
868   PetscFunctionBegin;
869   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
870   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
871   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
872   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
873   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
874   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
875   PetscFunctionReturn(0);
876 }
877 
878 static PetscErrorCode TSDestroy_Theta(TS ts) {
879   PetscFunctionBegin;
880   PetscCall(TSReset_Theta(ts));
881   if (ts->dm) {
882     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
883     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
884   }
885   PetscCall(PetscFree(ts->data));
886   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
887   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
888   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
889   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
890   PetscFunctionReturn(0);
891 }
892 
893 /*
894   This defines the nonlinear equation that is to be solved with SNES
895   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
896 
897   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
898   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
899   U = (U_{n+1} + U0)/2
900 */
901 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts) {
902   TS_Theta *th = (TS_Theta *)ts->data;
903   Vec       X0, Xdot;
904   DM        dm, dmsave;
905   PetscReal shift = th->shift;
906 
907   PetscFunctionBegin;
908   PetscCall(SNESGetDM(snes, &dm));
909   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
910   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
911   if (x != X0) {
912     PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
913   } else {
914     PetscCall(VecZeroEntries(Xdot));
915   }
916   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
917   dmsave = ts->dm;
918   ts->dm = dm;
919   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
920   ts->dm = dmsave;
921   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
922   PetscFunctionReturn(0);
923 }
924 
925 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts) {
926   TS_Theta *th = (TS_Theta *)ts->data;
927   Vec       Xdot;
928   DM        dm, dmsave;
929   PetscReal shift = th->shift;
930 
931   PetscFunctionBegin;
932   PetscCall(SNESGetDM(snes, &dm));
933   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
934   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));
935 
936   dmsave = ts->dm;
937   ts->dm = dm;
938   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
939   ts->dm = dmsave;
940   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
941   PetscFunctionReturn(0);
942 }
943 
944 static PetscErrorCode TSForwardSetUp_Theta(TS ts) {
945   TS_Theta *th     = (TS_Theta *)ts->data;
946   TS        quadts = ts->quadraturets;
947 
948   PetscFunctionBegin;
949   /* combine sensitivities to parameters and sensitivities to initial values into one array */
950   th->num_tlm = ts->num_parameters;
951   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
952   if (quadts && quadts->mat_sensip) {
953     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
954     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
955   }
956   /* backup sensitivity results for roll-backs */
957   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));
958 
959   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
960   PetscFunctionReturn(0);
961 }
962 
963 static PetscErrorCode TSForwardReset_Theta(TS ts) {
964   TS_Theta *th     = (TS_Theta *)ts->data;
965   TS        quadts = ts->quadraturets;
966 
967   PetscFunctionBegin;
968   if (quadts && quadts->mat_sensip) {
969     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
970     PetscCall(MatDestroy(&th->MatIntegralSensip0));
971   }
972   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
973   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
974   PetscCall(MatDestroy(&th->MatFwdSensip0));
975   PetscFunctionReturn(0);
976 }
977 
978 static PetscErrorCode TSSetUp_Theta(TS ts) {
979   TS_Theta *th     = (TS_Theta *)ts->data;
980   TS        quadts = ts->quadraturets;
981   PetscBool match;
982 
983   PetscFunctionBegin;
984   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
985     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
986   }
987   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
988   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
989   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
990   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
991 
992   th->order = (th->Theta == 0.5) ? 2 : 1;
993   th->shift = 1 / (th->Theta * ts->time_step);
994 
995   PetscCall(TSGetDM(ts, &ts->dm));
996   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
997   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
998 
999   PetscCall(TSGetAdapt(ts, &ts->adapt));
1000   PetscCall(TSAdaptCandidatesClear(ts->adapt));
1001   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
1002   if (!match) {
1003     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
1004     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
1005   }
1006   PetscCall(TSGetSNES(ts, &ts->snes));
1007 
1008   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 /*------------------------------------------------------------*/
1013 
1014 static PetscErrorCode TSAdjointSetUp_Theta(TS ts) {
1015   TS_Theta *th = (TS_Theta *)ts->data;
1016 
1017   PetscFunctionBegin;
1018   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
1019   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
1020   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1021   if (ts->vecs_sensi2) {
1022     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
1023     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
1024     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1025     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
1026     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1027   }
1028   if (ts->vecs_sensi2p) {
1029     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
1030     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1031     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
1032     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1033   }
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems *PetscOptionsObject) {
1038   TS_Theta *th = (TS_Theta *)ts->data;
1039 
1040   PetscFunctionBegin;
1041   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1042   {
1043     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
1044     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
1045     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1046   }
1047   PetscOptionsHeadEnd();
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer) {
1052   TS_Theta *th = (TS_Theta *)ts->data;
1053   PetscBool iascii;
1054 
1055   PetscFunctionBegin;
1056   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1057   if (iascii) {
1058     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
1059     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1060   }
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta) {
1065   TS_Theta *th = (TS_Theta *)ts->data;
1066 
1067   PetscFunctionBegin;
1068   *theta = th->Theta;
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta) {
1073   TS_Theta *th = (TS_Theta *)ts->data;
1074 
1075   PetscFunctionBegin;
1076   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
1077   th->Theta = theta;
1078   th->order = (th->Theta == 0.5) ? 2 : 1;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint) {
1083   TS_Theta *th = (TS_Theta *)ts->data;
1084 
1085   PetscFunctionBegin;
1086   *endpoint = th->endpoint;
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg) {
1091   TS_Theta *th = (TS_Theta *)ts->data;
1092 
1093   PetscFunctionBegin;
1094   th->endpoint = flg;
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #if defined(PETSC_HAVE_COMPLEX)
1099 static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi) {
1100   PetscComplex    z   = xr + xi * PETSC_i, f;
1101   TS_Theta       *th  = (TS_Theta *)ts->data;
1102   const PetscReal one = 1.0;
1103 
1104   PetscFunctionBegin;
1105   f   = (one + (one - th->Theta) * z) / (one - th->Theta * z);
1106   *yr = PetscRealPartComplex(f);
1107   *yi = PetscImaginaryPartComplex(f);
1108   PetscFunctionReturn(0);
1109 }
1110 #endif
1111 
1112 static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[]) {
1113   TS_Theta *th = (TS_Theta *)ts->data;
1114 
1115   PetscFunctionBegin;
1116   if (ns) {
1117     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
1118     else *ns = 2;                                   /* endpoint form */
1119   }
1120   if (Y) {
1121     if (!th->endpoint && th->Theta != 1.0) {
1122       th->Stages[0] = th->X;
1123     } else {
1124       th->Stages[0] = th->X0;
1125       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1126     }
1127     *Y = th->Stages;
1128   }
1129   PetscFunctionReturn(0);
1130 }
1131 
1132 /* ------------------------------------------------------------ */
1133 /*MC
1134       TSTHETA - DAE solver using the implicit Theta method
1135 
1136    Level: beginner
1137 
1138    Options Database:
1139 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1140 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
1141 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
1142 
1143    Notes:
1144 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
1145 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
1146 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
1147 
1148    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.
1149 
1150    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1151 
1152 .vb
1153   Theta | Theta
1154   -------------
1155         |  1
1156 .ve
1157 
1158    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1159 
1160    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1161 
1162 .vb
1163   0 | 0         0
1164   1 | 1-Theta   Theta
1165   -------------------
1166     | 1-Theta   Theta
1167 .ve
1168 
1169    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1170 
1171    To apply a diagonally implicit RK method to DAE, the stage formula
1172 
1173 $  Y_i = X + h sum_j a_ij Y'_j
1174 
1175    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1176 
1177 .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1178 
1179 M*/
1180 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) {
1181   TS_Theta *th;
1182 
1183   PetscFunctionBegin;
1184   ts->ops->reset          = TSReset_Theta;
1185   ts->ops->adjointreset   = TSAdjointReset_Theta;
1186   ts->ops->destroy        = TSDestroy_Theta;
1187   ts->ops->view           = TSView_Theta;
1188   ts->ops->setup          = TSSetUp_Theta;
1189   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1190   ts->ops->adjointreset   = TSAdjointReset_Theta;
1191   ts->ops->step           = TSStep_Theta;
1192   ts->ops->interpolate    = TSInterpolate_Theta;
1193   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
1194   ts->ops->rollback       = TSRollBack_Theta;
1195   ts->ops->setfromoptions = TSSetFromOptions_Theta;
1196   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
1197   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1198 #if defined(PETSC_HAVE_COMPLEX)
1199   ts->ops->linearstability = TSComputeLinearStability_Theta;
1200 #endif
1201   ts->ops->getstages       = TSGetStages_Theta;
1202   ts->ops->adjointstep     = TSAdjointStep_Theta;
1203   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1204   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
1205   ts->default_adapt_type   = TSADAPTNONE;
1206 
1207   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
1208   ts->ops->forwardreset     = TSForwardReset_Theta;
1209   ts->ops->forwardstep      = TSForwardStep_Theta;
1210   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1211 
1212   ts->usessnes = PETSC_TRUE;
1213 
1214   PetscCall(PetscNew(&th));
1215   ts->data = (void *)th;
1216 
1217   th->VecsDeltaLam   = NULL;
1218   th->VecsDeltaMu    = NULL;
1219   th->VecsSensiTemp  = NULL;
1220   th->VecsSensi2Temp = NULL;
1221 
1222   th->extrapolate = PETSC_FALSE;
1223   th->Theta       = 0.5;
1224   th->order       = 2;
1225   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
1226   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
1227   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
1228   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
1229   PetscFunctionReturn(0);
1230 }
1231 
1232 /*@
1233   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
1234 
1235   Not Collective
1236 
1237   Input Parameter:
1238 .  ts - timestepping context
1239 
1240   Output Parameter:
1241 .  theta - stage abscissa
1242 
1243   Note:
1244   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
1245 
1246   Level: Advanced
1247 
1248 .seealso: `TSThetaSetTheta()`
1249 @*/
1250 PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta) {
1251   PetscFunctionBegin;
1252   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1253   PetscValidRealPointer(theta, 2);
1254   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
1255   PetscFunctionReturn(0);
1256 }
1257 
1258 /*@
1259   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
1260 
1261   Not Collective
1262 
1263   Input Parameters:
1264 +  ts - timestepping context
1265 -  theta - stage abscissa
1266 
1267   Options Database:
1268 .  -ts_theta_theta <theta> - set theta
1269 
1270   Level: Intermediate
1271 
1272 .seealso: `TSThetaGetTheta()`
1273 @*/
1274 PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta) {
1275   PetscFunctionBegin;
1276   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1277   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
1278   PetscFunctionReturn(0);
1279 }
1280 
1281 /*@
1282   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1283 
1284   Not Collective
1285 
1286   Input Parameter:
1287 .  ts - timestepping context
1288 
1289   Output Parameter:
1290 .  endpoint - PETSC_TRUE when using the endpoint variant
1291 
1292   Level: Advanced
1293 
1294 .seealso: `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
1295 @*/
1296 PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint) {
1297   PetscFunctionBegin;
1298   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1299   PetscValidBoolPointer(endpoint, 2);
1300   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 /*@
1305   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1306 
1307   Not Collective
1308 
1309   Input Parameters:
1310 +  ts - timestepping context
1311 -  flg - PETSC_TRUE to use the endpoint variant
1312 
1313   Options Database:
1314 .  -ts_theta_endpoint <flg> - use the endpoint variant
1315 
1316   Level: Intermediate
1317 
1318 .seealso: `TSTHETA`, `TSCN`
1319 @*/
1320 PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg) {
1321   PetscFunctionBegin;
1322   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1323   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 /*
1328  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1329  * The creation functions for these specializations are below.
1330  */
1331 
1332 static PetscErrorCode TSSetUp_BEuler(TS ts) {
1333   TS_Theta *th = (TS_Theta *)ts->data;
1334 
1335   PetscFunctionBegin;
1336   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");
1337   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");
1338   PetscCall(TSSetUp_Theta(ts));
1339   PetscFunctionReturn(0);
1340 }
1341 
1342 static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer) {
1343   PetscFunctionBegin;
1344   PetscFunctionReturn(0);
1345 }
1346 
1347 /*MC
1348       TSBEULER - ODE solver using the implicit backward Euler method
1349 
1350   Level: beginner
1351 
1352   Notes:
1353   TSBEULER is equivalent to TSTHETA with Theta=1.0
1354 
1355 $  -ts_type theta -ts_theta_theta 1.0
1356 
1357 .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1358 
1359 M*/
1360 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) {
1361   PetscFunctionBegin;
1362   PetscCall(TSCreate_Theta(ts));
1363   PetscCall(TSThetaSetTheta(ts, 1.0));
1364   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
1365   ts->ops->setup = TSSetUp_BEuler;
1366   ts->ops->view  = TSView_BEuler;
1367   PetscFunctionReturn(0);
1368 }
1369 
1370 static PetscErrorCode TSSetUp_CN(TS ts) {
1371   TS_Theta *th = (TS_Theta *)ts->data;
1372 
1373   PetscFunctionBegin;
1374   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");
1375   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");
1376   PetscCall(TSSetUp_Theta(ts));
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer) {
1381   PetscFunctionBegin;
1382   PetscFunctionReturn(0);
1383 }
1384 
1385 /*MC
1386       TSCN - ODE solver using the implicit Crank-Nicolson method.
1387 
1388   Level: beginner
1389 
1390   Notes:
1391   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1392 
1393 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1394 
1395 .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`
1396 
1397 M*/
1398 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) {
1399   PetscFunctionBegin;
1400   PetscCall(TSCreate_Theta(ts));
1401   PetscCall(TSThetaSetTheta(ts, 0.5));
1402   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
1403   ts->ops->setup = TSSetUp_CN;
1404   ts->ops->view  = TSView_CN;
1405   PetscFunctionReturn(0);
1406 }
1407