xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 1a352fa80c5fd62651f9d25423fa30f3a9234bec)
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;
15*1a352fa8SHong Zhang   PetscReal    shift;                    /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
161566a47fSLisandro Dalcin   PetscInt     order;
171566a47fSLisandro Dalcin   PetscBool    endpoint;
181566a47fSLisandro Dalcin   PetscBool    extrapolate;
19715f1b00SHong Zhang   TSStepStatus status;
20*1a352fa8SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events, used by cost integral */
21*1a352fa8SHong Zhang   PetscReal    ptime0;                   /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
22*1a352fa8SHong Zhang   PetscReal    time_step0;               /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/
231566a47fSLisandro Dalcin 
24715f1b00SHong Zhang   /* context for sensitivity analysis */
25022c081aSHong Zhang   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
26b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
27b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
2813b1b0ffSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
2913b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
30b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
3113b1b0ffSHong Zhang   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
327207e0fdSHong Zhang   Mat          MatIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
337207e0fdSHong Zhang   Mat          MatIntegralSensip0;       /* backup for roll-backs due to events */
34b5b4017aSHong Zhang   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
35b5b4017aSHong Zhang   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
36b5b4017aSHong Zhang   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
37b5b4017aSHong Zhang   Vec          *VecsAffine;              /* Working vectors to store residuals */
38715f1b00SHong Zhang   /* context for error estimation */
391566a47fSLisandro Dalcin   Vec          vec_sol_prev;
401566a47fSLisandro Dalcin   Vec          vec_lte_work;
41316643e7SJed Brown } TS_Theta;
42316643e7SJed Brown 
437445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
447445fe48SJed Brown {
457445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
467445fe48SJed Brown   PetscErrorCode ierr;
477445fe48SJed Brown 
487445fe48SJed Brown   PetscFunctionBegin;
497445fe48SJed Brown   if (X0) {
507445fe48SJed Brown     if (dm && dm != ts->dm) {
510d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
527445fe48SJed Brown     } else *X0 = ts->vec_sol;
537445fe48SJed Brown   }
547445fe48SJed Brown   if (Xdot) {
557445fe48SJed Brown     if (dm && dm != ts->dm) {
560d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
577445fe48SJed Brown     } else *Xdot = th->Xdot;
587445fe48SJed Brown   }
597445fe48SJed Brown   PetscFunctionReturn(0);
607445fe48SJed Brown }
617445fe48SJed Brown 
620d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
630d0b770aSPeter Brune {
640d0b770aSPeter Brune   PetscErrorCode ierr;
650d0b770aSPeter Brune 
660d0b770aSPeter Brune   PetscFunctionBegin;
670d0b770aSPeter Brune   if (X0) {
680d0b770aSPeter Brune     if (dm && dm != ts->dm) {
690d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
700d0b770aSPeter Brune     }
710d0b770aSPeter Brune   }
720d0b770aSPeter Brune   if (Xdot) {
730d0b770aSPeter Brune     if (dm && dm != ts->dm) {
740d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
750d0b770aSPeter Brune     }
760d0b770aSPeter Brune   }
770d0b770aSPeter Brune   PetscFunctionReturn(0);
780d0b770aSPeter Brune }
790d0b770aSPeter Brune 
807445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
817445fe48SJed Brown {
827445fe48SJed Brown   PetscFunctionBegin;
837445fe48SJed Brown   PetscFunctionReturn(0);
847445fe48SJed Brown }
857445fe48SJed Brown 
867445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
877445fe48SJed Brown {
887445fe48SJed Brown   TS             ts = (TS)ctx;
897445fe48SJed Brown   PetscErrorCode ierr;
907445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
917445fe48SJed Brown 
927445fe48SJed Brown   PetscFunctionBegin;
937445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
977445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
987445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
990d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
1000d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1017445fe48SJed Brown   PetscFunctionReturn(0);
1027445fe48SJed Brown }
1037445fe48SJed Brown 
104258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
105258e1594SPeter Brune {
106258e1594SPeter Brune   PetscFunctionBegin;
107258e1594SPeter Brune   PetscFunctionReturn(0);
108258e1594SPeter Brune }
109258e1594SPeter Brune 
110258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
111258e1594SPeter Brune {
112258e1594SPeter Brune   TS             ts = (TS)ctx;
113258e1594SPeter Brune   PetscErrorCode ierr;
114258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
115258e1594SPeter Brune 
116258e1594SPeter Brune   PetscFunctionBegin;
117258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
118258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
119258e1594SPeter Brune 
120258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122258e1594SPeter Brune 
123258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125258e1594SPeter Brune 
126258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
127258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
128258e1594SPeter Brune   PetscFunctionReturn(0);
129258e1594SPeter Brune }
130258e1594SPeter Brune 
131022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
132022c081aSHong Zhang {
133022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
134cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
135022c081aSHong Zhang   PetscErrorCode ierr;
136022c081aSHong Zhang 
137022c081aSHong Zhang   PetscFunctionBegin;
138022c081aSHong Zhang   if (th->endpoint) {
139022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
140022c081aSHong Zhang     if (th->Theta!=1.0) {
141*1a352fa8SHong Zhang       ierr = TSComputeRHSFunction(quadts,th->ptime0,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
142*1a352fa8SHong Zhang       ierr = VecAXPY(quadts->vec_sol,th->time_step0*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
143022c081aSHong Zhang     }
144cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
145*1a352fa8SHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step0*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
146022c081aSHong Zhang   } else {
147cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
148*1a352fa8SHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step0,ts->vec_costintegrand);CHKERRQ(ierr);
149022c081aSHong Zhang   }
150022c081aSHong Zhang   PetscFunctionReturn(0);
151022c081aSHong Zhang }
152022c081aSHong Zhang 
153b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
154b1cb36f3SHong Zhang {
155b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
156cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
157b1cb36f3SHong Zhang   PetscErrorCode ierr;
158b1cb36f3SHong Zhang 
159b1cb36f3SHong Zhang   PetscFunctionBegin;
160b1cb36f3SHong Zhang   /* backup cost integral */
161cd4cee2dSHong Zhang   ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr);
162022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
163b1cb36f3SHong Zhang   PetscFunctionReturn(0);
164b1cb36f3SHong Zhang }
165b1cb36f3SHong Zhang 
166b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
167b1cb36f3SHong Zhang {
168*1a352fa8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
169b1cb36f3SHong Zhang   PetscErrorCode ierr;
170b1cb36f3SHong Zhang 
171b1cb36f3SHong Zhang   PetscFunctionBegin;
172*1a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
173*1a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
174*1a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
175022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
17624655328SShri   PetscFunctionReturn(0);
17724655328SShri }
17824655328SShri 
179470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1801566a47fSLisandro Dalcin {
1811566a47fSLisandro Dalcin   PetscInt       nits,lits;
1821566a47fSLisandro Dalcin   PetscErrorCode ierr;
1831566a47fSLisandro Dalcin 
1841566a47fSLisandro Dalcin   PetscFunctionBegin;
1851566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1861566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1871566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1881566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1891566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1901566a47fSLisandro Dalcin }
1911566a47fSLisandro Dalcin 
192193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
193316643e7SJed Brown {
194316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1951566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1964957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1971566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
198051f2191SLisandro Dalcin   PetscErrorCode ierr;
199316643e7SJed Brown 
200316643e7SJed Brown   PetscFunctionBegin;
2011566a47fSLisandro Dalcin   if (!ts->steprollback) {
2022ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
2033b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
2041566a47fSLisandro Dalcin   }
2051566a47fSLisandro Dalcin 
2061566a47fSLisandro Dalcin   th->status     = TS_STEP_INCOMPLETE;
2073e05ceb1SHong Zhang   th->shift      = 1/(th->Theta*ts->time_step);
2081566a47fSLisandro Dalcin   th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
209*1a352fa8SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
210be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
211fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
2123e05ceb1SHong Zhang       ierr = VecAXPY(th->X,1/th->shift,th->Xdot);CHKERRQ(ierr);
213be5899b3SLisandro Dalcin     }
214eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
215eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2161566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2171566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2181566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2191566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2201566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
221eb284becSJed Brown     }
222be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
223470880abSPatrick Sanan     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
224be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2251566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
226fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
227051f2191SLisandro Dalcin 
228051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2291566a47fSLisandro Dalcin     if (th->endpoint) {
2301566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2311566a47fSLisandro Dalcin     } else {
2323e05ceb1SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr);
2331566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2341566a47fSLisandro Dalcin     }
2351566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2361566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2371566a47fSLisandro Dalcin     if (!accept) {
2381566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
239be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
240be5899b3SLisandro Dalcin       goto reject_step;
241051f2191SLisandro Dalcin     }
2423b1890cdSShri Abhyankar 
243715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
244*1a352fa8SHong Zhang       th->ptime0     = ts->ptime;
245*1a352fa8SHong Zhang       th->time_step0 = ts->time_step;
24617145e75SHong Zhang     }
2472b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
248cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
249051f2191SLisandro Dalcin     break;
250051f2191SLisandro Dalcin 
251051f2191SLisandro Dalcin   reject_step:
252fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2531566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
254051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2551566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2563b1890cdSShri Abhyankar     }
2573b1890cdSShri Abhyankar   }
258316643e7SJed Brown   PetscFunctionReturn(0);
259316643e7SJed Brown }
260316643e7SJed Brown 
261c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
262c9681e5dSHong Zhang {
263c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
264cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
265c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
266c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
267c9681e5dSHong Zhang   PetscInt       nadj;
2687207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
269c9681e5dSHong Zhang   KSP            ksp;
270c9681e5dSHong Zhang   PetscScalar    *xarr;
271077c4feaSHong Zhang   TSEquationType eqtype;
272077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
273*1a352fa8SHong Zhang   PetscReal      adjoint_time_step;
274c9681e5dSHong Zhang   PetscErrorCode ierr;
275c9681e5dSHong Zhang 
276c9681e5dSHong Zhang   PetscFunctionBegin;
277077c4feaSHong Zhang   ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr);
278077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
279077c4feaSHong Zhang     isexplicitode  = PETSC_TRUE;
280077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
281077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
282077c4feaSHong Zhang   }
283c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
284c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
285cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
2867207e0fdSHong Zhang   if (quadts) {
287cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
288cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
2897207e0fdSHong Zhang   }
290c9681e5dSHong Zhang 
291*1a352fa8SHong Zhang   th->stage_time    = ts->ptime;
292*1a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
293c9681e5dSHong Zhang 
29433f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2957207e0fdSHong Zhang   if (quadts) {
296cd4cee2dSHong Zhang     ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
2977207e0fdSHong Zhang   }
298cd4cee2dSHong Zhang 
299c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
300c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
301*1a352fa8SHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./adjoint_time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
302cd4cee2dSHong Zhang     if (quadJ) {
303cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
304cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
305cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
306cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
307cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
308c9681e5dSHong Zhang     }
309c9681e5dSHong Zhang   }
310c9681e5dSHong Zhang 
311c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
31204f7482eSHong Zhang   ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr);
31304f7482eSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
314c9681e5dSHong Zhang 
315c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
316c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
317c9681e5dSHong Zhang     KSPConvergedReason kspreason;
318c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
319c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
320c9681e5dSHong Zhang     if (kspreason < 0) {
321c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32229b3c8a1SHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
323c9681e5dSHong Zhang     }
324c9681e5dSHong Zhang   }
325c9681e5dSHong Zhang 
326c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
327c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
328c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
329c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
330c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
331b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
332c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
333b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
334c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
335c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
336*1a352fa8SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],1./adjoint_time_step);CHKERRQ(ierr);
33733f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
338c9681e5dSHong Zhang       if (ts->vecs_fup) {
33933f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
340c9681e5dSHong Zhang       }
341c9681e5dSHong Zhang     }
342c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
343c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
34405c9950eSHong Zhang       KSPConvergedReason kspreason;
345c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
34605c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
34705c9950eSHong Zhang       if (kspreason < 0) {
34805c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
34929b3c8a1SHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
35005c9950eSHong Zhang       }
351c9681e5dSHong Zhang     }
352c9681e5dSHong Zhang   }
353c9681e5dSHong Zhang 
354c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
355077c4feaSHong Zhang   if (!isexplicitode) {
3563e05ceb1SHong Zhang     th->shift = 0.0;
35704f7482eSHong Zhang     ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr);
35804f7482eSHong Zhang     ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
359c9681e5dSHong Zhang     ierr = MatScale(J,-1.);CHKERRQ(ierr);
360c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
36133f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
362c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
363*1a352fa8SHong Zhang       ierr = VecScale(VecsSensiTemp[nadj],adjoint_time_step);CHKERRQ(ierr);
364c9681e5dSHong Zhang       ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
365c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
366c9681e5dSHong Zhang         ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
367*1a352fa8SHong Zhang         ierr = VecScale(VecsSensi2Temp[nadj],adjoint_time_step);CHKERRQ(ierr);
368c9681e5dSHong Zhang         ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
369c9681e5dSHong Zhang       }
370c9681e5dSHong Zhang     }
371077c4feaSHong Zhang   }
372c9681e5dSHong Zhang   if (ts->vecs_sensip) {
373*1a352fa8SHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,1./adjoint_time_step,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
3747207e0fdSHong Zhang     if (quadts) {
375cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
3767207e0fdSHong Zhang     }
377c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
378c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
379b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
38033f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
381b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
382c9681e5dSHong Zhang     }
383cd4cee2dSHong Zhang 
384c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
385c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
386*1a352fa8SHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
387cd4cee2dSHong Zhang       if (quadJp) {
388cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
389cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
390*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr);
391cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
392cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
39333f72c5dSHong Zhang       }
394c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
395c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
396*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
397c9681e5dSHong Zhang         if (ts->vecs_fpu) {
398*1a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
399c9681e5dSHong Zhang         }
400c9681e5dSHong Zhang         if (ts->vecs_fpp) {
401*1a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
402c9681e5dSHong Zhang         }
403c9681e5dSHong Zhang       }
404c9681e5dSHong Zhang     }
405c9681e5dSHong Zhang   }
406c9681e5dSHong Zhang 
407c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
408c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
409c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
410c9681e5dSHong Zhang   }
411c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
412c9681e5dSHong Zhang   PetscFunctionReturn(0);
413c9681e5dSHong Zhang }
414c9681e5dSHong Zhang 
41542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4162ca6e920SHong Zhang {
4172ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
418cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
419b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
420b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4212ca6e920SHong Zhang   PetscInt       nadj;
4227207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4232ca6e920SHong Zhang   KSP            ksp;
424b5b4017aSHong Zhang   PetscScalar    *xarr;
425*1a352fa8SHong Zhang   PetscReal      adjoint_time_step;
426*1a352fa8SHong Zhang   PetscReal      adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, ususally ts->ptime is larger than adjoint_ptime) */
427b5b4017aSHong Zhang   PetscErrorCode ierr;
4282ca6e920SHong Zhang 
4292ca6e920SHong Zhang   PetscFunctionBegin;
430c9681e5dSHong Zhang   if (th->Theta == 1.) {
431c9681e5dSHong Zhang     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
432c9681e5dSHong Zhang     PetscFunctionReturn(0);
433c9681e5dSHong Zhang   }
4342ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
435302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
436cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
4377207e0fdSHong Zhang   if (quadts) {
438cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
439cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
4407207e0fdSHong Zhang   }
441b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
442*1a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step);
443*1a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
444*1a352fa8SHong Zhang   adjoint_time_step = -ts->time_step;  /* always positive since time_step is negative */
4452ca6e920SHong Zhang 
446b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
44733f72c5dSHong Zhang   /* Cost function has an integral term */
4487207e0fdSHong Zhang   if (quadts) {
44905755b9cSHong Zhang     if (th->endpoint) {
450cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
45105755b9cSHong Zhang     } else {
452cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
45305755b9cSHong Zhang     }
4547207e0fdSHong Zhang   }
45533f72c5dSHong Zhang 
456abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
457b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
458*1a352fa8SHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*adjoint_time_step));CHKERRQ(ierr);
459cd4cee2dSHong Zhang     if (quadJ) {
460cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
461cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
462cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
463cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
464cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
46536eaed60SHong Zhang     }
4662ca6e920SHong Zhang   }
4673c54f92cSHong Zhang 
468b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
469*1a352fa8SHong Zhang   th->shift = 1./(th->Theta*adjoint_time_step);
4703c54f92cSHong Zhang   if (th->endpoint) {
47104f7482eSHong Zhang     ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr);
4723c54f92cSHong Zhang   } else {
47304f7482eSHong Zhang     ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr);
4743c54f92cSHong Zhang   }
475cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
4762ca6e920SHong Zhang 
477b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
478abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4791d14d03bSHong Zhang     KSPConvergedReason kspreason;
480b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
4811d14d03bSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
4821d14d03bSHong Zhang     if (kspreason < 0) {
4831d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
48429b3c8a1SHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
4851d14d03bSHong Zhang     }
4862ca6e920SHong Zhang   }
4873c54f92cSHong Zhang 
4881d14d03bSHong Zhang   /* Second-order adjoint */
489b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
490b5b4017aSHong Zhang     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
491b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
492b5b4017aSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
493b5b4017aSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
494b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
495b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
496b5b4017aSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
497b5b4017aSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
498b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
499b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
500b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
501b5b4017aSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
5023e05ceb1SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],th->shift);CHKERRQ(ierr);
50333f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
504b5b4017aSHong Zhang       if (ts->vecs_fup) {
50533f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
506b5b4017aSHong Zhang       }
507b5b4017aSHong Zhang     }
508b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
509b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
5101d14d03bSHong Zhang       KSPConvergedReason kspreason;
5111d14d03bSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
5121d14d03bSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
5131d14d03bSHong Zhang       if (kspreason < 0) {
5141d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
51529b3c8a1SHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
5161d14d03bSHong Zhang       }
517b5b4017aSHong Zhang     }
518b5b4017aSHong Zhang   }
519b5b4017aSHong Zhang 
52036eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5211d14d03bSHong Zhang   if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
522*1a352fa8SHong Zhang     th->shift      = 1./((th->Theta-1.)*adjoint_time_step);
523*1a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
52404f7482eSHong Zhang     ierr           = TSComputeSNESJacobian(ts,th->X0,J,Jpre);CHKERRQ(ierr);
52504f7482eSHong Zhang     ierr           = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
52633f72c5dSHong Zhang     /* R_U at t_n */
5277207e0fdSHong Zhang     if (quadts) {
528*1a352fa8SHong Zhang       ierr = TSComputeRHSJacobian(quadts,adjoint_ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
5297207e0fdSHong Zhang     }
530abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
531b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
532cd4cee2dSHong Zhang       if (quadJ) {
533cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
534cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
535cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr);
536cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
537cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
538b5b4017aSHong Zhang       }
5393e05ceb1SHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./th->shift);CHKERRQ(ierr);
540b5b4017aSHong Zhang     }
5411d14d03bSHong Zhang 
5421d14d03bSHong Zhang     /* Second-order adjoint */
5431d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
544b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
545b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
546b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
547b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
548*1a352fa8SHong Zhang       ierr = TSComputeIHessianProductFunctionUU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
549b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
55033f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
551b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
552*1a352fa8SHong Zhang       ierr = TSComputeIHessianProductFunctionUP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
553b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
55433f72c5dSHong 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  */
555b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
55633f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
557b5b4017aSHong Zhang         if (ts->vecs_fup) {
55833f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
559b5b4017aSHong Zhang         }
5603e05ceb1SHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./th->shift);CHKERRQ(ierr);
561b5b4017aSHong Zhang       }
562b5e0532dSHong Zhang     }
5633fd52205SHong Zhang 
564c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
565c586f2e8SHong Zhang 
5663fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
567b5b4017aSHong Zhang       /* U_{n+1} */
568*1a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,-1./(th->Theta*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5697207e0fdSHong Zhang       if (quadts) {
570cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
5717207e0fdSHong Zhang       }
572abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
573b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
574*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5753fd52205SHong Zhang       }
57633f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
57733f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
57833f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
57933f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
580b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
581b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
58233f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
58333f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
58433f72c5dSHong Zhang 
585b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
586b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
587b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
58833f72c5dSHong 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)  */
589b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
590*1a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
591b5b4017aSHong Zhang           if (ts->vecs_fpu) {
592*1a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
593b5b4017aSHong Zhang           }
594b5b4017aSHong Zhang           if (ts->vecs_fpp) {
595*1a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
596b5b4017aSHong Zhang           }
597b5b4017aSHong Zhang         }
598b5b4017aSHong Zhang       }
599b5b4017aSHong Zhang 
600b5b4017aSHong Zhang       /* U_s */
601*1a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,adjoint_ptime,th->X0,th->Xdot,1./((th->Theta-1.0)*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6027207e0fdSHong Zhang       if (quadts) {
603*1a352fa8SHong Zhang         ierr = TSComputeRHSJacobianP(quadts,adjoint_ptime,th->X0,quadJp);CHKERRQ(ierr);
6047207e0fdSHong Zhang       }
605abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
606b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
607*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
60833f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
60933f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
61033f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
61133f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
612b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
613*1a352fa8SHong Zhang           ierr = TSComputeIHessianProductFunctionPU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
61433f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
61533f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
616b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
617*1a352fa8SHong Zhang           ierr = TSComputeIHessianProductFunctionPP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
618b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
61933f72c5dSHong 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) */
620b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
621*1a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
622b5b4017aSHong Zhang             if (ts->vecs_fpu) {
623*1a352fa8SHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
624b5e0532dSHong Zhang             }
625b5b4017aSHong Zhang             if (ts->vecs_fpp) {
626*1a352fa8SHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
627b5b4017aSHong Zhang             }
62836eaed60SHong Zhang           }
62936eaed60SHong Zhang         }
6303fd52205SHong Zhang       }
631b5e0532dSHong Zhang     }
6323fd52205SHong Zhang   } else { /* one-stage case */
6333e05ceb1SHong Zhang     th->shift = 0.0;
63404f7482eSHong Zhang     ierr      = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr); /* get -f_y */
63504f7482eSHong Zhang     ierr      = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
6367207e0fdSHong Zhang     if (quadts) {
637cd4cee2dSHong Zhang       ierr  = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
6387207e0fdSHong Zhang     }
639abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
640b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
641*1a352fa8SHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-adjoint_time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
642cd4cee2dSHong Zhang       if (quadJ) {
643cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
644cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
645*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],adjoint_time_step,ts->vec_drdu_col);CHKERRQ(ierr);
646cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
647cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
64836eaed60SHong Zhang       }
6492ca6e920SHong Zhang     }
6503fd52205SHong Zhang     if (ts->vecs_sensip) {
6513e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6527207e0fdSHong Zhang       if (quadts) {
653cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
6547207e0fdSHong Zhang       }
655abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
65633f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
657*1a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
658cd4cee2dSHong Zhang         if (quadJp) {
659cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
660cd4cee2dSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
661*1a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr);
662cd4cee2dSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
663cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
66436eaed60SHong Zhang         }
66536eaed60SHong Zhang       }
6663fd52205SHong Zhang     }
6672ca6e920SHong Zhang   }
6682ca6e920SHong Zhang 
6692ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6702ca6e920SHong Zhang   PetscFunctionReturn(0);
6712ca6e920SHong Zhang }
6722ca6e920SHong Zhang 
673cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
674cd652676SJed Brown {
675cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
676be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
677cd652676SJed Brown   PetscErrorCode ierr;
678cd652676SJed Brown 
679cd652676SJed Brown   PetscFunctionBegin;
680a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
681be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
682be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
683cd652676SJed Brown   PetscFunctionReturn(0);
684cd652676SJed Brown }
685cd652676SJed Brown 
6861566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
6871566a47fSLisandro Dalcin {
6881566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
6891566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
6901566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
6917453f775SEmil Constantinescu   PetscReal      wltea,wlter;
6921566a47fSLisandro Dalcin   PetscErrorCode ierr;
6931566a47fSLisandro Dalcin 
6941566a47fSLisandro Dalcin   PetscFunctionBegin;
6952ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
6961566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
697fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
6981566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
699fecfb714SLisandro Dalcin   {
700be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
701be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
7021566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
7031566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
7041566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
7051566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
7061566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
7077453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
7081566a47fSLisandro Dalcin   }
7091566a47fSLisandro Dalcin   if (order) *order = 2;
7101566a47fSLisandro Dalcin   PetscFunctionReturn(0);
7111566a47fSLisandro Dalcin }
7121566a47fSLisandro Dalcin 
7131566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7141566a47fSLisandro Dalcin {
7151566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7167207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7171566a47fSLisandro Dalcin   PetscErrorCode ierr;
7181566a47fSLisandro Dalcin 
7191566a47fSLisandro Dalcin   PetscFunctionBegin;
7201566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
721cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
722cd4cee2dSHong Zhang     ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr);
7231566a47fSLisandro Dalcin   }
724715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
72513b1b0ffSHong Zhang   if (ts->mat_sensip) {
72613b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7276f92202bSHong Zhang   }
7287207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7297207e0fdSHong Zhang     ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7306f92202bSHong Zhang   }
731715f1b00SHong Zhang   PetscFunctionReturn(0);
732715f1b00SHong Zhang }
733715f1b00SHong Zhang 
734715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
735715f1b00SHong Zhang {
736715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
737cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
73813b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
739b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7407207e0fdSHong Zhang   PetscInt       ntlm;
741715f1b00SHong Zhang   KSP            ksp;
7427207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
74313b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
744*1a352fa8SHong Zhang   PetscReal      previous_shift;
745715f1b00SHong Zhang   PetscErrorCode ierr;
746715f1b00SHong Zhang 
747715f1b00SHong Zhang   PetscFunctionBegin;
748*1a352fa8SHong Zhang   previous_shift = th->shift;
74913b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
75013b1b0ffSHong Zhang 
7517207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7527207e0fdSHong Zhang     ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7536f92202bSHong Zhang   }
754715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
755cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
7567207e0fdSHong Zhang   if (quadts) {
757cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
758cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
7597207e0fdSHong Zhang   }
760715f1b00SHong Zhang 
761715f1b00SHong Zhang   /* Build RHS */
762715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
763*1a352fa8SHong Zhang     th->shift = 1./((th->Theta-1.)*th->time_step0);
764*1a352fa8SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime0,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
76513b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
76613b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
767715f1b00SHong Zhang 
768022c081aSHong Zhang     /* Add the f_p forcing terms */
7690e11c21fSHong Zhang     if (ts->Jacp) {
770*1a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime0,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
77133f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7723e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
77333f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7740e11c21fSHong Zhang     }
775715f1b00SHong Zhang   } else { /* 1-stage method */
7763e05ceb1SHong Zhang     th->shift = 0.0;
7773e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
77813b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
77913b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
78013b1b0ffSHong Zhang 
781022c081aSHong Zhang     /* Add the f_p forcing terms */
7820e11c21fSHong Zhang     if (ts->Jacp) {
7833e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
78433f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
785715f1b00SHong Zhang     }
7860e11c21fSHong Zhang   }
787715f1b00SHong Zhang 
788715f1b00SHong Zhang   /* Build LHS */
789*1a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
790715f1b00SHong Zhang   if (th->endpoint) {
7913e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
792715f1b00SHong Zhang   } else {
7933e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
794715f1b00SHong Zhang   }
795cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
796715f1b00SHong Zhang 
797715f1b00SHong Zhang   /*
798715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
799c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
800715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
801715f1b00SHong Zhang   */
802715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8037207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
804*1a352fa8SHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime0,th->X0,quadJ,NULL);CHKERRQ(ierr);
805*1a352fa8SHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->ptime0,th->X0,quadJp);CHKERRQ(ierr);
8067207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8077207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
808*1a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
809715f1b00SHong Zhang     }
810715f1b00SHong Zhang   }
811715f1b00SHong Zhang 
812715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
813022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
814b5b4017aSHong Zhang     KSPConvergedReason kspreason;
81513b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
816b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
817715f1b00SHong Zhang     if (th->endpoint) {
81813b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
819b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
820b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
821b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
82213b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
823715f1b00SHong Zhang     } else {
824b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
825715f1b00SHong Zhang     }
826b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
827b5b4017aSHong Zhang     if (kspreason < 0) {
828b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
829b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
830b5b4017aSHong Zhang     }
831b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
83213b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
833715f1b00SHong Zhang   }
834715f1b00SHong Zhang 
83513b1b0ffSHong Zhang   /*
83613b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
837c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
83813b1b0ffSHong Zhang   */
8397207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
84013b1b0ffSHong Zhang     if (!th->endpoint) {
8417207e0fdSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */
842cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
843cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
8447207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8457207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
846*1a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
84713b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
84813b1b0ffSHong Zhang     } else {
849cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
850cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
8517207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8527207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
853*1a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
85413b1b0ffSHong Zhang     }
85513b1b0ffSHong Zhang   } else {
85613b1b0ffSHong Zhang     if (!th->endpoint) {
85713b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
858715f1b00SHong Zhang     }
859715f1b00SHong Zhang   }
8601566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8611566a47fSLisandro Dalcin }
8621566a47fSLisandro Dalcin 
8636affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
8646affb6f8SHong Zhang {
8656affb6f8SHong Zhang   TS_Theta *th = (TS_Theta*)ts->data;
8666affb6f8SHong Zhang 
8676affb6f8SHong Zhang   PetscFunctionBegin;
8686affb6f8SHong Zhang   if (ns) *ns = 1;
8696affb6f8SHong Zhang   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
8706affb6f8SHong Zhang   PetscFunctionReturn(0);
8716affb6f8SHong Zhang }
8726affb6f8SHong Zhang 
873316643e7SJed Brown /*------------------------------------------------------------*/
874277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
875316643e7SJed Brown {
876316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
877316643e7SJed Brown   PetscErrorCode ierr;
878316643e7SJed Brown 
879316643e7SJed Brown   PetscFunctionBegin;
8806bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
8816bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
8823b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
883eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
8841566a47fSLisandro Dalcin 
8851566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
8861566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
8871566a47fSLisandro Dalcin 
888b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
889277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
890277b19d0SLisandro Dalcin }
891277b19d0SLisandro Dalcin 
892ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
893ece46509SHong Zhang {
894ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
895ece46509SHong Zhang   PetscErrorCode ierr;
896ece46509SHong Zhang 
897ece46509SHong Zhang   PetscFunctionBegin;
898ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
899ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
900ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
901ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
902ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
903ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
904ece46509SHong Zhang   PetscFunctionReturn(0);
905ece46509SHong Zhang }
906ece46509SHong Zhang 
907277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
908277b19d0SLisandro Dalcin {
909277b19d0SLisandro Dalcin   PetscErrorCode ierr;
910277b19d0SLisandro Dalcin 
911277b19d0SLisandro Dalcin   PetscFunctionBegin;
912277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
913b3a6b972SJed Brown   if (ts->dm) {
914b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
915b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
916b3a6b972SJed Brown   }
917277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
918bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
919bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
920bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
921bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
922316643e7SJed Brown   PetscFunctionReturn(0);
923316643e7SJed Brown }
924316643e7SJed Brown 
925316643e7SJed Brown /*
926316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9272b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
928316643e7SJed Brown */
9290f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
930316643e7SJed Brown {
931316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
932316643e7SJed Brown   PetscErrorCode ierr;
9337445fe48SJed Brown   Vec            X0,Xdot;
9347445fe48SJed Brown   DM             dm,dmsave;
9353e05ceb1SHong Zhang   PetscReal      shift = th->shift;
936316643e7SJed Brown 
937316643e7SJed Brown   PetscFunctionBegin;
9387445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9395a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9407445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
941494ce9fcSHong Zhang   if (x != X0) {
942b296d7d5SJed Brown     ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
943494ce9fcSHong Zhang   } else {
944494ce9fcSHong Zhang     ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
945494ce9fcSHong Zhang   }
9467445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9477445fe48SJed Brown   dmsave = ts->dm;
9487445fe48SJed Brown   ts->dm = dm;
9497445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
9507445fe48SJed Brown   ts->dm = dmsave;
9510d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
952316643e7SJed Brown   PetscFunctionReturn(0);
953316643e7SJed Brown }
954316643e7SJed Brown 
955d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
956316643e7SJed Brown {
957316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
958316643e7SJed Brown   PetscErrorCode ierr;
9597445fe48SJed Brown   Vec            Xdot;
9607445fe48SJed Brown   DM             dm,dmsave;
9613e05ceb1SHong Zhang   PetscReal      shift = th->shift;
962316643e7SJed Brown 
963316643e7SJed Brown   PetscFunctionBegin;
9647445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
965be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9660298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
9677445fe48SJed Brown 
9687445fe48SJed Brown   dmsave = ts->dm;
9697445fe48SJed Brown   ts->dm = dm;
970d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
9717445fe48SJed Brown   ts->dm = dmsave;
9720298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
973316643e7SJed Brown   PetscFunctionReturn(0);
974316643e7SJed Brown }
975316643e7SJed Brown 
976715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
977715f1b00SHong Zhang {
978715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9797207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
980715f1b00SHong Zhang   PetscErrorCode ierr;
981715f1b00SHong Zhang 
982715f1b00SHong Zhang   PetscFunctionBegin;
983022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
98413b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
98513b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
9867207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9877207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
9887207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
989715f1b00SHong Zhang   }
9906f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
99113b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
99213b1b0ffSHong Zhang 
993b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
994715f1b00SHong Zhang   PetscFunctionReturn(0);
995715f1b00SHong Zhang }
996715f1b00SHong Zhang 
9977adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
9987adebddeSHong Zhang {
9997adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10007207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10017adebddeSHong Zhang   PetscErrorCode ierr;
10027adebddeSHong Zhang 
10037adebddeSHong Zhang   PetscFunctionBegin;
10047207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10057207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10067207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
10077adebddeSHong Zhang   }
10087adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
10097adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10107adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
10117adebddeSHong Zhang   PetscFunctionReturn(0);
10127adebddeSHong Zhang }
10137adebddeSHong Zhang 
1014316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1015316643e7SJed Brown {
1016316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1017cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10182ffb9264SLisandro Dalcin   PetscBool      match;
1019316643e7SJed Brown   PetscErrorCode ierr;
1020316643e7SJed Brown 
1021316643e7SJed Brown   PetscFunctionBegin;
1022cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1023cd4cee2dSHong Zhang     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1024b1cb36f3SHong Zhang   }
102539d13666SHong Zhang   if (!th->X) {
1026316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
102739d13666SHong Zhang   }
102839d13666SHong Zhang   if (!th->Xdot) {
1029316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
103039d13666SHong Zhang   }
103139d13666SHong Zhang   if (!th->X0) {
10323b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
103339d13666SHong Zhang   }
10341566a47fSLisandro Dalcin   if (th->endpoint) {
10351566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
10367445fe48SJed Brown   }
10373b1890cdSShri Abhyankar 
10381566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
103920d49056SMatthew G. Knepley   th->shift = 1/(th->Theta*ts->time_step);
10401566a47fSLisandro Dalcin 
10411566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
10421566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10431566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10441566a47fSLisandro Dalcin 
10451566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
10461566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
10472ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
10482ffb9264SLisandro Dalcin   if (!match) {
10491566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
10501566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
10513b1890cdSShri Abhyankar   }
10521566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1053b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1054b4dd3bf9SBarry Smith }
10550c7376e5SHong Zhang 
1056b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1057b4dd3bf9SBarry Smith 
1058b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1059b4dd3bf9SBarry Smith {
1060b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1061b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1062b4dd3bf9SBarry Smith 
1063b4dd3bf9SBarry Smith   PetscFunctionBegin;
1064b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1065b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
10662ca6e920SHong Zhang   if (ts->vecs_sensip) {
1067b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
10682ca6e920SHong Zhang   }
1069b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1070b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1071b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
107267633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
107367633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
107467633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1075b5b4017aSHong Zhang   }
1076c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1077c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
107867633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
107967633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
108067633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1081b5b4017aSHong Zhang   }
1082316643e7SJed Brown   PetscFunctionReturn(0);
1083316643e7SJed Brown }
1084316643e7SJed Brown 
10854416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1086316643e7SJed Brown {
1087316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1088316643e7SJed Brown   PetscErrorCode ierr;
1089316643e7SJed Brown 
1090316643e7SJed Brown   PetscFunctionBegin;
1091e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1092316643e7SJed Brown   {
10930298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
10940298fd71SBarry 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);
109503842d09SLisandro 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);
1096316643e7SJed Brown   }
1097316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1098316643e7SJed Brown   PetscFunctionReturn(0);
1099316643e7SJed Brown }
1100316643e7SJed Brown 
1101316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1102316643e7SJed Brown {
1103316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1104ace3abfcSBarry Smith   PetscBool      iascii;
1105316643e7SJed Brown   PetscErrorCode ierr;
1106316643e7SJed Brown 
1107316643e7SJed Brown   PetscFunctionBegin;
1108251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1109316643e7SJed Brown   if (iascii) {
11107c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1111ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1112316643e7SJed Brown   }
1113316643e7SJed Brown   PetscFunctionReturn(0);
1114316643e7SJed Brown }
1115316643e7SJed Brown 
1116560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11170de4c49aSJed Brown {
11180de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11190de4c49aSJed Brown 
11200de4c49aSJed Brown   PetscFunctionBegin;
11210de4c49aSJed Brown   *theta = th->Theta;
11220de4c49aSJed Brown   PetscFunctionReturn(0);
11230de4c49aSJed Brown }
11240de4c49aSJed Brown 
1125560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11260de4c49aSJed Brown {
11270de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11280de4c49aSJed Brown 
11290de4c49aSJed Brown   PetscFunctionBegin;
11307c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11310de4c49aSJed Brown   th->Theta = theta;
11321566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11330de4c49aSJed Brown   PetscFunctionReturn(0);
11340de4c49aSJed Brown }
1135eb284becSJed Brown 
1136560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
113726f2ff8fSLisandro Dalcin {
113826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
113926f2ff8fSLisandro Dalcin 
114026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
114126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
114226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
114326f2ff8fSLisandro Dalcin }
114426f2ff8fSLisandro Dalcin 
1145560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1146eb284becSJed Brown {
1147eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1148eb284becSJed Brown 
1149eb284becSJed Brown   PetscFunctionBegin;
1150eb284becSJed Brown   th->endpoint = flg;
1151eb284becSJed Brown   PetscFunctionReturn(0);
1152eb284becSJed Brown }
11530de4c49aSJed Brown 
1154f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1155f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1156f9c1d6abSBarry Smith {
1157f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1158f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11593fd8ae06SJed Brown   const PetscReal one = 1.0;
1160f9c1d6abSBarry Smith 
1161f9c1d6abSBarry Smith   PetscFunctionBegin;
11623fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1163f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1164f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1165f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1166f9c1d6abSBarry Smith }
1167f9c1d6abSBarry Smith #endif
1168f9c1d6abSBarry Smith 
116942682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
117042682096SHong Zhang {
117142682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
117242682096SHong Zhang 
117342682096SHong Zhang   PetscFunctionBegin;
11741566a47fSLisandro Dalcin   if (ns) *ns = 1;
11751566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
117642682096SHong Zhang   PetscFunctionReturn(0);
117742682096SHong Zhang }
1178f9c1d6abSBarry Smith 
1179316643e7SJed Brown /* ------------------------------------------------------------ */
1180316643e7SJed Brown /*MC
118196f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1182316643e7SJed Brown 
1183316643e7SJed Brown    Level: beginner
1184316643e7SJed Brown 
11854eb428fdSBarry Smith    Options Database:
1186c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1187c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
118803842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11894eb428fdSBarry Smith 
1190eb284becSJed Brown    Notes:
11910c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
11920c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
11934eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
11944eb428fdSBarry Smith 
1195eb284becSJed Brown    This method can be applied to DAE.
1196eb284becSJed Brown 
1197eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
1198eb284becSJed Brown 
1199eb284becSJed Brown .vb
1200eb284becSJed Brown   Theta | Theta
1201eb284becSJed Brown   -------------
1202eb284becSJed Brown         |  1
1203eb284becSJed Brown .ve
1204eb284becSJed Brown 
1205eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1206eb284becSJed Brown 
1207eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1208eb284becSJed Brown 
1209eb284becSJed Brown .vb
1210eb284becSJed Brown   0 | 0         0
1211eb284becSJed Brown   1 | 1-Theta   Theta
1212eb284becSJed Brown   -------------------
1213eb284becSJed Brown     | 1-Theta   Theta
1214eb284becSJed Brown .ve
1215eb284becSJed Brown 
1216eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1217eb284becSJed Brown 
1218eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1219eb284becSJed Brown 
1220eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1221eb284becSJed Brown 
12224eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1223eb284becSJed Brown 
1224eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1225316643e7SJed Brown 
1226316643e7SJed Brown M*/
12278cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1228316643e7SJed Brown {
1229316643e7SJed Brown   TS_Theta       *th;
1230316643e7SJed Brown   PetscErrorCode ierr;
1231316643e7SJed Brown 
1232316643e7SJed Brown   PetscFunctionBegin;
1233277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1234ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1235316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1236316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1237316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
123842f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1239ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1240316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1241cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12421566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
124324655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1244316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12450f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12460f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1247f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1248f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1249f9c1d6abSBarry Smith #endif
125042682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
125142f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1252b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1253b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12542ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12556affb6f8SHong Zhang 
1256715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12577adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1258715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12596affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1260316643e7SJed Brown 
1261efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1262efd4aadfSBarry Smith 
1263b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1264316643e7SJed Brown   ts->data = (void*)th;
1265316643e7SJed Brown 
1266715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1267715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1268715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1269b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1270715f1b00SHong Zhang 
12716f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1272316643e7SJed Brown   th->Theta       = 0.5;
12731566a47fSLisandro Dalcin   th->order       = 2;
1274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1277bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1278316643e7SJed Brown   PetscFunctionReturn(0);
1279316643e7SJed Brown }
12800de4c49aSJed Brown 
12810de4c49aSJed Brown /*@
12820de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12830de4c49aSJed Brown 
12840de4c49aSJed Brown   Not Collective
12850de4c49aSJed Brown 
12860de4c49aSJed Brown   Input Parameter:
12870de4c49aSJed Brown .  ts - timestepping context
12880de4c49aSJed Brown 
12890de4c49aSJed Brown   Output Parameter:
12900de4c49aSJed Brown .  theta - stage abscissa
12910de4c49aSJed Brown 
12920de4c49aSJed Brown   Note:
12930de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
12940de4c49aSJed Brown 
12950de4c49aSJed Brown   Level: Advanced
12960de4c49aSJed Brown 
12970de4c49aSJed Brown .seealso: TSThetaSetTheta()
12980de4c49aSJed Brown @*/
12997087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13000de4c49aSJed Brown {
13014ac538c5SBarry Smith   PetscErrorCode ierr;
13020de4c49aSJed Brown 
13030de4c49aSJed Brown   PetscFunctionBegin;
1304afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13050de4c49aSJed Brown   PetscValidPointer(theta,2);
13064ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
13070de4c49aSJed Brown   PetscFunctionReturn(0);
13080de4c49aSJed Brown }
13090de4c49aSJed Brown 
13100de4c49aSJed Brown /*@
13110de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13120de4c49aSJed Brown 
13130de4c49aSJed Brown   Not Collective
13140de4c49aSJed Brown 
13150de4c49aSJed Brown   Input Parameter:
13160de4c49aSJed Brown +  ts - timestepping context
13170de4c49aSJed Brown -  theta - stage abscissa
13180de4c49aSJed Brown 
13190de4c49aSJed Brown   Options Database:
13200de4c49aSJed Brown .  -ts_theta_theta <theta>
13210de4c49aSJed Brown 
13220de4c49aSJed Brown   Level: Intermediate
13230de4c49aSJed Brown 
13240de4c49aSJed Brown .seealso: TSThetaGetTheta()
13250de4c49aSJed Brown @*/
13267087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13270de4c49aSJed Brown {
13284ac538c5SBarry Smith   PetscErrorCode ierr;
13290de4c49aSJed Brown 
13300de4c49aSJed Brown   PetscFunctionBegin;
1331afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13324ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
13330de4c49aSJed Brown   PetscFunctionReturn(0);
13340de4c49aSJed Brown }
1335f33bbcb6SJed Brown 
133626f2ff8fSLisandro Dalcin /*@
133726f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
133826f2ff8fSLisandro Dalcin 
133926f2ff8fSLisandro Dalcin   Not Collective
134026f2ff8fSLisandro Dalcin 
134126f2ff8fSLisandro Dalcin   Input Parameter:
134226f2ff8fSLisandro Dalcin .  ts - timestepping context
134326f2ff8fSLisandro Dalcin 
134426f2ff8fSLisandro Dalcin   Output Parameter:
134526f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
134626f2ff8fSLisandro Dalcin 
134726f2ff8fSLisandro Dalcin   Level: Advanced
134826f2ff8fSLisandro Dalcin 
134926f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
135026f2ff8fSLisandro Dalcin @*/
135126f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
135226f2ff8fSLisandro Dalcin {
135326f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
135426f2ff8fSLisandro Dalcin 
135526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
135626f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
135726f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1358163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
135926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
136026f2ff8fSLisandro Dalcin }
136126f2ff8fSLisandro Dalcin 
1362eb284becSJed Brown /*@
1363eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1364eb284becSJed Brown 
1365eb284becSJed Brown   Not Collective
1366eb284becSJed Brown 
1367eb284becSJed Brown   Input Parameter:
1368eb284becSJed Brown +  ts - timestepping context
1369eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1370eb284becSJed Brown 
1371eb284becSJed Brown   Options Database:
1372eb284becSJed Brown .  -ts_theta_endpoint <flg>
1373eb284becSJed Brown 
1374eb284becSJed Brown   Level: Intermediate
1375eb284becSJed Brown 
1376eb284becSJed Brown .seealso: TSTHETA, TSCN
1377eb284becSJed Brown @*/
1378eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1379eb284becSJed Brown {
1380eb284becSJed Brown   PetscErrorCode ierr;
1381eb284becSJed Brown 
1382eb284becSJed Brown   PetscFunctionBegin;
1383eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1384eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1385eb284becSJed Brown   PetscFunctionReturn(0);
1386eb284becSJed Brown }
1387eb284becSJed Brown 
1388f33bbcb6SJed Brown /*
1389f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1390f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1391f33bbcb6SJed Brown  */
1392f33bbcb6SJed Brown 
13931566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
13941566a47fSLisandro Dalcin {
13951566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
13961566a47fSLisandro Dalcin   PetscErrorCode ierr;
13971566a47fSLisandro Dalcin 
13981566a47fSLisandro Dalcin   PetscFunctionBegin;
13991566a47fSLisandro 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");
14001566a47fSLisandro 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");
14011566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14021566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14031566a47fSLisandro Dalcin }
14041566a47fSLisandro Dalcin 
1405f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1406f33bbcb6SJed Brown {
1407f33bbcb6SJed Brown   PetscFunctionBegin;
1408f33bbcb6SJed Brown   PetscFunctionReturn(0);
1409f33bbcb6SJed Brown }
1410f33bbcb6SJed Brown 
1411f33bbcb6SJed Brown /*MC
1412f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1413f33bbcb6SJed Brown 
1414f33bbcb6SJed Brown   Level: beginner
1415f33bbcb6SJed Brown 
14164eb428fdSBarry Smith   Notes:
1417c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14184eb428fdSBarry Smith 
14191566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14204eb428fdSBarry Smith 
1421f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1422f33bbcb6SJed Brown 
1423f33bbcb6SJed Brown M*/
14248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1425f33bbcb6SJed Brown {
1426f33bbcb6SJed Brown   PetscErrorCode ierr;
1427f33bbcb6SJed Brown 
1428f33bbcb6SJed Brown   PetscFunctionBegin;
1429f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1430f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
14311566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
14320c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1433f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1434f33bbcb6SJed Brown   PetscFunctionReturn(0);
1435f33bbcb6SJed Brown }
1436f33bbcb6SJed Brown 
14371566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14381566a47fSLisandro Dalcin {
14391566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14401566a47fSLisandro Dalcin   PetscErrorCode ierr;
14411566a47fSLisandro Dalcin 
14421566a47fSLisandro Dalcin   PetscFunctionBegin;
14431566a47fSLisandro 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");
14441566a47fSLisandro 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");
14451566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14461566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14471566a47fSLisandro Dalcin }
14481566a47fSLisandro Dalcin 
1449f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1450f33bbcb6SJed Brown {
1451f33bbcb6SJed Brown   PetscFunctionBegin;
1452f33bbcb6SJed Brown   PetscFunctionReturn(0);
1453f33bbcb6SJed Brown }
1454f33bbcb6SJed Brown 
1455f33bbcb6SJed Brown /*MC
1456f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1457f33bbcb6SJed Brown 
1458f33bbcb6SJed Brown   Level: beginner
1459f33bbcb6SJed Brown 
1460f33bbcb6SJed Brown   Notes:
14617cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14627cf5af47SJed Brown 
14637cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1464f33bbcb6SJed Brown 
1465f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1466f33bbcb6SJed Brown 
1467f33bbcb6SJed Brown M*/
14688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1469f33bbcb6SJed Brown {
1470f33bbcb6SJed Brown   PetscErrorCode ierr;
1471f33bbcb6SJed Brown 
1472f33bbcb6SJed Brown   PetscFunctionBegin;
1473f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1474f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1475eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
14760c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1477f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1478f33bbcb6SJed Brown   PetscFunctionReturn(0);
1479f33bbcb6SJed Brown }
1480