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