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