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