xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 6f92202b937c8db5b63afc621286d06689bf89d4)
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 {
10715f1b00SHong Zhang   /* context for time stepping */
111566a47fSLisandro Dalcin   PetscReal    stage_time;
121566a47fSLisandro Dalcin   Vec          X0,X,Xdot;                /* Storage for stages and time derivative */
131566a47fSLisandro Dalcin   Vec          affine;                   /* Affine vector needed for residual at beginning of step in endpoint formulation */
141566a47fSLisandro Dalcin   PetscReal    Theta;
15715f1b00SHong Zhang   PetscReal    ptime;
16715f1b00SHong Zhang   PetscReal    time_step;
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20715f1b00SHong Zhang   TSStepStatus status;
21715f1b00SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
221566a47fSLisandro Dalcin 
23715f1b00SHong Zhang   /* context for sensitivity analysis */
24b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
25b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
26b5e0532dSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be timed with Jacobian transpose */
27715f1b00SHong Zhang   Vec          *VecsDeltaFwdSensi;      /* Increment of the forward sensitivity at stage */
28715f1b00SHong Zhang   Vec          *VecsDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
29715f1b00SHong Zhang   Vec          VecIntegralSensiTemp;     /* Working vector for forward integral sensitivity */
30715f1b00SHong Zhang   Vec          VecIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
31*6f92202bSHong Zhang   Vec          *VecsFwdSensi0;          /* backup for roll-backs due to events */
32*6f92202bSHong Zhang   Vec          *VecsFwdSensip0;         /* backup for roll-backs due to events */
33*6f92202bSHong Zhang   Vec          *VecsIntegralSensi0;      /* backup for roll-backs due to events */
34*6f92202bSHong Zhang   Vec          *VecsIntegralSensip0;     /* backup for roll-backs due to events */
351566a47fSLisandro Dalcin 
36715f1b00SHong Zhang   /* context for error estimation */
371566a47fSLisandro Dalcin   Vec          vec_sol_prev;
381566a47fSLisandro Dalcin   Vec          vec_lte_work;
39316643e7SJed Brown } TS_Theta;
40316643e7SJed Brown 
417445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
427445fe48SJed Brown {
437445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
447445fe48SJed Brown   PetscErrorCode ierr;
457445fe48SJed Brown 
467445fe48SJed Brown   PetscFunctionBegin;
477445fe48SJed Brown   if (X0) {
487445fe48SJed Brown     if (dm && dm != ts->dm) {
490d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
507445fe48SJed Brown     } else *X0 = ts->vec_sol;
517445fe48SJed Brown   }
527445fe48SJed Brown   if (Xdot) {
537445fe48SJed Brown     if (dm && dm != ts->dm) {
540d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
557445fe48SJed Brown     } else *Xdot = th->Xdot;
567445fe48SJed Brown   }
577445fe48SJed Brown   PetscFunctionReturn(0);
587445fe48SJed Brown }
597445fe48SJed Brown 
600d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
610d0b770aSPeter Brune {
620d0b770aSPeter Brune   PetscErrorCode ierr;
630d0b770aSPeter Brune 
640d0b770aSPeter Brune   PetscFunctionBegin;
650d0b770aSPeter Brune   if (X0) {
660d0b770aSPeter Brune     if (dm && dm != ts->dm) {
670d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
680d0b770aSPeter Brune     }
690d0b770aSPeter Brune   }
700d0b770aSPeter Brune   if (Xdot) {
710d0b770aSPeter Brune     if (dm && dm != ts->dm) {
720d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
730d0b770aSPeter Brune     }
740d0b770aSPeter Brune   }
750d0b770aSPeter Brune   PetscFunctionReturn(0);
760d0b770aSPeter Brune }
770d0b770aSPeter Brune 
787445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
797445fe48SJed Brown {
807445fe48SJed Brown   PetscFunctionBegin;
817445fe48SJed Brown   PetscFunctionReturn(0);
827445fe48SJed Brown }
837445fe48SJed Brown 
847445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
857445fe48SJed Brown {
867445fe48SJed Brown   TS             ts = (TS)ctx;
877445fe48SJed Brown   PetscErrorCode ierr;
887445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
897445fe48SJed Brown 
907445fe48SJed Brown   PetscFunctionBegin;
917445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
927445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
937445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
970d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
980d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
997445fe48SJed Brown   PetscFunctionReturn(0);
1007445fe48SJed Brown }
1017445fe48SJed Brown 
102258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
103258e1594SPeter Brune {
104258e1594SPeter Brune   PetscFunctionBegin;
105258e1594SPeter Brune   PetscFunctionReturn(0);
106258e1594SPeter Brune }
107258e1594SPeter Brune 
108258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
109258e1594SPeter Brune {
110258e1594SPeter Brune   TS             ts = (TS)ctx;
111258e1594SPeter Brune   PetscErrorCode ierr;
112258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
113258e1594SPeter Brune 
114258e1594SPeter Brune   PetscFunctionBegin;
115258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
116258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
117258e1594SPeter Brune 
118258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune 
121258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune 
124258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
125258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
126258e1594SPeter Brune   PetscFunctionReturn(0);
127258e1594SPeter Brune }
128258e1594SPeter Brune 
129b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
130b1cb36f3SHong Zhang {
131b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
132b1cb36f3SHong Zhang   PetscErrorCode ierr;
133b1cb36f3SHong Zhang 
134b1cb36f3SHong Zhang   PetscFunctionBegin;
135b1cb36f3SHong Zhang   /* backup cost integral */
136b1cb36f3SHong Zhang   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
137b1cb36f3SHong Zhang   if (th->endpoint) {
138b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1391566a47fSLisandro Dalcin     ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
140b1cb36f3SHong Zhang   }
141b1cb36f3SHong Zhang   ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
142b1cb36f3SHong Zhang   if (th->endpoint) {
143b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
144b1cb36f3SHong Zhang   } else {
145b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
146b1cb36f3SHong Zhang   }
147b1cb36f3SHong Zhang   PetscFunctionReturn(0);
148b1cb36f3SHong Zhang }
149b1cb36f3SHong Zhang 
150b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
151b1cb36f3SHong Zhang {
152b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
153b1cb36f3SHong Zhang   PetscErrorCode ierr;
154b1cb36f3SHong Zhang 
155b1cb36f3SHong Zhang   PetscFunctionBegin;
156b1cb36f3SHong Zhang   if (th->endpoint) {
157b1cb36f3SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
158b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
159b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
160b1cb36f3SHong Zhang     if (th->Theta!=1) {
161b1cb36f3SHong Zhang       ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1621566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr);
163b1cb36f3SHong Zhang     }
164b1cb36f3SHong Zhang   } else {
165b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
166b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
167b1cb36f3SHong Zhang   }
16824655328SShri   PetscFunctionReturn(0);
16924655328SShri }
17024655328SShri 
1711566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
1721566a47fSLisandro Dalcin {
1731566a47fSLisandro Dalcin   PetscInt       nits,lits;
1741566a47fSLisandro Dalcin   PetscErrorCode ierr;
1751566a47fSLisandro Dalcin 
1761566a47fSLisandro Dalcin   PetscFunctionBegin;
1771566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1781566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1791566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1801566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1811566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1821566a47fSLisandro Dalcin }
1831566a47fSLisandro Dalcin 
184193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
185316643e7SJed Brown {
186316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1871566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1884957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1891566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
190051f2191SLisandro Dalcin   PetscErrorCode ierr;
191316643e7SJed Brown 
192316643e7SJed Brown   PetscFunctionBegin;
1931566a47fSLisandro Dalcin   if (!ts->steprollback) {
1942ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
1953b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
1961566a47fSLisandro Dalcin   }
1971566a47fSLisandro Dalcin 
1981566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
1991566a47fSLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
200715f1b00SHong Zhang 
2011566a47fSLisandro Dalcin     PetscReal shift = 1/(th->Theta*ts->time_step);
2021566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
203316643e7SJed Brown 
204be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
205fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
206be5899b3SLisandro Dalcin       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
207be5899b3SLisandro Dalcin     }
208eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
209eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2101566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2111566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2121566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2131566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2141566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
215eb284becSJed Brown     }
216be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
2171566a47fSLisandro Dalcin     ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
218be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2191566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
220fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
221051f2191SLisandro Dalcin 
222051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2231566a47fSLisandro Dalcin     if (th->endpoint) {
2241566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2251566a47fSLisandro Dalcin     } else {
226cb3a5882SLisandro Dalcin       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
2271566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2281566a47fSLisandro Dalcin     }
2291566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2301566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2311566a47fSLisandro Dalcin     if (!accept) {
2321566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
233be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
234be5899b3SLisandro Dalcin       goto reject_step;
235051f2191SLisandro Dalcin     }
2363b1890cdSShri Abhyankar 
237715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
238b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
239b1cb36f3SHong Zhang       th->time_step = ts->time_step;
24017145e75SHong Zhang     }
2411566a47fSLisandro Dalcin 
2422b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
243cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
244051f2191SLisandro Dalcin     break;
245051f2191SLisandro Dalcin 
246051f2191SLisandro Dalcin   reject_step:
247fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2481566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
249051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2501566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2513b1890cdSShri Abhyankar     }
2523b1890cdSShri Abhyankar   }
253316643e7SJed Brown   PetscFunctionReturn(0);
254316643e7SJed Brown }
255316643e7SJed Brown 
25642f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
2572ca6e920SHong Zhang {
2582ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
259b5e0532dSHong Zhang   Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
2602ca6e920SHong Zhang   PetscInt            nadj;
2612ca6e920SHong Zhang   PetscErrorCode      ierr;
2622ca6e920SHong Zhang   Mat                 J,Jp;
2632ca6e920SHong Zhang   KSP                 ksp;
2642ca6e920SHong Zhang   PetscReal           shift;
2652ca6e920SHong Zhang 
2662ca6e920SHong Zhang   PetscFunctionBegin;
2672ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
268302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
2692ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
270b5e0532dSHong Zhang 
271b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
2723fd52205SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
273b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
2742ca6e920SHong Zhang 
275a4cab896SHong Zhang   /* Build RHS */
2762c39e106SBarry Smith   if (ts->vec_costintegral) { /* Cost function has an integral term */
27705755b9cSHong Zhang     if (th->endpoint) {
278d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
27905755b9cSHong Zhang     } else {
280d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
28105755b9cSHong Zhang     }
28205755b9cSHong Zhang   }
283abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
284b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
285b5e0532dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
2862c39e106SBarry Smith     if (ts->vec_costintegral) {
287b5e0532dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
28836eaed60SHong Zhang     }
2892ca6e920SHong Zhang   }
2903c54f92cSHong Zhang 
2912ca6e920SHong Zhang   /* Build LHS */
2922ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
2933c54f92cSHong Zhang   if (th->endpoint) {
2943c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2953c54f92cSHong Zhang   } else {
2963c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2973c54f92cSHong Zhang   }
2982ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
2992ca6e920SHong Zhang 
3002ca6e920SHong Zhang   /* Solve LHS X = RHS */
301abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
302b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
3032ca6e920SHong Zhang   }
3043c54f92cSHong Zhang 
30536eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
306b5e0532dSHong Zhang   if(th->endpoint) { /* two-stage case */
307b5e0532dSHong Zhang     if (th->Theta!=1.) {
3083fd52205SHong Zhang       shift = -1./((th->Theta-1.)*ts->time_step);
309b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3102c39e106SBarry Smith       if (ts->vec_costintegral) {
311b5e0532dSHong Zhang         ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
3128263dbbfSHong Zhang       }
313abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
314b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
3152c39e106SBarry Smith         if (ts->vec_costintegral) {
31636eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
31736eaed60SHong Zhang         }
3183c54f92cSHong Zhang         ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3192ca6e920SHong Zhang       }
320b5e0532dSHong Zhang     } else { /* backward Euler */
321b5e0532dSHong Zhang       shift = 0.0;
322b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
323b5e0532dSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
324b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
325b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
326b5e0532dSHong Zhang         if (ts->vec_costintegral) {
327b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
328b5e0532dSHong Zhang         }
329b5e0532dSHong Zhang       }
330b5e0532dSHong Zhang     }
3313fd52205SHong Zhang 
3323fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
3335bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
334abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
335b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
336b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3373fd52205SHong Zhang       }
338b5e0532dSHong Zhang       if (th->Theta!=1.) {
339b5e0532dSHong Zhang         ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
340abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
341b5e0532dSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
342b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
343b5e0532dSHong Zhang         }
3443fd52205SHong Zhang       }
3452c39e106SBarry Smith       if (ts->vec_costintegral) {
346d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
347abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
34836eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
34936eaed60SHong Zhang         }
350b5e0532dSHong Zhang         if (th->Theta!=1.) {
351b5e0532dSHong Zhang           ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
352abc2977eSBarry Smith           for (nadj=0; nadj<ts->numcost; nadj++) {
35336eaed60SHong Zhang             ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
35436eaed60SHong Zhang           }
35536eaed60SHong Zhang         }
3563fd52205SHong Zhang       }
357b5e0532dSHong Zhang     }
3583fd52205SHong Zhang   } else { /* one-stage case */
3593c54f92cSHong Zhang     shift = 0.0;
360a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
3612c39e106SBarry Smith     if (ts->vec_costintegral) {
362d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
3633c54f92cSHong Zhang     }
364abc2977eSBarry Smith     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);
3672c39e106SBarry Smith       if (ts->vec_costintegral) {
36836eaed60SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
36936eaed60SHong Zhang       }
3702ca6e920SHong Zhang     }
3713fd52205SHong Zhang     if (ts->vecs_sensip) {
3725bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
373abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
374b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
375b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3763fd52205SHong Zhang       }
3772c39e106SBarry Smith       if (ts->vec_costintegral) {
378d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
379abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
38036eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
38136eaed60SHong Zhang         }
38236eaed60SHong Zhang       }
3833fd52205SHong Zhang     }
3842ca6e920SHong Zhang   }
3852ca6e920SHong Zhang 
3862ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
3872ca6e920SHong Zhang   PetscFunctionReturn(0);
3882ca6e920SHong Zhang }
3892ca6e920SHong Zhang 
390cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
391cd652676SJed Brown {
392cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
393be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
394cd652676SJed Brown   PetscErrorCode ierr;
395cd652676SJed Brown 
396cd652676SJed Brown   PetscFunctionBegin;
397a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
398be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
399be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
400cd652676SJed Brown   PetscFunctionReturn(0);
401cd652676SJed Brown }
402cd652676SJed Brown 
4031566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
4041566a47fSLisandro Dalcin {
4051566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
4061566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
4071566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
4087453f775SEmil Constantinescu   PetscReal      wltea,wlter;
4091566a47fSLisandro Dalcin   PetscErrorCode ierr;
4101566a47fSLisandro Dalcin 
4111566a47fSLisandro Dalcin   PetscFunctionBegin;
4122ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
4131566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
414fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
4151566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
416fecfb714SLisandro Dalcin   {
417be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
418be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
4191566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
4201566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
4211566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
4221566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
4231566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
4247453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
4251566a47fSLisandro Dalcin   }
4261566a47fSLisandro Dalcin   if (order) *order = 2;
4271566a47fSLisandro Dalcin   PetscFunctionReturn(0);
4281566a47fSLisandro Dalcin }
4291566a47fSLisandro Dalcin 
4301566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
4311566a47fSLisandro Dalcin {
4321566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
433*6f92202bSHong Zhang   PetscInt       ndir,ncost;
4341566a47fSLisandro Dalcin   PetscErrorCode ierr;
4351566a47fSLisandro Dalcin 
4361566a47fSLisandro Dalcin   PetscFunctionBegin;
4371566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
4381566a47fSLisandro Dalcin   if (ts->vec_costintegral && ts->costintegralfwd) {
4391566a47fSLisandro Dalcin     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4401566a47fSLisandro Dalcin   }
441715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
442*6f92202bSHong Zhang 
443*6f92202bSHong Zhang   if (ts->vecs_fwdsensi) {
444*6f92202bSHong Zhang     for (ndir=0;ndir<ts->num_initialvalues;ndir++) {
445*6f92202bSHong Zhang       ierr = VecCopy(th->VecsFwdSensi0[ndir],ts->vecs_fwdsensi[ndir]);CHKERRQ(ierr);
446*6f92202bSHong Zhang     }
447*6f92202bSHong Zhang   }
448*6f92202bSHong Zhang   if (ts->vecs_fwdsensip) {
449*6f92202bSHong Zhang     for (ndir=0;ndir<ts->num_parameters;ndir++) {
450*6f92202bSHong Zhang       ierr = VecCopy(th->VecsFwdSensip0[ndir],ts->vecs_fwdsensip[ndir]);CHKERRQ(ierr);
451*6f92202bSHong Zhang     }
452*6f92202bSHong Zhang   }
453*6f92202bSHong Zhang   if (ts->vecs_integral_sensi) {
454*6f92202bSHong Zhang     for (ncost=0;ncost<ts->numcost;ncost++) {
455*6f92202bSHong Zhang       ierr = VecCopy(th->VecsIntegralSensi0[ncost],ts->vecs_integral_sensi[ncost]);CHKERRQ(ierr);
456*6f92202bSHong Zhang     }
457*6f92202bSHong Zhang   }
458*6f92202bSHong Zhang   if (ts->vecs_integral_sensip) {
459*6f92202bSHong Zhang     for (ncost=0;ncost<ts->numcost;ncost++) {
460*6f92202bSHong Zhang       ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr);
461*6f92202bSHong Zhang     }
462*6f92202bSHong Zhang   }
463715f1b00SHong Zhang   PetscFunctionReturn(0);
464715f1b00SHong Zhang }
465715f1b00SHong Zhang 
466715f1b00SHong Zhang static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt nump,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative)
467715f1b00SHong Zhang {
468715f1b00SHong Zhang   PetscInt    ndir;
469715f1b00SHong Zhang   PetscScalar *v;
470715f1b00SHong Zhang   PetscErrorCode ierr;
471715f1b00SHong Zhang 
472715f1b00SHong Zhang   PetscFunctionBegin;
473715f1b00SHong Zhang 
474715f1b00SHong Zhang   ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
475715f1b00SHong Zhang   for (ndir=0; ndir<nump; ndir++) {
476715f1b00SHong Zhang     ierr = VecDot(DRncostDY,s[ndir],&v[ndir]);CHKERRQ(ierr);
477715f1b00SHong Zhang   }
478715f1b00SHong Zhang   ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
479715f1b00SHong Zhang   if (DRncostDP) {
480715f1b00SHong Zhang     ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr);
481715f1b00SHong Zhang   }
482715f1b00SHong Zhang   PetscFunctionReturn(0);
483715f1b00SHong Zhang }
484715f1b00SHong Zhang 
485715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
486715f1b00SHong Zhang {
487715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
488715f1b00SHong Zhang   Vec            *VecsDeltaFwdSensi  = th->VecsDeltaFwdSensi;
489715f1b00SHong Zhang   Vec            *VecsDeltaFwdSensip = th->VecsDeltaFwdSensip;
490715f1b00SHong Zhang   PetscInt       ndir,ncost;
491715f1b00SHong Zhang   KSP            ksp;
492715f1b00SHong Zhang   Mat            J,Jp;
493715f1b00SHong Zhang   PetscReal      shift;
494715f1b00SHong Zhang   PetscErrorCode ierr;
495715f1b00SHong Zhang 
496715f1b00SHong Zhang   PetscFunctionBegin;
497*6f92202bSHong Zhang   for (ndir=0;ndir<ts->num_initialvalues;ndir++) {
498*6f92202bSHong Zhang     ierr = VecCopy(ts->vecs_fwdsensi[ndir],th->VecsFwdSensi0[ndir]);CHKERRQ(ierr);
499*6f92202bSHong Zhang   }
500*6f92202bSHong Zhang   for (ndir=0;ndir<ts->num_parameters;ndir++) {
501*6f92202bSHong Zhang     ierr = VecCopy(ts->vecs_fwdsensip[ndir],th->VecsFwdSensip0[ndir]);CHKERRQ(ierr);
502*6f92202bSHong Zhang   }
503*6f92202bSHong Zhang   if (ts->vecs_integral_sensi) {
504*6f92202bSHong Zhang     for (ncost=0; ncost<ts->numcost;ncost++) {
505*6f92202bSHong Zhang       ierr = VecCopy(ts->vecs_integral_sensi[ncost],th->VecsIntegralSensi0[ncost]);CHKERRQ(ierr);
506*6f92202bSHong Zhang     }
507*6f92202bSHong Zhang   }
508*6f92202bSHong Zhang   if (ts->vecs_integral_sensip) {
509*6f92202bSHong Zhang     for (ncost=0;ncost<ts->numcost;ncost++) {
510*6f92202bSHong Zhang       ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr);
511*6f92202bSHong Zhang     }
512*6f92202bSHong Zhang   }
513*6f92202bSHong Zhang 
514715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
515715f1b00SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
516715f1b00SHong Zhang 
517715f1b00SHong Zhang   /* Build RHS */
518715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
519715f1b00SHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
520715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
521715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
522715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
523715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensip[ndir],(th->Theta-1.)/th->Theta);
524715f1b00SHong Zhang     }
525715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
526715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
527715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensi[ndir],(th->Theta-1.)/th->Theta);
528715f1b00SHong Zhang     }
529715f1b00SHong Zhang 
530715f1b00SHong Zhang     if (ts->num_parameters) { /* Add the f_p forcing terms */
531715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr);
532715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
533715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],(1.-th->Theta)/th->Theta,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
534715f1b00SHong Zhang       }
535715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr);
536715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
537715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
538715f1b00SHong Zhang       }
539715f1b00SHong Zhang     }
540715f1b00SHong Zhang   } else { /* 1-stage method */
541715f1b00SHong Zhang     shift = 0.0;
542715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
543715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
544715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
545715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensip[ndir],-1);
546715f1b00SHong Zhang     }
547715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
548715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
549715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensi[ndir],-1);
550715f1b00SHong Zhang     }
551715f1b00SHong Zhang     if (ts->num_parameters) { /* Add the f_p forcing terms */
552715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr);
553715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
554715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
555715f1b00SHong Zhang       }
556715f1b00SHong Zhang     }
557715f1b00SHong Zhang   }
558715f1b00SHong Zhang 
559715f1b00SHong Zhang   /* Build LHS */
560715f1b00SHong Zhang   shift = 1/(th->Theta*th->time_step);
561715f1b00SHong Zhang   if (th->endpoint) {
562715f1b00SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
563715f1b00SHong Zhang   } else {
564715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
565715f1b00SHong Zhang   }
566715f1b00SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
567715f1b00SHong Zhang 
568715f1b00SHong Zhang   /*
569715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
570715f1b00SHong Zhang     drdy|t_n*S(t_n) + drdp|t_n
571715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
572715f1b00SHong Zhang     It is important to note that sensitivitesi to parameters (sensip) and sensitivities to initial values (sensi) are independent of each other, but integral sensip relies on sensip and integral sensi requires sensi
573715f1b00SHong Zhang    */
574715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
575715f1b00SHong Zhang     if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
576715f1b00SHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
577715f1b00SHong Zhang     }
578715f1b00SHong Zhang     if (ts->vecs_integral_sensip) {
579715f1b00SHong Zhang       ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
580715f1b00SHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
581715f1b00SHong Zhang         ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
582715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
583715f1b00SHong Zhang       }
584715f1b00SHong Zhang     }
585715f1b00SHong Zhang     if (ts->vecs_integral_sensi) {
586715f1b00SHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
587715f1b00SHong Zhang         ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
588715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr);
589715f1b00SHong Zhang       }
590715f1b00SHong Zhang     }
591715f1b00SHong Zhang   }
592715f1b00SHong Zhang 
593715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
594715f1b00SHong Zhang   for (ndir=0; ndir<ts->num_parameters; ndir++) {
595715f1b00SHong Zhang     if (th->endpoint) {
596715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],ts->vecs_fwdsensip[ndir]);
597715f1b00SHong Zhang     } else {
598715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],VecsDeltaFwdSensip[ndir]);
599715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],1.,VecsDeltaFwdSensip[ndir]);
600715f1b00SHong Zhang     }
601715f1b00SHong Zhang   }
602715f1b00SHong Zhang   /* Solve the linear system for trajectory sensitivities to initial values */
603715f1b00SHong Zhang   for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
604715f1b00SHong Zhang     if (th->endpoint) {
605715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],ts->vecs_fwdsensi[ndir]);
606715f1b00SHong Zhang     } else {
607715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],VecsDeltaFwdSensi[ndir]);
608715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],1.,VecsDeltaFwdSensi[ndir]);
609715f1b00SHong Zhang     }
610715f1b00SHong Zhang   }
611715f1b00SHong Zhang 
612715f1b00SHong Zhang   /*Evaluate the second stage of integral gradients with the 2-stage method:
613715f1b00SHong Zhang     drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
614715f1b00SHong Zhang   */
615715f1b00SHong Zhang   if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
616715f1b00SHong Zhang     ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
617715f1b00SHong Zhang   }
618715f1b00SHong Zhang   if (ts->vecs_integral_sensip) {
619715f1b00SHong Zhang     ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
620715f1b00SHong Zhang     for (ncost=0; ncost<ts->numcost; ncost++) {
621715f1b00SHong Zhang       ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
622715f1b00SHong Zhang       if (th->endpoint) {
623715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
624715f1b00SHong Zhang       } else {
625715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
626715f1b00SHong Zhang       }
627715f1b00SHong Zhang     }
628715f1b00SHong Zhang   }
629715f1b00SHong Zhang   if (ts->vecs_integral_sensi) {
630715f1b00SHong Zhang     for (ncost=0; ncost<ts->numcost; ncost++) {
631715f1b00SHong Zhang       ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
632715f1b00SHong Zhang       if (th->endpoint) {
633715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr);
634715f1b00SHong Zhang       } else {
635715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr);
636715f1b00SHong Zhang       }
637715f1b00SHong Zhang     }
638715f1b00SHong Zhang   }
639715f1b00SHong Zhang 
640715f1b00SHong Zhang   if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */
641715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
642715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensip[ndir]);
643715f1b00SHong Zhang     }
644715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
645715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ndir]);
646715f1b00SHong Zhang     }
647715f1b00SHong Zhang   }
648715f1b00SHong Zhang 
6491566a47fSLisandro Dalcin   PetscFunctionReturn(0);
6501566a47fSLisandro Dalcin }
6511566a47fSLisandro Dalcin 
652316643e7SJed Brown /*------------------------------------------------------------*/
653277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
654316643e7SJed Brown {
655316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
656316643e7SJed Brown   PetscErrorCode ierr;
657316643e7SJed Brown 
658316643e7SJed Brown   PetscFunctionBegin;
6596bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
6606bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
6613b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
662eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
6631566a47fSLisandro Dalcin 
6641566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
6651566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
6661566a47fSLisandro Dalcin 
667b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
668715f1b00SHong Zhang   if (ts->forward_solve) {
669715f1b00SHong Zhang     if (ts->vecs_fwdsensi) {
670*6f92202bSHong Zhang       ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
671*6f92202bSHong Zhang       ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsFwdSensi0);CHKERRQ(ierr);
672*6f92202bSHong Zhang     }
673715f1b00SHong Zhang     if (ts->vecs_fwdsensip) {
674715f1b00SHong Zhang       ierr = VecDestroyVecs(ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
675*6f92202bSHong Zhang       ierr = VecDestroyVecs(ts->num_parameters,&th->VecsFwdSensip0);CHKERRQ(ierr);
676715f1b00SHong Zhang     }
677715f1b00SHong Zhang     if (ts->vecs_integral_sensi) {
678715f1b00SHong Zhang       ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr);
679*6f92202bSHong Zhang       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
680715f1b00SHong Zhang     }
681715f1b00SHong Zhang     if (ts->vecs_integral_sensip) {
682715f1b00SHong Zhang       ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
683*6f92202bSHong Zhang       ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
684715f1b00SHong Zhang     }
685715f1b00SHong Zhang   }
686b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
687b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
688b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
689715f1b00SHong Zhang 
690277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
691277b19d0SLisandro Dalcin }
692277b19d0SLisandro Dalcin 
693277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
694277b19d0SLisandro Dalcin {
695277b19d0SLisandro Dalcin   PetscErrorCode ierr;
696277b19d0SLisandro Dalcin 
697277b19d0SLisandro Dalcin   PetscFunctionBegin;
698277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
699277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
700bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
701bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
702bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
704316643e7SJed Brown   PetscFunctionReturn(0);
705316643e7SJed Brown }
706316643e7SJed Brown 
707316643e7SJed Brown /*
708316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
7092b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
710316643e7SJed Brown */
7110f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
712316643e7SJed Brown {
713316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
714316643e7SJed Brown   PetscErrorCode ierr;
7157445fe48SJed Brown   Vec            X0,Xdot;
7167445fe48SJed Brown   DM             dm,dmsave;
7171566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
718316643e7SJed Brown 
719316643e7SJed Brown   PetscFunctionBegin;
7207445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
7215a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
7227445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
723b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
7247445fe48SJed Brown 
7257445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
7267445fe48SJed Brown   dmsave = ts->dm;
7277445fe48SJed Brown   ts->dm = dm;
7287445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
7297445fe48SJed Brown   ts->dm = dmsave;
7300d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
731316643e7SJed Brown   PetscFunctionReturn(0);
732316643e7SJed Brown }
733316643e7SJed Brown 
734d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
735316643e7SJed Brown {
736316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
737316643e7SJed Brown   PetscErrorCode ierr;
7387445fe48SJed Brown   Vec            Xdot;
7397445fe48SJed Brown   DM             dm,dmsave;
7401566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
741316643e7SJed Brown 
742316643e7SJed Brown   PetscFunctionBegin;
7437445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
744be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
7450298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
7467445fe48SJed Brown 
7477445fe48SJed Brown   dmsave = ts->dm;
7487445fe48SJed Brown   ts->dm = dm;
749d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
7507445fe48SJed Brown   ts->dm = dmsave;
7510298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
752316643e7SJed Brown   PetscFunctionReturn(0);
753316643e7SJed Brown }
754316643e7SJed Brown 
755715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
756715f1b00SHong Zhang {
757715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
758715f1b00SHong Zhang   PetscErrorCode ierr;
759715f1b00SHong Zhang 
760715f1b00SHong Zhang   PetscFunctionBegin;
761715f1b00SHong Zhang   if (ts->vecs_fwdsensi) {
762715f1b00SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
763715f1b00SHong Zhang   }
764715f1b00SHong Zhang   if (ts->vecs_fwdsensip) {
765715f1b00SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
766715f1b00SHong Zhang   }
767715f1b00SHong Zhang   if (ts->vecs_integral_sensi) {
768715f1b00SHong Zhang     ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr);
769715f1b00SHong Zhang   }
770715f1b00SHong Zhang   if (ts->vecs_integral_sensip) {
771715f1b00SHong Zhang     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
772715f1b00SHong Zhang   }
773*6f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
774*6f92202bSHong Zhang   if (ts->vecs_fwdsensi) {
775*6f92202bSHong Zhang       ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsFwdSensi0);CHKERRQ(ierr);
776*6f92202bSHong Zhang   }
777*6f92202bSHong Zhang   if (ts->vecs_fwdsensip) {
778*6f92202bSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsFwdSensip0);CHKERRQ(ierr);
779*6f92202bSHong Zhang   }
780*6f92202bSHong Zhang   if (ts->vecs_integral_sensi) {
781*6f92202bSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_integral_sensi[0],ts->numcost,&th->VecsIntegralSensi0);CHKERRQ(ierr);
782*6f92202bSHong Zhang   }
783*6f92202bSHong Zhang   if (ts->vecs_integral_sensip) {
784*6f92202bSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
785*6f92202bSHong Zhang   }
786715f1b00SHong Zhang   PetscFunctionReturn(0);
787715f1b00SHong Zhang }
788715f1b00SHong Zhang 
789316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
790316643e7SJed Brown {
791316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
7922ffb9264SLisandro Dalcin   PetscBool      match;
793316643e7SJed Brown   PetscErrorCode ierr;
794316643e7SJed Brown 
795316643e7SJed Brown   PetscFunctionBegin;
796b1cb36f3SHong Zhang   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
797b1cb36f3SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
798b1cb36f3SHong Zhang   }
79939d13666SHong Zhang   if (!th->X) {
800316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
80139d13666SHong Zhang   }
80239d13666SHong Zhang   if (!th->Xdot) {
803316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
80439d13666SHong Zhang   }
80539d13666SHong Zhang   if (!th->X0) {
8063b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
80739d13666SHong Zhang   }
8081566a47fSLisandro Dalcin   if (th->endpoint) {
8091566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
8107445fe48SJed Brown   }
8113b1890cdSShri Abhyankar 
8121566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
8131566a47fSLisandro Dalcin 
8141566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
8151566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
8161566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
8171566a47fSLisandro Dalcin 
8181566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
8191566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
8202ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
8212ffb9264SLisandro Dalcin   if (!match) {
8221566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
8231566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
8243b1890cdSShri Abhyankar   }
8251566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
826b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
827b4dd3bf9SBarry Smith }
8280c7376e5SHong Zhang 
829b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
830b4dd3bf9SBarry Smith 
831b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
832b4dd3bf9SBarry Smith {
833b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
834b4dd3bf9SBarry Smith   PetscErrorCode ierr;
835b4dd3bf9SBarry Smith 
836b4dd3bf9SBarry Smith   PetscFunctionBegin;
837b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
8382ca6e920SHong Zhang   if(ts->vecs_sensip) {
839b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
8402ca6e920SHong Zhang   }
841b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
842316643e7SJed Brown   PetscFunctionReturn(0);
843316643e7SJed Brown }
844316643e7SJed Brown 
8454416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
846316643e7SJed Brown {
847316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
848316643e7SJed Brown   PetscErrorCode ierr;
849316643e7SJed Brown 
850316643e7SJed Brown   PetscFunctionBegin;
851e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
852316643e7SJed Brown   {
8530298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
8540298fd71SBarry 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);
85503842d09SLisandro Dalcin     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
856316643e7SJed Brown   }
857316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
858316643e7SJed Brown   PetscFunctionReturn(0);
859316643e7SJed Brown }
860316643e7SJed Brown 
861316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
862316643e7SJed Brown {
863316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
864ace3abfcSBarry Smith   PetscBool      iascii;
865316643e7SJed Brown   PetscErrorCode ierr;
866316643e7SJed Brown 
867316643e7SJed Brown   PetscFunctionBegin;
868251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
869316643e7SJed Brown   if (iascii) {
8707c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
871ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
872316643e7SJed Brown   }
873316643e7SJed Brown   PetscFunctionReturn(0);
874316643e7SJed Brown }
875316643e7SJed Brown 
876560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
8770de4c49aSJed Brown {
8780de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
8790de4c49aSJed Brown 
8800de4c49aSJed Brown   PetscFunctionBegin;
8810de4c49aSJed Brown   *theta = th->Theta;
8820de4c49aSJed Brown   PetscFunctionReturn(0);
8830de4c49aSJed Brown }
8840de4c49aSJed Brown 
885560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
8860de4c49aSJed Brown {
8870de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
8880de4c49aSJed Brown 
8890de4c49aSJed Brown   PetscFunctionBegin;
8907c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
8910de4c49aSJed Brown   th->Theta = theta;
8921566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
8930de4c49aSJed Brown   PetscFunctionReturn(0);
8940de4c49aSJed Brown }
895eb284becSJed Brown 
896560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
89726f2ff8fSLisandro Dalcin {
89826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
89926f2ff8fSLisandro Dalcin 
90026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
90126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
90226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
90326f2ff8fSLisandro Dalcin }
90426f2ff8fSLisandro Dalcin 
905560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
906eb284becSJed Brown {
907eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
908eb284becSJed Brown 
909eb284becSJed Brown   PetscFunctionBegin;
910eb284becSJed Brown   th->endpoint = flg;
911eb284becSJed Brown   PetscFunctionReturn(0);
912eb284becSJed Brown }
9130de4c49aSJed Brown 
914f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
915f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
916f9c1d6abSBarry Smith {
917f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
918f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
9193fd8ae06SJed Brown   const PetscReal one = 1.0;
920f9c1d6abSBarry Smith 
921f9c1d6abSBarry Smith   PetscFunctionBegin;
9223fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
923f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
924f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
925f9c1d6abSBarry Smith   PetscFunctionReturn(0);
926f9c1d6abSBarry Smith }
927f9c1d6abSBarry Smith #endif
928f9c1d6abSBarry Smith 
92942682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
93042682096SHong Zhang {
93142682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
93242682096SHong Zhang 
93342682096SHong Zhang   PetscFunctionBegin;
9341566a47fSLisandro Dalcin   if (ns) *ns = 1;
9351566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
93642682096SHong Zhang   PetscFunctionReturn(0);
93742682096SHong Zhang }
938f9c1d6abSBarry Smith 
939316643e7SJed Brown /* ------------------------------------------------------------ */
940316643e7SJed Brown /*MC
94196f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
942316643e7SJed Brown 
943316643e7SJed Brown    Level: beginner
944316643e7SJed Brown 
9454eb428fdSBarry Smith    Options Database:
946c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
947c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
94803842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
9494eb428fdSBarry Smith 
950eb284becSJed Brown    Notes:
9510c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
9520c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
9534eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
9544eb428fdSBarry Smith 
955eb284becSJed Brown    This method can be applied to DAE.
956eb284becSJed Brown 
957eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
958eb284becSJed Brown 
959eb284becSJed Brown .vb
960eb284becSJed Brown   Theta | Theta
961eb284becSJed Brown   -------------
962eb284becSJed Brown         |  1
963eb284becSJed Brown .ve
964eb284becSJed Brown 
965eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
966eb284becSJed Brown 
967eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
968eb284becSJed Brown 
969eb284becSJed Brown .vb
970eb284becSJed Brown   0 | 0         0
971eb284becSJed Brown   1 | 1-Theta   Theta
972eb284becSJed Brown   -------------------
973eb284becSJed Brown     | 1-Theta   Theta
974eb284becSJed Brown .ve
975eb284becSJed Brown 
976eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
977eb284becSJed Brown 
978eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
979eb284becSJed Brown 
980eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
981eb284becSJed Brown 
9824eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
983eb284becSJed Brown 
984eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
985316643e7SJed Brown 
986316643e7SJed Brown M*/
9878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
988316643e7SJed Brown {
989316643e7SJed Brown   TS_Theta       *th;
990316643e7SJed Brown   PetscErrorCode ierr;
991316643e7SJed Brown 
992316643e7SJed Brown   PetscFunctionBegin;
993277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
994316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
995316643e7SJed Brown   ts->ops->view            = TSView_Theta;
996316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
99742f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
998316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
999cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
10001566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
100124655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1002316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
10030f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
10040f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1005f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1006f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1007f9c1d6abSBarry Smith #endif
100842682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
100942f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1010b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1011b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
10122ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
1013715f1b00SHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_Theta;
1014715f1b00SHong Zhang   ts->ops->forwardstep     = TSForwardStep_Theta;
1015316643e7SJed Brown 
1016efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1017efd4aadfSBarry Smith 
1018b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1019316643e7SJed Brown   ts->data = (void*)th;
1020316643e7SJed Brown 
1021715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
1022715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
1023715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
1024715f1b00SHong Zhang 
10256f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1026316643e7SJed Brown   th->Theta       = 0.5;
10271566a47fSLisandro Dalcin   th->order       = 2;
1028bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1029bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1030bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1031bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1032316643e7SJed Brown   PetscFunctionReturn(0);
1033316643e7SJed Brown }
10340de4c49aSJed Brown 
10350de4c49aSJed Brown /*@
10360de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
10370de4c49aSJed Brown 
10380de4c49aSJed Brown   Not Collective
10390de4c49aSJed Brown 
10400de4c49aSJed Brown   Input Parameter:
10410de4c49aSJed Brown .  ts - timestepping context
10420de4c49aSJed Brown 
10430de4c49aSJed Brown   Output Parameter:
10440de4c49aSJed Brown .  theta - stage abscissa
10450de4c49aSJed Brown 
10460de4c49aSJed Brown   Note:
10470de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
10480de4c49aSJed Brown 
10490de4c49aSJed Brown   Level: Advanced
10500de4c49aSJed Brown 
10510de4c49aSJed Brown .seealso: TSThetaSetTheta()
10520de4c49aSJed Brown @*/
10537087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
10540de4c49aSJed Brown {
10554ac538c5SBarry Smith   PetscErrorCode ierr;
10560de4c49aSJed Brown 
10570de4c49aSJed Brown   PetscFunctionBegin;
1058afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10590de4c49aSJed Brown   PetscValidPointer(theta,2);
10604ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
10610de4c49aSJed Brown   PetscFunctionReturn(0);
10620de4c49aSJed Brown }
10630de4c49aSJed Brown 
10640de4c49aSJed Brown /*@
10650de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
10660de4c49aSJed Brown 
10670de4c49aSJed Brown   Not Collective
10680de4c49aSJed Brown 
10690de4c49aSJed Brown   Input Parameter:
10700de4c49aSJed Brown +  ts - timestepping context
10710de4c49aSJed Brown -  theta - stage abscissa
10720de4c49aSJed Brown 
10730de4c49aSJed Brown   Options Database:
10740de4c49aSJed Brown .  -ts_theta_theta <theta>
10750de4c49aSJed Brown 
10760de4c49aSJed Brown   Level: Intermediate
10770de4c49aSJed Brown 
10780de4c49aSJed Brown .seealso: TSThetaGetTheta()
10790de4c49aSJed Brown @*/
10807087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
10810de4c49aSJed Brown {
10824ac538c5SBarry Smith   PetscErrorCode ierr;
10830de4c49aSJed Brown 
10840de4c49aSJed Brown   PetscFunctionBegin;
1085afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10864ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
10870de4c49aSJed Brown   PetscFunctionReturn(0);
10880de4c49aSJed Brown }
1089f33bbcb6SJed Brown 
109026f2ff8fSLisandro Dalcin /*@
109126f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
109226f2ff8fSLisandro Dalcin 
109326f2ff8fSLisandro Dalcin   Not Collective
109426f2ff8fSLisandro Dalcin 
109526f2ff8fSLisandro Dalcin   Input Parameter:
109626f2ff8fSLisandro Dalcin .  ts - timestepping context
109726f2ff8fSLisandro Dalcin 
109826f2ff8fSLisandro Dalcin   Output Parameter:
109926f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
110026f2ff8fSLisandro Dalcin 
110126f2ff8fSLisandro Dalcin   Level: Advanced
110226f2ff8fSLisandro Dalcin 
110326f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
110426f2ff8fSLisandro Dalcin @*/
110526f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
110626f2ff8fSLisandro Dalcin {
110726f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
110826f2ff8fSLisandro Dalcin 
110926f2ff8fSLisandro Dalcin   PetscFunctionBegin;
111026f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
111126f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1112163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
111326f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
111426f2ff8fSLisandro Dalcin }
111526f2ff8fSLisandro Dalcin 
1116eb284becSJed Brown /*@
1117eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1118eb284becSJed Brown 
1119eb284becSJed Brown   Not Collective
1120eb284becSJed Brown 
1121eb284becSJed Brown   Input Parameter:
1122eb284becSJed Brown +  ts - timestepping context
1123eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1124eb284becSJed Brown 
1125eb284becSJed Brown   Options Database:
1126eb284becSJed Brown .  -ts_theta_endpoint <flg>
1127eb284becSJed Brown 
1128eb284becSJed Brown   Level: Intermediate
1129eb284becSJed Brown 
1130eb284becSJed Brown .seealso: TSTHETA, TSCN
1131eb284becSJed Brown @*/
1132eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1133eb284becSJed Brown {
1134eb284becSJed Brown   PetscErrorCode ierr;
1135eb284becSJed Brown 
1136eb284becSJed Brown   PetscFunctionBegin;
1137eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1138eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1139eb284becSJed Brown   PetscFunctionReturn(0);
1140eb284becSJed Brown }
1141eb284becSJed Brown 
1142f33bbcb6SJed Brown /*
1143f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1144f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1145f33bbcb6SJed Brown  */
1146f33bbcb6SJed Brown 
11471566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
11481566a47fSLisandro Dalcin {
11491566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
11501566a47fSLisandro Dalcin   PetscErrorCode ierr;
11511566a47fSLisandro Dalcin 
11521566a47fSLisandro Dalcin   PetscFunctionBegin;
11531566a47fSLisandro Dalcin   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");
11541566a47fSLisandro Dalcin   if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n");
11551566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
11561566a47fSLisandro Dalcin   PetscFunctionReturn(0);
11571566a47fSLisandro Dalcin }
11581566a47fSLisandro Dalcin 
1159f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1160f33bbcb6SJed Brown {
1161f33bbcb6SJed Brown   PetscFunctionBegin;
1162f33bbcb6SJed Brown   PetscFunctionReturn(0);
1163f33bbcb6SJed Brown }
1164f33bbcb6SJed Brown 
1165f33bbcb6SJed Brown /*MC
1166f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1167f33bbcb6SJed Brown 
1168f33bbcb6SJed Brown   Level: beginner
1169f33bbcb6SJed Brown 
11704eb428fdSBarry Smith   Notes:
1171c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
11724eb428fdSBarry Smith 
11731566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
11744eb428fdSBarry Smith 
1175f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1176f33bbcb6SJed Brown 
1177f33bbcb6SJed Brown M*/
11788cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1179f33bbcb6SJed Brown {
1180f33bbcb6SJed Brown   PetscErrorCode ierr;
1181f33bbcb6SJed Brown 
1182f33bbcb6SJed Brown   PetscFunctionBegin;
1183f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1184f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
11851566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
11860c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1187f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1188f33bbcb6SJed Brown   PetscFunctionReturn(0);
1189f33bbcb6SJed Brown }
1190f33bbcb6SJed Brown 
11911566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
11921566a47fSLisandro Dalcin {
11931566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
11941566a47fSLisandro Dalcin   PetscErrorCode ierr;
11951566a47fSLisandro Dalcin 
11961566a47fSLisandro Dalcin   PetscFunctionBegin;
11971566a47fSLisandro Dalcin   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");
11981566a47fSLisandro Dalcin   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");
11991566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
12001566a47fSLisandro Dalcin   PetscFunctionReturn(0);
12011566a47fSLisandro Dalcin }
12021566a47fSLisandro Dalcin 
1203f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1204f33bbcb6SJed Brown {
1205f33bbcb6SJed Brown   PetscFunctionBegin;
1206f33bbcb6SJed Brown   PetscFunctionReturn(0);
1207f33bbcb6SJed Brown }
1208f33bbcb6SJed Brown 
1209f33bbcb6SJed Brown /*MC
1210f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1211f33bbcb6SJed Brown 
1212f33bbcb6SJed Brown   Level: beginner
1213f33bbcb6SJed Brown 
1214f33bbcb6SJed Brown   Notes:
12157cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
12167cf5af47SJed Brown 
12177cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1218f33bbcb6SJed Brown 
1219f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1220f33bbcb6SJed Brown 
1221f33bbcb6SJed Brown M*/
12228cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1223f33bbcb6SJed Brown {
1224f33bbcb6SJed Brown   PetscErrorCode ierr;
1225f33bbcb6SJed Brown 
1226f33bbcb6SJed Brown   PetscFunctionBegin;
1227f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1228f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1229eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
12300c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1231f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1232f33bbcb6SJed Brown   PetscFunctionReturn(0);
1233f33bbcb6SJed Brown }
1234