xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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__ "TSStep_Theta"
152 static PetscErrorCode TSStep_Theta(TS ts)
153 {
154   TS_Theta            *th = (TS_Theta*)ts->data;
155   PetscInt            its,lits,reject,next_scheme;
156   PetscReal           next_time_step;
157   SNESConvergedReason snesreason;
158   PetscErrorCode      ierr;
159   TSAdapt             adapt;
160   PetscBool           accept;
161 
162   PetscFunctionBegin;
163   th->status = TS_STEP_INCOMPLETE;
164   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
165   for (reject=0; reject<ts->max_reject && !ts->reason && th->status != TS_STEP_COMPLETE; reject++,ts->reject++) {
166     PetscReal shift = 1./(th->Theta*ts->time_step);
167     next_time_step = ts->time_step;
168     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
169     ierr = TSPreStep(ts);CHKERRQ(ierr);
170     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
171 
172     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
173       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
174       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
175       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
176       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
177     }
178     if (th->extrapolate) {
179       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
180     } else {
181       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
182     }
183     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
184     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
185     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
186     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
187     ts->snes_its += its; ts->ksp_its += lits;
188     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
189     ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
190     if (!accept) continue;
191     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
192     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
193     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
194     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
195     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
196     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
197 
198     if (accept) {
199       /* ignore next_scheme for now */
200       ts->ptime    += ts->time_step;
201       ts->time_step = next_time_step;
202       ts->steps++;
203       th->status = TS_STEP_COMPLETE;
204     } else {                    /* Roll back the current step */
205       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
206       ts->time_step = next_time_step;
207       th->status    = TS_STEP_INCOMPLETE;
208     }
209   }
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "TSInterpolate_Theta"
215 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
216 {
217   TS_Theta       *th   = (TS_Theta*)ts->data;
218   PetscReal      alpha = t - ts->ptime;
219   PetscErrorCode ierr;
220 
221   PetscFunctionBegin;
222   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
223   if (th->endpoint) alpha *= th->Theta;
224   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
225   PetscFunctionReturn(0);
226 }
227 
228 /*------------------------------------------------------------*/
229 #undef __FUNCT__
230 #define __FUNCT__ "TSReset_Theta"
231 static PetscErrorCode TSReset_Theta(TS ts)
232 {
233   TS_Theta       *th = (TS_Theta*)ts->data;
234   PetscErrorCode ierr;
235 
236   PetscFunctionBegin;
237   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
238   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
239   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
240   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "TSDestroy_Theta"
246 static PetscErrorCode TSDestroy_Theta(TS ts)
247 {
248   PetscErrorCode ierr;
249 
250   PetscFunctionBegin;
251   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
252   ierr = PetscFree(ts->data);CHKERRQ(ierr);
253   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C","",NULL);CHKERRQ(ierr);
254   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C","",NULL);CHKERRQ(ierr);
255   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C","",NULL);CHKERRQ(ierr);
256   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C","",NULL);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 /*
261   This defines the nonlinear equation that is to be solved with SNES
262   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
263 */
264 #undef __FUNCT__
265 #define __FUNCT__ "SNESTSFormFunction_Theta"
266 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
267 {
268   TS_Theta       *th = (TS_Theta*)ts->data;
269   PetscErrorCode ierr;
270   Vec            X0,Xdot;
271   DM             dm,dmsave;
272   PetscReal      shift = 1./(th->Theta*ts->time_step);
273 
274   PetscFunctionBegin;
275   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
276   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
277   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
278   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
279 
280   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
281   dmsave = ts->dm;
282   ts->dm = dm;
283   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
284   ts->dm = dmsave;
285   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "SNESTSFormJacobian_Theta"
291 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
292 {
293   TS_Theta       *th = (TS_Theta*)ts->data;
294   PetscErrorCode ierr;
295   Vec            Xdot;
296   DM             dm,dmsave;
297   PetscReal      shift = 1./(th->Theta*ts->time_step);
298 
299   PetscFunctionBegin;
300   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
301 
302   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
303   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
304 
305   dmsave = ts->dm;
306   ts->dm = dm;
307   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
308   ts->dm = dmsave;
309   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "TSSetUp_Theta"
315 static PetscErrorCode TSSetUp_Theta(TS ts)
316 {
317   TS_Theta       *th = (TS_Theta*)ts->data;
318   PetscErrorCode ierr;
319   SNES           snes;
320   DM             dm;
321 
322   PetscFunctionBegin;
323   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
324   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
325   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
326   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
327   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
328   if (dm) {
329     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
330     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
331   }
332   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
333   else th->order = 1;
334 
335   if (!th->adapt) {
336     TSAdapt adapt;
337     ierr = TSAdaptDestroy(&ts->adapt);CHKERRQ(ierr);
338     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
339     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
340   }
341   PetscFunctionReturn(0);
342 }
343 /*------------------------------------------------------------*/
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "TSSetFromOptions_Theta"
347 static PetscErrorCode TSSetFromOptions_Theta(TS ts)
348 {
349   TS_Theta       *th = (TS_Theta*)ts->data;
350   PetscErrorCode ierr;
351 
352   PetscFunctionBegin;
353   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
354   {
355     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
356     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
357     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
358     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
359     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
360   }
361   ierr = PetscOptionsTail();CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "TSView_Theta"
367 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
368 {
369   TS_Theta       *th = (TS_Theta*)ts->data;
370   PetscBool      iascii;
371   PetscErrorCode ierr;
372 
373   PetscFunctionBegin;
374   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
375   if (iascii) {
376     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
377     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
378   }
379   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
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(PetscObjectComm((PetscObject)ts),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 
428 #if defined(PETSC_HAVE_COMPLEX)
429 #undef __FUNCT__
430 #define __FUNCT__ "TSComputeLinearStability_Theta"
431 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
432 {
433   PetscComplex z   = xr + xi*PETSC_i,f;
434   TS_Theta     *th = (TS_Theta*)ts->data;
435   const PetscReal one = 1.0;
436 
437   PetscFunctionBegin;
438   f   = (one + (one - th->Theta)*z)/(one - 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 #undef __FUNCT__
497 #define __FUNCT__ "TSCreate_Theta"
498 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
499 {
500   TS_Theta       *th;
501   PetscErrorCode ierr;
502 
503   PetscFunctionBegin;
504   ts->ops->reset          = TSReset_Theta;
505   ts->ops->destroy        = TSDestroy_Theta;
506   ts->ops->view           = TSView_Theta;
507   ts->ops->setup          = TSSetUp_Theta;
508   ts->ops->step           = TSStep_Theta;
509   ts->ops->interpolate    = TSInterpolate_Theta;
510   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
511   ts->ops->setfromoptions = TSSetFromOptions_Theta;
512   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
513   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
514 #if defined(PETSC_HAVE_COMPLEX)
515   ts->ops->linearstability = TSComputeLinearStability_Theta;
516 #endif
517 
518   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
519   ts->data = (void*)th;
520 
521   th->extrapolate = PETSC_FALSE;
522   th->Theta       = 0.5;
523   th->ccfl        = 1.0;
524   th->adapt       = PETSC_FALSE;
525   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
526   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
527   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
528   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "TSThetaGetTheta"
534 /*@
535   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
536 
537   Not Collective
538 
539   Input Parameter:
540 .  ts - timestepping context
541 
542   Output Parameter:
543 .  theta - stage abscissa
544 
545   Note:
546   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
547 
548   Level: Advanced
549 
550 .seealso: TSThetaSetTheta()
551 @*/
552 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
553 {
554   PetscErrorCode ierr;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
558   PetscValidPointer(theta,2);
559   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "TSThetaSetTheta"
565 /*@
566   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
567 
568   Not Collective
569 
570   Input Parameter:
571 +  ts - timestepping context
572 -  theta - stage abscissa
573 
574   Options Database:
575 .  -ts_theta_theta <theta>
576 
577   Level: Intermediate
578 
579 .seealso: TSThetaGetTheta()
580 @*/
581 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
582 {
583   PetscErrorCode ierr;
584 
585   PetscFunctionBegin;
586   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
587   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
588   PetscFunctionReturn(0);
589 }
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "TSThetaGetEndpoint"
593 /*@
594   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
595 
596   Not Collective
597 
598   Input Parameter:
599 .  ts - timestepping context
600 
601   Output Parameter:
602 .  endpoint - PETSC_TRUE when using the endpoint variant
603 
604   Level: Advanced
605 
606 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
607 @*/
608 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
609 {
610   PetscErrorCode ierr;
611 
612   PetscFunctionBegin;
613   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
614   PetscValidPointer(endpoint,2);
615   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "TSThetaSetEndpoint"
621 /*@
622   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
623 
624   Not Collective
625 
626   Input Parameter:
627 +  ts - timestepping context
628 -  flg - PETSC_TRUE to use the endpoint variant
629 
630   Options Database:
631 .  -ts_theta_endpoint <flg>
632 
633   Level: Intermediate
634 
635 .seealso: TSTHETA, TSCN
636 @*/
637 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
638 {
639   PetscErrorCode ierr;
640 
641   PetscFunctionBegin;
642   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
643   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
644   PetscFunctionReturn(0);
645 }
646 
647 /*
648  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
649  * The creation functions for these specializations are below.
650  */
651 
652 #undef __FUNCT__
653 #define __FUNCT__ "TSView_BEuler"
654 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
655 {
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
660   PetscFunctionReturn(0);
661 }
662 
663 /*MC
664       TSBEULER - ODE solver using the implicit backward Euler method
665 
666   Level: beginner
667 
668   Notes:
669   TSBEULER is equivalent to TSTHETA with Theta=1.0
670 
671 $  -ts_type theta -ts_theta_theta 1.
672 
673 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
674 
675 M*/
676 #undef __FUNCT__
677 #define __FUNCT__ "TSCreate_BEuler"
678 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
679 {
680   PetscErrorCode ierr;
681 
682   PetscFunctionBegin;
683   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
684   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
685   ts->ops->view = TSView_BEuler;
686   PetscFunctionReturn(0);
687 }
688 
689 #undef __FUNCT__
690 #define __FUNCT__ "TSView_CN"
691 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
692 {
693   PetscErrorCode ierr;
694 
695   PetscFunctionBegin;
696   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
697   PetscFunctionReturn(0);
698 }
699 
700 /*MC
701       TSCN - ODE solver using the implicit Crank-Nicolson method.
702 
703   Level: beginner
704 
705   Notes:
706   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
707 
708 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
709 
710 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
711 
712 M*/
713 #undef __FUNCT__
714 #define __FUNCT__ "TSCreate_CN"
715 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
716 {
717   PetscErrorCode ierr;
718 
719   PetscFunctionBegin;
720   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
721   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
722   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
723   ts->ops->view = TSView_CN;
724   PetscFunctionReturn(0);
725 }
726