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