xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 715f1b00af42c8aa385eef4857042688c0e6692c)
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 {
10*715f1b00SHong 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;
15*715f1b00SHong Zhang   PetscReal    ptime;
16*715f1b00SHong Zhang   PetscReal    time_step;
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20*715f1b00SHong Zhang   TSStepStatus status;
21*715f1b00SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
221566a47fSLisandro Dalcin 
23*715f1b00SHong 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 */
27*715f1b00SHong Zhang   Vec          *VecsDeltaFwdSensi;      /* Increment of the forward sensitivity at stage */
28*715f1b00SHong Zhang   Vec          *VecsDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
29*715f1b00SHong Zhang   Vec          VecIntegralSensiTemp;     /* Working vector for forward integral sensitivity */
30*715f1b00SHong Zhang   Vec          VecIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
311566a47fSLisandro Dalcin 
32*715f1b00SHong Zhang   /* context for error estimation */
331566a47fSLisandro Dalcin   Vec          vec_sol_prev;
341566a47fSLisandro Dalcin   Vec          vec_lte_work;
35316643e7SJed Brown } TS_Theta;
36316643e7SJed Brown 
377445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
387445fe48SJed Brown {
397445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
407445fe48SJed Brown   PetscErrorCode ierr;
417445fe48SJed Brown 
427445fe48SJed Brown   PetscFunctionBegin;
437445fe48SJed Brown   if (X0) {
447445fe48SJed Brown     if (dm && dm != ts->dm) {
450d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
467445fe48SJed Brown     } else *X0 = ts->vec_sol;
477445fe48SJed Brown   }
487445fe48SJed Brown   if (Xdot) {
497445fe48SJed Brown     if (dm && dm != ts->dm) {
500d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
517445fe48SJed Brown     } else *Xdot = th->Xdot;
527445fe48SJed Brown   }
537445fe48SJed Brown   PetscFunctionReturn(0);
547445fe48SJed Brown }
557445fe48SJed Brown 
560d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
570d0b770aSPeter Brune {
580d0b770aSPeter Brune   PetscErrorCode ierr;
590d0b770aSPeter Brune 
600d0b770aSPeter Brune   PetscFunctionBegin;
610d0b770aSPeter Brune   if (X0) {
620d0b770aSPeter Brune     if (dm && dm != ts->dm) {
630d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
640d0b770aSPeter Brune     }
650d0b770aSPeter Brune   }
660d0b770aSPeter Brune   if (Xdot) {
670d0b770aSPeter Brune     if (dm && dm != ts->dm) {
680d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
690d0b770aSPeter Brune     }
700d0b770aSPeter Brune   }
710d0b770aSPeter Brune   PetscFunctionReturn(0);
720d0b770aSPeter Brune }
730d0b770aSPeter Brune 
747445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
757445fe48SJed Brown {
767445fe48SJed Brown   PetscFunctionBegin;
777445fe48SJed Brown   PetscFunctionReturn(0);
787445fe48SJed Brown }
797445fe48SJed Brown 
807445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
817445fe48SJed Brown {
827445fe48SJed Brown   TS             ts = (TS)ctx;
837445fe48SJed Brown   PetscErrorCode ierr;
847445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
857445fe48SJed Brown 
867445fe48SJed Brown   PetscFunctionBegin;
877445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
897445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
907445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
917445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
927445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
930d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
940d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   PetscFunctionReturn(0);
967445fe48SJed Brown }
977445fe48SJed Brown 
98258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
99258e1594SPeter Brune {
100258e1594SPeter Brune   PetscFunctionBegin;
101258e1594SPeter Brune   PetscFunctionReturn(0);
102258e1594SPeter Brune }
103258e1594SPeter Brune 
104258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
105258e1594SPeter Brune {
106258e1594SPeter Brune   TS             ts = (TS)ctx;
107258e1594SPeter Brune   PetscErrorCode ierr;
108258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
109258e1594SPeter Brune 
110258e1594SPeter Brune   PetscFunctionBegin;
111258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
112258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
113258e1594SPeter Brune 
114258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
115258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116258e1594SPeter Brune 
117258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119258e1594SPeter Brune 
120258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
121258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
122258e1594SPeter Brune   PetscFunctionReturn(0);
123258e1594SPeter Brune }
124258e1594SPeter Brune 
125b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
126b1cb36f3SHong Zhang {
127b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
128b1cb36f3SHong Zhang   PetscErrorCode ierr;
129b1cb36f3SHong Zhang 
130b1cb36f3SHong Zhang   PetscFunctionBegin;
131b1cb36f3SHong Zhang   /* backup cost integral */
132b1cb36f3SHong Zhang   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
133b1cb36f3SHong Zhang   if (th->endpoint) {
134b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1351566a47fSLisandro Dalcin     ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
136b1cb36f3SHong Zhang   }
137b1cb36f3SHong Zhang   ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
138b1cb36f3SHong Zhang   if (th->endpoint) {
139b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
140b1cb36f3SHong Zhang   } else {
141b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
142b1cb36f3SHong Zhang   }
143b1cb36f3SHong Zhang   PetscFunctionReturn(0);
144b1cb36f3SHong Zhang }
145b1cb36f3SHong Zhang 
146b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
147b1cb36f3SHong Zhang {
148b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
149b1cb36f3SHong Zhang   PetscErrorCode ierr;
150b1cb36f3SHong Zhang 
151b1cb36f3SHong Zhang   PetscFunctionBegin;
152b1cb36f3SHong Zhang   if (th->endpoint) {
153b1cb36f3SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
154b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
155b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
156b1cb36f3SHong Zhang     if (th->Theta!=1) {
157b1cb36f3SHong Zhang       ierr = TSAdjointComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1581566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_costintegral,ts->time_step*(th->Theta-1),ts->vec_costintegrand);CHKERRQ(ierr);
159b1cb36f3SHong Zhang     }
160b1cb36f3SHong Zhang   } else {
161b1cb36f3SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
162b1cb36f3SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
163b1cb36f3SHong Zhang   }
16424655328SShri   PetscFunctionReturn(0);
16524655328SShri }
16624655328SShri 
1671566a47fSLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
1681566a47fSLisandro Dalcin {
1691566a47fSLisandro Dalcin   PetscInt       nits,lits;
1701566a47fSLisandro Dalcin   PetscErrorCode ierr;
1711566a47fSLisandro Dalcin 
1721566a47fSLisandro Dalcin   PetscFunctionBegin;
1731566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1741566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1751566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1761566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1771566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1781566a47fSLisandro Dalcin }
1791566a47fSLisandro Dalcin 
180193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
181316643e7SJed Brown {
182316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1831566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1844957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1851566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
186051f2191SLisandro Dalcin   PetscErrorCode ierr;
187316643e7SJed Brown 
188316643e7SJed Brown   PetscFunctionBegin;
1891566a47fSLisandro Dalcin   if (!ts->steprollback) {
1902ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
1913b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
1921566a47fSLisandro Dalcin   }
1931566a47fSLisandro Dalcin 
1941566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
1951566a47fSLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
196*715f1b00SHong Zhang 
1971566a47fSLisandro Dalcin     PetscReal shift = 1/(th->Theta*ts->time_step);
1981566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
199316643e7SJed Brown 
200be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
201fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
202be5899b3SLisandro Dalcin       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
203be5899b3SLisandro Dalcin     }
204eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
205eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2061566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2071566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2081566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2091566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2101566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
211eb284becSJed Brown     }
212be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
2131566a47fSLisandro Dalcin     ierr = TS_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
214be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2151566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
216fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
217051f2191SLisandro Dalcin 
218051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2191566a47fSLisandro Dalcin     if (th->endpoint) {
2201566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2211566a47fSLisandro Dalcin     } else {
222cb3a5882SLisandro Dalcin       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
2231566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2241566a47fSLisandro Dalcin     }
2251566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2261566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2271566a47fSLisandro Dalcin     if (!accept) {
2281566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
229be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
230be5899b3SLisandro Dalcin       goto reject_step;
231051f2191SLisandro Dalcin     }
2323b1890cdSShri Abhyankar 
233*715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
234b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
235b1cb36f3SHong Zhang       th->time_step = ts->time_step;
23617145e75SHong Zhang     }
2371566a47fSLisandro Dalcin 
2382b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
239cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
240051f2191SLisandro Dalcin     break;
241051f2191SLisandro Dalcin 
242051f2191SLisandro Dalcin   reject_step:
243fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2441566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
245051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2461566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2473b1890cdSShri Abhyankar     }
2483b1890cdSShri Abhyankar   }
249316643e7SJed Brown   PetscFunctionReturn(0);
250316643e7SJed Brown }
251316643e7SJed Brown 
25242f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
2532ca6e920SHong Zhang {
2542ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
255b5e0532dSHong Zhang   Vec                 *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
2562ca6e920SHong Zhang   PetscInt            nadj;
2572ca6e920SHong Zhang   PetscErrorCode      ierr;
2582ca6e920SHong Zhang   Mat                 J,Jp;
2592ca6e920SHong Zhang   KSP                 ksp;
2602ca6e920SHong Zhang   PetscReal           shift;
2612ca6e920SHong Zhang 
2622ca6e920SHong Zhang   PetscFunctionBegin;
2632ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
264302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
2652ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
266b5e0532dSHong Zhang 
267b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
2683fd52205SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
269b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
2702ca6e920SHong Zhang 
271a4cab896SHong Zhang   /* Build RHS */
2722c39e106SBarry Smith   if (ts->vec_costintegral) { /* Cost function has an integral term */
27305755b9cSHong Zhang     if (th->endpoint) {
274d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
27505755b9cSHong Zhang     } else {
276d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
27705755b9cSHong Zhang     }
27805755b9cSHong Zhang   }
279abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
280b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
281b5e0532dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
2822c39e106SBarry Smith     if (ts->vec_costintegral) {
283b5e0532dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
28436eaed60SHong Zhang     }
2852ca6e920SHong Zhang   }
2863c54f92cSHong Zhang 
2872ca6e920SHong Zhang   /* Build LHS */
2882ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
2893c54f92cSHong Zhang   if (th->endpoint) {
2903c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2913c54f92cSHong Zhang   } else {
2923c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2933c54f92cSHong Zhang   }
2942ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
2952ca6e920SHong Zhang 
2962ca6e920SHong Zhang   /* Solve LHS X = RHS */
297abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
298b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
2992ca6e920SHong Zhang   }
3003c54f92cSHong Zhang 
30136eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
302b5e0532dSHong Zhang   if(th->endpoint) { /* two-stage case */
303b5e0532dSHong Zhang     if (th->Theta!=1.) {
3043fd52205SHong Zhang       shift = -1./((th->Theta-1.)*ts->time_step);
305b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3062c39e106SBarry Smith       if (ts->vec_costintegral) {
307b5e0532dSHong Zhang         ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
3088263dbbfSHong Zhang       }
309abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
310b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
3112c39e106SBarry Smith         if (ts->vec_costintegral) {
31236eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
31336eaed60SHong Zhang         }
3143c54f92cSHong Zhang         ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3152ca6e920SHong Zhang       }
316b5e0532dSHong Zhang     } else { /* backward Euler */
317b5e0532dSHong Zhang       shift = 0.0;
318b5e0532dSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
319b5e0532dSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
320b5e0532dSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
321b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
322b5e0532dSHong Zhang         if (ts->vec_costintegral) {
323b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
324b5e0532dSHong Zhang         }
325b5e0532dSHong Zhang       }
326b5e0532dSHong Zhang     }
3273fd52205SHong Zhang 
3283fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
3295bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
330abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
331b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
332b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3333fd52205SHong Zhang       }
334b5e0532dSHong Zhang       if (th->Theta!=1.) {
335b5e0532dSHong Zhang         ierr = TSAdjointComputeRHSJacobian(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr);
336abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
337b5e0532dSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
338b5e0532dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
339b5e0532dSHong Zhang         }
3403fd52205SHong Zhang       }
3412c39e106SBarry Smith       if (ts->vec_costintegral) {
342d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
343abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
34436eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
34536eaed60SHong Zhang         }
346b5e0532dSHong Zhang         if (th->Theta!=1.) {
347b5e0532dSHong Zhang           ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
348abc2977eSBarry Smith           for (nadj=0; nadj<ts->numcost; nadj++) {
34936eaed60SHong Zhang             ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
35036eaed60SHong Zhang           }
35136eaed60SHong Zhang         }
3523fd52205SHong Zhang       }
353b5e0532dSHong Zhang     }
3543fd52205SHong Zhang   } else { /* one-stage case */
3553c54f92cSHong Zhang     shift = 0.0;
356a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
3572c39e106SBarry Smith     if (ts->vec_costintegral) {
358d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
3593c54f92cSHong Zhang     }
360abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
361b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
362b5e0532dSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
3632c39e106SBarry Smith       if (ts->vec_costintegral) {
36436eaed60SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
36536eaed60SHong Zhang       }
3662ca6e920SHong Zhang     }
3673fd52205SHong Zhang     if (ts->vecs_sensip) {
3685bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
369abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
370b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
371b5e0532dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
3723fd52205SHong Zhang       }
3732c39e106SBarry Smith       if (ts->vec_costintegral) {
374d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
375abc2977eSBarry Smith         for (nadj=0; nadj<ts->numcost; nadj++) {
37636eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
37736eaed60SHong Zhang         }
37836eaed60SHong Zhang       }
3793fd52205SHong Zhang     }
3802ca6e920SHong Zhang   }
3812ca6e920SHong Zhang 
3822ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
3832ca6e920SHong Zhang   PetscFunctionReturn(0);
3842ca6e920SHong Zhang }
3852ca6e920SHong Zhang 
386cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
387cd652676SJed Brown {
388cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
389be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
390cd652676SJed Brown   PetscErrorCode ierr;
391cd652676SJed Brown 
392cd652676SJed Brown   PetscFunctionBegin;
393a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
394be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
395be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
396cd652676SJed Brown   PetscFunctionReturn(0);
397cd652676SJed Brown }
398cd652676SJed Brown 
3991566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
4001566a47fSLisandro Dalcin {
4011566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
4021566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
4031566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
4047453f775SEmil Constantinescu   PetscReal      wltea,wlter;
4051566a47fSLisandro Dalcin   PetscErrorCode ierr;
4061566a47fSLisandro Dalcin 
4071566a47fSLisandro Dalcin   PetscFunctionBegin;
4082ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
4091566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
410fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
4111566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
412fecfb714SLisandro Dalcin   {
413be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
414be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
4151566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
4161566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
4171566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
4181566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
4191566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
4207453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
4211566a47fSLisandro Dalcin   }
4221566a47fSLisandro Dalcin   if (order) *order = 2;
4231566a47fSLisandro Dalcin   PetscFunctionReturn(0);
4241566a47fSLisandro Dalcin }
4251566a47fSLisandro Dalcin 
4261566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
4271566a47fSLisandro Dalcin {
4281566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
4291566a47fSLisandro Dalcin   PetscErrorCode ierr;
4301566a47fSLisandro Dalcin 
4311566a47fSLisandro Dalcin   PetscFunctionBegin;
4321566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
4331566a47fSLisandro Dalcin   if (ts->vec_costintegral && ts->costintegralfwd) {
4341566a47fSLisandro Dalcin     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4351566a47fSLisandro Dalcin   }
436*715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
437*715f1b00SHong Zhang   PetscFunctionReturn(0);
438*715f1b00SHong Zhang }
439*715f1b00SHong Zhang 
440*715f1b00SHong Zhang static PetscErrorCode TSThetaIntegrandDerivative(TS ts,PetscInt nump,Vec DRncostDY,Vec* s,Vec DRncostDP,Vec VecIntegrandDerivative)
441*715f1b00SHong Zhang {
442*715f1b00SHong Zhang   PetscInt    ndir;
443*715f1b00SHong Zhang   PetscScalar *v;
444*715f1b00SHong Zhang   PetscErrorCode ierr;
445*715f1b00SHong Zhang 
446*715f1b00SHong Zhang   PetscFunctionBegin;
447*715f1b00SHong Zhang 
448*715f1b00SHong Zhang   ierr = VecGetArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
449*715f1b00SHong Zhang   for (ndir=0; ndir<nump; ndir++) {
450*715f1b00SHong Zhang     ierr = VecDot(DRncostDY,s[ndir],&v[ndir]);CHKERRQ(ierr);
451*715f1b00SHong Zhang   }
452*715f1b00SHong Zhang   ierr = VecRestoreArray(VecIntegrandDerivative,&v);CHKERRQ(ierr);
453*715f1b00SHong Zhang   if (DRncostDP) {
454*715f1b00SHong Zhang     ierr = VecAXPY(VecIntegrandDerivative,1.,DRncostDP);CHKERRQ(ierr);
455*715f1b00SHong Zhang   }
456*715f1b00SHong Zhang   PetscFunctionReturn(0);
457*715f1b00SHong Zhang }
458*715f1b00SHong Zhang 
459*715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
460*715f1b00SHong Zhang {
461*715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
462*715f1b00SHong Zhang   Vec            *VecsDeltaFwdSensi  = th->VecsDeltaFwdSensi;
463*715f1b00SHong Zhang   Vec            *VecsDeltaFwdSensip = th->VecsDeltaFwdSensip;
464*715f1b00SHong Zhang   PetscInt       ndir,ncost;
465*715f1b00SHong Zhang   KSP            ksp;
466*715f1b00SHong Zhang   Mat            J,Jp;
467*715f1b00SHong Zhang   PetscReal      shift;
468*715f1b00SHong Zhang   PetscErrorCode ierr;
469*715f1b00SHong Zhang 
470*715f1b00SHong Zhang   PetscFunctionBegin;
471*715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
472*715f1b00SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
473*715f1b00SHong Zhang 
474*715f1b00SHong Zhang   /* Build RHS */
475*715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
476*715f1b00SHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
477*715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
478*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
479*715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
480*715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensip[ndir],(th->Theta-1.)/th->Theta);
481*715f1b00SHong Zhang     }
482*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
483*715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
484*715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensi[ndir],(th->Theta-1.)/th->Theta);
485*715f1b00SHong Zhang     }
486*715f1b00SHong Zhang 
487*715f1b00SHong Zhang     if (ts->num_parameters) { /* Add the f_p forcing terms */
488*715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->ptime,th->X0,ts->vecs_jacp);CHKERRQ(ierr);
489*715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
490*715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],(1.-th->Theta)/th->Theta,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
491*715f1b00SHong Zhang       }
492*715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->vecs_jacp);CHKERRQ(ierr);
493*715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
494*715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
495*715f1b00SHong Zhang       }
496*715f1b00SHong Zhang     }
497*715f1b00SHong Zhang   } else { /* 1-stage method */
498*715f1b00SHong Zhang     shift = 0.0;
499*715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
500*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
501*715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensip[ndir],VecsDeltaFwdSensip[ndir]);
502*715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensip[ndir],-1);
503*715f1b00SHong Zhang     }
504*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
505*715f1b00SHong Zhang       ierr = MatMult(J,ts->vecs_fwdsensi[ndir],VecsDeltaFwdSensi[ndir]);
506*715f1b00SHong Zhang       ierr = VecScale(VecsDeltaFwdSensi[ndir],-1);
507*715f1b00SHong Zhang     }
508*715f1b00SHong Zhang     if (ts->num_parameters) { /* Add the f_p forcing terms */
509*715f1b00SHong Zhang       ierr = TSForwardComputeRHSJacobianP(ts,th->stage_time,th->X,ts->vecs_jacp);CHKERRQ(ierr);
510*715f1b00SHong Zhang       for (ndir=0; ndir<ts->num_parameters; ndir++) {
511*715f1b00SHong Zhang         ierr = VecAXPY(VecsDeltaFwdSensip[ndir],1,ts->vecs_jacp[ndir]);CHKERRQ(ierr);
512*715f1b00SHong Zhang       }
513*715f1b00SHong Zhang     }
514*715f1b00SHong Zhang   }
515*715f1b00SHong Zhang 
516*715f1b00SHong Zhang   /* Build LHS */
517*715f1b00SHong Zhang   shift = 1/(th->Theta*th->time_step);
518*715f1b00SHong Zhang   if (th->endpoint) {
519*715f1b00SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
520*715f1b00SHong Zhang   } else {
521*715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
522*715f1b00SHong Zhang   }
523*715f1b00SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
524*715f1b00SHong Zhang 
525*715f1b00SHong Zhang   /*
526*715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
527*715f1b00SHong Zhang     drdy|t_n*S(t_n) + drdp|t_n
528*715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
529*715f1b00SHong 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
530*715f1b00SHong Zhang    */
531*715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
532*715f1b00SHong Zhang     if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
533*715f1b00SHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,th->ptime,th->X0,ts->vecs_drdy);CHKERRQ(ierr);
534*715f1b00SHong Zhang     }
535*715f1b00SHong Zhang     if (ts->vecs_integral_sensip) {
536*715f1b00SHong Zhang       ierr = TSAdjointComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
537*715f1b00SHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
538*715f1b00SHong Zhang         ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
539*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
540*715f1b00SHong Zhang       }
541*715f1b00SHong Zhang     }
542*715f1b00SHong Zhang     if (ts->vecs_integral_sensi) {
543*715f1b00SHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
544*715f1b00SHong Zhang         ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
545*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensiTemp);CHKERRQ(ierr);
546*715f1b00SHong Zhang       }
547*715f1b00SHong Zhang     }
548*715f1b00SHong Zhang   }
549*715f1b00SHong Zhang 
550*715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
551*715f1b00SHong Zhang   for (ndir=0; ndir<ts->num_parameters; ndir++) {
552*715f1b00SHong Zhang     if (th->endpoint) {
553*715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],ts->vecs_fwdsensip[ndir]);
554*715f1b00SHong Zhang     } else {
555*715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensip[ndir],VecsDeltaFwdSensip[ndir]);
556*715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],1.,VecsDeltaFwdSensip[ndir]);
557*715f1b00SHong Zhang     }
558*715f1b00SHong Zhang   }
559*715f1b00SHong Zhang   /* Solve the linear system for trajectory sensitivities to initial values */
560*715f1b00SHong Zhang   for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
561*715f1b00SHong Zhang     if (th->endpoint) {
562*715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],ts->vecs_fwdsensi[ndir]);
563*715f1b00SHong Zhang     } else {
564*715f1b00SHong Zhang       ierr = KSPSolve(ksp,VecsDeltaFwdSensi[ndir],VecsDeltaFwdSensi[ndir]);
565*715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],1.,VecsDeltaFwdSensi[ndir]);
566*715f1b00SHong Zhang     }
567*715f1b00SHong Zhang   }
568*715f1b00SHong Zhang 
569*715f1b00SHong Zhang   /*Evaluate the second stage of integral gradients with the 2-stage method:
570*715f1b00SHong Zhang     drdy|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
571*715f1b00SHong Zhang   */
572*715f1b00SHong Zhang   if (ts->vecs_integral_sensi || ts->vecs_integral_sensip) {
573*715f1b00SHong Zhang     ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
574*715f1b00SHong Zhang   }
575*715f1b00SHong Zhang   if (ts->vecs_integral_sensip) {
576*715f1b00SHong Zhang     ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
577*715f1b00SHong Zhang     for (ncost=0; ncost<ts->numcost; ncost++) {
578*715f1b00SHong Zhang       ierr = TSThetaIntegrandDerivative(ts,ts->num_parameters,ts->vecs_drdy[ncost],ts->vecs_fwdsensip,ts->vecs_drdp[ncost],th->VecIntegralSensipTemp);
579*715f1b00SHong Zhang       if (th->endpoint) {
580*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
581*715f1b00SHong Zhang       } else {
582*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
583*715f1b00SHong Zhang       }
584*715f1b00SHong Zhang     }
585*715f1b00SHong Zhang   }
586*715f1b00SHong Zhang   if (ts->vecs_integral_sensi) {
587*715f1b00SHong Zhang     for (ncost=0; ncost<ts->numcost; ncost++) {
588*715f1b00SHong Zhang       ierr = TSThetaIntegrandDerivative(ts,ts->num_initialvalues,ts->vecs_drdy[ncost],ts->vecs_fwdsensi,NULL,th->VecIntegralSensiTemp);
589*715f1b00SHong Zhang       if (th->endpoint) {
590*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step*th->Theta,th->VecIntegralSensiTemp);CHKERRQ(ierr);
591*715f1b00SHong Zhang       } else {
592*715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensi[ncost],th->time_step,th->VecIntegralSensiTemp);CHKERRQ(ierr);
593*715f1b00SHong Zhang       }
594*715f1b00SHong Zhang     }
595*715f1b00SHong Zhang   }
596*715f1b00SHong Zhang 
597*715f1b00SHong Zhang   if (!th->endpoint) { /* VecsDeltaFwdSensip are the increment for the stage values for the 1-stage method */
598*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_parameters; ndir++) {
599*715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensip[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensip[ndir]);
600*715f1b00SHong Zhang     }
601*715f1b00SHong Zhang     for (ndir=0; ndir<ts->num_initialvalues; ndir++) {
602*715f1b00SHong Zhang       ierr = VecAXPY(ts->vecs_fwdsensi[ndir],(1.-th->Theta)/th->Theta,VecsDeltaFwdSensi[ndir]);
603*715f1b00SHong Zhang     }
604*715f1b00SHong Zhang   }
605*715f1b00SHong Zhang 
6061566a47fSLisandro Dalcin   PetscFunctionReturn(0);
6071566a47fSLisandro Dalcin }
6081566a47fSLisandro Dalcin 
609316643e7SJed Brown /*------------------------------------------------------------*/
610277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
611316643e7SJed Brown {
612316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
613316643e7SJed Brown   PetscErrorCode ierr;
614316643e7SJed Brown 
615316643e7SJed Brown   PetscFunctionBegin;
6166bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
6176bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
6183b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
619eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
6201566a47fSLisandro Dalcin 
6211566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
6221566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
6231566a47fSLisandro Dalcin 
624b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
625*715f1b00SHong Zhang   if (ts->forward_solve) {
626*715f1b00SHong Zhang     if (ts->vecs_fwdsensi) {
627*715f1b00SHong Zhang       ierr = VecDestroyVecs(ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);    }
628*715f1b00SHong Zhang     if (ts->vecs_fwdsensip) {
629*715f1b00SHong Zhang       ierr = VecDestroyVecs(ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
630*715f1b00SHong Zhang     }
631*715f1b00SHong Zhang     if (ts->vecs_integral_sensi) {
632*715f1b00SHong Zhang       ierr = VecDestroy(&th->VecIntegralSensiTemp);CHKERRQ(ierr);
633*715f1b00SHong Zhang     }
634*715f1b00SHong Zhang     if (ts->vecs_integral_sensip) {
635*715f1b00SHong Zhang       ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
636*715f1b00SHong Zhang     }
637*715f1b00SHong Zhang   }
638b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
639b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
640b5e0532dSHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
641*715f1b00SHong Zhang 
642277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
643277b19d0SLisandro Dalcin }
644277b19d0SLisandro Dalcin 
645277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
646277b19d0SLisandro Dalcin {
647277b19d0SLisandro Dalcin   PetscErrorCode ierr;
648277b19d0SLisandro Dalcin 
649277b19d0SLisandro Dalcin   PetscFunctionBegin;
650277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
651277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
652bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
653bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
654bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
655bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
656316643e7SJed Brown   PetscFunctionReturn(0);
657316643e7SJed Brown }
658316643e7SJed Brown 
659316643e7SJed Brown /*
660316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
6612b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
662316643e7SJed Brown */
6630f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
664316643e7SJed Brown {
665316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
666316643e7SJed Brown   PetscErrorCode ierr;
6677445fe48SJed Brown   Vec            X0,Xdot;
6687445fe48SJed Brown   DM             dm,dmsave;
6691566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
670316643e7SJed Brown 
671316643e7SJed Brown   PetscFunctionBegin;
6727445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
6735a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
6747445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
675b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
6767445fe48SJed Brown 
6777445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
6787445fe48SJed Brown   dmsave = ts->dm;
6797445fe48SJed Brown   ts->dm = dm;
6807445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
6817445fe48SJed Brown   ts->dm = dmsave;
6820d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
683316643e7SJed Brown   PetscFunctionReturn(0);
684316643e7SJed Brown }
685316643e7SJed Brown 
686d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
687316643e7SJed Brown {
688316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
689316643e7SJed Brown   PetscErrorCode ierr;
6907445fe48SJed Brown   Vec            Xdot;
6917445fe48SJed Brown   DM             dm,dmsave;
6921566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
693316643e7SJed Brown 
694316643e7SJed Brown   PetscFunctionBegin;
6957445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
696be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
6970298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
6987445fe48SJed Brown 
6997445fe48SJed Brown   dmsave = ts->dm;
7007445fe48SJed Brown   ts->dm = dm;
701d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
7027445fe48SJed Brown   ts->dm = dmsave;
7030298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
704316643e7SJed Brown   PetscFunctionReturn(0);
705316643e7SJed Brown }
706316643e7SJed Brown 
707*715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
708*715f1b00SHong Zhang {
709*715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
710*715f1b00SHong Zhang   PetscErrorCode ierr;
711*715f1b00SHong Zhang 
712*715f1b00SHong Zhang   PetscFunctionBegin;
713*715f1b00SHong Zhang   if (ts->vecs_fwdsensi) {
714*715f1b00SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_fwdsensi[0],ts->num_initialvalues,&th->VecsDeltaFwdSensi);CHKERRQ(ierr);
715*715f1b00SHong Zhang   }
716*715f1b00SHong Zhang   if (ts->vecs_fwdsensip) {
717*715f1b00SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_fwdsensip[0],ts->num_parameters,&th->VecsDeltaFwdSensip);CHKERRQ(ierr);
718*715f1b00SHong Zhang   }
719*715f1b00SHong Zhang   if (ts->vecs_integral_sensi) {
720*715f1b00SHong Zhang     ierr = VecDuplicate(ts->vecs_integral_sensi[0],&th->VecIntegralSensiTemp);CHKERRQ(ierr);
721*715f1b00SHong Zhang   }
722*715f1b00SHong Zhang   if (ts->vecs_integral_sensip) {
723*715f1b00SHong Zhang     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
724*715f1b00SHong Zhang   }
725*715f1b00SHong Zhang   PetscFunctionReturn(0);
726*715f1b00SHong Zhang }
727*715f1b00SHong Zhang 
728316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
729316643e7SJed Brown {
730316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
7312ffb9264SLisandro Dalcin   PetscBool      match;
732316643e7SJed Brown   PetscErrorCode ierr;
733316643e7SJed Brown 
734316643e7SJed Brown   PetscFunctionBegin;
735b1cb36f3SHong Zhang   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
736b1cb36f3SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
737b1cb36f3SHong Zhang   }
73839d13666SHong Zhang   if (!th->X) {
739316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
74039d13666SHong Zhang   }
74139d13666SHong Zhang   if (!th->Xdot) {
742316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
74339d13666SHong Zhang   }
74439d13666SHong Zhang   if (!th->X0) {
7453b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
74639d13666SHong Zhang   }
7471566a47fSLisandro Dalcin   if (th->endpoint) {
7481566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
7497445fe48SJed Brown   }
7503b1890cdSShri Abhyankar 
7511566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
7521566a47fSLisandro Dalcin 
7531566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
7541566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
7551566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
7561566a47fSLisandro Dalcin 
7571566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
7581566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
7592ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
7602ffb9264SLisandro Dalcin   if (!match) {
7611566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
7621566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
7633b1890cdSShri Abhyankar   }
7641566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
765b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
766b4dd3bf9SBarry Smith }
7670c7376e5SHong Zhang 
768b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
769b4dd3bf9SBarry Smith 
770b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
771b4dd3bf9SBarry Smith {
772b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
773b4dd3bf9SBarry Smith   PetscErrorCode ierr;
774b4dd3bf9SBarry Smith 
775b4dd3bf9SBarry Smith   PetscFunctionBegin;
776b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
7772ca6e920SHong Zhang   if(ts->vecs_sensip) {
778b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
7792ca6e920SHong Zhang   }
780b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
781316643e7SJed Brown   PetscFunctionReturn(0);
782316643e7SJed Brown }
783316643e7SJed Brown 
7844416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
785316643e7SJed Brown {
786316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
787316643e7SJed Brown   PetscErrorCode ierr;
788316643e7SJed Brown 
789316643e7SJed Brown   PetscFunctionBegin;
790e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
791316643e7SJed Brown   {
7920298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
7930298fd71SBarry 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);
79403842d09SLisandro 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);
795316643e7SJed Brown   }
796316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
797316643e7SJed Brown   PetscFunctionReturn(0);
798316643e7SJed Brown }
799316643e7SJed Brown 
800316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
801316643e7SJed Brown {
802316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
803ace3abfcSBarry Smith   PetscBool      iascii;
804316643e7SJed Brown   PetscErrorCode ierr;
805316643e7SJed Brown 
806316643e7SJed Brown   PetscFunctionBegin;
807251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
808316643e7SJed Brown   if (iascii) {
8097c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
810ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
811316643e7SJed Brown   }
812316643e7SJed Brown   PetscFunctionReturn(0);
813316643e7SJed Brown }
814316643e7SJed Brown 
815560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
8160de4c49aSJed Brown {
8170de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
8180de4c49aSJed Brown 
8190de4c49aSJed Brown   PetscFunctionBegin;
8200de4c49aSJed Brown   *theta = th->Theta;
8210de4c49aSJed Brown   PetscFunctionReturn(0);
8220de4c49aSJed Brown }
8230de4c49aSJed Brown 
824560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
8250de4c49aSJed Brown {
8260de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
8270de4c49aSJed Brown 
8280de4c49aSJed Brown   PetscFunctionBegin;
8297c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
8300de4c49aSJed Brown   th->Theta = theta;
8311566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
8320de4c49aSJed Brown   PetscFunctionReturn(0);
8330de4c49aSJed Brown }
834eb284becSJed Brown 
835560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
83626f2ff8fSLisandro Dalcin {
83726f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
83826f2ff8fSLisandro Dalcin 
83926f2ff8fSLisandro Dalcin   PetscFunctionBegin;
84026f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
84126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
84226f2ff8fSLisandro Dalcin }
84326f2ff8fSLisandro Dalcin 
844560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
845eb284becSJed Brown {
846eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
847eb284becSJed Brown 
848eb284becSJed Brown   PetscFunctionBegin;
849eb284becSJed Brown   th->endpoint = flg;
850eb284becSJed Brown   PetscFunctionReturn(0);
851eb284becSJed Brown }
8520de4c49aSJed Brown 
853f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
854f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
855f9c1d6abSBarry Smith {
856f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
857f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
8583fd8ae06SJed Brown   const PetscReal one = 1.0;
859f9c1d6abSBarry Smith 
860f9c1d6abSBarry Smith   PetscFunctionBegin;
8613fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
862f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
863f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
864f9c1d6abSBarry Smith   PetscFunctionReturn(0);
865f9c1d6abSBarry Smith }
866f9c1d6abSBarry Smith #endif
867f9c1d6abSBarry Smith 
86842682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
86942682096SHong Zhang {
87042682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
87142682096SHong Zhang 
87242682096SHong Zhang   PetscFunctionBegin;
8731566a47fSLisandro Dalcin   if (ns) *ns = 1;
8741566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
87542682096SHong Zhang   PetscFunctionReturn(0);
87642682096SHong Zhang }
877f9c1d6abSBarry Smith 
878316643e7SJed Brown /* ------------------------------------------------------------ */
879316643e7SJed Brown /*MC
88096f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
881316643e7SJed Brown 
882316643e7SJed Brown    Level: beginner
883316643e7SJed Brown 
8844eb428fdSBarry Smith    Options Database:
885c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
886c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
88703842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
8884eb428fdSBarry Smith 
889eb284becSJed Brown    Notes:
8900c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
8910c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
8924eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
8934eb428fdSBarry Smith 
894eb284becSJed Brown    This method can be applied to DAE.
895eb284becSJed Brown 
896eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
897eb284becSJed Brown 
898eb284becSJed Brown .vb
899eb284becSJed Brown   Theta | Theta
900eb284becSJed Brown   -------------
901eb284becSJed Brown         |  1
902eb284becSJed Brown .ve
903eb284becSJed Brown 
904eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
905eb284becSJed Brown 
906eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
907eb284becSJed Brown 
908eb284becSJed Brown .vb
909eb284becSJed Brown   0 | 0         0
910eb284becSJed Brown   1 | 1-Theta   Theta
911eb284becSJed Brown   -------------------
912eb284becSJed Brown     | 1-Theta   Theta
913eb284becSJed Brown .ve
914eb284becSJed Brown 
915eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
916eb284becSJed Brown 
917eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
918eb284becSJed Brown 
919eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
920eb284becSJed Brown 
9214eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
922eb284becSJed Brown 
923eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
924316643e7SJed Brown 
925316643e7SJed Brown M*/
9268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
927316643e7SJed Brown {
928316643e7SJed Brown   TS_Theta       *th;
929316643e7SJed Brown   PetscErrorCode ierr;
930316643e7SJed Brown 
931316643e7SJed Brown   PetscFunctionBegin;
932277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
933316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
934316643e7SJed Brown   ts->ops->view            = TSView_Theta;
935316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
93642f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
937316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
938cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
9391566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
94024655328SShri   ts->ops->rollback        = TSRollBack_Theta;
941316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
9420f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
9430f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
944f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
945f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
946f9c1d6abSBarry Smith #endif
94742682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
94842f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
949b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
950b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
9512ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
952*715f1b00SHong Zhang   ts->ops->forwardsetup    = TSForwardSetUp_Theta;
953*715f1b00SHong Zhang   ts->ops->forwardstep     = TSForwardStep_Theta;
954316643e7SJed Brown 
955efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
956efd4aadfSBarry Smith 
957b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
958316643e7SJed Brown   ts->data = (void*)th;
959316643e7SJed Brown 
960*715f1b00SHong Zhang   th->VecsDeltaLam   = NULL;
961*715f1b00SHong Zhang   th->VecsDeltaMu    = NULL;
962*715f1b00SHong Zhang   th->VecsSensiTemp  = NULL;
963*715f1b00SHong Zhang 
9646f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
965316643e7SJed Brown   th->Theta       = 0.5;
9661566a47fSLisandro Dalcin   th->order       = 2;
967bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
968bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
969bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
970bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
971316643e7SJed Brown   PetscFunctionReturn(0);
972316643e7SJed Brown }
9730de4c49aSJed Brown 
9740de4c49aSJed Brown /*@
9750de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
9760de4c49aSJed Brown 
9770de4c49aSJed Brown   Not Collective
9780de4c49aSJed Brown 
9790de4c49aSJed Brown   Input Parameter:
9800de4c49aSJed Brown .  ts - timestepping context
9810de4c49aSJed Brown 
9820de4c49aSJed Brown   Output Parameter:
9830de4c49aSJed Brown .  theta - stage abscissa
9840de4c49aSJed Brown 
9850de4c49aSJed Brown   Note:
9860de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
9870de4c49aSJed Brown 
9880de4c49aSJed Brown   Level: Advanced
9890de4c49aSJed Brown 
9900de4c49aSJed Brown .seealso: TSThetaSetTheta()
9910de4c49aSJed Brown @*/
9927087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
9930de4c49aSJed Brown {
9944ac538c5SBarry Smith   PetscErrorCode ierr;
9950de4c49aSJed Brown 
9960de4c49aSJed Brown   PetscFunctionBegin;
997afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
9980de4c49aSJed Brown   PetscValidPointer(theta,2);
9994ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
10000de4c49aSJed Brown   PetscFunctionReturn(0);
10010de4c49aSJed Brown }
10020de4c49aSJed Brown 
10030de4c49aSJed Brown /*@
10040de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
10050de4c49aSJed Brown 
10060de4c49aSJed Brown   Not Collective
10070de4c49aSJed Brown 
10080de4c49aSJed Brown   Input Parameter:
10090de4c49aSJed Brown +  ts - timestepping context
10100de4c49aSJed Brown -  theta - stage abscissa
10110de4c49aSJed Brown 
10120de4c49aSJed Brown   Options Database:
10130de4c49aSJed Brown .  -ts_theta_theta <theta>
10140de4c49aSJed Brown 
10150de4c49aSJed Brown   Level: Intermediate
10160de4c49aSJed Brown 
10170de4c49aSJed Brown .seealso: TSThetaGetTheta()
10180de4c49aSJed Brown @*/
10197087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
10200de4c49aSJed Brown {
10214ac538c5SBarry Smith   PetscErrorCode ierr;
10220de4c49aSJed Brown 
10230de4c49aSJed Brown   PetscFunctionBegin;
1024afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
10254ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
10260de4c49aSJed Brown   PetscFunctionReturn(0);
10270de4c49aSJed Brown }
1028f33bbcb6SJed Brown 
102926f2ff8fSLisandro Dalcin /*@
103026f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
103126f2ff8fSLisandro Dalcin 
103226f2ff8fSLisandro Dalcin   Not Collective
103326f2ff8fSLisandro Dalcin 
103426f2ff8fSLisandro Dalcin   Input Parameter:
103526f2ff8fSLisandro Dalcin .  ts - timestepping context
103626f2ff8fSLisandro Dalcin 
103726f2ff8fSLisandro Dalcin   Output Parameter:
103826f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
103926f2ff8fSLisandro Dalcin 
104026f2ff8fSLisandro Dalcin   Level: Advanced
104126f2ff8fSLisandro Dalcin 
104226f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
104326f2ff8fSLisandro Dalcin @*/
104426f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
104526f2ff8fSLisandro Dalcin {
104626f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
104726f2ff8fSLisandro Dalcin 
104826f2ff8fSLisandro Dalcin   PetscFunctionBegin;
104926f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
105026f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1051163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
105226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
105326f2ff8fSLisandro Dalcin }
105426f2ff8fSLisandro Dalcin 
1055eb284becSJed Brown /*@
1056eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1057eb284becSJed Brown 
1058eb284becSJed Brown   Not Collective
1059eb284becSJed Brown 
1060eb284becSJed Brown   Input Parameter:
1061eb284becSJed Brown +  ts - timestepping context
1062eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1063eb284becSJed Brown 
1064eb284becSJed Brown   Options Database:
1065eb284becSJed Brown .  -ts_theta_endpoint <flg>
1066eb284becSJed Brown 
1067eb284becSJed Brown   Level: Intermediate
1068eb284becSJed Brown 
1069eb284becSJed Brown .seealso: TSTHETA, TSCN
1070eb284becSJed Brown @*/
1071eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1072eb284becSJed Brown {
1073eb284becSJed Brown   PetscErrorCode ierr;
1074eb284becSJed Brown 
1075eb284becSJed Brown   PetscFunctionBegin;
1076eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1077eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1078eb284becSJed Brown   PetscFunctionReturn(0);
1079eb284becSJed Brown }
1080eb284becSJed Brown 
1081f33bbcb6SJed Brown /*
1082f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1083f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1084f33bbcb6SJed Brown  */
1085f33bbcb6SJed Brown 
10861566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
10871566a47fSLisandro Dalcin {
10881566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
10891566a47fSLisandro Dalcin   PetscErrorCode ierr;
10901566a47fSLisandro Dalcin 
10911566a47fSLisandro Dalcin   PetscFunctionBegin;
10921566a47fSLisandro 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");
10931566a47fSLisandro 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");
10941566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
10951566a47fSLisandro Dalcin   PetscFunctionReturn(0);
10961566a47fSLisandro Dalcin }
10971566a47fSLisandro Dalcin 
1098f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1099f33bbcb6SJed Brown {
1100f33bbcb6SJed Brown   PetscFunctionBegin;
1101f33bbcb6SJed Brown   PetscFunctionReturn(0);
1102f33bbcb6SJed Brown }
1103f33bbcb6SJed Brown 
1104f33bbcb6SJed Brown /*MC
1105f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1106f33bbcb6SJed Brown 
1107f33bbcb6SJed Brown   Level: beginner
1108f33bbcb6SJed Brown 
11094eb428fdSBarry Smith   Notes:
1110c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
11114eb428fdSBarry Smith 
11121566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
11134eb428fdSBarry Smith 
1114f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1115f33bbcb6SJed Brown 
1116f33bbcb6SJed Brown M*/
11178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1118f33bbcb6SJed Brown {
1119f33bbcb6SJed Brown   PetscErrorCode ierr;
1120f33bbcb6SJed Brown 
1121f33bbcb6SJed Brown   PetscFunctionBegin;
1122f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1123f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
11241566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
11250c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1126f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1127f33bbcb6SJed Brown   PetscFunctionReturn(0);
1128f33bbcb6SJed Brown }
1129f33bbcb6SJed Brown 
11301566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
11311566a47fSLisandro Dalcin {
11321566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
11331566a47fSLisandro Dalcin   PetscErrorCode ierr;
11341566a47fSLisandro Dalcin 
11351566a47fSLisandro Dalcin   PetscFunctionBegin;
11361566a47fSLisandro 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");
11371566a47fSLisandro 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");
11381566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
11391566a47fSLisandro Dalcin   PetscFunctionReturn(0);
11401566a47fSLisandro Dalcin }
11411566a47fSLisandro Dalcin 
1142f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1143f33bbcb6SJed Brown {
1144f33bbcb6SJed Brown   PetscFunctionBegin;
1145f33bbcb6SJed Brown   PetscFunctionReturn(0);
1146f33bbcb6SJed Brown }
1147f33bbcb6SJed Brown 
1148f33bbcb6SJed Brown /*MC
1149f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1150f33bbcb6SJed Brown 
1151f33bbcb6SJed Brown   Level: beginner
1152f33bbcb6SJed Brown 
1153f33bbcb6SJed Brown   Notes:
11547cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
11557cf5af47SJed Brown 
11567cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1157f33bbcb6SJed Brown 
1158f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1159f33bbcb6SJed Brown 
1160f33bbcb6SJed Brown M*/
11618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1162f33bbcb6SJed Brown {
1163f33bbcb6SJed Brown   PetscErrorCode ierr;
1164f33bbcb6SJed Brown 
1165f33bbcb6SJed Brown   PetscFunctionBegin;
1166f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1167f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1168eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
11690c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1170f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1171f33bbcb6SJed Brown   PetscFunctionReturn(0);
1172f33bbcb6SJed Brown }
1173