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