xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 73a71a0f875e342bef2af9dbfc05b9368f8a1bc4)
1 /*
2   Code for timestepping with implicit Theta method
3 */
4 #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
5 #include <petscsnesfas.h>
6 
7 typedef struct {
8   Vec       X,Xdot;                   /* Storage for one stage */
9   Vec       affine;                   /* Affine vector needed for residual at beginning of step */
10   PetscBool extrapolate;
11   PetscBool endpoint;
12   PetscReal Theta;
13   PetscReal shift;
14   PetscReal stage_time;
15 } TS_Theta;
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "TSThetaGetX0AndXdot"
19 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
20 {
21   TS_Theta       *th = (TS_Theta*)ts->data;
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   if (X0) {
26     if (dm && dm != ts->dm) {
27       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
28     } else *X0 = ts->vec_sol;
29   }
30   if (Xdot) {
31     if (dm && dm != ts->dm) {
32       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
33     } else *Xdot = th->Xdot;
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "TSThetaRestoreX0AndXdot"
41 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
42 {
43   PetscErrorCode ierr;
44 
45   PetscFunctionBegin;
46   if (X0) {
47     if (dm && dm != ts->dm) {
48       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
49     }
50   }
51   if (Xdot) {
52     if (dm && dm != ts->dm) {
53       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
54     }
55   }
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "DMCoarsenHook_TSTheta"
61 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
62 {
63 
64   PetscFunctionBegin;
65   PetscFunctionReturn(0);
66 }
67 
68 #undef __FUNCT__
69 #define __FUNCT__ "DMRestrictHook_TSTheta"
70 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
71 {
72   TS ts = (TS)ctx;
73   PetscErrorCode ierr;
74   Vec X0,Xdot,X0_c,Xdot_c;
75 
76   PetscFunctionBegin;
77   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
78   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
79   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
80   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
81   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
82   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
83   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
84   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "TSStep_Theta"
90 static PetscErrorCode TSStep_Theta(TS ts)
91 {
92   TS_Theta            *th = (TS_Theta*)ts->data;
93   PetscInt            its,lits;
94   PetscReal           next_time_step;
95   SNESConvergedReason snesreason;
96   PetscErrorCode      ierr;
97 
98   PetscFunctionBegin;
99   next_time_step = ts->time_step;
100   th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
101   th->shift = 1./(th->Theta*ts->time_step);
102   ierr = TSPreStep(ts);CHKERRQ(ierr);
103   ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
104 
105   if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
106     ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
107     if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
108     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
109     ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
110   }
111   if (th->extrapolate) {
112     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
113   } else {
114     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
115   }
116   ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
117   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
118   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
119   ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
120   ts->snes_its += its; ts->ksp_its += lits;
121   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
122     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
123     ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
124     PetscFunctionReturn(0);
125   }
126   if (th->endpoint) {
127     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
128   } else {
129     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
130     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
131   }
132   ts->ptime += ts->time_step;
133   ts->time_step = next_time_step;
134   ts->steps++;
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "TSInterpolate_Theta"
140 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
141 {
142   TS_Theta       *th = (TS_Theta*)ts->data;
143   PetscReal      alpha = t - ts->ptime;
144   PetscErrorCode ierr;
145 
146   PetscFunctionBegin;
147   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
148   if (th->endpoint) alpha *= th->Theta;
149   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
150   PetscFunctionReturn(0);
151 }
152 
153 /*------------------------------------------------------------*/
154 #undef __FUNCT__
155 #define __FUNCT__ "TSReset_Theta"
156 static PetscErrorCode TSReset_Theta(TS ts)
157 {
158   TS_Theta       *th = (TS_Theta*)ts->data;
159   PetscErrorCode  ierr;
160 
161   PetscFunctionBegin;
162   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
163   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
164   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNCT__
169 #define __FUNCT__ "TSDestroy_Theta"
170 static PetscErrorCode TSDestroy_Theta(TS ts)
171 {
172   PetscErrorCode  ierr;
173 
174   PetscFunctionBegin;
175   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
176   ierr = PetscFree(ts->data);CHKERRQ(ierr);
177   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
178   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
179   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
180   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 /*
185   This defines the nonlinear equation that is to be solved with SNES
186   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
187 */
188 #undef __FUNCT__
189 #define __FUNCT__ "SNESTSFormFunction_Theta"
190 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
191 {
192   TS_Theta       *th = (TS_Theta*)ts->data;
193   PetscErrorCode ierr;
194   Vec            X0,Xdot;
195   DM             dm,dmsave;
196 
197   PetscFunctionBegin;
198   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
199   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
200   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
201   ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr);
202 
203   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
204   dmsave = ts->dm;
205   ts->dm = dm;
206   ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
207   ts->dm = dmsave;
208   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "SNESTSFormJacobian_Theta"
214 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
215 {
216   TS_Theta       *th = (TS_Theta*)ts->data;
217   PetscErrorCode ierr;
218   Vec            Xdot;
219   DM             dm,dmsave;
220 
221   PetscFunctionBegin;
222   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
223 
224   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
225   ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
226 
227   dmsave = ts->dm;
228   ts->dm = dm;
229   ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
230   ts->dm = dmsave;
231   ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "TSSetUp_Theta"
237 static PetscErrorCode TSSetUp_Theta(TS ts)
238 {
239   TS_Theta       *th = (TS_Theta*)ts->data;
240   PetscErrorCode ierr;
241   SNES           snes;
242   DM             dm;
243 
244   PetscFunctionBegin;
245   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
246   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
247   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
248   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
249   if (dm) {
250     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
251   }
252   PetscFunctionReturn(0);
253 }
254 /*------------------------------------------------------------*/
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "TSSetFromOptions_Theta"
258 static PetscErrorCode TSSetFromOptions_Theta(TS ts)
259 {
260   TS_Theta       *th = (TS_Theta*)ts->data;
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
265   {
266     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
267     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
268     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr);
269     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
270   }
271   ierr = PetscOptionsTail();CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "TSView_Theta"
277 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
278 {
279   TS_Theta       *th = (TS_Theta*)ts->data;
280   PetscBool       iascii;
281   PetscErrorCode  ierr;
282 
283   PetscFunctionBegin;
284   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
285   if (iascii) {
286     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
287     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
288   }
289   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293 EXTERN_C_BEGIN
294 #undef __FUNCT__
295 #define __FUNCT__ "TSThetaGetTheta_Theta"
296 PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
297 {
298   TS_Theta *th = (TS_Theta*)ts->data;
299 
300   PetscFunctionBegin;
301   *theta = th->Theta;
302   PetscFunctionReturn(0);
303 }
304 
305 #undef __FUNCT__
306 #define __FUNCT__ "TSThetaSetTheta_Theta"
307 PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
308 {
309   TS_Theta *th = (TS_Theta*)ts->data;
310 
311   PetscFunctionBegin;
312   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
313   th->Theta = theta;
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNCT__
318 #define __FUNCT__ "TSThetaGetEndpoint_Theta"
319 PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
320 {
321   TS_Theta *th = (TS_Theta*)ts->data;
322 
323   PetscFunctionBegin;
324   *endpoint = th->endpoint;
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "TSThetaSetEndpoint_Theta"
330 PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
331 {
332   TS_Theta *th = (TS_Theta*)ts->data;
333 
334   PetscFunctionBegin;
335   th->endpoint = flg;
336   PetscFunctionReturn(0);
337 }
338 EXTERN_C_END
339 
340 /* ------------------------------------------------------------ */
341 /*MC
342       TSTHETA - DAE solver using the implicit Theta method
343 
344    Level: beginner
345 
346    Options Database:
347       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1);  Theta = 1.0 (
348       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
349       -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method
350 
351    Notes:
352 $  -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER)
353 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
354 
355 
356 
357    This method can be applied to DAE.
358 
359    This method is cast as a 1-stage implicit Runge-Kutta method.
360 
361 .vb
362   Theta | Theta
363   -------------
364         |  1
365 .ve
366 
367    For the default Theta=0.5, this is also known as the implicit midpoint rule.
368 
369    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
370 
371 .vb
372   0 | 0         0
373   1 | 1-Theta   Theta
374   -------------------
375     | 1-Theta   Theta
376 .ve
377 
378    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
379 
380    To apply a diagonally implicit RK method to DAE, the stage formula
381 
382 $  Y_i = X + h sum_j a_ij Y'_j
383 
384    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
385 
386 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
387 
388 M*/
389 EXTERN_C_BEGIN
390 #undef __FUNCT__
391 #define __FUNCT__ "TSCreate_Theta"
392 PetscErrorCode  TSCreate_Theta(TS ts)
393 {
394   TS_Theta       *th;
395   PetscErrorCode ierr;
396 
397   PetscFunctionBegin;
398   ts->ops->reset          = TSReset_Theta;
399   ts->ops->destroy        = TSDestroy_Theta;
400   ts->ops->view           = TSView_Theta;
401   ts->ops->setup          = TSSetUp_Theta;
402   ts->ops->step           = TSStep_Theta;
403   ts->ops->interpolate    = TSInterpolate_Theta;
404   ts->ops->setfromoptions = TSSetFromOptions_Theta;
405   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
406   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
407 
408   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
409   ts->data = (void*)th;
410 
411   th->extrapolate = PETSC_FALSE;
412   th->Theta       = 0.5;
413 
414   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
415   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
416   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
417   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420 EXTERN_C_END
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "TSThetaGetTheta"
424 /*@
425   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
426 
427   Not Collective
428 
429   Input Parameter:
430 .  ts - timestepping context
431 
432   Output Parameter:
433 .  theta - stage abscissa
434 
435   Note:
436   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
437 
438   Level: Advanced
439 
440 .seealso: TSThetaSetTheta()
441 @*/
442 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
443 {
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
448   PetscValidPointer(theta,2);
449   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
450   PetscFunctionReturn(0);
451 }
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "TSThetaSetTheta"
455 /*@
456   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
457 
458   Not Collective
459 
460   Input Parameter:
461 +  ts - timestepping context
462 -  theta - stage abscissa
463 
464   Options Database:
465 .  -ts_theta_theta <theta>
466 
467   Level: Intermediate
468 
469 .seealso: TSThetaGetTheta()
470 @*/
471 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
472 {
473   PetscErrorCode ierr;
474 
475   PetscFunctionBegin;
476   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
477   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 #undef __FUNCT__
482 #define __FUNCT__ "TSThetaGetEndpoint"
483 /*@
484   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
485 
486   Not Collective
487 
488   Input Parameter:
489 .  ts - timestepping context
490 
491   Output Parameter:
492 .  endpoint - PETSC_TRUE when using the endpoint variant
493 
494   Level: Advanced
495 
496 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
497 @*/
498 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
499 {
500   PetscErrorCode ierr;
501 
502   PetscFunctionBegin;
503   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
504   PetscValidPointer(endpoint,2);
505   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
506   PetscFunctionReturn(0);
507 }
508 
509 #undef __FUNCT__
510 #define __FUNCT__ "TSThetaSetEndpoint"
511 /*@
512   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
513 
514   Not Collective
515 
516   Input Parameter:
517 +  ts - timestepping context
518 -  flg - PETSC_TRUE to use the endpoint variant
519 
520   Options Database:
521 .  -ts_theta_endpoint <flg>
522 
523   Level: Intermediate
524 
525 .seealso: TSTHETA, TSCN
526 @*/
527 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
528 {
529   PetscErrorCode ierr;
530 
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
533   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
534   PetscFunctionReturn(0);
535 }
536 
537 /*
538  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
539  * The creation functions for these specializations are below.
540  */
541 
542 #undef __FUNCT__
543 #define __FUNCT__ "TSView_BEuler"
544 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
545 {
546   PetscErrorCode ierr;
547 
548   PetscFunctionBegin;
549   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
550   PetscFunctionReturn(0);
551 }
552 
553 /*MC
554       TSBEULER - ODE solver using the implicit backward Euler method
555 
556   Level: beginner
557 
558   Notes:
559   TSCN is equivalent to TSTHETA with Theta=1.0
560 
561 $  -ts_type theta -ts_theta_theta 1.
562 
563 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
564 
565 M*/
566 EXTERN_C_BEGIN
567 #undef __FUNCT__
568 #define __FUNCT__ "TSCreate_BEuler"
569 PetscErrorCode  TSCreate_BEuler(TS ts)
570 {
571   PetscErrorCode ierr;
572 
573   PetscFunctionBegin;
574   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
575   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
576   ts->ops->view = TSView_BEuler;
577   PetscFunctionReturn(0);
578 }
579 EXTERN_C_END
580 
581 #undef __FUNCT__
582 #define __FUNCT__ "TSView_CN"
583 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
584 {
585   PetscErrorCode ierr;
586 
587   PetscFunctionBegin;
588   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 /*MC
593       TSCN - ODE solver using the implicit Crank-Nicolson method.
594 
595   Level: beginner
596 
597   Notes:
598   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
599 
600 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
601 
602 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
603 
604 M*/
605 EXTERN_C_BEGIN
606 #undef __FUNCT__
607 #define __FUNCT__ "TSCreate_CN"
608 PetscErrorCode  TSCreate_CN(TS ts)
609 {
610   PetscErrorCode ierr;
611 
612   PetscFunctionBegin;
613   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
614   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
615   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
616   ts->ops->view = TSView_CN;
617   PetscFunctionReturn(0);
618 }
619 EXTERN_C_END
620