xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 4fc747eaadbeca11629f314a99edccbc2ed7b3d3)
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   PetscErrorCode ierr;
430 
431   PetscFunctionBegin;
432   /* Cannot compute LTE in first step or in restart after event */
433   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
434   /* Compute LTE using backward differences with non-constant time step */
435   {
436     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
437     PetscReal   a = 1 + h_prev/h;
438     PetscScalar scal[3]; Vec vecs[3];
439     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
440     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
441     ierr = VecCopy(X,Y);CHKERRQ(ierr);
442     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
443     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);CHKERRQ(ierr);
444   }
445   if (order) *order = 2;
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "TSRollBack_Theta"
451 static PetscErrorCode TSRollBack_Theta(TS ts)
452 {
453   TS_Theta       *th = (TS_Theta*)ts->data;
454   PetscErrorCode ierr;
455 
456   PetscFunctionBegin;
457   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
458   if (ts->vec_costintegral && ts->costintegralfwd) {
459     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 /*------------------------------------------------------------*/
465 #undef __FUNCT__
466 #define __FUNCT__ "TSReset_Theta"
467 static PetscErrorCode TSReset_Theta(TS ts)
468 {
469   TS_Theta       *th = (TS_Theta*)ts->data;
470   PetscErrorCode ierr;
471 
472   PetscFunctionBegin;
473   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
474   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
475   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
476   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
477 
478   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
479   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
480 
481   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
482   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
483   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
484   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "TSDestroy_Theta"
490 static PetscErrorCode TSDestroy_Theta(TS ts)
491 {
492   PetscErrorCode ierr;
493 
494   PetscFunctionBegin;
495   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
496   ierr = PetscFree(ts->data);CHKERRQ(ierr);
497   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
498   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
499   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
500   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
501   PetscFunctionReturn(0);
502 }
503 
504 /*
505   This defines the nonlinear equation that is to be solved with SNES
506   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
507 */
508 #undef __FUNCT__
509 #define __FUNCT__ "SNESTSFormFunction_Theta"
510 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
511 {
512   TS_Theta       *th = (TS_Theta*)ts->data;
513   PetscErrorCode ierr;
514   Vec            X0,Xdot;
515   DM             dm,dmsave;
516   PetscReal      shift = 1/(th->Theta*ts->time_step);
517 
518   PetscFunctionBegin;
519   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
520   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
521   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
522   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
523 
524   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
525   dmsave = ts->dm;
526   ts->dm = dm;
527   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
528   ts->dm = dmsave;
529   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "SNESTSFormJacobian_Theta"
535 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
536 {
537   TS_Theta       *th = (TS_Theta*)ts->data;
538   PetscErrorCode ierr;
539   Vec            Xdot;
540   DM             dm,dmsave;
541   PetscReal      shift = 1/(th->Theta*ts->time_step);
542 
543   PetscFunctionBegin;
544   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
545   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
546   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
547 
548   dmsave = ts->dm;
549   ts->dm = dm;
550   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
551   ts->dm = dmsave;
552   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
553   PetscFunctionReturn(0);
554 }
555 
556 #undef __FUNCT__
557 #define __FUNCT__ "TSSetUp_Theta"
558 static PetscErrorCode TSSetUp_Theta(TS ts)
559 {
560   TS_Theta       *th = (TS_Theta*)ts->data;
561   PetscErrorCode ierr;
562 
563   PetscFunctionBegin;
564   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
565     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
566   }
567   if (!th->X) {
568     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
569   }
570   if (!th->Xdot) {
571     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
572   }
573   if (!th->X0) {
574     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
575   }
576   if (th->endpoint) {
577     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
578   }
579 
580   th->order = (th->Theta == 0.5) ? 2 : 1;
581 
582   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
583   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
584   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
585 
586   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
587   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
588   if (!th->adapt) {
589     ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr);
590   } else {
591     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
592     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
593     if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
594       ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
595   }
596 
597   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
598   PetscFunctionReturn(0);
599 }
600 
601 /*------------------------------------------------------------*/
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "TSAdjointSetUp_Theta"
605 static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
606 {
607   TS_Theta       *th = (TS_Theta*)ts->data;
608   PetscErrorCode ierr;
609 
610   PetscFunctionBegin;
611   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
612   if(ts->vecs_sensip) {
613     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
614   }
615   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
616   PetscFunctionReturn(0);
617 }
618 /*------------------------------------------------------------*/
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "TSSetFromOptions_Theta"
622 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
623 {
624   TS_Theta       *th = (TS_Theta*)ts->data;
625   PetscErrorCode ierr;
626 
627   PetscFunctionBegin;
628   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
629   {
630     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
631     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
632     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
633     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
634   }
635   ierr = PetscOptionsTail();CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 #undef __FUNCT__
640 #define __FUNCT__ "TSView_Theta"
641 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
642 {
643   TS_Theta       *th = (TS_Theta*)ts->data;
644   PetscBool      iascii;
645   PetscErrorCode ierr;
646 
647   PetscFunctionBegin;
648   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
649   if (iascii) {
650     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
651     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
652   }
653   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
654   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "TSThetaGetTheta_Theta"
660 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
661 {
662   TS_Theta *th = (TS_Theta*)ts->data;
663 
664   PetscFunctionBegin;
665   *theta = th->Theta;
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "TSThetaSetTheta_Theta"
671 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
672 {
673   TS_Theta *th = (TS_Theta*)ts->data;
674 
675   PetscFunctionBegin;
676   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
677   th->Theta = theta;
678   th->order = (th->Theta == 0.5) ? 2 : 1;
679   PetscFunctionReturn(0);
680 }
681 
682 #undef __FUNCT__
683 #define __FUNCT__ "TSThetaGetEndpoint_Theta"
684 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
685 {
686   TS_Theta *th = (TS_Theta*)ts->data;
687 
688   PetscFunctionBegin;
689   *endpoint = th->endpoint;
690   PetscFunctionReturn(0);
691 }
692 
693 #undef __FUNCT__
694 #define __FUNCT__ "TSThetaSetEndpoint_Theta"
695 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
696 {
697   TS_Theta *th = (TS_Theta*)ts->data;
698 
699   PetscFunctionBegin;
700   th->endpoint = flg;
701   PetscFunctionReturn(0);
702 }
703 
704 #if defined(PETSC_HAVE_COMPLEX)
705 #undef __FUNCT__
706 #define __FUNCT__ "TSComputeLinearStability_Theta"
707 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
708 {
709   PetscComplex   z   = xr + xi*PETSC_i,f;
710   TS_Theta       *th = (TS_Theta*)ts->data;
711   const PetscReal one = 1.0;
712 
713   PetscFunctionBegin;
714   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
715   *yr = PetscRealPartComplex(f);
716   *yi = PetscImaginaryPartComplex(f);
717   PetscFunctionReturn(0);
718 }
719 #endif
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "TSGetStages_Theta"
723 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
724 {
725   TS_Theta     *th = (TS_Theta*)ts->data;
726 
727   PetscFunctionBegin;
728   if (ns) *ns = 1;
729   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
730   PetscFunctionReturn(0);
731 }
732 
733 /* ------------------------------------------------------------ */
734 /*MC
735       TSTHETA - DAE solver using the implicit Theta method
736 
737    Level: beginner
738 
739    Options Database:
740 +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
741 .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
742 .  -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method
743 -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
744 
745    Notes:
746 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
747 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
748 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
749 
750    This method can be applied to DAE.
751 
752    This method is cast as a 1-stage implicit Runge-Kutta method.
753 
754 .vb
755   Theta | Theta
756   -------------
757         |  1
758 .ve
759 
760    For the default Theta=0.5, this is also known as the implicit midpoint rule.
761 
762    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
763 
764 .vb
765   0 | 0         0
766   1 | 1-Theta   Theta
767   -------------------
768     | 1-Theta   Theta
769 .ve
770 
771    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
772 
773    To apply a diagonally implicit RK method to DAE, the stage formula
774 
775 $  Y_i = X + h sum_j a_ij Y'_j
776 
777    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
778 
779 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
780 
781 M*/
782 #undef __FUNCT__
783 #define __FUNCT__ "TSCreate_Theta"
784 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
785 {
786   TS_Theta       *th;
787   PetscErrorCode ierr;
788 
789   PetscFunctionBegin;
790   ts->ops->reset           = TSReset_Theta;
791   ts->ops->destroy         = TSDestroy_Theta;
792   ts->ops->view            = TSView_Theta;
793   ts->ops->setup           = TSSetUp_Theta;
794   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
795   ts->ops->step            = TSStep_Theta;
796   ts->ops->interpolate     = TSInterpolate_Theta;
797   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
798   ts->ops->rollback        = TSRollBack_Theta;
799   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
800   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
801   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
802 #if defined(PETSC_HAVE_COMPLEX)
803   ts->ops->linearstability = TSComputeLinearStability_Theta;
804 #endif
805   ts->ops->getstages       = TSGetStages_Theta;
806   ts->ops->adjointstep     = TSAdjointStep_Theta;
807   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
808   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
809 
810   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
811   ts->data = (void*)th;
812 
813   th->extrapolate = PETSC_FALSE;
814   th->Theta       = 0.5;
815   th->order       = 2;
816   th->adapt       = PETSC_FALSE;
817   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
818   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
819   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
820   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "TSThetaGetTheta"
826 /*@
827   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
828 
829   Not Collective
830 
831   Input Parameter:
832 .  ts - timestepping context
833 
834   Output Parameter:
835 .  theta - stage abscissa
836 
837   Note:
838   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
839 
840   Level: Advanced
841 
842 .seealso: TSThetaSetTheta()
843 @*/
844 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
845 {
846   PetscErrorCode ierr;
847 
848   PetscFunctionBegin;
849   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
850   PetscValidPointer(theta,2);
851   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
852   PetscFunctionReturn(0);
853 }
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "TSThetaSetTheta"
857 /*@
858   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
859 
860   Not Collective
861 
862   Input Parameter:
863 +  ts - timestepping context
864 -  theta - stage abscissa
865 
866   Options Database:
867 .  -ts_theta_theta <theta>
868 
869   Level: Intermediate
870 
871 .seealso: TSThetaGetTheta()
872 @*/
873 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
874 {
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
879   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
880   PetscFunctionReturn(0);
881 }
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "TSThetaGetEndpoint"
885 /*@
886   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
887 
888   Not Collective
889 
890   Input Parameter:
891 .  ts - timestepping context
892 
893   Output Parameter:
894 .  endpoint - PETSC_TRUE when using the endpoint variant
895 
896   Level: Advanced
897 
898 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
899 @*/
900 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
901 {
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906   PetscValidPointer(endpoint,2);
907   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "TSThetaSetEndpoint"
913 /*@
914   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
915 
916   Not Collective
917 
918   Input Parameter:
919 +  ts - timestepping context
920 -  flg - PETSC_TRUE to use the endpoint variant
921 
922   Options Database:
923 .  -ts_theta_endpoint <flg>
924 
925   Level: Intermediate
926 
927 .seealso: TSTHETA, TSCN
928 @*/
929 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
930 {
931   PetscErrorCode ierr;
932 
933   PetscFunctionBegin;
934   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
935   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
936   PetscFunctionReturn(0);
937 }
938 
939 /*
940  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
941  * The creation functions for these specializations are below.
942  */
943 
944 #undef __FUNCT__
945 #define __FUNCT__ "TSSetUp_BEuler"
946 static PetscErrorCode TSSetUp_BEuler(TS ts)
947 {
948   TS_Theta       *th = (TS_Theta*)ts->data;
949   PetscErrorCode ierr;
950 
951   PetscFunctionBegin;
952   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");
953   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");
954   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
955   PetscFunctionReturn(0);
956 }
957 
958 #undef __FUNCT__
959 #define __FUNCT__ "TSView_BEuler"
960 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
961 {
962   PetscErrorCode ierr;
963 
964   PetscFunctionBegin;
965   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
966   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
967   PetscFunctionReturn(0);
968 }
969 
970 /*MC
971       TSBEULER - ODE solver using the implicit backward Euler method
972 
973   Level: beginner
974 
975   Notes:
976   TSBEULER is equivalent to TSTHETA with Theta=1.0
977 
978 $  -ts_type theta -ts_theta_theta 1.0
979 
980 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
981 
982 M*/
983 #undef __FUNCT__
984 #define __FUNCT__ "TSCreate_BEuler"
985 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
986 {
987   PetscErrorCode ierr;
988 
989   PetscFunctionBegin;
990   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
991   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
992   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
993   ts->ops->setup = TSSetUp_BEuler;
994   ts->ops->view  = TSView_BEuler;
995   PetscFunctionReturn(0);
996 }
997 
998 #undef __FUNCT__
999 #define __FUNCT__ "TSSetUp_CN"
1000 static PetscErrorCode TSSetUp_CN(TS ts)
1001 {
1002   TS_Theta       *th = (TS_Theta*)ts->data;
1003   PetscErrorCode ierr;
1004 
1005   PetscFunctionBegin;
1006   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");
1007   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");
1008   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 #undef __FUNCT__
1013 #define __FUNCT__ "TSView_CN"
1014 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1015 {
1016   PetscErrorCode ierr;
1017 
1018   PetscFunctionBegin;
1019   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
1020   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 /*MC
1025       TSCN - ODE solver using the implicit Crank-Nicolson method.
1026 
1027   Level: beginner
1028 
1029   Notes:
1030   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
1031 
1032 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1033 
1034 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1035 
1036 M*/
1037 #undef __FUNCT__
1038 #define __FUNCT__ "TSCreate_CN"
1039 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1040 {
1041   PetscErrorCode ierr;
1042 
1043   PetscFunctionBegin;
1044   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1045   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1046   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
1047   ts->ops->setup = TSSetUp_CN;
1048   ts->ops->view  = TSView_CN;
1049   PetscFunctionReturn(0);
1050 }
1051