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