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