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