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