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