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