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