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