xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 9f4d3c52fa2fe0bb72fec4f4e85d8e495867af35)
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 timed with Jacobian transpose */
28   Vec          *VecsDeltaFwdSensi;       /* 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          *VecsIntegralSensi0;      /* backup for roll-backs due to events */
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 TS_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 = TS_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       ntlm,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   for (ntlm=0;ntlm<th->num_tlm;ntlm++) {
444     ierr = VecCopy(th->VecsFwdSensi0[ntlm],ts->vecs_fwdsensipacked[ntlm]);CHKERRQ(ierr);
445   }
446   if (ts->vecs_integral_sensi) {
447     for (ncost=0;ncost<ts->numcost;ncost++) {
448       ierr = VecCopy(th->VecsIntegralSensi0[ncost],ts->vecs_integral_sensi[ncost]);CHKERRQ(ierr);
449     }
450   }
451   if (ts->vecs_integral_sensip) {
452     for (ncost=0;ncost<ts->numcost;ncost++) {
453       ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr);
454     }
455   }
456   PetscFunctionReturn(0);
457 }
458 
459 static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt numtlm,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative)
460 {
461   PetscInt    ntlm,low,high;
462   PetscScalar *v;
463   PetscErrorCode ierr;
464 
465   PetscFunctionBegin;
466 
467   ierr = VecGetOwnershipRange(VecIntegrandDerivative,&low,&high);CHKERRQ(ierr);
468   ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
469   for (ntlm=low; ntlm<high; ntlm++) {
470     ierr = VecDot(DRncostDY,s[ntlm],&v[ntlm-low]);CHKERRQ(ierr);
471   }
472   ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
473   if (DRncostDP) {
474     ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr);
475   }
476   PetscFunctionReturn(0);
477 }
478 
479 static PetscErrorCode TSForwardStep_Theta(TS ts)
480 {
481   TS_Theta       *th = (TS_Theta*)ts->data;
482   Vec            *VecsDeltaFwdSensi  = th->VecsDeltaFwdSensi;
483   PetscInt       ntlm,ncost,np;
484   KSP            ksp;
485   Mat            J,Jp;
486   PetscReal      shift;
487   PetscErrorCode ierr;
488 
489   PetscFunctionBegin;
490   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
491     ierr = VecCopy(ts->vecs_fwdsensipacked[ntlm],th->VecsFwdSensi0[ntlm]);CHKERRQ(ierr);
492   }
493   for (ncost=0; ncost<ts->numcost; ncost++) {
494     if (ts->vecs_integral_sensi) {
495       ierr = VecCopy(ts->vecs_integral_sensi[ncost],th->VecsIntegralSensi0[ncost]);CHKERRQ(ierr);
496     }
497     if (ts->vecs_integral_sensip) {
498       ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr);
499     }
500   }
501 
502   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
503   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
504 
505   /* Build RHS */
506   if (th->endpoint) { /* 2-stage method*/
507     shift = 1./((th->Theta-1.)*th->time_step);
508     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
509 
510     for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
511       ierr = MatMult(J,ts->vecs_fwdsensipacked[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr);
512       ierr = VecScale(VecsDeltaFwdSensi[ntlm],(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
513     }
514     /* Add the f_p forcing terms */
515     ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr);
516     for (np=0; np<ts->num_parameters; np++) {
517       ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],(1.-th->Theta)/th->Theta,ts->vecs_jacp[np]);CHKERRQ(ierr);
518     }
519     ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr);
520     for (np=0; np<ts->num_parameters; np++) {
521       ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],1,ts->vecs_jacp[np]);CHKERRQ(ierr);
522     }
523   } else { /* 1-stage method */
524     shift = 0.0;
525     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
526     for (ntlm=0; ntlm<ts->num_parameters; ntlm++) {
527       ierr = MatMult(J,ts->vecs_fwdsensipacked[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr);
528       ierr = VecScale(VecsDeltaFwdSensi[ntlm],-1);CHKERRQ(ierr);
529     }
530     /* Add the f_p forcing terms */
531     ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr);
532     for (np=0; np<ts->num_parameters; np++) {
533       ierr = VecAXPY(VecsDeltaFwdSensi[np+ts->num_initialvalues],1,ts->vecs_jacp[np]);CHKERRQ(ierr);
534     }
535   }
536 
537   /* Build LHS */
538   shift = 1/(th->Theta*th->time_step);
539   if (th->endpoint) {
540     ierr  = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
541   } else {
542     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
543   }
544   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
545 
546   /*
547     Evaluate the first stage of integral gradients with the 2-stage method:
548     drdy|t_n*S(t_n) + drdp|t_n
549     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
550     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
551    */
552   if (th->endpoint) { /* 2-stage method only */
553     if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
554       ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
555     }
556     if (ts->vecs_integral_sensip) {
557       ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
558       for (ncost=0; ncost<ts->numcost; ncost++) {
559         ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],&ts->vecs_fwdsensipacked[ts->num_initialvalues],ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
560         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
561       }
562     }
563     if (ts->vecs_integral_sensi) {
564       for (ncost=0; ncost<ts->numcost; ncost++) {
565         ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensipacked,NULL,th->VecIntegralSensiTemp);CHKERRQ(ierr);
566         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr);
567       }
568     }
569   }
570 
571   /* Solve the tangent linear equation for forward sensitivities to parameters */
572   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
573     if (th->endpoint) {
574       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ntlm],ts->vecs_fwdsensipacked[ntlm]);CHKERRQ(ierr);
575     } else {
576       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ntlm],VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr);
577       ierr = VecAXPY(ts->vecs_fwdsensipacked[ntlm],1.,VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr);
578     }
579   }
580   /* Evaluate the second stage of integral gradients with the 2-stage method:
581     drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
582   */
583   if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
584     ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
585   }
586   if (ts->vecs_integral_sensip) {
587     ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
588     for (ncost=0; ncost<ts->numcost; ncost++) {
589       ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],&ts->vecs_fwdsensipacked[ts->num_initialvalues],ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
590       if (th->endpoint) {
591         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
592       } else {
593         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
594       }
595     }
596   }
597   if (ts->vecs_integral_sensi) {
598     for (ncost=0; ncost<ts->numcost; ncost++) {
599       ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensipacked,NULL,th->VecIntegralSensiTemp);CHKERRQ(ierr);
600       if (th->endpoint) {
601         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr);
602       } else {
603         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr);
604       }
605     }
606   }
607 
608   if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */
609     for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
610       ierr = VecAXPY(ts->vecs_fwdsensipacked[ntlm],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ntlm]);CHKERRQ(ierr);
611     }
612   }
613   PetscFunctionReturn(0);
614 }
615 
616 /*------------------------------------------------------------*/
617 static PetscErrorCode TSReset_Theta(TS ts)
618 {
619   TS_Theta       *th = (TS_Theta*)ts->data;
620   PetscErrorCode ierr;
621 
622   PetscFunctionBegin;
623   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
624   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
625   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
626   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
627 
628   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
629   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
630 
631   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
632   if (ts->forward_solve) {
633     ierr = VecDestroyVecs(th->num_tlm,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
634     ierr = VecDestroyVecs(th->num_tlm,&th->VecsFwdSensi0);CHKERRQ(ierr);
635     if (ts->vecs_integral_sensi) {
636       ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr);
637       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
638     }
639     if (ts->vecs_integral_sensip) {
640       ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
641       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
642     }
643   }
644   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
645   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
646   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
647 
648   PetscFunctionReturn(0);
649 }
650 
651 static PetscErrorCode TSDestroy_Theta(TS ts)
652 {
653   PetscErrorCode ierr;
654 
655   PetscFunctionBegin;
656   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
657   if (ts->dm) {
658     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
659     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
660   }
661   ierr = PetscFree(ts->data);CHKERRQ(ierr);
662   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
663   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
664   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
665   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
666   PetscFunctionReturn(0);
667 }
668 
669 /*
670   This defines the nonlinear equation that is to be solved with SNES
671   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
672 */
673 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
674 {
675   TS_Theta       *th = (TS_Theta*)ts->data;
676   PetscErrorCode ierr;
677   Vec            X0,Xdot;
678   DM             dm,dmsave;
679   PetscReal      shift = 1/(th->Theta*ts->time_step);
680 
681   PetscFunctionBegin;
682   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
683   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
684   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
685   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
686 
687   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
688   dmsave = ts->dm;
689   ts->dm = dm;
690   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
691   ts->dm = dmsave;
692   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
693   PetscFunctionReturn(0);
694 }
695 
696 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
697 {
698   TS_Theta       *th = (TS_Theta*)ts->data;
699   PetscErrorCode ierr;
700   Vec            Xdot;
701   DM             dm,dmsave;
702   PetscReal      shift = 1/(th->Theta*ts->time_step);
703 
704   PetscFunctionBegin;
705   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
706   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
707   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
708 
709   dmsave = ts->dm;
710   ts->dm = dm;
711   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
712   ts->dm = dmsave;
713   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
714   PetscFunctionReturn(0);
715 }
716 
717 static PetscErrorCode TSForwardSetUp_Theta(TS ts)
718 {
719   TS_Theta       *th = (TS_Theta*)ts->data;
720   PetscErrorCode ierr;
721 
722   PetscFunctionBegin;
723   /* combine sensitivities to parameters and sensitivities to initial values into one array */
724   th->num_tlm = ts->num_parameters + ts->num_initialvalues;
725   ierr = VecDuplicateVecs(ts->vecs_fwdsensipacked[0],th->num_tlm,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
726   if (ts->vecs_integral_sensi) {
727     ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr);
728   }
729   if (ts->vecs_integral_sensip) {
730     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
731   }
732   /* backup sensitivity results for roll-backs */
733   ierr = VecDuplicateVecs(ts->vecs_fwdsensipacked[0],th->num_tlm,&th->VecsFwdSensi0);CHKERRQ(ierr);
734   if (ts->vecs_integral_sensi) {
735     ierr = VecDuplicateVecs(ts->vecs_integral_sensi[0],ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
736   }
737   if (ts->vecs_integral_sensip) {
738     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
739   }
740   PetscFunctionReturn(0);
741 }
742 
743 static PetscErrorCode TSSetUp_Theta(TS ts)
744 {
745   TS_Theta       *th = (TS_Theta*)ts->data;
746   PetscBool      match;
747   PetscErrorCode ierr;
748 
749   PetscFunctionBegin;
750   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
751     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
752   }
753   if (!th->X) {
754     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
755   }
756   if (!th->Xdot) {
757     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
758   }
759   if (!th->X0) {
760     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
761   }
762   if (th->endpoint) {
763     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
764   }
765 
766   th->order = (th->Theta == 0.5) ? 2 : 1;
767 
768   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
769   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
770   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
771 
772   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
773   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
774   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
775   if (!match) {
776     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
777     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
778   }
779   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
780   PetscFunctionReturn(0);
781 }
782 
783 /*------------------------------------------------------------*/
784 
785 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
786 {
787   TS_Theta       *th = (TS_Theta*)ts->data;
788   PetscErrorCode ierr;
789 
790   PetscFunctionBegin;
791   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
792   if(ts->vecs_sensip) {
793     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
794   }
795   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
796   PetscFunctionReturn(0);
797 }
798 
799 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
800 {
801   TS_Theta       *th = (TS_Theta*)ts->data;
802   PetscErrorCode ierr;
803 
804   PetscFunctionBegin;
805   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
806   {
807     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
808     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
809     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
810   }
811   ierr = PetscOptionsTail();CHKERRQ(ierr);
812   PetscFunctionReturn(0);
813 }
814 
815 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
816 {
817   TS_Theta       *th = (TS_Theta*)ts->data;
818   PetscBool      iascii;
819   PetscErrorCode ierr;
820 
821   PetscFunctionBegin;
822   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
823   if (iascii) {
824     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
825     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
826   }
827   PetscFunctionReturn(0);
828 }
829 
830 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
831 {
832   TS_Theta *th = (TS_Theta*)ts->data;
833 
834   PetscFunctionBegin;
835   *theta = th->Theta;
836   PetscFunctionReturn(0);
837 }
838 
839 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
840 {
841   TS_Theta *th = (TS_Theta*)ts->data;
842 
843   PetscFunctionBegin;
844   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
845   th->Theta = theta;
846   th->order = (th->Theta == 0.5) ? 2 : 1;
847   PetscFunctionReturn(0);
848 }
849 
850 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
851 {
852   TS_Theta *th = (TS_Theta*)ts->data;
853 
854   PetscFunctionBegin;
855   *endpoint = th->endpoint;
856   PetscFunctionReturn(0);
857 }
858 
859 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
860 {
861   TS_Theta *th = (TS_Theta*)ts->data;
862 
863   PetscFunctionBegin;
864   th->endpoint = flg;
865   PetscFunctionReturn(0);
866 }
867 
868 #if defined(PETSC_HAVE_COMPLEX)
869 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
870 {
871   PetscComplex   z   = xr + xi*PETSC_i,f;
872   TS_Theta       *th = (TS_Theta*)ts->data;
873   const PetscReal one = 1.0;
874 
875   PetscFunctionBegin;
876   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
877   *yr = PetscRealPartComplex(f);
878   *yi = PetscImaginaryPartComplex(f);
879   PetscFunctionReturn(0);
880 }
881 #endif
882 
883 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
884 {
885   TS_Theta     *th = (TS_Theta*)ts->data;
886 
887   PetscFunctionBegin;
888   if (ns) *ns = 1;
889   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
890   PetscFunctionReturn(0);
891 }
892 
893 /* ------------------------------------------------------------ */
894 /*MC
895       TSTHETA - DAE solver using the implicit Theta method
896 
897    Level: beginner
898 
899    Options Database:
900 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
901 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
902 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
903 
904    Notes:
905 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
906 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
907 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
908 
909    This method can be applied to DAE.
910 
911    This method is cast as a 1-stage implicit Runge-Kutta method.
912 
913 .vb
914   Theta | Theta
915   -------------
916         |  1
917 .ve
918 
919    For the default Theta=0.5, this is also known as the implicit midpoint rule.
920 
921    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
922 
923 .vb
924   0 | 0         0
925   1 | 1-Theta   Theta
926   -------------------
927     | 1-Theta   Theta
928 .ve
929 
930    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
931 
932    To apply a diagonally implicit RK method to DAE, the stage formula
933 
934 $  Y_i = X + h sum_j a_ij Y'_j
935 
936    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
937 
938 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
939 
940 M*/
941 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
942 {
943   TS_Theta       *th;
944   PetscErrorCode ierr;
945 
946   PetscFunctionBegin;
947   ts->ops->reset           = TSReset_Theta;
948   ts->ops->destroy         = TSDestroy_Theta;
949   ts->ops->view            = TSView_Theta;
950   ts->ops->setup           = TSSetUp_Theta;
951   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
952   ts->ops->step            = TSStep_Theta;
953   ts->ops->interpolate     = TSInterpolate_Theta;
954   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
955   ts->ops->rollback        = TSRollBack_Theta;
956   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
957   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
958   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
959 #if defined(PETSC_HAVE_COMPLEX)
960   ts->ops->linearstability = TSComputeLinearStability_Theta;
961 #endif
962   ts->ops->getstages       = TSGetStages_Theta;
963   ts->ops->adjointstep     = TSAdjointStep_Theta;
964   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
965   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
966   ts->default_adapt_type   = TSADAPTNONE;
967   ts->ops->forwardsetup    = TSForwardSetUp_Theta;
968   ts->ops->forwardstep     = TSForwardStep_Theta;
969 
970   ts->usessnes = PETSC_TRUE;
971 
972   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
973   ts->data = (void*)th;
974 
975   th->VecsDeltaLam   = NULL;
976   th->VecsDeltaMu    = NULL;
977   th->VecsSensiTemp  = NULL;
978 
979   th->extrapolate = PETSC_FALSE;
980   th->Theta       = 0.5;
981   th->order       = 2;
982   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
983   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
984   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
985   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
986   PetscFunctionReturn(0);
987 }
988 
989 /*@
990   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
991 
992   Not Collective
993 
994   Input Parameter:
995 .  ts - timestepping context
996 
997   Output Parameter:
998 .  theta - stage abscissa
999 
1000   Note:
1001   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
1002 
1003   Level: Advanced
1004 
1005 .seealso: TSThetaSetTheta()
1006 @*/
1007 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
1008 {
1009   PetscErrorCode ierr;
1010 
1011   PetscFunctionBegin;
1012   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1013   PetscValidPointer(theta,2);
1014   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 /*@
1019   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
1020 
1021   Not Collective
1022 
1023   Input Parameter:
1024 +  ts - timestepping context
1025 -  theta - stage abscissa
1026 
1027   Options Database:
1028 .  -ts_theta_theta <theta>
1029 
1030   Level: Intermediate
1031 
1032 .seealso: TSThetaGetTheta()
1033 @*/
1034 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
1035 {
1036   PetscErrorCode ierr;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1040   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 /*@
1045   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1046 
1047   Not Collective
1048 
1049   Input Parameter:
1050 .  ts - timestepping context
1051 
1052   Output Parameter:
1053 .  endpoint - PETSC_TRUE when using the endpoint variant
1054 
1055   Level: Advanced
1056 
1057 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
1058 @*/
1059 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
1060 {
1061   PetscErrorCode ierr;
1062 
1063   PetscFunctionBegin;
1064   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1065   PetscValidPointer(endpoint,2);
1066   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
1067   PetscFunctionReturn(0);
1068 }
1069 
1070 /*@
1071   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1072 
1073   Not Collective
1074 
1075   Input Parameter:
1076 +  ts - timestepping context
1077 -  flg - PETSC_TRUE to use the endpoint variant
1078 
1079   Options Database:
1080 .  -ts_theta_endpoint <flg>
1081 
1082   Level: Intermediate
1083 
1084 .seealso: TSTHETA, TSCN
1085 @*/
1086 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1087 {
1088   PetscErrorCode ierr;
1089 
1090   PetscFunctionBegin;
1091   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1092   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1093   PetscFunctionReturn(0);
1094 }
1095 
1096 /*
1097  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1098  * The creation functions for these specializations are below.
1099  */
1100 
1101 static PetscErrorCode TSSetUp_BEuler(TS ts)
1102 {
1103   TS_Theta       *th = (TS_Theta*)ts->data;
1104   PetscErrorCode ierr;
1105 
1106   PetscFunctionBegin;
1107   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");
1108   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");
1109   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1114 {
1115   PetscFunctionBegin;
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /*MC
1120       TSBEULER - ODE solver using the implicit backward Euler method
1121 
1122   Level: beginner
1123 
1124   Notes:
1125   TSBEULER is equivalent to TSTHETA with Theta=1.0
1126 
1127 $  -ts_type theta -ts_theta_theta 1.0
1128 
1129 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1130 
1131 M*/
1132 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1133 {
1134   PetscErrorCode ierr;
1135 
1136   PetscFunctionBegin;
1137   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1138   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
1139   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
1140   ts->ops->setup = TSSetUp_BEuler;
1141   ts->ops->view  = TSView_BEuler;
1142   PetscFunctionReturn(0);
1143 }
1144 
1145 static PetscErrorCode TSSetUp_CN(TS ts)
1146 {
1147   TS_Theta       *th = (TS_Theta*)ts->data;
1148   PetscErrorCode ierr;
1149 
1150   PetscFunctionBegin;
1151   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");
1152   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");
1153   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1154   PetscFunctionReturn(0);
1155 }
1156 
1157 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1158 {
1159   PetscFunctionBegin;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 /*MC
1164       TSCN - ODE solver using the implicit Crank-Nicolson method.
1165 
1166   Level: beginner
1167 
1168   Notes:
1169   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1170 
1171 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1172 
1173 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1174 
1175 M*/
1176 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1177 {
1178   PetscErrorCode ierr;
1179 
1180   PetscFunctionBegin;
1181   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1182   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1183   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
1184   ts->ops->setup = TSSetUp_CN;
1185   ts->ops->view  = TSView_CN;
1186   PetscFunctionReturn(0);
1187 }
1188