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