xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision b1cb36f36a6b7dfc48c2470e77f0d724e5408ff8)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5a055c8caSBarry Smith #include <petscsnes.h>
61e25c274SJed Brown #include <petscdm.h>
73c54f92cSHong Zhang #include <petscmat.h>
8316643e7SJed Brown 
9316643e7SJed Brown typedef struct {
10316643e7SJed Brown   Vec          X,Xdot;                   /* Storage for one stage */
113b1890cdSShri Abhyankar   Vec          X0;                       /* work vector to store X0 */
12eb284becSJed Brown   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
13b5e0532dSHong Zhang   Vec          *VecsDeltaLam;             /* Increment of the adjoint sensitivity w.r.t IC at stage*/
14b5e0532dSHong Zhang   Vec          *VecsDeltaMu;              /* Increment of the adjoint sensitivity w.r.t P at stage*/
15b5e0532dSHong Zhang   Vec          *VecsSensiTemp;            /* Vector to be timed with Jacobian transpose*/
16*b1cb36f3SHong Zhang   Vec          VecCostIntegral0;          /* backup for roll-backs due to events */
17ace3abfcSBarry Smith   PetscBool    extrapolate;
18eb284becSJed Brown   PetscBool    endpoint;
19316643e7SJed Brown   PetscReal    Theta;
20316643e7SJed Brown   PetscReal    stage_time;
213b1890cdSShri Abhyankar   TSStepStatus status;
223b1890cdSShri Abhyankar   char         *name;
233b1890cdSShri Abhyankar   PetscInt     order;
243b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
253b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
26b5e0532dSHong Zhang   PetscReal    ptime;
27*b1cb36f3SHong Zhang   PetscReal    time_step;
28316643e7SJed Brown } TS_Theta;
29316643e7SJed Brown 
30316643e7SJed Brown #undef __FUNCT__
317445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
327445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
337445fe48SJed Brown {
347445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
357445fe48SJed Brown   PetscErrorCode ierr;
367445fe48SJed Brown 
377445fe48SJed Brown   PetscFunctionBegin;
387445fe48SJed Brown   if (X0) {
397445fe48SJed Brown     if (dm && dm != ts->dm) {
400d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
417445fe48SJed Brown     } else *X0 = ts->vec_sol;
427445fe48SJed Brown   }
437445fe48SJed Brown   if (Xdot) {
447445fe48SJed Brown     if (dm && dm != ts->dm) {
450d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
467445fe48SJed Brown     } else *Xdot = th->Xdot;
477445fe48SJed Brown   }
487445fe48SJed Brown   PetscFunctionReturn(0);
497445fe48SJed Brown }
507445fe48SJed Brown 
510d0b770aSPeter Brune #undef __FUNCT__
520d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
530d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
540d0b770aSPeter Brune {
550d0b770aSPeter Brune   PetscErrorCode ierr;
560d0b770aSPeter Brune 
570d0b770aSPeter Brune   PetscFunctionBegin;
580d0b770aSPeter Brune   if (X0) {
590d0b770aSPeter Brune     if (dm && dm != ts->dm) {
600d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
610d0b770aSPeter Brune     }
620d0b770aSPeter Brune   }
630d0b770aSPeter Brune   if (Xdot) {
640d0b770aSPeter Brune     if (dm && dm != ts->dm) {
650d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
660d0b770aSPeter Brune     }
670d0b770aSPeter Brune   }
680d0b770aSPeter Brune   PetscFunctionReturn(0);
690d0b770aSPeter Brune }
700d0b770aSPeter Brune 
717445fe48SJed Brown #undef __FUNCT__
727445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
737445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
747445fe48SJed Brown {
757445fe48SJed Brown 
767445fe48SJed Brown   PetscFunctionBegin;
777445fe48SJed Brown   PetscFunctionReturn(0);
787445fe48SJed Brown }
797445fe48SJed Brown 
807445fe48SJed Brown #undef __FUNCT__
817445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
827445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
837445fe48SJed Brown {
847445fe48SJed Brown   TS             ts = (TS)ctx;
857445fe48SJed Brown   PetscErrorCode ierr;
867445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
877445fe48SJed Brown 
887445fe48SJed Brown   PetscFunctionBegin;
897445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
907445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
917445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
927445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
937445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
950d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
960d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
977445fe48SJed Brown   PetscFunctionReturn(0);
987445fe48SJed Brown }
997445fe48SJed Brown 
1007445fe48SJed Brown #undef __FUNCT__
101258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
102258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
103258e1594SPeter Brune {
104258e1594SPeter Brune 
105258e1594SPeter Brune   PetscFunctionBegin;
106258e1594SPeter Brune   PetscFunctionReturn(0);
107258e1594SPeter Brune }
108258e1594SPeter Brune 
109258e1594SPeter Brune #undef __FUNCT__
110258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
111258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
112258e1594SPeter Brune {
113258e1594SPeter Brune   TS             ts = (TS)ctx;
114258e1594SPeter Brune   PetscErrorCode ierr;
115258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
116258e1594SPeter Brune 
117258e1594SPeter Brune   PetscFunctionBegin;
118258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
119258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
120258e1594SPeter Brune 
121258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune 
124258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126258e1594SPeter Brune 
127258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
128258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
129258e1594SPeter Brune   PetscFunctionReturn(0);
130258e1594SPeter Brune }
131258e1594SPeter Brune 
1323b1890cdSShri Abhyankar #undef __FUNCT__
1333b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1343b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1353b1890cdSShri Abhyankar {
1363b1890cdSShri Abhyankar   PetscErrorCode ierr;
1373b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1383b1890cdSShri Abhyankar 
1393b1890cdSShri Abhyankar   PetscFunctionBegin;
140ce94432eSBarry Smith   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");
1413b1890cdSShri Abhyankar   if (order == th->order) {
1423b1890cdSShri Abhyankar     if (th->endpoint) {
1433b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1443b1890cdSShri Abhyankar     } else {
1453b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1463b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1473b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1483b1890cdSShri Abhyankar     }
1493b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1503b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1513b1890cdSShri Abhyankar   }
1523b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1533b1890cdSShri Abhyankar }
154258e1594SPeter Brune 
155258e1594SPeter Brune #undef __FUNCT__
15624655328SShri #define __FUNCT__ "TSRollBack_Theta"
15724655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15824655328SShri {
15924655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
16024655328SShri   PetscErrorCode ierr;
16124655328SShri 
16224655328SShri   PetscFunctionBegin;
16324655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
16424655328SShri   th->status    = TS_STEP_INCOMPLETE;
165*b1cb36f3SHong Zhang   if (ts->vec_costintegral && ts->costintegralfwd) {
166*b1cb36f3SHong Zhang     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
167*b1cb36f3SHong Zhang   }
168*b1cb36f3SHong Zhang   PetscFunctionReturn(0);
169*b1cb36f3SHong Zhang }
170*b1cb36f3SHong Zhang 
171*b1cb36f3SHong Zhang #undef __FUNCT__
172*b1cb36f3SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_Theta"
173*b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
174*b1cb36f3SHong Zhang {
175*b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
176*b1cb36f3SHong Zhang   PetscErrorCode ierr;
177*b1cb36f3SHong Zhang 
178*b1cb36f3SHong Zhang   PetscFunctionBegin;
179*b1cb36f3SHong Zhang   /* backup cost integral */
180*b1cb36f3SHong Zhang   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
181*b1cb36f3SHong Zhang   if (th->endpoint) {
182*b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
183*b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
184*b1cb36f3SHong Zhang   }
185*b1cb36f3SHong Zhang   ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
186*b1cb36f3SHong Zhang   if (th->endpoint) {
187*b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
188*b1cb36f3SHong Zhang   } else {
189*b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
190*b1cb36f3SHong Zhang   }
191*b1cb36f3SHong Zhang   PetscFunctionReturn(0);
192*b1cb36f3SHong Zhang }
193*b1cb36f3SHong Zhang 
194*b1cb36f3SHong Zhang #undef __FUNCT__
195*b1cb36f3SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_Theta"
196*b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
197*b1cb36f3SHong Zhang {
198*b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
199*b1cb36f3SHong Zhang   PetscErrorCode ierr;
200*b1cb36f3SHong Zhang 
201*b1cb36f3SHong Zhang   PetscFunctionBegin;
202*b1cb36f3SHong Zhang   if (th->endpoint) {
203*b1cb36f3SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
204*b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
205*b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
206*b1cb36f3SHong Zhang     if (th->Theta!=1) {
207*b1cb36f3SHong Zhang       ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
208*b1cb36f3SHong Zhang       ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr);
209*b1cb36f3SHong Zhang     }
210*b1cb36f3SHong Zhang   }else {
211*b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
212*b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
213*b1cb36f3SHong Zhang   }
21424655328SShri   PetscFunctionReturn(0);
21524655328SShri }
21624655328SShri 
21724655328SShri #undef __FUNCT__
218316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
219193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
220316643e7SJed Brown {
221316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
2223b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
2233b1890cdSShri Abhyankar   PetscReal      next_time_step;
2243b1890cdSShri Abhyankar   TSAdapt        adapt;
2254957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
226051f2191SLisandro Dalcin   PetscErrorCode ierr;
227316643e7SJed Brown 
228316643e7SJed Brown   PetscFunctionBegin;
2293b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
2303b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
231051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
232b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
233eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
234b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
235b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
236316643e7SJed Brown 
237eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
238eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
239eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
240eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
241eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
242eb284becSJed Brown     }
243ace68cafSJed Brown     if (th->extrapolate) {
244b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
245ace68cafSJed Brown     } else {
2462b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
247ace68cafSJed Brown     }
248eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
249316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
250316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
2515ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
252051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
253552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
254b295832fSPierre Barbier de Reuille     ierr = TSAdaptCheckStage(adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
255051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
256051f2191SLisandro Dalcin 
2570298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
258051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2593b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
260552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2613b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2620298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2633b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
264051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
265051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
266051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
267051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
268051f2191SLisandro Dalcin       goto reject_step;
269051f2191SLisandro Dalcin     }
2703b1890cdSShri Abhyankar 
271*b1cb36f3SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
272*b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
273*b1cb36f3SHong Zhang       th->time_step = ts->time_step;
27417145e75SHong Zhang     }
2753b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2762b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
277cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
278316643e7SJed Brown     ts->steps++;
2793b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
280051f2191SLisandro Dalcin     break;
281051f2191SLisandro Dalcin 
282051f2191SLisandro Dalcin reject_step:
283051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
284051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
285051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2863b1890cdSShri Abhyankar     }
287051f2191SLisandro Dalcin     continue;
2883b1890cdSShri Abhyankar   }
289316643e7SJed Brown   PetscFunctionReturn(0);
290316643e7SJed Brown }
291316643e7SJed Brown 
292cd652676SJed Brown #undef __FUNCT__
29342f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_Theta"
29442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
2952ca6e920SHong Zhang {
2962ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
297b5e0532dSHong Zhang   Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
2982ca6e920SHong Zhang   PetscInt            nadj;
2992ca6e920SHong Zhang   PetscErrorCode      ierr;
3002ca6e920SHong Zhang   Mat                 J,Jp;
3012ca6e920SHong Zhang   KSP                 ksp;
3022ca6e920SHong Zhang   PetscReal           shift;
3032ca6e920SHong Zhang 
3042ca6e920SHong Zhang   PetscFunctionBegin;
3052ca6e920SHong Zhang 
3062ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
307302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
3082ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
309b5e0532dSHong Zhang 
310b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
3113fd52205SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
312b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
3132ca6e920SHong Zhang 
3142ca6e920SHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
3152ca6e920SHong Zhang 
316a4cab896SHong Zhang   /* Build RHS */
3172c39e106SBarry Smith   if (ts->vec_costintegral) { /* Cost function has an integral term */
31805755b9cSHong Zhang     if (th->endpoint) {
319d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
32005755b9cSHong Zhang     }else {
321d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
32205755b9cSHong Zhang     }
32305755b9cSHong Zhang   }
324abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
325b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
326b5e0532dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
3272c39e106SBarry Smith     if (ts->vec_costintegral) {
328b5e0532dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
32936eaed60SHong Zhang     }
3302ca6e920SHong Zhang   }
3313c54f92cSHong Zhang 
3322ca6e920SHong Zhang   /* Build LHS */
3332ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
3343c54f92cSHong Zhang   if (th->endpoint) {
3353c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3363c54f92cSHong Zhang   }else {
3373c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3383c54f92cSHong Zhang   }
3392ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
3402ca6e920SHong Zhang 
3412ca6e920SHong Zhang   /* Solve LHS X = RHS */
342abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
343b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
3442ca6e920SHong Zhang   }
3453c54f92cSHong Zhang 
34636eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
347b5e0532dSHong Zhang   if(th->endpoint) { /* two-stage case */
348b5e0532dSHong Zhang     if (th->Theta!=1.) {
3493fd52205SHong Zhang       shift = -1./((th->Theta-1.)*ts->time_step);
350b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3512c39e106SBarry Smith       if (ts->vec_costintegral) {
352b5e0532dSHong Zhang         ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
3538263dbbfSHong Zhang       }
354abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
355b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
3562c39e106SBarry Smith         if (ts->vec_costintegral) {
35736eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
35836eaed60SHong Zhang         }
3593c54f92cSHong Zhang         ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3602ca6e920SHong Zhang       }
361b5e0532dSHong Zhang     }else { /* backward Euler */
362b5e0532dSHong Zhang       shift = 0.0;
363b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
364b5e0532dSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
365b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
366b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
367b5e0532dSHong Zhang         if (ts->vec_costintegral) {
368b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
369b5e0532dSHong Zhang         }
370b5e0532dSHong Zhang       }
371b5e0532dSHong Zhang     }
3723fd52205SHong Zhang 
3733fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
3745bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
375abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
376b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
377b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3783fd52205SHong Zhang       }
379b5e0532dSHong Zhang       if (th->Theta!=1.) {
380b5e0532dSHong Zhang         ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
381abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
382b5e0532dSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
383b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
384b5e0532dSHong Zhang         }
3853fd52205SHong Zhang       }
3862c39e106SBarry Smith       if (ts->vec_costintegral) {
387d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
388abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
38936eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
39036eaed60SHong Zhang         }
391b5e0532dSHong Zhang         if (th->Theta!=1.) {
392b5e0532dSHong Zhang           ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
393abc2977eSBarry Smith           for (nadj=0; nadj<ts->numcost; nadj++) {
39436eaed60SHong Zhang             ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
39536eaed60SHong Zhang           }
39636eaed60SHong Zhang         }
3973fd52205SHong Zhang       }
398b5e0532dSHong Zhang     }
3993fd52205SHong Zhang   }else { /* one-stage case */
4003c54f92cSHong Zhang     shift = 0.0;
401a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
4022c39e106SBarry Smith     if (ts->vec_costintegral) {
403d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
4043c54f92cSHong Zhang     }
405abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
406b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
407b5e0532dSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
4082c39e106SBarry Smith       if (ts->vec_costintegral) {
40936eaed60SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
41036eaed60SHong Zhang       }
4112ca6e920SHong Zhang     }
4123fd52205SHong Zhang     if (ts->vecs_sensip) {
4135bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
414abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
415b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
416b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
4173fd52205SHong Zhang       }
4182c39e106SBarry Smith       if (ts->vec_costintegral) {
419d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
420abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
42136eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
42236eaed60SHong Zhang         }
42336eaed60SHong Zhang       }
4243fd52205SHong Zhang     }
4252ca6e920SHong Zhang   }
4262ca6e920SHong Zhang 
4272ca6e920SHong Zhang   ts->steps++;
4282ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
4292ca6e920SHong Zhang   PetscFunctionReturn(0);
4302ca6e920SHong Zhang }
4312ca6e920SHong Zhang 
4322ca6e920SHong Zhang #undef __FUNCT__
433cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
434cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
435cd652676SJed Brown {
436cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
4375a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
438cd652676SJed Brown   PetscErrorCode ierr;
439cd652676SJed Brown 
440cd652676SJed Brown   PetscFunctionBegin;
441a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
4425a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
4435a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
444cd652676SJed Brown   PetscFunctionReturn(0);
445cd652676SJed Brown }
446cd652676SJed Brown 
447316643e7SJed Brown /*------------------------------------------------------------*/
448316643e7SJed Brown #undef __FUNCT__
449277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
450277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
451316643e7SJed Brown {
452316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
453316643e7SJed Brown   PetscErrorCode ierr;
454316643e7SJed Brown 
455316643e7SJed Brown   PetscFunctionBegin;
4566bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
4576bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
4583b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
459eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
460*b1cb36f3SHong Zhang   if (th->VecCostIntegral0) {
461*b1cb36f3SHong Zhang     ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
462*b1cb36f3SHong Zhang   }
463b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
464b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
465b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
466277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
467277b19d0SLisandro Dalcin }
468277b19d0SLisandro Dalcin 
469277b19d0SLisandro Dalcin #undef __FUNCT__
470277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
471277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
472277b19d0SLisandro Dalcin {
473277b19d0SLisandro Dalcin   PetscErrorCode ierr;
474277b19d0SLisandro Dalcin 
475277b19d0SLisandro Dalcin   PetscFunctionBegin;
476277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
477277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
478bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
479bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
480bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
481bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
482316643e7SJed Brown   PetscFunctionReturn(0);
483316643e7SJed Brown }
484316643e7SJed Brown 
485316643e7SJed Brown /*
486316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
4872b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
488316643e7SJed Brown */
489316643e7SJed Brown #undef __FUNCT__
4900f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
4910f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
492316643e7SJed Brown {
493316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
494316643e7SJed Brown   PetscErrorCode ierr;
4957445fe48SJed Brown   Vec            X0,Xdot;
4967445fe48SJed Brown   DM             dm,dmsave;
497b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
498316643e7SJed Brown 
499316643e7SJed Brown   PetscFunctionBegin;
5007445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
5015a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
5027445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
503b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
5047445fe48SJed Brown 
5057445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
5067445fe48SJed Brown   dmsave = ts->dm;
5077445fe48SJed Brown   ts->dm = dm;
5087445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
5097445fe48SJed Brown   ts->dm = dmsave;
5100d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
511316643e7SJed Brown   PetscFunctionReturn(0);
512316643e7SJed Brown }
513316643e7SJed Brown 
514316643e7SJed Brown #undef __FUNCT__
5150f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
516d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
517316643e7SJed Brown {
518316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
519316643e7SJed Brown   PetscErrorCode ierr;
5207445fe48SJed Brown   Vec            Xdot;
5217445fe48SJed Brown   DM             dm,dmsave;
522b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
523316643e7SJed Brown 
524316643e7SJed Brown   PetscFunctionBegin;
5257445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
5267445fe48SJed Brown 
5270f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
5280298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
5297445fe48SJed Brown 
5307445fe48SJed Brown   dmsave = ts->dm;
5317445fe48SJed Brown   ts->dm = dm;
532d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
5337445fe48SJed Brown   ts->dm = dmsave;
5340298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
535316643e7SJed Brown   PetscFunctionReturn(0);
536316643e7SJed Brown }
537316643e7SJed Brown 
538316643e7SJed Brown #undef __FUNCT__
539316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
540316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
541316643e7SJed Brown {
542316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
543316643e7SJed Brown   PetscErrorCode ierr;
5447445fe48SJed Brown   SNES           snes;
545ef749922SLisandro Dalcin   TSAdapt        adapt;
5467445fe48SJed Brown   DM             dm;
547316643e7SJed Brown 
548316643e7SJed Brown   PetscFunctionBegin;
549*b1cb36f3SHong Zhang   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
550*b1cb36f3SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
551*b1cb36f3SHong Zhang   }
55239d13666SHong Zhang   if (!th->X) {
553316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
55439d13666SHong Zhang   }
55539d13666SHong Zhang   if (!th->Xdot) {
556316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
55739d13666SHong Zhang   }
55839d13666SHong Zhang   if (!th->X0) {
5593b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
56039d13666SHong Zhang   }
5617445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
5627445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
5637445fe48SJed Brown   if (dm) {
5647445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
565258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
5667445fe48SJed Brown   }
5673b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
5683b1890cdSShri Abhyankar   else th->order = 1;
5693b1890cdSShri Abhyankar 
570552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
571ef749922SLisandro Dalcin   if (!th->adapt) {
5723b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
5733b1890cdSShri Abhyankar   }
574b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
575b4dd3bf9SBarry Smith }
5760c7376e5SHong Zhang 
5770c7376e5SHong Zhang #undef __FUNCT__
5780c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_BEuler"
5790c7376e5SHong Zhang static PetscErrorCode TSSetUp_BEuler(TS ts)
5800c7376e5SHong Zhang {
5810c7376e5SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
5820c7376e5SHong Zhang   PetscErrorCode ierr;
5830c7376e5SHong Zhang 
5840c7376e5SHong Zhang   PetscFunctionBegin;
585b9eb5ee8SHong Zhang   if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n");
586421bc905SSatish Balay   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
5870c7376e5SHong Zhang   PetscFunctionReturn(0);
5880c7376e5SHong Zhang }
5890c7376e5SHong Zhang 
5900c7376e5SHong Zhang #undef __FUNCT__
5910c7376e5SHong Zhang #define __FUNCT__ "TSSetUp_CN"
5920c7376e5SHong Zhang static PetscErrorCode TSSetUp_CN(TS ts)
5930c7376e5SHong Zhang {
5940c7376e5SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
5950c7376e5SHong Zhang   PetscErrorCode ierr;
5960c7376e5SHong Zhang 
5970c7376e5SHong Zhang   PetscFunctionBegin;
598b9eb5ee8SHong Zhang   if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n");
599b9eb5ee8SHong Zhang   if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n");
600421bc905SSatish Balay   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
6010c7376e5SHong Zhang   PetscFunctionReturn(0);
6020c7376e5SHong Zhang }
603b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
604b4dd3bf9SBarry Smith 
605b4dd3bf9SBarry Smith #undef __FUNCT__
606b4dd3bf9SBarry Smith #define __FUNCT__ "TSAdjointSetUp_Theta"
607b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
608b4dd3bf9SBarry Smith {
609b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
610b4dd3bf9SBarry Smith   PetscErrorCode ierr;
611b4dd3bf9SBarry Smith 
612b4dd3bf9SBarry Smith   PetscFunctionBegin;
613b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
6142ca6e920SHong Zhang   if(ts->vecs_sensip) {
615b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
6162ca6e920SHong Zhang   }
617b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
618316643e7SJed Brown   PetscFunctionReturn(0);
619316643e7SJed Brown }
620316643e7SJed Brown /*------------------------------------------------------------*/
621316643e7SJed Brown 
622316643e7SJed Brown #undef __FUNCT__
623316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
6244416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
625316643e7SJed Brown {
626316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
627316643e7SJed Brown   PetscErrorCode ierr;
628316643e7SJed Brown 
629316643e7SJed Brown   PetscFunctionBegin;
630e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
631316643e7SJed Brown   {
6320298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
6330298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
6340298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
6350298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
636316643e7SJed Brown   }
637316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
638316643e7SJed Brown   PetscFunctionReturn(0);
639316643e7SJed Brown }
640316643e7SJed Brown 
641316643e7SJed Brown #undef __FUNCT__
642316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
643316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
644316643e7SJed Brown {
645316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
646ace3abfcSBarry Smith   PetscBool      iascii;
647316643e7SJed Brown   PetscErrorCode ierr;
648316643e7SJed Brown 
649316643e7SJed Brown   PetscFunctionBegin;
650251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
651316643e7SJed Brown   if (iascii) {
6527c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
653ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
654316643e7SJed Brown   }
655ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
656316643e7SJed Brown   PetscFunctionReturn(0);
657316643e7SJed Brown }
658316643e7SJed Brown 
6590de4c49aSJed Brown #undef __FUNCT__
6600de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
6617087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
6620de4c49aSJed Brown {
6630de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6640de4c49aSJed Brown 
6650de4c49aSJed Brown   PetscFunctionBegin;
6660de4c49aSJed Brown   *theta = th->Theta;
6670de4c49aSJed Brown   PetscFunctionReturn(0);
6680de4c49aSJed Brown }
6690de4c49aSJed Brown 
6700de4c49aSJed Brown #undef __FUNCT__
6710de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
6727087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
6730de4c49aSJed Brown {
6740de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
6750de4c49aSJed Brown 
6760de4c49aSJed Brown   PetscFunctionBegin;
6777c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
6780de4c49aSJed Brown   th->Theta = theta;
6790de4c49aSJed Brown   PetscFunctionReturn(0);
6800de4c49aSJed Brown }
681eb284becSJed Brown 
682eb284becSJed Brown #undef __FUNCT__
68378e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
68426f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
68526f2ff8fSLisandro Dalcin {
68626f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
68726f2ff8fSLisandro Dalcin 
68826f2ff8fSLisandro Dalcin   PetscFunctionBegin;
68926f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
69026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
69126f2ff8fSLisandro Dalcin }
69226f2ff8fSLisandro Dalcin 
69326f2ff8fSLisandro Dalcin #undef __FUNCT__
69426f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
695eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
696eb284becSJed Brown {
697eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
698eb284becSJed Brown 
699eb284becSJed Brown   PetscFunctionBegin;
700eb284becSJed Brown   th->endpoint = flg;
701eb284becSJed Brown   PetscFunctionReturn(0);
702eb284becSJed Brown }
7030de4c49aSJed Brown 
704f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
705f9c1d6abSBarry Smith #undef __FUNCT__
706f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
707f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
708f9c1d6abSBarry Smith {
709f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
710f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
7113fd8ae06SJed Brown   const PetscReal one = 1.0;
712f9c1d6abSBarry Smith 
713f9c1d6abSBarry Smith   PetscFunctionBegin;
7143fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
715f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
716f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
717f9c1d6abSBarry Smith   PetscFunctionReturn(0);
718f9c1d6abSBarry Smith }
719f9c1d6abSBarry Smith #endif
720f9c1d6abSBarry Smith 
72142682096SHong Zhang #undef __FUNCT__
72242682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
72342682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
72442682096SHong Zhang {
72542682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
72642682096SHong Zhang 
72742682096SHong Zhang   PetscFunctionBegin;
7282ca6e920SHong Zhang   *ns = 1;
7292ca6e920SHong Zhang   if(Y) {
730b5e0532dSHong Zhang     *Y  = (th->endpoint)?&(th->X0):&(th->X);
7312ca6e920SHong Zhang   }
73242682096SHong Zhang   PetscFunctionReturn(0);
73342682096SHong Zhang }
734f9c1d6abSBarry Smith 
735316643e7SJed Brown /* ------------------------------------------------------------ */
736316643e7SJed Brown /*MC
73796f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
738316643e7SJed Brown 
739316643e7SJed Brown    Level: beginner
740316643e7SJed Brown 
7414eb428fdSBarry Smith    Options Database:
742c82ed962SBarry Smith +      -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
743c82ed962SBarry Smith .      -ts_theta_extrapolate <flg> - Extrapolate stage solution from previous solution (sometimes unstable)
744c82ed962SBarry Smith .      -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
745c82ed962SBarry Smith -     -ts_theta_adapt <flg> - Use time-step adaptivity with the Theta method
7464eb428fdSBarry Smith 
747eb284becSJed Brown    Notes:
7480c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
7490c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
7504eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
7514eb428fdSBarry Smith 
752eb284becSJed Brown    This method can be applied to DAE.
753eb284becSJed Brown 
754eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
755eb284becSJed Brown 
756eb284becSJed Brown .vb
757eb284becSJed Brown   Theta | Theta
758eb284becSJed Brown   -------------
759eb284becSJed Brown         |  1
760eb284becSJed Brown .ve
761eb284becSJed Brown 
762eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
763eb284becSJed Brown 
764eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
765eb284becSJed Brown 
766eb284becSJed Brown .vb
767eb284becSJed Brown   0 | 0         0
768eb284becSJed Brown   1 | 1-Theta   Theta
769eb284becSJed Brown   -------------------
770eb284becSJed Brown     | 1-Theta   Theta
771eb284becSJed Brown .ve
772eb284becSJed Brown 
773eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
774eb284becSJed Brown 
775eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
776eb284becSJed Brown 
777eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
778eb284becSJed Brown 
7794eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
780eb284becSJed Brown 
781eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
782316643e7SJed Brown 
783316643e7SJed Brown M*/
784316643e7SJed Brown #undef __FUNCT__
785316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
7868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
787316643e7SJed Brown {
788316643e7SJed Brown   TS_Theta       *th;
789316643e7SJed Brown   PetscErrorCode ierr;
790316643e7SJed Brown 
791316643e7SJed Brown   PetscFunctionBegin;
792277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
793316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
794316643e7SJed Brown   ts->ops->view            = TSView_Theta;
795316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
79642f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
797316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
798cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
7993b1890cdSShri Abhyankar   ts->ops->evaluatestep    = TSEvaluateStep_Theta;
80024655328SShri   ts->ops->rollback        = TSRollBack_Theta;
801316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
8020f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
8030f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
804f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
805f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
806f9c1d6abSBarry Smith #endif
80742682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
80842f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
809*b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
810*b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
811316643e7SJed Brown 
812b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
813316643e7SJed Brown   ts->data = (void*)th;
814316643e7SJed Brown 
8156f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
816316643e7SJed Brown   th->Theta       = 0.5;
8173b1890cdSShri Abhyankar   th->ccfl        = 1.0;
8183b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
819bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
820bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
821bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
822bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
823316643e7SJed Brown   PetscFunctionReturn(0);
824316643e7SJed Brown }
8250de4c49aSJed Brown 
8260de4c49aSJed Brown #undef __FUNCT__
8270de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
8280de4c49aSJed Brown /*@
8290de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
8300de4c49aSJed Brown 
8310de4c49aSJed Brown   Not Collective
8320de4c49aSJed Brown 
8330de4c49aSJed Brown   Input Parameter:
8340de4c49aSJed Brown .  ts - timestepping context
8350de4c49aSJed Brown 
8360de4c49aSJed Brown   Output Parameter:
8370de4c49aSJed Brown .  theta - stage abscissa
8380de4c49aSJed Brown 
8390de4c49aSJed Brown   Note:
8400de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
8410de4c49aSJed Brown 
8420de4c49aSJed Brown   Level: Advanced
8430de4c49aSJed Brown 
8440de4c49aSJed Brown .seealso: TSThetaSetTheta()
8450de4c49aSJed Brown @*/
8467087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
8470de4c49aSJed Brown {
8484ac538c5SBarry Smith   PetscErrorCode ierr;
8490de4c49aSJed Brown 
8500de4c49aSJed Brown   PetscFunctionBegin;
851afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8520de4c49aSJed Brown   PetscValidPointer(theta,2);
8534ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
8540de4c49aSJed Brown   PetscFunctionReturn(0);
8550de4c49aSJed Brown }
8560de4c49aSJed Brown 
8570de4c49aSJed Brown #undef __FUNCT__
8580de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
8590de4c49aSJed Brown /*@
8600de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
8610de4c49aSJed Brown 
8620de4c49aSJed Brown   Not Collective
8630de4c49aSJed Brown 
8640de4c49aSJed Brown   Input Parameter:
8650de4c49aSJed Brown +  ts - timestepping context
8660de4c49aSJed Brown -  theta - stage abscissa
8670de4c49aSJed Brown 
8680de4c49aSJed Brown   Options Database:
8690de4c49aSJed Brown .  -ts_theta_theta <theta>
8700de4c49aSJed Brown 
8710de4c49aSJed Brown   Level: Intermediate
8720de4c49aSJed Brown 
8730de4c49aSJed Brown .seealso: TSThetaGetTheta()
8740de4c49aSJed Brown @*/
8757087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
8760de4c49aSJed Brown {
8774ac538c5SBarry Smith   PetscErrorCode ierr;
8780de4c49aSJed Brown 
8790de4c49aSJed Brown   PetscFunctionBegin;
880afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
8814ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
8820de4c49aSJed Brown   PetscFunctionReturn(0);
8830de4c49aSJed Brown }
884f33bbcb6SJed Brown 
885eb284becSJed Brown #undef __FUNCT__
88626f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
88726f2ff8fSLisandro Dalcin /*@
88826f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
88926f2ff8fSLisandro Dalcin 
89026f2ff8fSLisandro Dalcin   Not Collective
89126f2ff8fSLisandro Dalcin 
89226f2ff8fSLisandro Dalcin   Input Parameter:
89326f2ff8fSLisandro Dalcin .  ts - timestepping context
89426f2ff8fSLisandro Dalcin 
89526f2ff8fSLisandro Dalcin   Output Parameter:
89626f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
89726f2ff8fSLisandro Dalcin 
89826f2ff8fSLisandro Dalcin   Level: Advanced
89926f2ff8fSLisandro Dalcin 
90026f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
90126f2ff8fSLisandro Dalcin @*/
90226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
90326f2ff8fSLisandro Dalcin {
90426f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
90526f2ff8fSLisandro Dalcin 
90626f2ff8fSLisandro Dalcin   PetscFunctionBegin;
90726f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
90826f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
90926f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
91026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
91126f2ff8fSLisandro Dalcin }
91226f2ff8fSLisandro Dalcin 
91326f2ff8fSLisandro Dalcin #undef __FUNCT__
914eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
915eb284becSJed Brown /*@
916eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
917eb284becSJed Brown 
918eb284becSJed Brown   Not Collective
919eb284becSJed Brown 
920eb284becSJed Brown   Input Parameter:
921eb284becSJed Brown +  ts - timestepping context
922eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
923eb284becSJed Brown 
924eb284becSJed Brown   Options Database:
925eb284becSJed Brown .  -ts_theta_endpoint <flg>
926eb284becSJed Brown 
927eb284becSJed Brown   Level: Intermediate
928eb284becSJed Brown 
929eb284becSJed Brown .seealso: TSTHETA, TSCN
930eb284becSJed Brown @*/
931eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
932eb284becSJed Brown {
933eb284becSJed Brown   PetscErrorCode ierr;
934eb284becSJed Brown 
935eb284becSJed Brown   PetscFunctionBegin;
936eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
937eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
938eb284becSJed Brown   PetscFunctionReturn(0);
939eb284becSJed Brown }
940eb284becSJed Brown 
941f33bbcb6SJed Brown /*
942f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
943f33bbcb6SJed Brown  * The creation functions for these specializations are below.
944f33bbcb6SJed Brown  */
945f33bbcb6SJed Brown 
946f33bbcb6SJed Brown #undef __FUNCT__
947f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
948f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
949f33bbcb6SJed Brown {
950d52bd9f3SBarry Smith   PetscErrorCode ierr;
951d52bd9f3SBarry Smith 
952f33bbcb6SJed Brown   PetscFunctionBegin;
953d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
954f33bbcb6SJed Brown   PetscFunctionReturn(0);
955f33bbcb6SJed Brown }
956f33bbcb6SJed Brown 
957f33bbcb6SJed Brown /*MC
958f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
959f33bbcb6SJed Brown 
960f33bbcb6SJed Brown   Level: beginner
961f33bbcb6SJed Brown 
9624eb428fdSBarry Smith   Notes:
963c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
9644eb428fdSBarry Smith 
9654eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
9664eb428fdSBarry Smith 
967f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
968f33bbcb6SJed Brown 
969f33bbcb6SJed Brown M*/
970f33bbcb6SJed Brown #undef __FUNCT__
971f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
9728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
973f33bbcb6SJed Brown {
974f33bbcb6SJed Brown   PetscErrorCode ierr;
975f33bbcb6SJed Brown 
976f33bbcb6SJed Brown   PetscFunctionBegin;
977f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
978f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
9790c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
980f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
981f33bbcb6SJed Brown   PetscFunctionReturn(0);
982f33bbcb6SJed Brown }
983f33bbcb6SJed Brown 
984f33bbcb6SJed Brown #undef __FUNCT__
985f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
986f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
987f33bbcb6SJed Brown {
988d52bd9f3SBarry Smith   PetscErrorCode ierr;
989d52bd9f3SBarry Smith 
990f33bbcb6SJed Brown   PetscFunctionBegin;
991d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
992f33bbcb6SJed Brown   PetscFunctionReturn(0);
993f33bbcb6SJed Brown }
994f33bbcb6SJed Brown 
995f33bbcb6SJed Brown /*MC
996f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
997f33bbcb6SJed Brown 
998f33bbcb6SJed Brown   Level: beginner
999f33bbcb6SJed Brown 
1000f33bbcb6SJed Brown   Notes:
10017cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
10027cf5af47SJed Brown 
10037cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1004f33bbcb6SJed Brown 
1005f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1006f33bbcb6SJed Brown 
1007f33bbcb6SJed Brown M*/
1008f33bbcb6SJed Brown #undef __FUNCT__
1009f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
10108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1011f33bbcb6SJed Brown {
1012f33bbcb6SJed Brown   PetscErrorCode ierr;
1013f33bbcb6SJed Brown 
1014f33bbcb6SJed Brown   PetscFunctionBegin;
1015f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1016f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1017eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
10180c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1019f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
1020f33bbcb6SJed Brown   PetscFunctionReturn(0);
1021f33bbcb6SJed Brown }
1022