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