xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 6f92202b937c8db5b63afc621286d06689bf89d4)
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   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
25   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
26   Vec          *VecsSensiTemp;           /* Vector to be timed with Jacobian transpose */
27   Vec          *VecsDeltaFwdSensi;      /* Increment of the forward sensitivity at stage */
28   Vec          *VecsDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
29   Vec          VecIntegralSensiTemp;     /* Working vector for forward integral sensitivity */
30   Vec          VecIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
31   Vec          *VecsFwdSensi0;          /* backup for roll-backs due to events */
32   Vec          *VecsFwdSensip0;         /* backup for roll-backs due to events */
33   Vec          *VecsIntegralSensi0;      /* backup for roll-backs due to events */
34   Vec          *VecsIntegralSensip0;     /* backup for roll-backs due to events */
35 
36   /* context for error estimation */
37   Vec          vec_sol_prev;
38   Vec          vec_lte_work;
39 } TS_Theta;
40 
41 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
42 {
43   TS_Theta       *th = (TS_Theta*)ts->data;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   if (X0) {
48     if (dm && dm != ts->dm) {
49       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
50     } else *X0 = ts->vec_sol;
51   }
52   if (Xdot) {
53     if (dm && dm != ts->dm) {
54       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
55     } else *Xdot = th->Xdot;
56   }
57   PetscFunctionReturn(0);
58 }
59 
60 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
61 {
62   PetscErrorCode ierr;
63 
64   PetscFunctionBegin;
65   if (X0) {
66     if (dm && dm != ts->dm) {
67       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
68     }
69   }
70   if (Xdot) {
71     if (dm && dm != ts->dm) {
72       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
73     }
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
79 {
80   PetscFunctionBegin;
81   PetscFunctionReturn(0);
82 }
83 
84 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
85 {
86   TS             ts = (TS)ctx;
87   PetscErrorCode ierr;
88   Vec            X0,Xdot,X0_c,Xdot_c;
89 
90   PetscFunctionBegin;
91   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
92   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
93   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
94   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
95   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
96   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
97   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
98   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
103 {
104   PetscFunctionBegin;
105   PetscFunctionReturn(0);
106 }
107 
108 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
109 {
110   TS             ts = (TS)ctx;
111   PetscErrorCode ierr;
112   Vec            X0,Xdot,X0_sub,Xdot_sub;
113 
114   PetscFunctionBegin;
115   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
116   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
117 
118   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120 
121   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123 
124   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
125   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
130 {
131   TS_Theta       *th = (TS_Theta*)ts->data;
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   /* backup cost integral */
136   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
137   if (th->endpoint) {
138     ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
139     ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
140   }
141   ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
142   if (th->endpoint) {
143     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
144   } else {
145     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
151 {
152   TS_Theta       *th = (TS_Theta*)ts->data;
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   if (th->endpoint) {
157     /* Evolve ts->vec_costintegral to compute integrals */
158     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
159     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
160     if (th->Theta!=1) {
161       ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
162       ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr);
163     }
164   } else {
165     ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
166     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
167   }
168   PetscFunctionReturn(0);
169 }
170 
171 static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
172 {
173   PetscInt       nits,lits;
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
178   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
179   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
180   ts->snes_its += nits; ts->ksp_its += lits;
181   PetscFunctionReturn(0);
182 }
183 
184 static PetscErrorCode TSStep_Theta(TS ts)
185 {
186   TS_Theta       *th = (TS_Theta*)ts->data;
187   PetscInt       rejections = 0;
188   PetscBool      stageok,accept = PETSC_TRUE;
189   PetscReal      next_time_step = ts->time_step;
190   PetscErrorCode ierr;
191 
192   PetscFunctionBegin;
193   if (!ts->steprollback) {
194     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
195     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
196   }
197 
198   th->status = TS_STEP_INCOMPLETE;
199   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
200 
201     PetscReal shift = 1/(th->Theta*ts->time_step);
202     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
203 
204     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
205     if (th->extrapolate && !ts->steprestart) {
206       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
207     }
208     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
209       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
210       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
211       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
212       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
213     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
214       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
215     }
216     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
217     ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
218     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
219     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
220     if (!stageok) goto reject_step;
221 
222     th->status = TS_STEP_PENDING;
223     if (th->endpoint) {
224       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
225     } else {
226       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
227       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
228     }
229     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
230     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
231     if (!accept) {
232       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
233       ts->time_step = next_time_step;
234       goto reject_step;
235     }
236 
237     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
238       th->ptime     = ts->ptime;
239       th->time_step = ts->time_step;
240     }
241 
242     ts->ptime += ts->time_step;
243     ts->time_step = next_time_step;
244     break;
245 
246   reject_step:
247     ts->reject++; accept = PETSC_FALSE;
248     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
249       ts->reason = TS_DIVERGED_STEP_REJECTED;
250       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
251     }
252   }
253   PetscFunctionReturn(0);
254 }
255 
256 static PetscErrorCode TSAdjointStep_Theta(TS ts)
257 {
258   TS_Theta            *th = (TS_Theta*)ts->data;
259   Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
260   PetscInt            nadj;
261   PetscErrorCode      ierr;
262   Mat                 J,Jp;
263   KSP                 ksp;
264   PetscReal           shift;
265 
266   PetscFunctionBegin;
267   th->status = TS_STEP_INCOMPLETE;
268   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
269   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
270 
271   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
272   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
273   th->ptime      = ts->ptime + ts->time_step;
274 
275   /* Build RHS */
276   if (ts->vec_costintegral) { /* Cost function has an integral term */
277     if (th->endpoint) {
278       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
279     } else {
280       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
281     }
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->Theta*ts->time_step));CHKERRQ(ierr);
286     if (ts->vec_costintegral) {
287       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
288     }
289   }
290 
291   /* Build LHS */
292   shift = -1./(th->Theta*ts->time_step);
293   if (th->endpoint) {
294     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
295   } else {
296     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
297   }
298   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
299 
300   /* Solve LHS X = RHS */
301   for (nadj=0; nadj<ts->numcost; nadj++) {
302     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
303   }
304 
305   /* Update sensitivities, and evaluate integrals if there is any */
306   if(th->endpoint) { /* two-stage case */
307     if (th->Theta!=1.) {
308       shift = -1./((th->Theta-1.)*ts->time_step);
309       ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
310       if (ts->vec_costintegral) {
311         ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
312       }
313       for (nadj=0; nadj<ts->numcost; nadj++) {
314         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
315         if (ts->vec_costintegral) {
316           ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
317         }
318         ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
319       }
320     } else { /* backward Euler */
321       shift = 0.0;
322       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
323       for (nadj=0; nadj<ts->numcost; nadj++) {
324         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
325         ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
326         if (ts->vec_costintegral) {
327           ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
328         }
329       }
330     }
331 
332     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
333       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
334       for (nadj=0; nadj<ts->numcost; nadj++) {
335         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
336         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
337       }
338       if (th->Theta!=1.) {
339         ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
340         for (nadj=0; nadj<ts->numcost; nadj++) {
341           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
342           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
343         }
344       }
345       if (ts->vec_costintegral) {
346         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
347         for (nadj=0; nadj<ts->numcost; nadj++) {
348           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
349         }
350         if (th->Theta!=1.) {
351           ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
352           for (nadj=0; nadj<ts->numcost; nadj++) {
353             ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
354           }
355         }
356       }
357     }
358   } else { /* one-stage case */
359     shift = 0.0;
360     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
361     if (ts->vec_costintegral) {
362       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
363     }
364     for (nadj=0; nadj<ts->numcost; nadj++) {
365       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
366       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
367       if (ts->vec_costintegral) {
368         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
369       }
370     }
371     if (ts->vecs_sensip) {
372       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
373       for (nadj=0; nadj<ts->numcost; nadj++) {
374         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
375         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
376       }
377       if (ts->vec_costintegral) {
378         ierr = TSAdjointComputeDRDPFunction(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],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
381         }
382       }
383     }
384   }
385 
386   th->status = TS_STEP_COMPLETE;
387   PetscFunctionReturn(0);
388 }
389 
390 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
391 {
392   TS_Theta       *th = (TS_Theta*)ts->data;
393   PetscReal      dt  = t - ts->ptime;
394   PetscErrorCode ierr;
395 
396   PetscFunctionBegin;
397   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
398   if (th->endpoint) dt *= th->Theta;
399   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
400   PetscFunctionReturn(0);
401 }
402 
403 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
404 {
405   TS_Theta       *th = (TS_Theta*)ts->data;
406   Vec            X = ts->vec_sol;      /* X = solution */
407   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
408   PetscReal      wltea,wlter;
409   PetscErrorCode ierr;
410 
411   PetscFunctionBegin;
412   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
413   /* Cannot compute LTE in first step or in restart after event */
414   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
415   /* Compute LTE using backward differences with non-constant time step */
416   {
417     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
418     PetscReal   a = 1 + h_prev/h;
419     PetscScalar scal[3]; Vec vecs[3];
420     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
421     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
422     ierr = VecCopy(X,Y);CHKERRQ(ierr);
423     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
424     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
425   }
426   if (order) *order = 2;
427   PetscFunctionReturn(0);
428 }
429 
430 static PetscErrorCode TSRollBack_Theta(TS ts)
431 {
432   TS_Theta       *th = (TS_Theta*)ts->data;
433   PetscInt       ndir,ncost;
434   PetscErrorCode ierr;
435 
436   PetscFunctionBegin;
437   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
438   if (ts->vec_costintegral && ts->costintegralfwd) {
439     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
440   }
441   th->status = TS_STEP_INCOMPLETE;
442 
443   if (ts->vecs_fwdsensi) {
444     for (ndir=0;ndir<ts->num_initialvalues;ndir++) {
445       ierr = VecCopy(th->VecsFwdSensi0[ndir],ts->vecs_fwdsensi[ndir]);CHKERRQ(ierr);
446     }
447   }
448   if (ts->vecs_fwdsensip) {
449     for (ndir=0;ndir<ts->num_parameters;ndir++) {
450       ierr = VecCopy(th->VecsFwdSensip0[ndir],ts->vecs_fwdsensip[ndir]);CHKERRQ(ierr);
451     }
452   }
453   if (ts->vecs_integral_sensi) {
454     for (ncost=0;ncost<ts->numcost;ncost++) {
455       ierr = VecCopy(th->VecsIntegralSensi0[ncost],ts->vecs_integral_sensi[ncost]);CHKERRQ(ierr);
456     }
457   }
458   if (ts->vecs_integral_sensip) {
459     for (ncost=0;ncost<ts->numcost;ncost++) {
460       ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr);
461     }
462   }
463   PetscFunctionReturn(0);
464 }
465 
466 static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt nump,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative)
467 {
468   PetscInt    ndir;
469   PetscScalar *v;
470   PetscErrorCode ierr;
471 
472   PetscFunctionBegin;
473 
474   ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
475   for (ndir=0; ndir<nump; ndir++) {
476     ierr = VecDot(DRncostDY,s[ndir],&v[ndir]);CHKERRQ(ierr);
477   }
478   ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
479   if (DRncostDP) {
480     ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr);
481   }
482   PetscFunctionReturn(0);
483 }
484 
485 static PetscErrorCode TSForwardStep_Theta(TS ts)
486 {
487   TS_Theta       *th = (TS_Theta*)ts->data;
488   Vec            *VecsDeltaFwdSensi  = th->VecsDeltaFwdSensi;
489   Vec            *VecsDeltaFwdSensip = th->VecsDeltaFwdSensip;
490   PetscInt       ndir,ncost;
491   KSP            ksp;
492   Mat            J,Jp;
493   PetscReal      shift;
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   for (ndir=0;ndir<ts->num_initialvalues;ndir++) {
498     ierr = VecCopy(ts->vecs_fwdsensi[ndir],th->VecsFwdSensi0[ndir]);CHKERRQ(ierr);
499   }
500   for (ndir=0;ndir<ts->num_parameters;ndir++) {
501     ierr = VecCopy(ts->vecs_fwdsensip[ndir],th->VecsFwdSensip0[ndir]);CHKERRQ(ierr);
502   }
503   if (ts->vecs_integral_sensi) {
504     for (ncost=0; ncost<ts->numcost;ncost++) {
505       ierr = VecCopy(ts->vecs_integral_sensi[ncost],th->VecsIntegralSensi0[ncost]);CHKERRQ(ierr);
506     }
507   }
508   if (ts->vecs_integral_sensip) {
509     for (ncost=0;ncost<ts->numcost;ncost++) {
510       ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr);
511     }
512   }
513 
514   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
515   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
516 
517   /* Build RHS */
518   if (th->endpoint) { /* 2-stage method*/
519     shift = 1./((th->Theta-1.)*th->time_step);
520     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
521     for (ndir=0; ndir<ts->num_parameters; ndir++) {
522       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
523       ierr = VecScale(VecsDeltaFwdSensip[ndir],(th->Theta-1.)/th->Theta);
524     }
525     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
526       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
527       ierr = VecScale(VecsDeltaFwdSensi[ndir],(th->Theta-1.)/th->Theta);
528     }
529 
530     if (ts->num_parameters) { /* Add the f_p forcing terms */
531       ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr);
532       for (ndir=0; ndir<ts->num_parameters; ndir++) {
533         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],(1.-th->Theta)/th->Theta,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
534       }
535       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr);
536       for (ndir=0; ndir<ts->num_parameters; ndir++) {
537         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
538       }
539     }
540   } else { /* 1-stage method */
541     shift = 0.0;
542     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
543     for (ndir=0; ndir<ts->num_parameters; ndir++) {
544       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
545       ierr = VecScale(VecsDeltaFwdSensip[ndir],-1);
546     }
547     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
548       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
549       ierr = VecScale(VecsDeltaFwdSensi[ndir],-1);
550     }
551     if (ts->num_parameters) { /* Add the f_p forcing terms */
552       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr);
553       for (ndir=0; ndir<ts->num_parameters; ndir++) {
554         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
555       }
556     }
557   }
558 
559   /* Build LHS */
560   shift = 1/(th->Theta*th->time_step);
561   if (th->endpoint) {
562     ierr  = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
563   } else {
564     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
565   }
566   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
567 
568   /*
569     Evaluate the first stage of integral gradients with the 2-stage method:
570     drdy|t_n*S(t_n) + drdp|t_n
571     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
572     It is important to note that sensitivitesi to parameters (sensip) and sensitivities to initial values (sensi) are independent of each other, but integral sensip relies on sensip and integral sensi requires sensi
573    */
574   if (th->endpoint) { /* 2-stage method only */
575     if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
576       ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
577     }
578     if (ts->vecs_integral_sensip) {
579       ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
580       for (ncost=0; ncost<ts->numcost; ncost++) {
581         ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
582         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
583       }
584     }
585     if (ts->vecs_integral_sensi) {
586       for (ncost=0; ncost<ts->numcost; ncost++) {
587         ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
588         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr);
589       }
590     }
591   }
592 
593   /* Solve the tangent linear equation for forward sensitivities to parameters */
594   for (ndir=0; ndir<ts->num_parameters; ndir++) {
595     if (th->endpoint) {
596       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],ts->vecs_fwdsensip[ndir]);
597     } else {
598       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],VecsDeltaFwdSensip[ndir]);
599       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],1.,VecsDeltaFwdSensip[ndir]);
600     }
601   }
602   /* Solve the linear system for trajectory sensitivities to initial values */
603   for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
604     if (th->endpoint) {
605       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],ts->vecs_fwdsensi[ndir]);
606     } else {
607       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],VecsDeltaFwdSensi[ndir]);
608       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],1.,VecsDeltaFwdSensi[ndir]);
609     }
610   }
611 
612   /*Evaluate the second stage of integral gradients with the 2-stage method:
613     drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
614   */
615   if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
616     ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
617   }
618   if (ts->vecs_integral_sensip) {
619     ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
620     for (ncost=0; ncost<ts->numcost; ncost++) {
621       ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
622       if (th->endpoint) {
623         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
624       } else {
625         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
626       }
627     }
628   }
629   if (ts->vecs_integral_sensi) {
630     for (ncost=0; ncost<ts->numcost; ncost++) {
631       ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
632       if (th->endpoint) {
633         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr);
634       } else {
635         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr);
636       }
637     }
638   }
639 
640   if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */
641     for (ndir=0; ndir<ts->num_parameters; ndir++) {
642       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensip[ndir]);
643     }
644     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
645       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ndir]);
646     }
647   }
648 
649   PetscFunctionReturn(0);
650 }
651 
652 /*------------------------------------------------------------*/
653 static PetscErrorCode TSReset_Theta(TS ts)
654 {
655   TS_Theta       *th = (TS_Theta*)ts->data;
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
660   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
661   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
662   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
663 
664   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
665   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
666 
667   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
668   if (ts->forward_solve) {
669     if (ts->vecs_fwdsensi) {
670       ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
671       ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsFwdSensi0);CHKERRQ(ierr);
672     }
673     if (ts->vecs_fwdsensip) {
674       ierr = VecDestroyVecs(ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
675       ierr = VecDestroyVecs(ts->num_parameters,&th->VecsFwdSensip0);CHKERRQ(ierr);
676     }
677     if (ts->vecs_integral_sensi) {
678       ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr);
679       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
680     }
681     if (ts->vecs_integral_sensip) {
682       ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
683       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
684     }
685   }
686   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
687   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
688   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
689 
690   PetscFunctionReturn(0);
691 }
692 
693 static PetscErrorCode TSDestroy_Theta(TS ts)
694 {
695   PetscErrorCode ierr;
696 
697   PetscFunctionBegin;
698   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
699   ierr = PetscFree(ts->data);CHKERRQ(ierr);
700   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
701   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
702   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
703   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 /*
708   This defines the nonlinear equation that is to be solved with SNES
709   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
710 */
711 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
712 {
713   TS_Theta       *th = (TS_Theta*)ts->data;
714   PetscErrorCode ierr;
715   Vec            X0,Xdot;
716   DM             dm,dmsave;
717   PetscReal      shift = 1/(th->Theta*ts->time_step);
718 
719   PetscFunctionBegin;
720   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
721   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
722   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
723   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
724 
725   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
726   dmsave = ts->dm;
727   ts->dm = dm;
728   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
729   ts->dm = dmsave;
730   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
731   PetscFunctionReturn(0);
732 }
733 
734 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
735 {
736   TS_Theta       *th = (TS_Theta*)ts->data;
737   PetscErrorCode ierr;
738   Vec            Xdot;
739   DM             dm,dmsave;
740   PetscReal      shift = 1/(th->Theta*ts->time_step);
741 
742   PetscFunctionBegin;
743   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
744   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
745   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
746 
747   dmsave = ts->dm;
748   ts->dm = dm;
749   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
750   ts->dm = dmsave;
751   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
752   PetscFunctionReturn(0);
753 }
754 
755 static PetscErrorCode TSForwardSetUp_Theta(TS ts)
756 {
757   TS_Theta       *th = (TS_Theta*)ts->data;
758   PetscErrorCode ierr;
759 
760   PetscFunctionBegin;
761   if (ts->vecs_fwdsensi) {
762     ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
763   }
764   if (ts->vecs_fwdsensip) {
765     ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
766   }
767   if (ts->vecs_integral_sensi) {
768     ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr);
769   }
770   if (ts->vecs_integral_sensip) {
771     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
772   }
773   /* backup sensitivity results for roll-backs */
774   if (ts->vecs_fwdsensi) {
775       ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsFwdSensi0);CHKERRQ(ierr);
776   }
777   if (ts->vecs_fwdsensip) {
778     ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsFwdSensip0);CHKERRQ(ierr);
779   }
780   if (ts->vecs_integral_sensi) {
781     ierr = VecDuplicateVecs(ts->vecs_integral_sensi[0],ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
782   }
783   if (ts->vecs_integral_sensip) {
784     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
785   }
786   PetscFunctionReturn(0);
787 }
788 
789 static PetscErrorCode TSSetUp_Theta(TS ts)
790 {
791   TS_Theta       *th = (TS_Theta*)ts->data;
792   PetscBool      match;
793   PetscErrorCode ierr;
794 
795   PetscFunctionBegin;
796   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
797     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
798   }
799   if (!th->X) {
800     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
801   }
802   if (!th->Xdot) {
803     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
804   }
805   if (!th->X0) {
806     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
807   }
808   if (th->endpoint) {
809     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
810   }
811 
812   th->order = (th->Theta == 0.5) ? 2 : 1;
813 
814   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
815   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
816   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
817 
818   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
819   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
820   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
821   if (!match) {
822     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
823     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
824   }
825   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
826   PetscFunctionReturn(0);
827 }
828 
829 /*------------------------------------------------------------*/
830 
831 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
832 {
833   TS_Theta       *th = (TS_Theta*)ts->data;
834   PetscErrorCode ierr;
835 
836   PetscFunctionBegin;
837   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
838   if(ts->vecs_sensip) {
839     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
840   }
841   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
842   PetscFunctionReturn(0);
843 }
844 
845 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
846 {
847   TS_Theta       *th = (TS_Theta*)ts->data;
848   PetscErrorCode ierr;
849 
850   PetscFunctionBegin;
851   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
852   {
853     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
854     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
855     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
856   }
857   ierr = PetscOptionsTail();CHKERRQ(ierr);
858   PetscFunctionReturn(0);
859 }
860 
861 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
862 {
863   TS_Theta       *th = (TS_Theta*)ts->data;
864   PetscBool      iascii;
865   PetscErrorCode ierr;
866 
867   PetscFunctionBegin;
868   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
869   if (iascii) {
870     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
871     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
872   }
873   PetscFunctionReturn(0);
874 }
875 
876 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
877 {
878   TS_Theta *th = (TS_Theta*)ts->data;
879 
880   PetscFunctionBegin;
881   *theta = th->Theta;
882   PetscFunctionReturn(0);
883 }
884 
885 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
886 {
887   TS_Theta *th = (TS_Theta*)ts->data;
888 
889   PetscFunctionBegin;
890   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
891   th->Theta = theta;
892   th->order = (th->Theta == 0.5) ? 2 : 1;
893   PetscFunctionReturn(0);
894 }
895 
896 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
897 {
898   TS_Theta *th = (TS_Theta*)ts->data;
899 
900   PetscFunctionBegin;
901   *endpoint = th->endpoint;
902   PetscFunctionReturn(0);
903 }
904 
905 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
906 {
907   TS_Theta *th = (TS_Theta*)ts->data;
908 
909   PetscFunctionBegin;
910   th->endpoint = flg;
911   PetscFunctionReturn(0);
912 }
913 
914 #if defined(PETSC_HAVE_COMPLEX)
915 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
916 {
917   PetscComplex   z   = xr + xi*PETSC_i,f;
918   TS_Theta       *th = (TS_Theta*)ts->data;
919   const PetscReal one = 1.0;
920 
921   PetscFunctionBegin;
922   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
923   *yr = PetscRealPartComplex(f);
924   *yi = PetscImaginaryPartComplex(f);
925   PetscFunctionReturn(0);
926 }
927 #endif
928 
929 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
930 {
931   TS_Theta     *th = (TS_Theta*)ts->data;
932 
933   PetscFunctionBegin;
934   if (ns) *ns = 1;
935   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
936   PetscFunctionReturn(0);
937 }
938 
939 /* ------------------------------------------------------------ */
940 /*MC
941       TSTHETA - DAE solver using the implicit Theta method
942 
943    Level: beginner
944 
945    Options Database:
946 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
947 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
948 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
949 
950    Notes:
951 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
952 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
953 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
954 
955    This method can be applied to DAE.
956 
957    This method is cast as a 1-stage implicit Runge-Kutta method.
958 
959 .vb
960   Theta | Theta
961   -------------
962         |  1
963 .ve
964 
965    For the default Theta=0.5, this is also known as the implicit midpoint rule.
966 
967    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
968 
969 .vb
970   0 | 0         0
971   1 | 1-Theta   Theta
972   -------------------
973     | 1-Theta   Theta
974 .ve
975 
976    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
977 
978    To apply a diagonally implicit RK method to DAE, the stage formula
979 
980 $  Y_i = X + h sum_j a_ij Y'_j
981 
982    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
983 
984 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
985 
986 M*/
987 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
988 {
989   TS_Theta       *th;
990   PetscErrorCode ierr;
991 
992   PetscFunctionBegin;
993   ts->ops->reset           = TSReset_Theta;
994   ts->ops->destroy         = TSDestroy_Theta;
995   ts->ops->view            = TSView_Theta;
996   ts->ops->setup           = TSSetUp_Theta;
997   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
998   ts->ops->step            = TSStep_Theta;
999   ts->ops->interpolate     = TSInterpolate_Theta;
1000   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
1001   ts->ops->rollback        = TSRollBack_Theta;
1002   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
1003   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
1004   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1005 #if defined(PETSC_HAVE_COMPLEX)
1006   ts->ops->linearstability = TSComputeLinearStability_Theta;
1007 #endif
1008   ts->ops->getstages       = TSGetStages_Theta;
1009   ts->ops->adjointstep     = TSAdjointStep_Theta;
1010   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1011   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
1012   ts->default_adapt_type   = TSADAPTNONE;
1013   ts->ops->forwardsetup    = TSForwardSetUp_Theta;
1014   ts->ops->forwardstep     = TSForwardStep_Theta;
1015 
1016   ts->usessnes = PETSC_TRUE;
1017 
1018   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1019   ts->data = (void*)th;
1020 
1021   th->VecsDeltaLam   = NULL;
1022   th->VecsDeltaMu    = NULL;
1023   th->VecsSensiTemp  = NULL;
1024 
1025   th->extrapolate = PETSC_FALSE;
1026   th->Theta       = 0.5;
1027   th->order       = 2;
1028   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1030   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1031   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 /*@
1036   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
1037 
1038   Not Collective
1039 
1040   Input Parameter:
1041 .  ts - timestepping context
1042 
1043   Output Parameter:
1044 .  theta - stage abscissa
1045 
1046   Note:
1047   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
1048 
1049   Level: Advanced
1050 
1051 .seealso: TSThetaSetTheta()
1052 @*/
1053 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
1054 {
1055   PetscErrorCode ierr;
1056 
1057   PetscFunctionBegin;
1058   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1059   PetscValidPointer(theta,2);
1060   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 /*@
1065   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
1066 
1067   Not Collective
1068 
1069   Input Parameter:
1070 +  ts - timestepping context
1071 -  theta - stage abscissa
1072 
1073   Options Database:
1074 .  -ts_theta_theta <theta>
1075 
1076   Level: Intermediate
1077 
1078 .seealso: TSThetaGetTheta()
1079 @*/
1080 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
1081 {
1082   PetscErrorCode ierr;
1083 
1084   PetscFunctionBegin;
1085   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1086   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 /*@
1091   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1092 
1093   Not Collective
1094 
1095   Input Parameter:
1096 .  ts - timestepping context
1097 
1098   Output Parameter:
1099 .  endpoint - PETSC_TRUE when using the endpoint variant
1100 
1101   Level: Advanced
1102 
1103 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
1104 @*/
1105 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
1106 {
1107   PetscErrorCode ierr;
1108 
1109   PetscFunctionBegin;
1110   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1111   PetscValidPointer(endpoint,2);
1112   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 /*@
1117   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1118 
1119   Not Collective
1120 
1121   Input Parameter:
1122 +  ts - timestepping context
1123 -  flg - PETSC_TRUE to use the endpoint variant
1124 
1125   Options Database:
1126 .  -ts_theta_endpoint <flg>
1127 
1128   Level: Intermediate
1129 
1130 .seealso: TSTHETA, TSCN
1131 @*/
1132 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1133 {
1134   PetscErrorCode ierr;
1135 
1136   PetscFunctionBegin;
1137   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1138   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 /*
1143  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1144  * The creation functions for these specializations are below.
1145  */
1146 
1147 static PetscErrorCode TSSetUp_BEuler(TS ts)
1148 {
1149   TS_Theta       *th = (TS_Theta*)ts->data;
1150   PetscErrorCode ierr;
1151 
1152   PetscFunctionBegin;
1153   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");
1154   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");
1155   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1156   PetscFunctionReturn(0);
1157 }
1158 
1159 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1160 {
1161   PetscFunctionBegin;
1162   PetscFunctionReturn(0);
1163 }
1164 
1165 /*MC
1166       TSBEULER - ODE solver using the implicit backward Euler method
1167 
1168   Level: beginner
1169 
1170   Notes:
1171   TSBEULER is equivalent to TSTHETA with Theta=1.0
1172 
1173 $  -ts_type theta -ts_theta_theta 1.0
1174 
1175 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1176 
1177 M*/
1178 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1179 {
1180   PetscErrorCode ierr;
1181 
1182   PetscFunctionBegin;
1183   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1184   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
1185   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
1186   ts->ops->setup = TSSetUp_BEuler;
1187   ts->ops->view  = TSView_BEuler;
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 static PetscErrorCode TSSetUp_CN(TS ts)
1192 {
1193   TS_Theta       *th = (TS_Theta*)ts->data;
1194   PetscErrorCode ierr;
1195 
1196   PetscFunctionBegin;
1197   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");
1198   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");
1199   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1200   PetscFunctionReturn(0);
1201 }
1202 
1203 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1204 {
1205   PetscFunctionBegin;
1206   PetscFunctionReturn(0);
1207 }
1208 
1209 /*MC
1210       TSCN - ODE solver using the implicit Crank-Nicolson method.
1211 
1212   Level: beginner
1213 
1214   Notes:
1215   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1216 
1217 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1218 
1219 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1220 
1221 M*/
1222 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1223 {
1224   PetscErrorCode ierr;
1225 
1226   PetscFunctionBegin;
1227   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1228   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1229   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
1230   ts->ops->setup = TSSetUp_CN;
1231   ts->ops->view  = TSView_CN;
1232   PetscFunctionReturn(0);
1233 }
1234