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