xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision a4cab8963d92c31bbf27668ad7bbb2395036c865)
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                 VecStage,*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 ? ts->time_step : (1-th->Theta)*ts->time_step); /* time_step is negative*/
254 
255   ierr = TSPreStep(ts);CHKERRQ(ierr);
256 
257   /* Build RHS */
258   if (th->endpoint && th->Theta != 1.0) { /* This formulation assumes linear time-independent mass matrix */
259     shift = -1./((th->Theta-1.0)*ts->time_step);
260     ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
261     for (nadj=0; nadj<ts->numberadjs; nadj++) {
262       ierr  = MatMultTranspose(J,ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
263       ierr = VecScale(VecSensiTemp[nadj],(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
264     }
265   }else { /* Assume mass matrix to be identity for now */
266     for (nadj=0; nadj<ts->numberadjs; nadj++) {
267       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
268       ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
269     }
270   }
271   /* Build LHS */
272   shift = -1./(th->Theta*ts->time_step);
273   VecStage = (th->endpoint) ? th->X0 : th->X;
274   ierr = TSComputeIJacobian(ts,th->stage_time,VecStage,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
275   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
276 
277   /* Solve LHS X = RHS */
278   for (nadj=0; nadj<ts->numberadjs; nadj++) {
279     ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr);
280   }
281   if(th->endpoint) {
282     for (nadj=0; nadj<ts->numberadjs; nadj++) {
283       ierr = VecCopy(VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
284     }
285   }else {
286     shift = -1./(th->Theta*ts->time_step);
287     for (nadj=0; nadj<ts->numberadjs; nadj++) {
288       ierr = VecAXPBYPCZ(VecSensiTemp[nadj],shift,-shift,0,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
289       ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr);
290     }
291   }
292 
293   ts->ptime += ts->time_step;
294   ts->steps++;
295   th->status = TS_STEP_COMPLETE;
296 
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "TSInterpolate_Theta"
302 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
303 {
304   TS_Theta       *th   = (TS_Theta*)ts->data;
305   PetscReal      alpha = t - ts->ptime;
306   PetscErrorCode ierr;
307 
308   PetscFunctionBegin;
309   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
310   if (th->endpoint) alpha *= th->Theta;
311   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
312   PetscFunctionReturn(0);
313 }
314 
315 /*------------------------------------------------------------*/
316 #undef __FUNCT__
317 #define __FUNCT__ "TSReset_Theta"
318 static PetscErrorCode TSReset_Theta(TS ts)
319 {
320   TS_Theta       *th = (TS_Theta*)ts->data;
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
325   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
326   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
327   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
328   if(ts->reverse_mode) {
329     ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
330     if(th->VecDeltaMu) {
331       ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
332     }
333     ierr = VecDestroyVecs(ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
334   }
335   PetscFunctionReturn(0);
336 }
337 
338 #undef __FUNCT__
339 #define __FUNCT__ "TSDestroy_Theta"
340 static PetscErrorCode TSDestroy_Theta(TS ts)
341 {
342   PetscErrorCode ierr;
343 
344   PetscFunctionBegin;
345   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
346   ierr = PetscFree(ts->data);CHKERRQ(ierr);
347   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
348   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
349   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
350   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
351   PetscFunctionReturn(0);
352 }
353 
354 /*
355   This defines the nonlinear equation that is to be solved with SNES
356   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
357 */
358 #undef __FUNCT__
359 #define __FUNCT__ "SNESTSFormFunction_Theta"
360 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
361 {
362   TS_Theta       *th = (TS_Theta*)ts->data;
363   PetscErrorCode ierr;
364   Vec            X0,Xdot;
365   DM             dm,dmsave;
366   PetscReal      shift = 1./(th->Theta*ts->time_step);
367 
368   PetscFunctionBegin;
369   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
370   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
371   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
372   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
373 
374   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
375   dmsave = ts->dm;
376   ts->dm = dm;
377   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
378   ts->dm = dmsave;
379   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "SNESTSFormJacobian_Theta"
385 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
386 {
387   TS_Theta       *th = (TS_Theta*)ts->data;
388   PetscErrorCode ierr;
389   Vec            Xdot;
390   DM             dm,dmsave;
391   PetscReal      shift = 1./(th->Theta*ts->time_step);
392 
393   PetscFunctionBegin;
394   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
395 
396   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
397   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
398 
399   dmsave = ts->dm;
400   ts->dm = dm;
401   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
402   ts->dm = dmsave;
403   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "TSSetUp_Theta"
409 static PetscErrorCode TSSetUp_Theta(TS ts)
410 {
411   TS_Theta       *th = (TS_Theta*)ts->data;
412   PetscErrorCode ierr;
413   SNES           snes;
414   TSAdapt        adapt;
415   DM             dm;
416 
417   PetscFunctionBegin;
418   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
419   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
420   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
421   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
422   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
423   if (dm) {
424     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
425     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
426   }
427   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
428   else th->order = 1;
429 
430   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
431   if (!th->adapt) {
432     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
433   }
434   if (ts->reverse_mode) {
435     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
436     if(ts->vecs_sensip) {
437       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
438     }
439     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
440   }
441   PetscFunctionReturn(0);
442 }
443 /*------------------------------------------------------------*/
444 
445 #undef __FUNCT__
446 #define __FUNCT__ "TSSetFromOptions_Theta"
447 static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
448 {
449   TS_Theta       *th = (TS_Theta*)ts->data;
450   PetscErrorCode ierr;
451 
452   PetscFunctionBegin;
453   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
454   {
455     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
456     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
457     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
458     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
459     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
460   }
461   ierr = PetscOptionsTail();CHKERRQ(ierr);
462   PetscFunctionReturn(0);
463 }
464 
465 #undef __FUNCT__
466 #define __FUNCT__ "TSView_Theta"
467 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
468 {
469   TS_Theta       *th = (TS_Theta*)ts->data;
470   PetscBool      iascii;
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
475   if (iascii) {
476     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
477     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
478   }
479   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "TSThetaGetTheta_Theta"
485 PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
486 {
487   TS_Theta *th = (TS_Theta*)ts->data;
488 
489   PetscFunctionBegin;
490   *theta = th->Theta;
491   PetscFunctionReturn(0);
492 }
493 
494 #undef __FUNCT__
495 #define __FUNCT__ "TSThetaSetTheta_Theta"
496 PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
497 {
498   TS_Theta *th = (TS_Theta*)ts->data;
499 
500   PetscFunctionBegin;
501   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
502   th->Theta = theta;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "TSThetaGetEndpoint_Theta"
508 PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
509 {
510   TS_Theta *th = (TS_Theta*)ts->data;
511 
512   PetscFunctionBegin;
513   *endpoint = th->endpoint;
514   PetscFunctionReturn(0);
515 }
516 
517 #undef __FUNCT__
518 #define __FUNCT__ "TSThetaSetEndpoint_Theta"
519 PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
520 {
521   TS_Theta *th = (TS_Theta*)ts->data;
522 
523   PetscFunctionBegin;
524   th->endpoint = flg;
525   PetscFunctionReturn(0);
526 }
527 
528 #if defined(PETSC_HAVE_COMPLEX)
529 #undef __FUNCT__
530 #define __FUNCT__ "TSComputeLinearStability_Theta"
531 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
532 {
533   PetscComplex z   = xr + xi*PETSC_i,f;
534   TS_Theta     *th = (TS_Theta*)ts->data;
535   const PetscReal one = 1.0;
536 
537   PetscFunctionBegin;
538   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
539   *yr = PetscRealPartComplex(f);
540   *yi = PetscImaginaryPartComplex(f);
541   PetscFunctionReturn(0);
542 }
543 #endif
544 
545 #undef __FUNCT__
546 #define __FUNCT__ "TSGetStages_Theta"
547 static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
548 {
549   TS_Theta     *th = (TS_Theta*)ts->data;
550 
551   PetscFunctionBegin;
552   *ns = 1;
553   if(Y) {
554     if(th->endpoint) { /* return the first (explicit) stage X0 for checkpointing */
555       *Y  = &(th->X0);
556     }else { /* return the stage value*/
557       *Y  = &(th->X);
558     }
559   }
560   PetscFunctionReturn(0);
561 }
562 
563 /* ------------------------------------------------------------ */
564 /*MC
565       TSTHETA - DAE solver using the implicit Theta method
566 
567    Level: beginner
568 
569    Options Database:
570       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
571       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
572       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
573 
574    Notes:
575 $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
576 $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
577 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
578 
579 
580 
581    This method can be applied to DAE.
582 
583    This method is cast as a 1-stage implicit Runge-Kutta method.
584 
585 .vb
586   Theta | Theta
587   -------------
588         |  1
589 .ve
590 
591    For the default Theta=0.5, this is also known as the implicit midpoint rule.
592 
593    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
594 
595 .vb
596   0 | 0         0
597   1 | 1-Theta   Theta
598   -------------------
599     | 1-Theta   Theta
600 .ve
601 
602    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
603 
604    To apply a diagonally implicit RK method to DAE, the stage formula
605 
606 $  Y_i = X + h sum_j a_ij Y'_j
607 
608    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
609 
610 .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
611 
612 M*/
613 #undef __FUNCT__
614 #define __FUNCT__ "TSCreate_Theta"
615 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
616 {
617   TS_Theta       *th;
618   PetscErrorCode ierr;
619 
620   PetscFunctionBegin;
621   ts->ops->reset          = TSReset_Theta;
622   ts->ops->destroy        = TSDestroy_Theta;
623   ts->ops->view           = TSView_Theta;
624   ts->ops->setup          = TSSetUp_Theta;
625   ts->ops->step           = TSStep_Theta;
626   ts->ops->interpolate    = TSInterpolate_Theta;
627   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
628   ts->ops->rollback       = TSRollBack_Theta;
629   ts->ops->setfromoptions = TSSetFromOptions_Theta;
630   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
631   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
632 #if defined(PETSC_HAVE_COMPLEX)
633   ts->ops->linearstability = TSComputeLinearStability_Theta;
634 #endif
635   ts->ops->getstages      = TSGetStages_Theta;
636   ts->ops->stepadj        = TSStepAdj_Theta;
637 
638   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
639   ts->data = (void*)th;
640 
641   th->extrapolate = PETSC_FALSE;
642   th->Theta       = 0.5;
643   th->ccfl        = 1.0;
644   th->adapt       = PETSC_FALSE;
645   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
646   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
647   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
648   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
649   PetscFunctionReturn(0);
650 }
651 
652 #undef __FUNCT__
653 #define __FUNCT__ "TSThetaGetTheta"
654 /*@
655   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
656 
657   Not Collective
658 
659   Input Parameter:
660 .  ts - timestepping context
661 
662   Output Parameter:
663 .  theta - stage abscissa
664 
665   Note:
666   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
667 
668   Level: Advanced
669 
670 .seealso: TSThetaSetTheta()
671 @*/
672 PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
673 {
674   PetscErrorCode ierr;
675 
676   PetscFunctionBegin;
677   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
678   PetscValidPointer(theta,2);
679   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "TSThetaSetTheta"
685 /*@
686   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
687 
688   Not Collective
689 
690   Input Parameter:
691 +  ts - timestepping context
692 -  theta - stage abscissa
693 
694   Options Database:
695 .  -ts_theta_theta <theta>
696 
697   Level: Intermediate
698 
699 .seealso: TSThetaGetTheta()
700 @*/
701 PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
707   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 #undef __FUNCT__
712 #define __FUNCT__ "TSThetaGetEndpoint"
713 /*@
714   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
715 
716   Not Collective
717 
718   Input Parameter:
719 .  ts - timestepping context
720 
721   Output Parameter:
722 .  endpoint - PETSC_TRUE when using the endpoint variant
723 
724   Level: Advanced
725 
726 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
727 @*/
728 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
729 {
730   PetscErrorCode ierr;
731 
732   PetscFunctionBegin;
733   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
734   PetscValidPointer(endpoint,2);
735   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
736   PetscFunctionReturn(0);
737 }
738 
739 #undef __FUNCT__
740 #define __FUNCT__ "TSThetaSetEndpoint"
741 /*@
742   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
743 
744   Not Collective
745 
746   Input Parameter:
747 +  ts - timestepping context
748 -  flg - PETSC_TRUE to use the endpoint variant
749 
750   Options Database:
751 .  -ts_theta_endpoint <flg>
752 
753   Level: Intermediate
754 
755 .seealso: TSTHETA, TSCN
756 @*/
757 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
758 {
759   PetscErrorCode ierr;
760 
761   PetscFunctionBegin;
762   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
763   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
764   PetscFunctionReturn(0);
765 }
766 
767 /*
768  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
769  * The creation functions for these specializations are below.
770  */
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "TSView_BEuler"
774 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
775 {
776   PetscErrorCode ierr;
777 
778   PetscFunctionBegin;
779   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
780   PetscFunctionReturn(0);
781 }
782 
783 /*MC
784       TSBEULER - ODE solver using the implicit backward Euler method
785 
786   Level: beginner
787 
788   Notes:
789   TSBEULER is equivalent to TSTHETA with Theta=1.0
790 
791 $  -ts_type theta -ts_theta_theta 1.
792 
793 .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
794 
795 M*/
796 #undef __FUNCT__
797 #define __FUNCT__ "TSCreate_BEuler"
798 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
799 {
800   PetscErrorCode ierr;
801 
802   PetscFunctionBegin;
803   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
804   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
805   ts->ops->view = TSView_BEuler;
806   PetscFunctionReturn(0);
807 }
808 
809 #undef __FUNCT__
810 #define __FUNCT__ "TSView_CN"
811 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
812 {
813   PetscErrorCode ierr;
814 
815   PetscFunctionBegin;
816   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
817   PetscFunctionReturn(0);
818 }
819 
820 /*MC
821       TSCN - ODE solver using the implicit Crank-Nicolson method.
822 
823   Level: beginner
824 
825   Notes:
826   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
827 
828 $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
829 
830 .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
831 
832 M*/
833 #undef __FUNCT__
834 #define __FUNCT__ "TSCreate_CN"
835 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
836 {
837   PetscErrorCode ierr;
838 
839   PetscFunctionBegin;
840   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
841   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
842   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
843   ts->ops->view = TSView_CN;
844   PetscFunctionReturn(0);
845 }
846