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