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