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