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