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