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