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