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