xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 33f72c5d48e626d09b7d7a978ff73ecaa099928a)
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 */
24022c081aSHong Zhang   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
25b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
26b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
2713b1b0ffSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
2813b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
29b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
3013b1b0ffSHong Zhang   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
31715f1b00SHong Zhang   Vec          VecIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
326f92202bSHong Zhang   Vec          *VecsIntegralSensip0;     /* backup for roll-backs due to events */
33b5b4017aSHong Zhang   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
34b5b4017aSHong Zhang   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
35b5b4017aSHong Zhang   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
36b5b4017aSHong Zhang   Vec          *VecsAffine;              /* Working vectors to store residuals */
37715f1b00SHong Zhang   /* context for error estimation */
381566a47fSLisandro Dalcin   Vec          vec_sol_prev;
391566a47fSLisandro Dalcin   Vec          vec_lte_work;
40316643e7SJed Brown } TS_Theta;
41316643e7SJed Brown 
427445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
437445fe48SJed Brown {
447445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
457445fe48SJed Brown   PetscErrorCode ierr;
467445fe48SJed Brown 
477445fe48SJed Brown   PetscFunctionBegin;
487445fe48SJed Brown   if (X0) {
497445fe48SJed Brown     if (dm && dm != ts->dm) {
500d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
517445fe48SJed Brown     } else *X0 = ts->vec_sol;
527445fe48SJed Brown   }
537445fe48SJed Brown   if (Xdot) {
547445fe48SJed Brown     if (dm && dm != ts->dm) {
550d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
567445fe48SJed Brown     } else *Xdot = th->Xdot;
577445fe48SJed Brown   }
587445fe48SJed Brown   PetscFunctionReturn(0);
597445fe48SJed Brown }
607445fe48SJed Brown 
610d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
620d0b770aSPeter Brune {
630d0b770aSPeter Brune   PetscErrorCode ierr;
640d0b770aSPeter Brune 
650d0b770aSPeter Brune   PetscFunctionBegin;
660d0b770aSPeter Brune   if (X0) {
670d0b770aSPeter Brune     if (dm && dm != ts->dm) {
680d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
690d0b770aSPeter Brune     }
700d0b770aSPeter Brune   }
710d0b770aSPeter Brune   if (Xdot) {
720d0b770aSPeter Brune     if (dm && dm != ts->dm) {
730d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
740d0b770aSPeter Brune     }
750d0b770aSPeter Brune   }
760d0b770aSPeter Brune   PetscFunctionReturn(0);
770d0b770aSPeter Brune }
780d0b770aSPeter Brune 
797445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
807445fe48SJed Brown {
817445fe48SJed Brown   PetscFunctionBegin;
827445fe48SJed Brown   PetscFunctionReturn(0);
837445fe48SJed Brown }
847445fe48SJed Brown 
857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
867445fe48SJed Brown {
877445fe48SJed Brown   TS             ts = (TS)ctx;
887445fe48SJed Brown   PetscErrorCode ierr;
897445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
907445fe48SJed Brown 
917445fe48SJed Brown   PetscFunctionBegin;
927445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
937445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
957445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
977445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
980d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
990d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1007445fe48SJed Brown   PetscFunctionReturn(0);
1017445fe48SJed Brown }
1027445fe48SJed Brown 
103258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
104258e1594SPeter Brune {
105258e1594SPeter Brune   PetscFunctionBegin;
106258e1594SPeter Brune   PetscFunctionReturn(0);
107258e1594SPeter Brune }
108258e1594SPeter Brune 
109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110258e1594SPeter Brune {
111258e1594SPeter Brune   TS             ts = (TS)ctx;
112258e1594SPeter Brune   PetscErrorCode ierr;
113258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
114258e1594SPeter Brune 
115258e1594SPeter Brune   PetscFunctionBegin;
116258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
117258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
118258e1594SPeter Brune 
119258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune 
122258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune 
125258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
126258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
127258e1594SPeter Brune   PetscFunctionReturn(0);
128258e1594SPeter Brune }
129258e1594SPeter Brune 
130022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
131022c081aSHong Zhang {
132022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
133022c081aSHong Zhang   PetscErrorCode ierr;
134022c081aSHong Zhang 
135022c081aSHong Zhang   PetscFunctionBegin;
136022c081aSHong Zhang   if (th->endpoint) {
137022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
138022c081aSHong Zhang     if (th->Theta!=1.0) {
139022c081aSHong Zhang       ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
140022c081aSHong Zhang       ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
141022c081aSHong Zhang     }
142022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
143022c081aSHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
144022c081aSHong Zhang   } else {
145022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
146022c081aSHong Zhang     ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
147022c081aSHong Zhang   }
148022c081aSHong Zhang   PetscFunctionReturn(0);
149022c081aSHong Zhang }
150022c081aSHong Zhang 
151b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
152b1cb36f3SHong Zhang {
153b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
154b1cb36f3SHong Zhang   PetscErrorCode ierr;
155b1cb36f3SHong Zhang 
156b1cb36f3SHong Zhang   PetscFunctionBegin;
157b1cb36f3SHong Zhang   /* backup cost integral */
158b1cb36f3SHong Zhang   ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr);
159022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
160b1cb36f3SHong Zhang   PetscFunctionReturn(0);
161b1cb36f3SHong Zhang }
162b1cb36f3SHong Zhang 
163b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
164b1cb36f3SHong Zhang {
165b1cb36f3SHong Zhang   PetscErrorCode ierr;
166b1cb36f3SHong Zhang 
167b1cb36f3SHong Zhang   PetscFunctionBegin;
168022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
16924655328SShri   PetscFunctionReturn(0);
17024655328SShri }
17124655328SShri 
172470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1731566a47fSLisandro Dalcin {
1741566a47fSLisandro Dalcin   PetscInt       nits,lits;
1751566a47fSLisandro Dalcin   PetscErrorCode ierr;
1761566a47fSLisandro Dalcin 
1771566a47fSLisandro Dalcin   PetscFunctionBegin;
1781566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1791566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1801566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1811566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1821566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1831566a47fSLisandro Dalcin }
1841566a47fSLisandro Dalcin 
185193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
186316643e7SJed Brown {
187316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1881566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1894957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1901566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
191051f2191SLisandro Dalcin   PetscErrorCode ierr;
192316643e7SJed Brown 
193316643e7SJed Brown   PetscFunctionBegin;
1941566a47fSLisandro Dalcin   if (!ts->steprollback) {
1952ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
1963b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
1971566a47fSLisandro Dalcin   }
1981566a47fSLisandro Dalcin 
1991566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
2001566a47fSLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
201715f1b00SHong Zhang 
2021566a47fSLisandro Dalcin     PetscReal shift = 1/(th->Theta*ts->time_step);
2031566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
204316643e7SJed Brown 
205be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
206fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
207be5899b3SLisandro Dalcin       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
208be5899b3SLisandro Dalcin     }
209eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
210eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2111566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2121566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2131566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2141566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2151566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
216eb284becSJed Brown     }
217be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
218470880abSPatrick Sanan     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
219be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2201566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
221fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
222051f2191SLisandro Dalcin 
223051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2241566a47fSLisandro Dalcin     if (th->endpoint) {
2251566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2261566a47fSLisandro Dalcin     } else {
227cb3a5882SLisandro Dalcin       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
2281566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2291566a47fSLisandro Dalcin     }
2301566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2311566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2321566a47fSLisandro Dalcin     if (!accept) {
2331566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
234be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
235be5899b3SLisandro Dalcin       goto reject_step;
236051f2191SLisandro Dalcin     }
2373b1890cdSShri Abhyankar 
238715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
239b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
240b1cb36f3SHong Zhang       th->time_step = ts->time_step;
24117145e75SHong Zhang     }
2421566a47fSLisandro Dalcin 
2432b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
244cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
245051f2191SLisandro Dalcin     break;
246051f2191SLisandro Dalcin 
247051f2191SLisandro Dalcin   reject_step:
248fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2491566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
250051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2511566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2523b1890cdSShri Abhyankar     }
2533b1890cdSShri Abhyankar   }
254316643e7SJed Brown   PetscFunctionReturn(0);
255316643e7SJed Brown }
256316643e7SJed Brown 
257c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
258c9681e5dSHong Zhang {
259c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
260c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
261c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
262c9681e5dSHong Zhang   PetscInt       nadj;
263c9681e5dSHong Zhang   Mat            J,Jp;
264c9681e5dSHong Zhang   KSP            ksp;
265c9681e5dSHong Zhang   PetscReal      shift;
266c9681e5dSHong Zhang   PetscScalar    *xarr;
267c9681e5dSHong Zhang   PetscErrorCode ierr;
268c9681e5dSHong Zhang 
269c9681e5dSHong Zhang   PetscFunctionBegin;
270c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
271c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
272c9681e5dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
273c9681e5dSHong Zhang 
274c9681e5dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
275c9681e5dSHong Zhang   th->stage_time = ts->ptime; /* time_step is negative*/
276c9681e5dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
277c9681e5dSHong Zhang   th->time_step  = -ts->time_step;
278c9681e5dSHong Zhang 
279*33f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
280c9681e5dSHong Zhang   ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
281c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
282c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
283c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
284*33f72c5dSHong Zhang     if (ts->vecs_drdu) {
285c9681e5dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
286c9681e5dSHong Zhang     }
287c9681e5dSHong Zhang   }
288c9681e5dSHong Zhang 
289c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
290c9681e5dSHong Zhang   shift = 1./th->time_step;
291c9681e5dSHong Zhang   ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
292c9681e5dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
293c9681e5dSHong Zhang 
294c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
295c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
296c9681e5dSHong Zhang     KSPConvergedReason kspreason;
297c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
298c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
299c9681e5dSHong Zhang     if (kspreason < 0) {
300c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
301c9681e5dSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
302c9681e5dSHong Zhang     }
303c9681e5dSHong Zhang   }
304c9681e5dSHong Zhang 
305c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
306c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
307c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
308c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
309c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
310c9681e5dSHong Zhang     ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
311c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
312c9681e5dSHong Zhang     ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
313c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
314c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
315c9681e5dSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
316*33f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
317c9681e5dSHong Zhang       if (ts->vecs_fup) {
318*33f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
319c9681e5dSHong Zhang       }
320c9681e5dSHong Zhang     }
321c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
322c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
32305c9950eSHong Zhang       KSPConvergedReason kspreason;
324c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
32505c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
32605c9950eSHong Zhang       if (kspreason < 0) {
32705c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32805c9950eSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
32905c9950eSHong Zhang       }
330c9681e5dSHong Zhang     }
331c9681e5dSHong Zhang   }
332c9681e5dSHong Zhang 
333c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
334c9681e5dSHong Zhang   shift = 0.0;
335*33f72c5dSHong Zhang   ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */
336c9681e5dSHong Zhang   ierr  = MatScale(J,-1.);CHKERRQ(ierr);
337c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
338*33f72c5dSHong Zhang     /* Add f_U \lambda_s to the original RHS */
339c9681e5dSHong Zhang     ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
340c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr);
341c9681e5dSHong Zhang     ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
342c9681e5dSHong Zhang     if (ts->vecs_sensi2) {
343c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
344c9681e5dSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr);
345c9681e5dSHong Zhang       ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
346c9681e5dSHong Zhang     }
347c9681e5dSHong Zhang   }
348c9681e5dSHong Zhang   if (ts->vecs_sensip) {
349c9681e5dSHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
350*33f72c5dSHong Zhang     ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
351c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
352c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
353c9681e5dSHong Zhang       ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
354*33f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
355c9681e5dSHong Zhang       ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
356c9681e5dSHong Zhang     }
357c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
358c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
359c9681e5dSHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
360*33f72c5dSHong Zhang       if (ts->vecs_drdp) {
361*33f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
362*33f72c5dSHong Zhang       }
363c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
364c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
365c9681e5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
366c9681e5dSHong Zhang         if (ts->vecs_fpu) {
367c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
368c9681e5dSHong Zhang         }
369c9681e5dSHong Zhang         if (ts->vecs_fpp) {
370c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
371c9681e5dSHong Zhang         }
372c9681e5dSHong Zhang       }
373c9681e5dSHong Zhang     }
374c9681e5dSHong Zhang   }
375c9681e5dSHong Zhang 
376c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
377c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
378c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
379c9681e5dSHong Zhang   }
380c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
381c9681e5dSHong Zhang   PetscFunctionReturn(0);
382c9681e5dSHong Zhang }
383c9681e5dSHong Zhang 
38442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
3852ca6e920SHong Zhang {
3862ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
387b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
388b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
3892ca6e920SHong Zhang   PetscInt       nadj;
3902ca6e920SHong Zhang   Mat            J,Jp;
3912ca6e920SHong Zhang   KSP            ksp;
3922ca6e920SHong Zhang   PetscReal      shift;
393b5b4017aSHong Zhang   PetscScalar    *xarr;
394b5b4017aSHong Zhang   PetscErrorCode ierr;
3952ca6e920SHong Zhang 
3962ca6e920SHong Zhang   PetscFunctionBegin;
397c9681e5dSHong Zhang   if (th->Theta == 1.) {
398c9681e5dSHong Zhang     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
399c9681e5dSHong Zhang     PetscFunctionReturn(0);
400c9681e5dSHong Zhang   }
4012ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
402302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
4032ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
404b5e0532dSHong Zhang 
405b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
406022c081aSHong Zhang   th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/
407b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
408022c081aSHong Zhang   th->time_step  = -ts->time_step;
4092ca6e920SHong Zhang 
410b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
411*33f72c5dSHong Zhang   /* Cost function has an integral term */
41205755b9cSHong Zhang   if (th->endpoint) {
413c9ad9de0SHong Zhang     ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr);
41405755b9cSHong Zhang   } else {
415c9ad9de0SHong Zhang     ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
41605755b9cSHong Zhang   }
417*33f72c5dSHong Zhang 
418abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
419b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
420022c081aSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr);
4212c39e106SBarry Smith     if (ts->vec_costintegral) {
422c9ad9de0SHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
42336eaed60SHong Zhang     }
4242ca6e920SHong Zhang   }
4253c54f92cSHong Zhang 
426b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
427022c081aSHong Zhang   shift = 1./(th->Theta*th->time_step);
4283c54f92cSHong Zhang   if (th->endpoint) {
429022c081aSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
4303c54f92cSHong Zhang   } else {
4313c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
4323c54f92cSHong Zhang   }
4332ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
4342ca6e920SHong Zhang 
435b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
436abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4371d14d03bSHong Zhang     KSPConvergedReason kspreason;
438b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
4391d14d03bSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
4401d14d03bSHong Zhang     if (kspreason < 0) {
4411d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
4421d14d03bSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
4431d14d03bSHong Zhang     }
4442ca6e920SHong Zhang   }
4453c54f92cSHong Zhang 
4461d14d03bSHong Zhang   /* Second-order adjoint */
447b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
448b5b4017aSHong Zhang     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
449b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
450b5b4017aSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
451b5b4017aSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
452b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
453b5b4017aSHong Zhang     ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
454b5b4017aSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
455b5b4017aSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
456b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
457b5b4017aSHong Zhang     ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
458b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
459b5b4017aSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
4601d14d03bSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
461*33f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
462b5b4017aSHong Zhang       if (ts->vecs_fup) {
463*33f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
464b5b4017aSHong Zhang       }
465b5b4017aSHong Zhang     }
466b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
467b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
4681d14d03bSHong Zhang       KSPConvergedReason kspreason;
4691d14d03bSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
4701d14d03bSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
4711d14d03bSHong Zhang       if (kspreason < 0) {
4721d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
4731d14d03bSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
4741d14d03bSHong Zhang       }
475b5b4017aSHong Zhang     }
476b5b4017aSHong Zhang   }
477b5b4017aSHong Zhang 
47836eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
4791d14d03bSHong Zhang   if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
480022c081aSHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
481b5e0532dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
482*33f72c5dSHong Zhang     /* R_U at t_n */
483c9ad9de0SHong Zhang     ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr);
484abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
485b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
486*33f72c5dSHong Zhang       if (ts->vecs_drdu) {
4871d14d03bSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
488b5b4017aSHong Zhang       }
4891d14d03bSHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
490b5b4017aSHong Zhang     }
4911d14d03bSHong Zhang 
4921d14d03bSHong Zhang     /* Second-order adjoint */
4931d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
494b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
495b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
496b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
497b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
498b5b4017aSHong Zhang       ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
499b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
500*33f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
501b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
502b5b4017aSHong Zhang       ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
503b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
504*33f72c5dSHong Zhang         /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2  */
505b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
506*33f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
507b5b4017aSHong Zhang         if (ts->vecs_fup) {
508*33f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
509b5b4017aSHong Zhang         }
5101d14d03bSHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr);
511b5b4017aSHong Zhang       }
512b5e0532dSHong Zhang     }
5133fd52205SHong Zhang 
5143fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
515b5b4017aSHong Zhang       /* U_{n+1} */
516*33f72c5dSHong Zhang       shift = -1./(th->Theta*th->time_step);
517*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
518b5b4017aSHong Zhang       ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
519abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
520b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
521*33f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5223fd52205SHong Zhang       }
523*33f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
524*33f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
525*33f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
526*33f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
527b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
528b5b4017aSHong Zhang         ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
529*33f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
530*33f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
531*33f72c5dSHong Zhang 
532b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
533b5b4017aSHong Zhang         ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
534b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
535*33f72c5dSHong Zhang           /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2)  */
536b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
537*33f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
538b5b4017aSHong Zhang           if (ts->vecs_fpu) {
539*33f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
540b5b4017aSHong Zhang           }
541b5b4017aSHong Zhang           if (ts->vecs_fpp) {
542*33f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
543b5b4017aSHong Zhang           }
544b5b4017aSHong Zhang         }
545b5b4017aSHong Zhang       }
546b5b4017aSHong Zhang 
547b5b4017aSHong Zhang       /* U_s */
548*33f72c5dSHong Zhang       shift = 1./((th->Theta-1.0)*th->time_step);
549*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
550b5b4017aSHong Zhang       ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
551abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
552b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
553*33f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
554*33f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
555*33f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
556*33f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
557*33f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
558b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
559*33f72c5dSHong Zhang           ierr = TSComputeIHessianProductFunction3(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
560*33f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
561*33f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
562b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
563*33f72c5dSHong Zhang           ierr = TSComputeIHessianProductFunction4(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
564b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
565*33f72c5dSHong Zhang             /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */
566b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
567*33f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
568b5b4017aSHong Zhang             if (ts->vecs_fpu) {
569*33f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
570b5e0532dSHong Zhang             }
571b5b4017aSHong Zhang             if (ts->vecs_fpp) {
572*33f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
573b5b4017aSHong Zhang             }
57436eaed60SHong Zhang           }
57536eaed60SHong Zhang         }
5763fd52205SHong Zhang       }
577b5e0532dSHong Zhang     }
5783fd52205SHong Zhang   } else { /* one-stage case */
5793c54f92cSHong Zhang     shift = 0.0;
580a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
581c9ad9de0SHong Zhang     ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
582abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
583b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
584022c081aSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
585*33f72c5dSHong Zhang       if (ts->vecs_drdu) {
586c9ad9de0SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr);
58736eaed60SHong Zhang       }
5882ca6e920SHong Zhang     }
5893fd52205SHong Zhang     if (ts->vecs_sensip) {
590*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
591c9ad9de0SHong Zhang       ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
592abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
593*33f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
594*33f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
595*33f72c5dSHong Zhang         if (ts->vecs_drdp) {
596022c081aSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
59736eaed60SHong Zhang         }
59836eaed60SHong Zhang       }
5993fd52205SHong Zhang     }
6002ca6e920SHong Zhang   }
6012ca6e920SHong Zhang 
6022ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6032ca6e920SHong Zhang   PetscFunctionReturn(0);
6042ca6e920SHong Zhang }
6052ca6e920SHong Zhang 
606cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
607cd652676SJed Brown {
608cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
609be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
610cd652676SJed Brown   PetscErrorCode ierr;
611cd652676SJed Brown 
612cd652676SJed Brown   PetscFunctionBegin;
613a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
614be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
615be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
616cd652676SJed Brown   PetscFunctionReturn(0);
617cd652676SJed Brown }
618cd652676SJed Brown 
6191566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
6201566a47fSLisandro Dalcin {
6211566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
6221566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
6231566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
6247453f775SEmil Constantinescu   PetscReal      wltea,wlter;
6251566a47fSLisandro Dalcin   PetscErrorCode ierr;
6261566a47fSLisandro Dalcin 
6271566a47fSLisandro Dalcin   PetscFunctionBegin;
6282ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
6291566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
630fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
6311566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
632fecfb714SLisandro Dalcin   {
633be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
634be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
6351566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
6361566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
6371566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
6381566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
6391566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
6407453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
6411566a47fSLisandro Dalcin   }
6421566a47fSLisandro Dalcin   if (order) *order = 2;
6431566a47fSLisandro Dalcin   PetscFunctionReturn(0);
6441566a47fSLisandro Dalcin }
6451566a47fSLisandro Dalcin 
6461566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
6471566a47fSLisandro Dalcin {
6481566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
64913b1b0ffSHong Zhang   PetscInt       ncost;
6501566a47fSLisandro Dalcin   PetscErrorCode ierr;
6511566a47fSLisandro Dalcin 
6521566a47fSLisandro Dalcin   PetscFunctionBegin;
6531566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
6541566a47fSLisandro Dalcin   if (ts->vec_costintegral && ts->costintegralfwd) {
6551566a47fSLisandro Dalcin     ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
6561566a47fSLisandro Dalcin   }
657715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
65813b1b0ffSHong Zhang   if (ts->mat_sensip) {
65913b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6606f92202bSHong Zhang   }
6616f92202bSHong Zhang   if (ts->vecs_integral_sensip) {
6626f92202bSHong Zhang     for (ncost=0;ncost<ts->numcost;ncost++) {
6636f92202bSHong Zhang       ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr);
6646f92202bSHong Zhang     }
6656f92202bSHong Zhang   }
666715f1b00SHong Zhang   PetscFunctionReturn(0);
667715f1b00SHong Zhang }
668715f1b00SHong Zhang 
669715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
670715f1b00SHong Zhang {
671715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
67213b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
673b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
67413b1b0ffSHong Zhang   PetscInt       ncost,ntlm;
675715f1b00SHong Zhang   KSP            ksp;
676715f1b00SHong Zhang   Mat            J,Jp;
677715f1b00SHong Zhang   PetscReal      shift;
67813b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
679715f1b00SHong Zhang   PetscErrorCode ierr;
680715f1b00SHong Zhang 
681715f1b00SHong Zhang   PetscFunctionBegin;
68213b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
68313b1b0ffSHong Zhang 
6846f92202bSHong Zhang   for (ncost=0; ncost<ts->numcost; ncost++) {
6856f92202bSHong Zhang     if (ts->vecs_integral_sensip) {
6866f92202bSHong Zhang       ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr);
6876f92202bSHong Zhang     }
6886f92202bSHong Zhang   }
6896f92202bSHong Zhang 
690715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
691715f1b00SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
692715f1b00SHong Zhang 
693715f1b00SHong Zhang   /* Build RHS */
694715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
695715f1b00SHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
696715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
69713b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
69813b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
699715f1b00SHong Zhang 
700022c081aSHong Zhang     /* Add the f_p forcing terms */
7010e11c21fSHong Zhang     if (ts->Jacp) {
702*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
703*33f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
704*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
705*33f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7060e11c21fSHong Zhang     }
707715f1b00SHong Zhang   } else { /* 1-stage method */
708715f1b00SHong Zhang     shift = 0.0;
709715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
71013b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
71113b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
71213b1b0ffSHong Zhang 
713022c081aSHong Zhang     /* Add the f_p forcing terms */
7140e11c21fSHong Zhang     if (ts->Jacp) {
715*33f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
716*33f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
717715f1b00SHong Zhang     }
7180e11c21fSHong Zhang   }
719715f1b00SHong Zhang 
720715f1b00SHong Zhang   /* Build LHS */
721715f1b00SHong Zhang   shift = 1/(th->Theta*th->time_step);
722715f1b00SHong Zhang   if (th->endpoint) {
723715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
724715f1b00SHong Zhang   } else {
725715f1b00SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
726715f1b00SHong Zhang   }
727715f1b00SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
728715f1b00SHong Zhang 
729715f1b00SHong Zhang   /*
730715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
731c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
732715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
733715f1b00SHong Zhang   */
734715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
735715f1b00SHong Zhang     if (ts->vecs_integral_sensip) {
736c9ad9de0SHong Zhang       ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr);
7370e11c21fSHong Zhang       if (ts->vecs_drdp) {
738c9ad9de0SHong Zhang         ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr);
7390e11c21fSHong Zhang       }
740715f1b00SHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
741c9ad9de0SHong Zhang         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
7420e11c21fSHong Zhang         if (ts->vecs_drdp) {
74313b1b0ffSHong Zhang           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
7440e11c21fSHong Zhang         }
745715f1b00SHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr);
746715f1b00SHong Zhang       }
747715f1b00SHong Zhang     }
748715f1b00SHong Zhang   }
749715f1b00SHong Zhang 
750715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
751022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
752b5b4017aSHong Zhang     KSPConvergedReason kspreason;
75313b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
754b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
755715f1b00SHong Zhang     if (th->endpoint) {
75613b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
757b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
758b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
759b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
76013b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
761715f1b00SHong Zhang     } else {
762b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
763715f1b00SHong Zhang     }
764b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
765b5b4017aSHong Zhang     if (kspreason < 0) {
766b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
767b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
768b5b4017aSHong Zhang     }
769b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
77013b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
771715f1b00SHong Zhang   }
772715f1b00SHong Zhang 
773b5b4017aSHong Zhang 
77413b1b0ffSHong Zhang   /*
77513b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
776c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
77713b1b0ffSHong Zhang   */
77813b1b0ffSHong Zhang   if (ts->vecs_integral_sensip) {
77913b1b0ffSHong Zhang     if (!th->endpoint) {
78013b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
781c9ad9de0SHong Zhang       ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr);
7820e11c21fSHong Zhang       if (ts->vecs_drdp) {
783c9ad9de0SHong Zhang         ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
7840e11c21fSHong Zhang       }
78513b1b0ffSHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
786c9ad9de0SHong Zhang         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
7870e11c21fSHong Zhang         if (ts->vecs_drdp) {
78813b1b0ffSHong Zhang           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
7890e11c21fSHong Zhang         }
79013b1b0ffSHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr);
79113b1b0ffSHong Zhang       }
79213b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
79313b1b0ffSHong Zhang     } else {
794c9ad9de0SHong Zhang       ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr);
7950e11c21fSHong Zhang       if (ts->vecs_drdp) {
796c9ad9de0SHong Zhang         ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
7970e11c21fSHong Zhang       }
79813b1b0ffSHong Zhang       for (ncost=0; ncost<ts->numcost; ncost++) {
799c9ad9de0SHong Zhang         ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr);
8000e11c21fSHong Zhang         if (ts->vecs_drdp) {
80113b1b0ffSHong Zhang           ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr);
8020e11c21fSHong Zhang         }
80313b1b0ffSHong Zhang         ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr);
80413b1b0ffSHong Zhang       }
80513b1b0ffSHong Zhang     }
80613b1b0ffSHong Zhang   } else {
80713b1b0ffSHong Zhang     if (!th->endpoint) {
80813b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
809715f1b00SHong Zhang     }
810715f1b00SHong Zhang   }
8111566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8121566a47fSLisandro Dalcin }
8131566a47fSLisandro Dalcin 
8146affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
8156affb6f8SHong Zhang {
8166affb6f8SHong Zhang   TS_Theta *th = (TS_Theta*)ts->data;
8176affb6f8SHong Zhang 
8186affb6f8SHong Zhang   PetscFunctionBegin;
8196affb6f8SHong Zhang   if (ns) *ns = 1;
8206affb6f8SHong Zhang   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
8216affb6f8SHong Zhang   PetscFunctionReturn(0);
8226affb6f8SHong Zhang }
8236affb6f8SHong Zhang 
824316643e7SJed Brown /*------------------------------------------------------------*/
825277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
826316643e7SJed Brown {
827316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
828316643e7SJed Brown   PetscErrorCode ierr;
829316643e7SJed Brown 
830316643e7SJed Brown   PetscFunctionBegin;
8316bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
8326bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
8333b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
834eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
8351566a47fSLisandro Dalcin 
8361566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
8371566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
8381566a47fSLisandro Dalcin 
839b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
840277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
841277b19d0SLisandro Dalcin }
842277b19d0SLisandro Dalcin 
843ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
844ece46509SHong Zhang {
845ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
846ece46509SHong Zhang   PetscErrorCode ierr;
847ece46509SHong Zhang 
848ece46509SHong Zhang   PetscFunctionBegin;
849ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
850ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
851ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
852ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
853ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
854ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
855ece46509SHong Zhang   PetscFunctionReturn(0);
856ece46509SHong Zhang }
857ece46509SHong Zhang 
858277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
859277b19d0SLisandro Dalcin {
860277b19d0SLisandro Dalcin   PetscErrorCode ierr;
861277b19d0SLisandro Dalcin 
862277b19d0SLisandro Dalcin   PetscFunctionBegin;
863277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
864b3a6b972SJed Brown   if (ts->dm) {
865b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
866b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
867b3a6b972SJed Brown   }
868277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
869bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
870bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
871bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
872bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
873316643e7SJed Brown   PetscFunctionReturn(0);
874316643e7SJed Brown }
875316643e7SJed Brown 
876316643e7SJed Brown /*
877316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
8782b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
879316643e7SJed Brown */
8800f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
881316643e7SJed Brown {
882316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
883316643e7SJed Brown   PetscErrorCode ierr;
8847445fe48SJed Brown   Vec            X0,Xdot;
8857445fe48SJed Brown   DM             dm,dmsave;
8861566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
887316643e7SJed Brown 
888316643e7SJed Brown   PetscFunctionBegin;
8897445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
8905a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
8917445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
892b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
8937445fe48SJed Brown 
8947445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
8957445fe48SJed Brown   dmsave = ts->dm;
8967445fe48SJed Brown   ts->dm = dm;
8977445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
8987445fe48SJed Brown   ts->dm = dmsave;
8990d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
900316643e7SJed Brown   PetscFunctionReturn(0);
901316643e7SJed Brown }
902316643e7SJed Brown 
903d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
904316643e7SJed Brown {
905316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
906316643e7SJed Brown   PetscErrorCode ierr;
9077445fe48SJed Brown   Vec            Xdot;
9087445fe48SJed Brown   DM             dm,dmsave;
9091566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
910316643e7SJed Brown 
911316643e7SJed Brown   PetscFunctionBegin;
9127445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
913be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9140298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
9157445fe48SJed Brown 
9167445fe48SJed Brown   dmsave = ts->dm;
9177445fe48SJed Brown   ts->dm = dm;
918d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
9197445fe48SJed Brown   ts->dm = dmsave;
9200298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
921316643e7SJed Brown   PetscFunctionReturn(0);
922316643e7SJed Brown }
923316643e7SJed Brown 
924715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
925715f1b00SHong Zhang {
926715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
927715f1b00SHong Zhang   PetscErrorCode ierr;
928715f1b00SHong Zhang 
929715f1b00SHong Zhang   PetscFunctionBegin;
930022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
93113b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
93213b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
933715f1b00SHong Zhang   if (ts->vecs_integral_sensip) {
934715f1b00SHong Zhang     ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr);
935715f1b00SHong Zhang   }
9366f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
93713b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
93813b1b0ffSHong Zhang 
9396f92202bSHong Zhang   if (ts->vecs_integral_sensip) {
9406f92202bSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
9416f92202bSHong Zhang   }
942b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
943715f1b00SHong Zhang   PetscFunctionReturn(0);
944715f1b00SHong Zhang }
945715f1b00SHong Zhang 
9467adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
9477adebddeSHong Zhang {
9487adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9497adebddeSHong Zhang   PetscErrorCode ierr;
9507adebddeSHong Zhang 
9517adebddeSHong Zhang   PetscFunctionBegin;
9527adebddeSHong Zhang   if (ts->vecs_integral_sensip) {
9537adebddeSHong Zhang     ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr);
9547adebddeSHong Zhang     ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr);
9557adebddeSHong Zhang   }
9567adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
9577adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
9587adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
9597adebddeSHong Zhang   PetscFunctionReturn(0);
9607adebddeSHong Zhang }
9617adebddeSHong Zhang 
962316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
963316643e7SJed Brown {
964316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
9652ffb9264SLisandro Dalcin   PetscBool      match;
966316643e7SJed Brown   PetscErrorCode ierr;
967316643e7SJed Brown 
968316643e7SJed Brown   PetscFunctionBegin;
969b1cb36f3SHong Zhang   if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
970b1cb36f3SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr);
971b1cb36f3SHong Zhang   }
97239d13666SHong Zhang   if (!th->X) {
973316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
97439d13666SHong Zhang   }
97539d13666SHong Zhang   if (!th->Xdot) {
976316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
97739d13666SHong Zhang   }
97839d13666SHong Zhang   if (!th->X0) {
9793b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
98039d13666SHong Zhang   }
9811566a47fSLisandro Dalcin   if (th->endpoint) {
9821566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
9837445fe48SJed Brown   }
9843b1890cdSShri Abhyankar 
9851566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
9861566a47fSLisandro Dalcin 
9871566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
9881566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
9891566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
9901566a47fSLisandro Dalcin 
9911566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
9921566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
9932ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
9942ffb9264SLisandro Dalcin   if (!match) {
9951566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
9961566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
9973b1890cdSShri Abhyankar   }
9981566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
999b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1000b4dd3bf9SBarry Smith }
10010c7376e5SHong Zhang 
1002b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1003b4dd3bf9SBarry Smith 
1004b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1005b4dd3bf9SBarry Smith {
1006b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1007b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1008b4dd3bf9SBarry Smith 
1009b4dd3bf9SBarry Smith   PetscFunctionBegin;
1010b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1011b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
10122ca6e920SHong Zhang   if (ts->vecs_sensip) {
1013b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
10142ca6e920SHong Zhang   }
1015b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1016b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1017b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
1018b5b4017aSHong Zhang   }
1019c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1020c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
1021b5b4017aSHong Zhang   }
1022316643e7SJed Brown   PetscFunctionReturn(0);
1023316643e7SJed Brown }
1024316643e7SJed Brown 
10254416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1026316643e7SJed Brown {
1027316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1028316643e7SJed Brown   PetscErrorCode ierr;
1029316643e7SJed Brown 
1030316643e7SJed Brown   PetscFunctionBegin;
1031e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1032316643e7SJed Brown   {
10330298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
10340298fd71SBarry 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);
103503842d09SLisandro 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);
1036316643e7SJed Brown   }
1037316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1038316643e7SJed Brown   PetscFunctionReturn(0);
1039316643e7SJed Brown }
1040316643e7SJed Brown 
1041316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1042316643e7SJed Brown {
1043316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1044ace3abfcSBarry Smith   PetscBool      iascii;
1045316643e7SJed Brown   PetscErrorCode ierr;
1046316643e7SJed Brown 
1047316643e7SJed Brown   PetscFunctionBegin;
1048251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1049316643e7SJed Brown   if (iascii) {
10507c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1051ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1052316643e7SJed Brown   }
1053316643e7SJed Brown   PetscFunctionReturn(0);
1054316643e7SJed Brown }
1055316643e7SJed Brown 
1056560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
10570de4c49aSJed Brown {
10580de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
10590de4c49aSJed Brown 
10600de4c49aSJed Brown   PetscFunctionBegin;
10610de4c49aSJed Brown   *theta = th->Theta;
10620de4c49aSJed Brown   PetscFunctionReturn(0);
10630de4c49aSJed Brown }
10640de4c49aSJed Brown 
1065560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
10660de4c49aSJed Brown {
10670de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
10680de4c49aSJed Brown 
10690de4c49aSJed Brown   PetscFunctionBegin;
10707c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
10710de4c49aSJed Brown   th->Theta = theta;
10721566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
10730de4c49aSJed Brown   PetscFunctionReturn(0);
10740de4c49aSJed Brown }
1075eb284becSJed Brown 
1076560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
107726f2ff8fSLisandro Dalcin {
107826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
107926f2ff8fSLisandro Dalcin 
108026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
108126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
108226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
108326f2ff8fSLisandro Dalcin }
108426f2ff8fSLisandro Dalcin 
1085560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1086eb284becSJed Brown {
1087eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1088eb284becSJed Brown 
1089eb284becSJed Brown   PetscFunctionBegin;
1090eb284becSJed Brown   th->endpoint = flg;
1091eb284becSJed Brown   PetscFunctionReturn(0);
1092eb284becSJed Brown }
10930de4c49aSJed Brown 
1094f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1095f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1096f9c1d6abSBarry Smith {
1097f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1098f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
10993fd8ae06SJed Brown   const PetscReal one = 1.0;
1100f9c1d6abSBarry Smith 
1101f9c1d6abSBarry Smith   PetscFunctionBegin;
11023fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1103f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1104f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1105f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1106f9c1d6abSBarry Smith }
1107f9c1d6abSBarry Smith #endif
1108f9c1d6abSBarry Smith 
110942682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
111042682096SHong Zhang {
111142682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
111242682096SHong Zhang 
111342682096SHong Zhang   PetscFunctionBegin;
11141566a47fSLisandro Dalcin   if (ns) *ns = 1;
11151566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
111642682096SHong Zhang   PetscFunctionReturn(0);
111742682096SHong Zhang }
1118f9c1d6abSBarry Smith 
1119316643e7SJed Brown /* ------------------------------------------------------------ */
1120316643e7SJed Brown /*MC
112196f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1122316643e7SJed Brown 
1123316643e7SJed Brown    Level: beginner
1124316643e7SJed Brown 
11254eb428fdSBarry Smith    Options Database:
1126c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1127c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
112803842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11294eb428fdSBarry Smith 
1130eb284becSJed Brown    Notes:
11310c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
11320c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
11334eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
11344eb428fdSBarry Smith 
1135eb284becSJed Brown    This method can be applied to DAE.
1136eb284becSJed Brown 
1137eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
1138eb284becSJed Brown 
1139eb284becSJed Brown .vb
1140eb284becSJed Brown   Theta | Theta
1141eb284becSJed Brown   -------------
1142eb284becSJed Brown         |  1
1143eb284becSJed Brown .ve
1144eb284becSJed Brown 
1145eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1146eb284becSJed Brown 
1147eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1148eb284becSJed Brown 
1149eb284becSJed Brown .vb
1150eb284becSJed Brown   0 | 0         0
1151eb284becSJed Brown   1 | 1-Theta   Theta
1152eb284becSJed Brown   -------------------
1153eb284becSJed Brown     | 1-Theta   Theta
1154eb284becSJed Brown .ve
1155eb284becSJed Brown 
1156eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1157eb284becSJed Brown 
1158eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1159eb284becSJed Brown 
1160eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1161eb284becSJed Brown 
11624eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1163eb284becSJed Brown 
1164eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1165316643e7SJed Brown 
1166316643e7SJed Brown M*/
11678cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1168316643e7SJed Brown {
1169316643e7SJed Brown   TS_Theta       *th;
1170316643e7SJed Brown   PetscErrorCode ierr;
1171316643e7SJed Brown 
1172316643e7SJed Brown   PetscFunctionBegin;
1173277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1174ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1175316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1176316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1177316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
117842f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1179ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1180316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1181cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
11821566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
118324655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1184316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
11850f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
11860f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1187f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1188f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1189f9c1d6abSBarry Smith #endif
119042682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
119142f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1192b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1193b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
11942ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
11956affb6f8SHong Zhang 
1196715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
11977adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1198715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
11996affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1200316643e7SJed Brown 
1201efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1202efd4aadfSBarry Smith 
1203b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1204316643e7SJed Brown   ts->data = (void*)th;
1205316643e7SJed Brown 
1206715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1207715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1208715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1209b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1210715f1b00SHong Zhang 
12116f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1212316643e7SJed Brown   th->Theta       = 0.5;
12131566a47fSLisandro Dalcin   th->order       = 2;
1214bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1215bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1216bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1217bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1218316643e7SJed Brown   PetscFunctionReturn(0);
1219316643e7SJed Brown }
12200de4c49aSJed Brown 
12210de4c49aSJed Brown /*@
12220de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12230de4c49aSJed Brown 
12240de4c49aSJed Brown   Not Collective
12250de4c49aSJed Brown 
12260de4c49aSJed Brown   Input Parameter:
12270de4c49aSJed Brown .  ts - timestepping context
12280de4c49aSJed Brown 
12290de4c49aSJed Brown   Output Parameter:
12300de4c49aSJed Brown .  theta - stage abscissa
12310de4c49aSJed Brown 
12320de4c49aSJed Brown   Note:
12330de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
12340de4c49aSJed Brown 
12350de4c49aSJed Brown   Level: Advanced
12360de4c49aSJed Brown 
12370de4c49aSJed Brown .seealso: TSThetaSetTheta()
12380de4c49aSJed Brown @*/
12397087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
12400de4c49aSJed Brown {
12414ac538c5SBarry Smith   PetscErrorCode ierr;
12420de4c49aSJed Brown 
12430de4c49aSJed Brown   PetscFunctionBegin;
1244afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12450de4c49aSJed Brown   PetscValidPointer(theta,2);
12464ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
12470de4c49aSJed Brown   PetscFunctionReturn(0);
12480de4c49aSJed Brown }
12490de4c49aSJed Brown 
12500de4c49aSJed Brown /*@
12510de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
12520de4c49aSJed Brown 
12530de4c49aSJed Brown   Not Collective
12540de4c49aSJed Brown 
12550de4c49aSJed Brown   Input Parameter:
12560de4c49aSJed Brown +  ts - timestepping context
12570de4c49aSJed Brown -  theta - stage abscissa
12580de4c49aSJed Brown 
12590de4c49aSJed Brown   Options Database:
12600de4c49aSJed Brown .  -ts_theta_theta <theta>
12610de4c49aSJed Brown 
12620de4c49aSJed Brown   Level: Intermediate
12630de4c49aSJed Brown 
12640de4c49aSJed Brown .seealso: TSThetaGetTheta()
12650de4c49aSJed Brown @*/
12667087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
12670de4c49aSJed Brown {
12684ac538c5SBarry Smith   PetscErrorCode ierr;
12690de4c49aSJed Brown 
12700de4c49aSJed Brown   PetscFunctionBegin;
1271afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12724ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
12730de4c49aSJed Brown   PetscFunctionReturn(0);
12740de4c49aSJed Brown }
1275f33bbcb6SJed Brown 
127626f2ff8fSLisandro Dalcin /*@
127726f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
127826f2ff8fSLisandro Dalcin 
127926f2ff8fSLisandro Dalcin   Not Collective
128026f2ff8fSLisandro Dalcin 
128126f2ff8fSLisandro Dalcin   Input Parameter:
128226f2ff8fSLisandro Dalcin .  ts - timestepping context
128326f2ff8fSLisandro Dalcin 
128426f2ff8fSLisandro Dalcin   Output Parameter:
128526f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
128626f2ff8fSLisandro Dalcin 
128726f2ff8fSLisandro Dalcin   Level: Advanced
128826f2ff8fSLisandro Dalcin 
128926f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
129026f2ff8fSLisandro Dalcin @*/
129126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
129226f2ff8fSLisandro Dalcin {
129326f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
129426f2ff8fSLisandro Dalcin 
129526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
129626f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
129726f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1298163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
129926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
130026f2ff8fSLisandro Dalcin }
130126f2ff8fSLisandro Dalcin 
1302eb284becSJed Brown /*@
1303eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1304eb284becSJed Brown 
1305eb284becSJed Brown   Not Collective
1306eb284becSJed Brown 
1307eb284becSJed Brown   Input Parameter:
1308eb284becSJed Brown +  ts - timestepping context
1309eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1310eb284becSJed Brown 
1311eb284becSJed Brown   Options Database:
1312eb284becSJed Brown .  -ts_theta_endpoint <flg>
1313eb284becSJed Brown 
1314eb284becSJed Brown   Level: Intermediate
1315eb284becSJed Brown 
1316eb284becSJed Brown .seealso: TSTHETA, TSCN
1317eb284becSJed Brown @*/
1318eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1319eb284becSJed Brown {
1320eb284becSJed Brown   PetscErrorCode ierr;
1321eb284becSJed Brown 
1322eb284becSJed Brown   PetscFunctionBegin;
1323eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1324eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1325eb284becSJed Brown   PetscFunctionReturn(0);
1326eb284becSJed Brown }
1327eb284becSJed Brown 
1328f33bbcb6SJed Brown /*
1329f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1330f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1331f33bbcb6SJed Brown  */
1332f33bbcb6SJed Brown 
13331566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
13341566a47fSLisandro Dalcin {
13351566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
13361566a47fSLisandro Dalcin   PetscErrorCode ierr;
13371566a47fSLisandro Dalcin 
13381566a47fSLisandro Dalcin   PetscFunctionBegin;
13391566a47fSLisandro 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");
13401566a47fSLisandro 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");
13411566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
13421566a47fSLisandro Dalcin   PetscFunctionReturn(0);
13431566a47fSLisandro Dalcin }
13441566a47fSLisandro Dalcin 
1345f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1346f33bbcb6SJed Brown {
1347f33bbcb6SJed Brown   PetscFunctionBegin;
1348f33bbcb6SJed Brown   PetscFunctionReturn(0);
1349f33bbcb6SJed Brown }
1350f33bbcb6SJed Brown 
1351f33bbcb6SJed Brown /*MC
1352f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1353f33bbcb6SJed Brown 
1354f33bbcb6SJed Brown   Level: beginner
1355f33bbcb6SJed Brown 
13564eb428fdSBarry Smith   Notes:
1357c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
13584eb428fdSBarry Smith 
13591566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
13604eb428fdSBarry Smith 
1361f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1362f33bbcb6SJed Brown 
1363f33bbcb6SJed Brown M*/
13648cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1365f33bbcb6SJed Brown {
1366f33bbcb6SJed Brown   PetscErrorCode ierr;
1367f33bbcb6SJed Brown 
1368f33bbcb6SJed Brown   PetscFunctionBegin;
1369f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1370f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
13711566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
13720c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1373f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1374f33bbcb6SJed Brown   PetscFunctionReturn(0);
1375f33bbcb6SJed Brown }
1376f33bbcb6SJed Brown 
13771566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
13781566a47fSLisandro Dalcin {
13791566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
13801566a47fSLisandro Dalcin   PetscErrorCode ierr;
13811566a47fSLisandro Dalcin 
13821566a47fSLisandro Dalcin   PetscFunctionBegin;
13831566a47fSLisandro 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");
13841566a47fSLisandro 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");
13851566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
13861566a47fSLisandro Dalcin   PetscFunctionReturn(0);
13871566a47fSLisandro Dalcin }
13881566a47fSLisandro Dalcin 
1389f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1390f33bbcb6SJed Brown {
1391f33bbcb6SJed Brown   PetscFunctionBegin;
1392f33bbcb6SJed Brown   PetscFunctionReturn(0);
1393f33bbcb6SJed Brown }
1394f33bbcb6SJed Brown 
1395f33bbcb6SJed Brown /*MC
1396f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1397f33bbcb6SJed Brown 
1398f33bbcb6SJed Brown   Level: beginner
1399f33bbcb6SJed Brown 
1400f33bbcb6SJed Brown   Notes:
14017cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14027cf5af47SJed Brown 
14037cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1404f33bbcb6SJed Brown 
1405f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1406f33bbcb6SJed Brown 
1407f33bbcb6SJed Brown M*/
14088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1409f33bbcb6SJed Brown {
1410f33bbcb6SJed Brown   PetscErrorCode ierr;
1411f33bbcb6SJed Brown 
1412f33bbcb6SJed Brown   PetscFunctionBegin;
1413f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1414f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1415eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
14160c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1417f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1418f33bbcb6SJed Brown   PetscFunctionReturn(0);
1419f33bbcb6SJed Brown }
1420