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