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