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