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