xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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;
12cc4f23bcSHong Zhang   Vec          Stages[2];                 /* Storage for stage solutions */
13cc4f23bcSHong Zhang   Vec          X0,X,Xdot;                /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */
141566a47fSLisandro Dalcin   Vec          affine;                   /* Affine vector needed for residual at beginning of step in endpoint formulation */
151566a47fSLisandro Dalcin   PetscReal    Theta;
161a352fa8SHong Zhang   PetscReal    shift;                    /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20715f1b00SHong Zhang   TSStepStatus status;
211a352fa8SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events, used by cost integral */
221a352fa8SHong Zhang   PetscReal    ptime0;                   /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
231a352fa8SHong Zhang   PetscReal    time_step0;               /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/
241566a47fSLisandro Dalcin 
25715f1b00SHong Zhang   /* context for sensitivity analysis */
26022c081aSHong Zhang   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
27b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
28b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
2913b1b0ffSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
30cc4f23bcSHong Zhang   Mat          MatFwdStages[2];          /* TLM Stages */
3113b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
32b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
3313b1b0ffSHong Zhang   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
347207e0fdSHong Zhang   Mat          MatIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
357207e0fdSHong Zhang   Mat          MatIntegralSensip0;       /* backup for roll-backs due to events */
36b5b4017aSHong Zhang   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
37b5b4017aSHong Zhang   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
38b5b4017aSHong Zhang   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
39b5b4017aSHong Zhang   Vec          *VecsAffine;              /* Working vectors to store residuals */
40715f1b00SHong Zhang   /* context for error estimation */
411566a47fSLisandro Dalcin   Vec          vec_sol_prev;
421566a47fSLisandro Dalcin   Vec          vec_lte_work;
43316643e7SJed Brown } TS_Theta;
44316643e7SJed Brown 
457445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
467445fe48SJed Brown {
477445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
487445fe48SJed Brown 
497445fe48SJed Brown   PetscFunctionBegin;
507445fe48SJed Brown   if (X0) {
517445fe48SJed Brown     if (dm && dm != ts->dm) {
529566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm,"TSTheta_X0",X0));
537445fe48SJed Brown     } else *X0 = ts->vec_sol;
547445fe48SJed Brown   }
557445fe48SJed Brown   if (Xdot) {
567445fe48SJed Brown     if (dm && dm != ts->dm) {
579566063dSJacob Faibussowitsch       PetscCall(DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot));
587445fe48SJed Brown     } else *Xdot = th->Xdot;
597445fe48SJed Brown   }
607445fe48SJed Brown   PetscFunctionReturn(0);
617445fe48SJed Brown }
627445fe48SJed Brown 
630d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
640d0b770aSPeter Brune {
650d0b770aSPeter Brune   PetscFunctionBegin;
660d0b770aSPeter Brune   if (X0) {
670d0b770aSPeter Brune     if (dm && dm != ts->dm) {
689566063dSJacob Faibussowitsch       PetscCall(DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0));
690d0b770aSPeter Brune     }
700d0b770aSPeter Brune   }
710d0b770aSPeter Brune   if (Xdot) {
720d0b770aSPeter Brune     if (dm && dm != ts->dm) {
739566063dSJacob Faibussowitsch       PetscCall(DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot));
740d0b770aSPeter Brune     }
750d0b770aSPeter Brune   }
760d0b770aSPeter Brune   PetscFunctionReturn(0);
770d0b770aSPeter Brune }
780d0b770aSPeter Brune 
797445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
807445fe48SJed Brown {
817445fe48SJed Brown   PetscFunctionBegin;
827445fe48SJed Brown   PetscFunctionReturn(0);
837445fe48SJed Brown }
847445fe48SJed Brown 
857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
867445fe48SJed Brown {
877445fe48SJed Brown   TS             ts = (TS)ctx;
887445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
897445fe48SJed Brown 
907445fe48SJed Brown   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot));
929566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c));
939566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct,X0,X0_c));
949566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct,Xdot,Xdot_c));
959566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(X0_c,rscale,X0_c));
969566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Xdot_c,rscale,Xdot_c));
979566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot));
989566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c));
997445fe48SJed Brown   PetscFunctionReturn(0);
1007445fe48SJed Brown }
1017445fe48SJed Brown 
102258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
103258e1594SPeter Brune {
104258e1594SPeter Brune   PetscFunctionBegin;
105258e1594SPeter Brune   PetscFunctionReturn(0);
106258e1594SPeter Brune }
107258e1594SPeter Brune 
108258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
109258e1594SPeter Brune {
110258e1594SPeter Brune   TS             ts = (TS)ctx;
111258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
112258e1594SPeter Brune 
113258e1594SPeter Brune   PetscFunctionBegin;
1149566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot));
1159566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub));
116258e1594SPeter Brune 
1179566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD));
1189566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD));
119258e1594SPeter Brune 
1209566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD));
1219566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD));
122258e1594SPeter Brune 
1239566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot));
1249566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub));
125258e1594SPeter Brune   PetscFunctionReturn(0);
126258e1594SPeter Brune }
127258e1594SPeter Brune 
128022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
129022c081aSHong Zhang {
130022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
131cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
132022c081aSHong Zhang 
133022c081aSHong Zhang   PetscFunctionBegin;
134022c081aSHong Zhang   if (th->endpoint) {
135022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
136022c081aSHong Zhang     if (th->Theta!=1.0) {
1379566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSFunction(quadts,th->ptime0,th->X0,ts->vec_costintegrand));
1389566063dSJacob Faibussowitsch       PetscCall(VecAXPY(quadts->vec_sol,th->time_step0*(1.0-th->Theta),ts->vec_costintegrand));
139022c081aSHong Zhang     }
1409566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand));
1419566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,th->time_step0*th->Theta,ts->vec_costintegrand));
142022c081aSHong Zhang   } else {
1439566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand));
1449566063dSJacob Faibussowitsch     PetscCall(VecAXPY(quadts->vec_sol,th->time_step0,ts->vec_costintegrand));
145022c081aSHong Zhang   }
146022c081aSHong Zhang   PetscFunctionReturn(0);
147022c081aSHong Zhang }
148022c081aSHong Zhang 
149b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
150b1cb36f3SHong Zhang {
151b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
152cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
153b1cb36f3SHong Zhang 
154b1cb36f3SHong Zhang   PetscFunctionBegin;
155b1cb36f3SHong Zhang   /* backup cost integral */
1569566063dSJacob Faibussowitsch   PetscCall(VecCopy(quadts->vec_sol,th->VecCostIntegral0));
1579566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
158b1cb36f3SHong Zhang   PetscFunctionReturn(0);
159b1cb36f3SHong Zhang }
160b1cb36f3SHong Zhang 
161b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
162b1cb36f3SHong Zhang {
1631a352fa8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
164b1cb36f3SHong Zhang 
165b1cb36f3SHong Zhang   PetscFunctionBegin;
1661a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
1671a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
1681a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
1699566063dSJacob Faibussowitsch   PetscCall(TSThetaEvaluateCostIntegral(ts));
17024655328SShri   PetscFunctionReturn(0);
17124655328SShri }
17224655328SShri 
173470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1741566a47fSLisandro Dalcin {
1751566a47fSLisandro Dalcin   PetscInt       nits,lits;
1761566a47fSLisandro Dalcin 
1771566a47fSLisandro Dalcin   PetscFunctionBegin;
1789566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes,b,x));
1799566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes,&nits));
1809566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits));
1811566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1821566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1831566a47fSLisandro Dalcin }
1841566a47fSLisandro Dalcin 
185193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
186316643e7SJed Brown {
187316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1881566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1894957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1901566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
191316643e7SJed Brown 
192316643e7SJed Brown   PetscFunctionBegin;
1931566a47fSLisandro Dalcin   if (!ts->steprollback) {
1949566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0,th->vec_sol_prev));
1959566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol,th->X0));
1961566a47fSLisandro Dalcin   }
1971566a47fSLisandro Dalcin 
1981566a47fSLisandro Dalcin   th->status     = TS_STEP_INCOMPLETE;
19999260352SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
2003e05ceb1SHong Zhang     th->shift      = 1/(th->Theta*ts->time_step);
2011566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
2029566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0,th->X));
203fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
2049566063dSJacob Faibussowitsch       PetscCall(VecAXPY(th->X,1/th->shift,th->Xdot));
205be5899b3SLisandro Dalcin     }
206eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
2079566063dSJacob Faibussowitsch       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol,&th->affine));
2089566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
2099566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE));
2109566063dSJacob Faibussowitsch       PetscCall(VecScale(th->affine,(th->Theta-1)/th->Theta));
2111566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2129566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->affine));
213eb284becSJed Brown     }
2149566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,th->stage_time));
2159566063dSJacob Faibussowitsch     PetscCall(TSTheta_SNESSolve(ts,th->affine,th->X));
2169566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts,th->stage_time,0,&th->X));
2179566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok));
218fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
219051f2191SLisandro Dalcin 
220051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2211566a47fSLisandro Dalcin     if (th->endpoint) {
2229566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X,ts->vec_sol));
2231566a47fSLisandro Dalcin     } else {
2249566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X));
2259566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vec_sol,ts->time_step,th->Xdot));
2261566a47fSLisandro Dalcin     }
2279566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
2281566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2291566a47fSLisandro Dalcin     if (!accept) {
2309566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0,ts->vec_sol));
231be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
232be5899b3SLisandro Dalcin       goto reject_step;
233051f2191SLisandro Dalcin     }
2343b1890cdSShri Abhyankar 
235715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2361a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2371a352fa8SHong Zhang       th->time_step0 = ts->time_step;
23817145e75SHong Zhang     }
2392b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
240cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
241051f2191SLisandro Dalcin     break;
242051f2191SLisandro Dalcin 
243051f2191SLisandro Dalcin   reject_step:
244fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2451566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
246051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
24763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,rejections));
2483b1890cdSShri Abhyankar     }
2493b1890cdSShri Abhyankar   }
250316643e7SJed Brown   PetscFunctionReturn(0);
251316643e7SJed Brown }
252316643e7SJed Brown 
253c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
254c9681e5dSHong Zhang {
255c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
256cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
257c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
258c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
259c9681e5dSHong Zhang   PetscInt       nadj;
2607207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
261c9681e5dSHong Zhang   KSP            ksp;
262c9681e5dSHong Zhang   PetscScalar    *xarr;
263077c4feaSHong Zhang   TSEquationType eqtype;
264077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2651a352fa8SHong Zhang   PetscReal      adjoint_time_step;
266c9681e5dSHong Zhang 
267c9681e5dSHong Zhang   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(TSGetEquationType(ts,&eqtype));
269077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
270077c4feaSHong Zhang     isexplicitode  = PETSC_TRUE;
271077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
272077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
273077c4feaSHong Zhang   }
274c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
2759566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes,&ksp));
2769566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,&Jpre,NULL,NULL));
2777207e0fdSHong Zhang   if (quadts) {
2789566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL));
2799566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL));
2807207e0fdSHong Zhang   }
281c9681e5dSHong Zhang 
2821a352fa8SHong Zhang   th->stage_time    = ts->ptime;
2831a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
284c9681e5dSHong Zhang 
28533f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
286*1baa6e33SBarry Smith   if (quadts) PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
287cd4cee2dSHong Zhang 
288c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
2899566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]));
2909566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj],1./adjoint_time_step)); /* lambda_{n+1}/h */
291cd4cee2dSHong Zhang     if (quadJ) {
2929566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
2939566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
2949566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col));
2959566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
2969566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
297c9681e5dSHong Zhang     }
298c9681e5dSHong Zhang   }
299c9681e5dSHong Zhang 
300c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
301c87ba875SIvan Yashchuk   th->shift = 1./adjoint_time_step;
3029566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
3039566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
304c9681e5dSHong Zhang 
305c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
306c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
307c9681e5dSHong Zhang     KSPConvergedReason kspreason;
3089566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]));
3099566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
310c9681e5dSHong Zhang     if (kspreason < 0) {
311c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
31263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj));
313c9681e5dSHong Zhang     }
314c9681e5dSHong Zhang   }
315c9681e5dSHong Zhang 
316c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
317c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3189566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
3199566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
320c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3219566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
322c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3239566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
324c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
3259566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
3269566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj],1./adjoint_time_step));
3279566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]));
328c9681e5dSHong Zhang       if (ts->vecs_fup) {
3299566063dSJacob Faibussowitsch         PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]));
330c9681e5dSHong Zhang       }
331c9681e5dSHong Zhang     }
332c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
333c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
33405c9950eSHong Zhang       KSPConvergedReason kspreason;
3359566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]));
3369566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp,&kspreason));
33705c9950eSHong Zhang       if (kspreason < 0) {
33805c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
33963a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj));
34005c9950eSHong Zhang       }
341c9681e5dSHong Zhang     }
342c9681e5dSHong Zhang   }
343c9681e5dSHong Zhang 
344c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
345077c4feaSHong Zhang   if (!isexplicitode) {
3463e05ceb1SHong Zhang     th->shift = 0.0;
3479566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
3489566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
349c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
35033f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3519566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj],-1.));
3529566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]));
3539566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj],-adjoint_time_step));
3549566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]));
355c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3569566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]));
3579566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj],-adjoint_time_step));
3589566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]));
359c9681e5dSHong Zhang       }
360c9681e5dSHong Zhang     }
361077c4feaSHong Zhang   }
362c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3639566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot,-1./adjoint_time_step,1.0/adjoint_time_step,0,th->X0,ts->vec_sol));
3649566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,1./adjoint_time_step,ts->Jacp,PETSC_FALSE)); /* get -f_p */
365*1baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
366c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
367c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3689566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
36933f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3709566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
371c9681e5dSHong Zhang     }
372cd4cee2dSHong Zhang 
373c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
3749566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
3759566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]));
376cd4cee2dSHong Zhang       if (quadJp) {
3779566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
3789566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
3799566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col));
3809566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3819566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
38233f72c5dSHong Zhang       }
383c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
3849566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
3859566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,VecsDeltaMu2[nadj]));
386c9681e5dSHong Zhang         if (ts->vecs_fpu) {
3879566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpu[nadj]));
388c9681e5dSHong Zhang         }
389c9681e5dSHong Zhang         if (ts->vecs_fpp) {
3909566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpp[nadj]));
391c9681e5dSHong Zhang         }
392c9681e5dSHong Zhang       }
393c9681e5dSHong Zhang     }
394c9681e5dSHong Zhang   }
395c9681e5dSHong Zhang 
396c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
3979566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
3989566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
399c9681e5dSHong Zhang   }
400c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
401c9681e5dSHong Zhang   PetscFunctionReturn(0);
402c9681e5dSHong Zhang }
403c9681e5dSHong Zhang 
40442f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4052ca6e920SHong Zhang {
4062ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
407cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
408b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
409b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4102ca6e920SHong Zhang   PetscInt       nadj;
4117207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4122ca6e920SHong Zhang   KSP            ksp;
413b5b4017aSHong Zhang   PetscScalar    *xarr;
4141a352fa8SHong Zhang   PetscReal      adjoint_time_step;
4151a352fa8SHong 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) */
4162ca6e920SHong Zhang 
4172ca6e920SHong Zhang   PetscFunctionBegin;
418c9681e5dSHong Zhang   if (th->Theta == 1.) {
4199566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
420c9681e5dSHong Zhang     PetscFunctionReturn(0);
421c9681e5dSHong Zhang   }
4222ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4239566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes,&ksp));
4249566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,&Jpre,NULL,NULL));
4257207e0fdSHong Zhang   if (quadts) {
4269566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL));
4279566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL));
4287207e0fdSHong Zhang   }
429b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4301a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step);
4311a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4321a352fa8SHong Zhang   adjoint_time_step = -ts->time_step;  /* always positive since time_step is negative */
4332ca6e920SHong Zhang 
43482ad9101SHong Zhang   if (!th->endpoint) {
4355072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4369566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0,1.0/(1.0-th->Theta),th->Theta/(th->Theta-1.0),0,th->X,ts->vec_sol));
43782ad9101SHong Zhang   }
43882ad9101SHong Zhang 
439b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
44033f72c5dSHong Zhang   /* Cost function has an integral term */
4417207e0fdSHong Zhang   if (quadts) {
44205755b9cSHong Zhang     if (th->endpoint) {
4439566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
44405755b9cSHong Zhang     } else {
4459566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
44605755b9cSHong Zhang     }
4477207e0fdSHong Zhang   }
44833f72c5dSHong Zhang 
449abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4509566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]));
4519566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj],1./(th->Theta*adjoint_time_step)));
452cd4cee2dSHong Zhang     if (quadJ) {
4539566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
4549566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
4559566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col));
4569566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4579566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
45836eaed60SHong Zhang     }
4592ca6e920SHong Zhang   }
4603c54f92cSHong Zhang 
461b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4621a352fa8SHong Zhang   th->shift = 1./(th->Theta*adjoint_time_step);
4633c54f92cSHong Zhang   if (th->endpoint) {
4649566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
4653c54f92cSHong Zhang   } else {
4669566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X,J,Jpre));
4673c54f92cSHong Zhang   }
4689566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
4692ca6e920SHong Zhang 
470b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
471abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4721d14d03bSHong Zhang     KSPConvergedReason kspreason;
4739566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]));
4749566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
4751d14d03bSHong Zhang     if (kspreason < 0) {
4761d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
47763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj));
4781d14d03bSHong Zhang     }
4792ca6e920SHong Zhang   }
4803c54f92cSHong Zhang 
4811d14d03bSHong Zhang   /* Second-order adjoint */
482b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4833c633725SBarry Smith     PetscCheck(th->endpoint,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
484b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4859566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
4869566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
487b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
4889566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
4899566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4909566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
491b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
4929566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
493b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
4949566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
4959566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj],th->shift));
4969566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]));
497b5b4017aSHong Zhang       if (ts->vecs_fup) {
4989566063dSJacob Faibussowitsch         PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]));
499b5b4017aSHong Zhang       }
500b5b4017aSHong Zhang     }
501b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
502b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
5031d14d03bSHong Zhang       KSPConvergedReason kspreason;
5049566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]));
5059566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp,&kspreason));
5061d14d03bSHong Zhang       if (kspreason < 0) {
5071d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
50863a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj));
5091d14d03bSHong Zhang       }
510b5b4017aSHong Zhang     }
511b5b4017aSHong Zhang   }
512b5b4017aSHong Zhang 
51336eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5141d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5151a352fa8SHong Zhang     th->shift      = 1./((th->Theta-1.)*adjoint_time_step);
5161a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
5179566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X0,J,Jpre));
5189566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
51933f72c5dSHong Zhang     /* R_U at t_n */
520*1baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts,adjoint_ptime,th->X0,quadJ,NULL));
521abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
5229566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]));
523cd4cee2dSHong Zhang       if (quadJ) {
5249566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
5259566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
5269566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col));
5279566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5289566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
529b5b4017aSHong Zhang       }
5309566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj],1./th->shift));
531b5b4017aSHong Zhang     }
5321d14d03bSHong Zhang 
5331d14d03bSHong Zhang     /* Second-order adjoint */
5341d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
535b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5369566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0,0,&xarr));
5379566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
538b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5399566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
5409566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5419566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0,&xarr));
542b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5439566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
544b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
54533f72c5dSHong 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  */
5469566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]));
5479566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]));
548b5b4017aSHong Zhang         if (ts->vecs_fup) {
5499566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]));
550b5b4017aSHong Zhang         }
5519566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj],1./th->shift));
552b5b4017aSHong Zhang       }
553b5e0532dSHong Zhang     }
5543fd52205SHong Zhang 
555c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
556c586f2e8SHong Zhang 
5573fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
558b5b4017aSHong Zhang       /* U_{n+1} */
55982ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
5609566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol));
5619566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,-1./(th->Theta*adjoint_time_step),ts->Jacp,PETSC_FALSE));
562*1baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
563abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
5649566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
5659566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu[nadj]));
5663b711c3fSHong Zhang         if (quadJp) {
5679566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
5689566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
5699566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*th->Theta,ts->vec_drdp_col));
5709566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5719566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
5723b711c3fSHong Zhang         }
5733fd52205SHong Zhang       }
57433f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
57533f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5769566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
5779566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
578b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5799566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
5809566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5819566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
58233f72c5dSHong Zhang 
583b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5849566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
585b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
58633f72c5dSHong 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)  */
5879566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
5889566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu2[nadj]));
589b5b4017aSHong Zhang           if (ts->vecs_fpu) {
5909566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpu[nadj]));
591b5b4017aSHong Zhang           }
592b5b4017aSHong Zhang           if (ts->vecs_fpp) {
5939566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpp[nadj]));
594b5b4017aSHong Zhang           }
595b5b4017aSHong Zhang         }
596b5b4017aSHong Zhang       }
597b5b4017aSHong Zhang 
598b5b4017aSHong Zhang       /* U_s */
5999566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
6009566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,adjoint_ptime,th->X0,th->Xdot,1./((th->Theta-1.0)*adjoint_time_step),ts->Jacp,PETSC_FALSE));
601*1baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts,adjoint_ptime,th->X0,quadJp));
602abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
6039566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
6049566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu[nadj]));
6053b711c3fSHong Zhang         if (quadJp) {
6069566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
6079566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
6089566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*(1.0-th->Theta),ts->vec_drdp_col));
6099566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6109566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
6113b711c3fSHong Zhang         }
61233f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
61333f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
6149566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0,0,&xarr));
6159566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
616b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6179566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
6189566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
6199566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0,&xarr));
620b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6219566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
622b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
62333f72c5dSHong 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) */
6249566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
6259566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]));
626b5b4017aSHong Zhang             if (ts->vecs_fpu) {
6279566063dSJacob Faibussowitsch               PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]));
628b5e0532dSHong Zhang             }
629b5b4017aSHong Zhang             if (ts->vecs_fpp) {
6309566063dSJacob Faibussowitsch               PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]));
631b5b4017aSHong Zhang             }
63236eaed60SHong Zhang           }
63336eaed60SHong Zhang         }
6343fd52205SHong Zhang       }
635b5e0532dSHong Zhang     }
6363fd52205SHong Zhang   } else { /* one-stage case */
6373e05ceb1SHong Zhang     th->shift = 0.0;
6389566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X,J,Jpre)); /* get -f_y */
6399566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
640*1baa6e33SBarry Smith     if (quadts) PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
641abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
6429566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]));
6439566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj],-adjoint_time_step,VecsSensiTemp[nadj]));
644cd4cee2dSHong Zhang       if (quadJ) {
6459566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
6469566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
6479566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj],adjoint_time_step,ts->vec_drdu_col));
6489566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6499566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
65036eaed60SHong Zhang       }
6512ca6e920SHong Zhang     }
6523fd52205SHong Zhang     if (ts->vecs_sensip) {
65382ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
6549566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X));
6559566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
656*1baa6e33SBarry Smith       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp));
657abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
6589566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
6599566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]));
660cd4cee2dSHong Zhang         if (quadJp) {
6619566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
6629566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
6639566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col));
6649566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6659566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
66636eaed60SHong Zhang         }
66736eaed60SHong Zhang       }
6683fd52205SHong Zhang     }
6692ca6e920SHong Zhang   }
6702ca6e920SHong Zhang 
6712ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6722ca6e920SHong Zhang   PetscFunctionReturn(0);
6732ca6e920SHong Zhang }
6742ca6e920SHong Zhang 
675cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
676cd652676SJed Brown {
677cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
678be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
679cd652676SJed Brown 
680cd652676SJed Brown   PetscFunctionBegin;
6819566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,th->X));
682be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6839566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X,dt,th->Xdot,th->X));
684cd652676SJed Brown   PetscFunctionReturn(0);
685cd652676SJed Brown }
686cd652676SJed Brown 
6871566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
6881566a47fSLisandro Dalcin {
6891566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
6901566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
6911566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
6927453f775SEmil Constantinescu   PetscReal      wltea,wlter;
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;
7059566063dSJacob Faibussowitsch     PetscCall(VecCopy(X,Y));
7069566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y,3,scal,vecs));
7079566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter));
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 
7181566a47fSLisandro Dalcin   PetscFunctionBegin;
7199566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0,ts->vec_sol));
720cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
7219566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->VecCostIntegral0,quadts->vec_sol));
7221566a47fSLisandro Dalcin   }
723715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
724*1baa6e33SBarry Smith   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN));
7257207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7269566063dSJacob Faibussowitsch     PetscCall(MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN));
7276f92202bSHong Zhang   }
728715f1b00SHong Zhang   PetscFunctionReturn(0);
729715f1b00SHong Zhang }
730715f1b00SHong Zhang 
731715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
732715f1b00SHong Zhang {
733715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
734cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
73513b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
736b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7377207e0fdSHong Zhang   PetscInt       ntlm;
738715f1b00SHong Zhang   KSP            ksp;
7397207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
74013b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
7411a352fa8SHong Zhang   PetscReal      previous_shift;
742715f1b00SHong Zhang 
743715f1b00SHong Zhang   PetscFunctionBegin;
7441a352fa8SHong Zhang   previous_shift = th->shift;
7459566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN));
74613b1b0ffSHong Zhang 
7477207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7489566063dSJacob Faibussowitsch     PetscCall(MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN));
7496f92202bSHong Zhang   }
7509566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes,&ksp));
7519566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,&Jpre,NULL,NULL));
7527207e0fdSHong Zhang   if (quadts) {
7539566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL));
7549566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL));
7557207e0fdSHong Zhang   }
756715f1b00SHong Zhang 
757715f1b00SHong Zhang   /* Build RHS */
758715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7591a352fa8SHong Zhang     th->shift = 1./((th->Theta-1.)*th->time_step0);
7609566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->ptime0,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
7619566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip));
7629566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta));
763715f1b00SHong Zhang 
764022c081aSHong Zhang     /* Add the f_p forcing terms */
7650e11c21fSHong Zhang     if (ts->Jacp) {
7669566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7679566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->ptime0,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
7689566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN));
76982ad9101SHong Zhang       th->shift = previous_shift;
7709566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol));
7719566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
7729566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN));
7730e11c21fSHong Zhang     }
774715f1b00SHong Zhang   } else { /* 1-stage method */
7753e05ceb1SHong Zhang     th->shift = 0.0;
7769566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
7779566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip));
7789566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip,-1.));
77913b1b0ffSHong Zhang 
780022c081aSHong Zhang     /* Add the f_p forcing terms */
7810e11c21fSHong Zhang     if (ts->Jacp) {
78282ad9101SHong Zhang       th->shift = previous_shift;
7839566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X));
7849566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
7859566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN));
786715f1b00SHong Zhang     }
7870e11c21fSHong Zhang   }
788715f1b00SHong Zhang 
789715f1b00SHong Zhang   /* Build LHS */
7901a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
791715f1b00SHong Zhang   if (th->endpoint) {
7929566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
793715f1b00SHong Zhang   } else {
7949566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
795715f1b00SHong Zhang   }
7969566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
797715f1b00SHong Zhang 
798715f1b00SHong Zhang   /*
799715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
800c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
801715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
802715f1b00SHong Zhang   */
803715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8047207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8059566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->ptime0,th->X0,quadJ,NULL));
8069566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->ptime0,th->X0,quadJp));
8079566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8089566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8099566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
810715f1b00SHong Zhang     }
811715f1b00SHong Zhang   }
812715f1b00SHong Zhang 
813715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
814022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
815b5b4017aSHong Zhang     KSPConvergedReason kspreason;
8169566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr));
8179566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol,barr));
818715f1b00SHong Zhang     if (th->endpoint) {
8199566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr));
8209566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
8219566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col));
8229566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
8239566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
824715f1b00SHong Zhang     } else {
8259566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol));
826715f1b00SHong Zhang     }
8279566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
828b5b4017aSHong Zhang     if (kspreason < 0) {
829b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
83063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm));
831b5b4017aSHong Zhang     }
8329566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
8339566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip,&barr));
834715f1b00SHong Zhang   }
835715f1b00SHong Zhang 
83613b1b0ffSHong Zhang   /*
83713b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
838c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
83913b1b0ffSHong Zhang   */
8407207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
84113b1b0ffSHong Zhang     if (!th->endpoint) {
8429566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN)); /* stage sensitivity */
8439566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
8449566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp));
8459566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8469566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8479566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
8489566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN));
84913b1b0ffSHong Zhang     } else {
8509566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
8519566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
8529566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8539566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8549566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
85513b1b0ffSHong Zhang     }
85613b1b0ffSHong Zhang   } else {
85713b1b0ffSHong Zhang     if (!th->endpoint) {
8589566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN));
859715f1b00SHong Zhang     }
860715f1b00SHong Zhang   }
8611566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8621566a47fSLisandro Dalcin }
8631566a47fSLisandro Dalcin 
864cc4f23bcSHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat *stagesensip[])
8656affb6f8SHong Zhang {
8666affb6f8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
8676affb6f8SHong Zhang 
8686affb6f8SHong Zhang   PetscFunctionBegin;
8697409eceaSHong Zhang   if (ns) {
8707409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8717409eceaSHong Zhang     else *ns = 2; /* endpoint form */
8727409eceaSHong Zhang   }
8735072d23cSHong Zhang   if (stagesensip) {
874cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8757409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
876cc4f23bcSHong Zhang     } else {
877cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
878cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
879cc4f23bcSHong Zhang     }
880cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8815072d23cSHong Zhang   }
8826affb6f8SHong Zhang   PetscFunctionReturn(0);
8836affb6f8SHong Zhang }
8846affb6f8SHong Zhang 
885316643e7SJed Brown /*------------------------------------------------------------*/
886277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
887316643e7SJed Brown {
888316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
889316643e7SJed Brown 
890316643e7SJed Brown   PetscFunctionBegin;
8919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
8929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
8939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
8949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
8951566a47fSLisandro Dalcin 
8969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
8979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
8981566a47fSLisandro Dalcin 
8999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
900277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
901277b19d0SLisandro Dalcin }
902277b19d0SLisandro Dalcin 
903ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
904ece46509SHong Zhang {
905ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
906ece46509SHong Zhang 
907ece46509SHong Zhang   PetscFunctionBegin;
9089566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaLam));
9099566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaMu));
9109566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2));
9119566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2));
9129566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsSensiTemp));
9139566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp));
914ece46509SHong Zhang   PetscFunctionReturn(0);
915ece46509SHong Zhang }
916ece46509SHong Zhang 
917277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
918277b19d0SLisandro Dalcin {
919277b19d0SLisandro Dalcin   PetscFunctionBegin;
9209566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
921b3a6b972SJed Brown   if (ts->dm) {
9229566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts));
9239566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts));
924b3a6b972SJed Brown   }
9259566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL));
9279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL));
9289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL));
9299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL));
930316643e7SJed Brown   PetscFunctionReturn(0);
931316643e7SJed Brown }
932316643e7SJed Brown 
933316643e7SJed Brown /*
934316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9352b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9360fd10508SMatthew Knepley 
9370fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9380fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9390fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
940316643e7SJed Brown */
9410f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
942316643e7SJed Brown {
943316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
9447445fe48SJed Brown   Vec            X0,Xdot;
9457445fe48SJed Brown   DM             dm,dmsave;
9463e05ceb1SHong Zhang   PetscReal      shift = th->shift;
947316643e7SJed Brown 
948316643e7SJed Brown   PetscFunctionBegin;
9499566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
9505a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9519566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot));
952494ce9fcSHong Zhang   if (x != X0) {
9539566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x));
954494ce9fcSHong Zhang   } else {
9559566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
956494ce9fcSHong Zhang   }
9577445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9587445fe48SJed Brown   dmsave = ts->dm;
9597445fe48SJed Brown   ts->dm = dm;
9609566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE));
9617445fe48SJed Brown   ts->dm = dmsave;
9629566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot));
963316643e7SJed Brown   PetscFunctionReturn(0);
964316643e7SJed Brown }
965316643e7SJed Brown 
966d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
967316643e7SJed Brown {
968316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
9697445fe48SJed Brown   Vec            Xdot;
9707445fe48SJed Brown   DM             dm,dmsave;
9713e05ceb1SHong Zhang   PetscReal      shift = th->shift;
972316643e7SJed Brown 
973316643e7SJed Brown   PetscFunctionBegin;
9749566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
975be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9769566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot));
9777445fe48SJed Brown 
9787445fe48SJed Brown   dmsave = ts->dm;
9797445fe48SJed Brown   ts->dm = dm;
9809566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE));
9817445fe48SJed Brown   ts->dm = dmsave;
9829566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot));
983316643e7SJed Brown   PetscFunctionReturn(0);
984316643e7SJed Brown }
985316643e7SJed Brown 
986715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
987715f1b00SHong Zhang {
988715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9897207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
990715f1b00SHong Zhang 
991715f1b00SHong Zhang   PetscFunctionBegin;
992022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
99313b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
9949566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip));
9957207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9969566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp));
9979566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0));
998715f1b00SHong Zhang   }
9996f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
10009566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0));
100113b1b0ffSHong Zhang 
10029566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol));
1003715f1b00SHong Zhang   PetscFunctionReturn(0);
1004715f1b00SHong Zhang }
1005715f1b00SHong Zhang 
10067adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
10077adebddeSHong Zhang {
10087adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10097207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10107adebddeSHong Zhang 
10117adebddeSHong Zhang   PetscFunctionBegin;
10127207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
10149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
10157adebddeSHong Zhang   }
10169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
10179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
10189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
10197adebddeSHong Zhang   PetscFunctionReturn(0);
10207adebddeSHong Zhang }
10217adebddeSHong Zhang 
1022316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1023316643e7SJed Brown {
1024316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1025cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10262ffb9264SLisandro Dalcin   PetscBool      match;
1027316643e7SJed Brown 
1028316643e7SJed Brown   PetscFunctionBegin;
1029cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
10309566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0));
1031b1cb36f3SHong Zhang   }
103239d13666SHong Zhang   if (!th->X) {
10339566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->X));
103439d13666SHong Zhang   }
103539d13666SHong Zhang   if (!th->Xdot) {
10369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->Xdot));
103739d13666SHong Zhang   }
103839d13666SHong Zhang   if (!th->X0) {
10399566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->X0));
104039d13666SHong Zhang   }
10411566a47fSLisandro Dalcin   if (th->endpoint) {
10429566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->affine));
10437445fe48SJed Brown   }
10443b1890cdSShri Abhyankar 
10451566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
104620d49056SMatthew G. Knepley   th->shift = 1/(th->Theta*ts->time_step);
10471566a47fSLisandro Dalcin 
10489566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&ts->dm));
10499566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts));
10509566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts));
10511566a47fSLisandro Dalcin 
10529566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&ts->adapt));
10539566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match));
10552ffb9264SLisandro Dalcin   if (!match) {
10569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->vec_sol_prev));
10579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->vec_lte_work));
10583b1890cdSShri Abhyankar   }
10599566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&ts->snes));
1060cc4f23bcSHong Zhang 
1061cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1062b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1063b4dd3bf9SBarry Smith }
10640c7376e5SHong Zhang 
1065b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1066b4dd3bf9SBarry Smith 
1067b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1068b4dd3bf9SBarry Smith {
1069b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1070b4dd3bf9SBarry Smith 
1071b4dd3bf9SBarry Smith   PetscFunctionBegin;
10729566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam));
10739566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp));
10742ca6e920SHong Zhang   if (ts->vecs_sensip) {
10759566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu));
10762ca6e920SHong Zhang   }
1077b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10789566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2));
10799566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp));
108067633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
108167633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
108267633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1083b5b4017aSHong Zhang   }
1084c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
10859566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2));
108667633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
108767633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
108867633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1089b5b4017aSHong Zhang   }
1090316643e7SJed Brown   PetscFunctionReturn(0);
1091316643e7SJed Brown }
1092316643e7SJed Brown 
10934416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1094316643e7SJed Brown {
1095316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1096316643e7SJed Brown 
1097316643e7SJed Brown   PetscFunctionBegin;
1098d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"Theta ODE solver options");
1099316643e7SJed Brown   {
11009566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL));
11019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL));
11029566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL));
1103316643e7SJed Brown   }
1104d0609cedSBarry Smith   PetscOptionsHeadEnd();
1105316643e7SJed Brown   PetscFunctionReturn(0);
1106316643e7SJed Brown }
1107316643e7SJed Brown 
1108316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1109316643e7SJed Brown {
1110316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1111ace3abfcSBarry Smith   PetscBool      iascii;
1112316643e7SJed Brown 
1113316643e7SJed Brown   PetscFunctionBegin;
11149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1115316643e7SJed Brown   if (iascii) {
11169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta));
11179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no"));
1118316643e7SJed Brown   }
1119316643e7SJed Brown   PetscFunctionReturn(0);
1120316643e7SJed Brown }
1121316643e7SJed Brown 
1122560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11230de4c49aSJed Brown {
11240de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11250de4c49aSJed Brown 
11260de4c49aSJed Brown   PetscFunctionBegin;
11270de4c49aSJed Brown   *theta = th->Theta;
11280de4c49aSJed Brown   PetscFunctionReturn(0);
11290de4c49aSJed Brown }
11300de4c49aSJed Brown 
1131560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11320de4c49aSJed Brown {
11330de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11340de4c49aSJed Brown 
11350de4c49aSJed Brown   PetscFunctionBegin;
1136cad9d221SBarry Smith   PetscCheck(theta > 0 && theta <= 1,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11370de4c49aSJed Brown   th->Theta = theta;
11381566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11390de4c49aSJed Brown   PetscFunctionReturn(0);
11400de4c49aSJed Brown }
1141eb284becSJed Brown 
1142560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
114326f2ff8fSLisandro Dalcin {
114426f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
114526f2ff8fSLisandro Dalcin 
114626f2ff8fSLisandro Dalcin   PetscFunctionBegin;
114726f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
114826f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
114926f2ff8fSLisandro Dalcin }
115026f2ff8fSLisandro Dalcin 
1151560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1152eb284becSJed Brown {
1153eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1154eb284becSJed Brown 
1155eb284becSJed Brown   PetscFunctionBegin;
1156eb284becSJed Brown   th->endpoint = flg;
1157eb284becSJed Brown   PetscFunctionReturn(0);
1158eb284becSJed Brown }
11590de4c49aSJed Brown 
1160f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1161f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1162f9c1d6abSBarry Smith {
1163f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1164f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11653fd8ae06SJed Brown   const PetscReal one = 1.0;
1166f9c1d6abSBarry Smith 
1167f9c1d6abSBarry Smith   PetscFunctionBegin;
11683fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1169f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1170f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1171f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1172f9c1d6abSBarry Smith }
1173f9c1d6abSBarry Smith #endif
1174f9c1d6abSBarry Smith 
1175cc4f23bcSHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec *Y[])
117642682096SHong Zhang {
117742682096SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
117842682096SHong Zhang 
117942682096SHong Zhang   PetscFunctionBegin;
11807409eceaSHong Zhang   if (ns) {
11817409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11827409eceaSHong Zhang     else *ns = 2; /* endpoint form */
11837409eceaSHong Zhang   }
11845072d23cSHong Zhang   if (Y) {
1185cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
11867409eceaSHong Zhang       th->Stages[0] = th->X;
1187cc4f23bcSHong Zhang     } else {
1188cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1189cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1190cc4f23bcSHong Zhang     }
1191cc4f23bcSHong Zhang     *Y = th->Stages;
11925072d23cSHong Zhang   }
119342682096SHong Zhang   PetscFunctionReturn(0);
119442682096SHong Zhang }
1195f9c1d6abSBarry Smith 
1196316643e7SJed Brown /* ------------------------------------------------------------ */
1197316643e7SJed Brown /*MC
119896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1199316643e7SJed Brown 
1200316643e7SJed Brown    Level: beginner
1201316643e7SJed Brown 
12024eb428fdSBarry Smith    Options Database:
1203c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1204c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
120503842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
12064eb428fdSBarry Smith 
1207eb284becSJed Brown    Notes:
12080c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
12090c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
12104eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
12114eb428fdSBarry Smith 
12127409eceaSHong Zhang    The endpoint variant of the Theta method and backward Euler can be applied to DAE. The midpoint variant is not suitable for DAEs because it is not stiffly accurate.
1213eb284becSJed Brown 
12147409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1215eb284becSJed Brown 
1216eb284becSJed Brown .vb
1217eb284becSJed Brown   Theta | Theta
1218eb284becSJed Brown   -------------
1219eb284becSJed Brown         |  1
1220eb284becSJed Brown .ve
1221eb284becSJed Brown 
1222eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1223eb284becSJed Brown 
1224eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1225eb284becSJed Brown 
1226eb284becSJed Brown .vb
1227eb284becSJed Brown   0 | 0         0
1228eb284becSJed Brown   1 | 1-Theta   Theta
1229eb284becSJed Brown   -------------------
1230eb284becSJed Brown     | 1-Theta   Theta
1231eb284becSJed Brown .ve
1232eb284becSJed Brown 
1233eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1234eb284becSJed Brown 
1235eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1236eb284becSJed Brown 
1237eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1238eb284becSJed Brown 
12394eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1240eb284becSJed Brown 
1241db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1242316643e7SJed Brown 
1243316643e7SJed Brown M*/
12448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1245316643e7SJed Brown {
1246316643e7SJed Brown   TS_Theta       *th;
1247316643e7SJed Brown 
1248316643e7SJed Brown   PetscFunctionBegin;
1249277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1250ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1251316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1252316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1253316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
125442f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1255ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1256316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1257cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12581566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
125924655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1260316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12610f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12620f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1263f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1264f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1265f9c1d6abSBarry Smith #endif
126642682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
126742f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1268b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1269b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12702ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12716affb6f8SHong Zhang 
1272715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12737adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1274715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12756affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1276316643e7SJed Brown 
1277efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1278efd4aadfSBarry Smith 
12799566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&th));
1280316643e7SJed Brown   ts->data = (void*)th;
1281316643e7SJed Brown 
1282715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1283715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1284715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1285b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1286715f1b00SHong Zhang 
12876f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1288316643e7SJed Brown   th->Theta       = 0.5;
12891566a47fSLisandro Dalcin   th->order       = 2;
12909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta));
12919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta));
12929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta));
12939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta));
1294316643e7SJed Brown   PetscFunctionReturn(0);
1295316643e7SJed Brown }
12960de4c49aSJed Brown 
12970de4c49aSJed Brown /*@
12980de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12990de4c49aSJed Brown 
13000de4c49aSJed Brown   Not Collective
13010de4c49aSJed Brown 
13020de4c49aSJed Brown   Input Parameter:
13030de4c49aSJed Brown .  ts - timestepping context
13040de4c49aSJed Brown 
13050de4c49aSJed Brown   Output Parameter:
13060de4c49aSJed Brown .  theta - stage abscissa
13070de4c49aSJed Brown 
13080de4c49aSJed Brown   Note:
13090de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
13100de4c49aSJed Brown 
13110de4c49aSJed Brown   Level: Advanced
13120de4c49aSJed Brown 
1313db781477SPatrick Sanan .seealso: `TSThetaSetTheta()`
13140de4c49aSJed Brown @*/
13157087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13160de4c49aSJed Brown {
13170de4c49aSJed Brown   PetscFunctionBegin;
1318afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1319dadcf809SJacob Faibussowitsch   PetscValidRealPointer(theta,2);
1320cac4c232SBarry Smith   PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
13210de4c49aSJed Brown   PetscFunctionReturn(0);
13220de4c49aSJed Brown }
13230de4c49aSJed Brown 
13240de4c49aSJed Brown /*@
13250de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13260de4c49aSJed Brown 
13270de4c49aSJed Brown   Not Collective
13280de4c49aSJed Brown 
1329d8d19677SJose E. Roman   Input Parameters:
13300de4c49aSJed Brown +  ts - timestepping context
13310de4c49aSJed Brown -  theta - stage abscissa
13320de4c49aSJed Brown 
13330de4c49aSJed Brown   Options Database:
133467b8a455SSatish Balay .  -ts_theta_theta <theta> - set theta
13350de4c49aSJed Brown 
13360de4c49aSJed Brown   Level: Intermediate
13370de4c49aSJed Brown 
1338db781477SPatrick Sanan .seealso: `TSThetaGetTheta()`
13390de4c49aSJed Brown @*/
13407087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13410de4c49aSJed Brown {
13420de4c49aSJed Brown   PetscFunctionBegin;
1343afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1344cac4c232SBarry Smith   PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
13450de4c49aSJed Brown   PetscFunctionReturn(0);
13460de4c49aSJed Brown }
1347f33bbcb6SJed Brown 
134826f2ff8fSLisandro Dalcin /*@
134926f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
135026f2ff8fSLisandro Dalcin 
135126f2ff8fSLisandro Dalcin   Not Collective
135226f2ff8fSLisandro Dalcin 
135326f2ff8fSLisandro Dalcin   Input Parameter:
135426f2ff8fSLisandro Dalcin .  ts - timestepping context
135526f2ff8fSLisandro Dalcin 
135626f2ff8fSLisandro Dalcin   Output Parameter:
135726f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
135826f2ff8fSLisandro Dalcin 
135926f2ff8fSLisandro Dalcin   Level: Advanced
136026f2ff8fSLisandro Dalcin 
1361db781477SPatrick Sanan .seealso: `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
136226f2ff8fSLisandro Dalcin @*/
136326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
136426f2ff8fSLisandro Dalcin {
136526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
136626f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1367dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(endpoint,2);
1368cac4c232SBarry Smith   PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
136926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
137026f2ff8fSLisandro Dalcin }
137126f2ff8fSLisandro Dalcin 
1372eb284becSJed Brown /*@
1373eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1374eb284becSJed Brown 
1375eb284becSJed Brown   Not Collective
1376eb284becSJed Brown 
1377d8d19677SJose E. Roman   Input Parameters:
1378eb284becSJed Brown +  ts - timestepping context
1379eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1380eb284becSJed Brown 
1381eb284becSJed Brown   Options Database:
138267b8a455SSatish Balay .  -ts_theta_endpoint <flg> - use the endpoint variant
1383eb284becSJed Brown 
1384eb284becSJed Brown   Level: Intermediate
1385eb284becSJed Brown 
1386db781477SPatrick Sanan .seealso: `TSTHETA`, `TSCN`
1387eb284becSJed Brown @*/
1388eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1389eb284becSJed Brown {
1390eb284becSJed Brown   PetscFunctionBegin;
1391eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1392cac4c232SBarry Smith   PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
1393eb284becSJed Brown   PetscFunctionReturn(0);
1394eb284becSJed Brown }
1395eb284becSJed Brown 
1396f33bbcb6SJed Brown /*
1397f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1398f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1399f33bbcb6SJed Brown  */
1400f33bbcb6SJed Brown 
14011566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
14021566a47fSLisandro Dalcin {
14031566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14041566a47fSLisandro Dalcin 
14051566a47fSLisandro Dalcin   PetscFunctionBegin;
140608401ef6SPierre Jolivet   PetscCheck(th->Theta == 1.0,PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler");
140728b400f6SJacob Faibussowitsch   PetscCheck(!th->endpoint,PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler");
14089566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14091566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14101566a47fSLisandro Dalcin }
14111566a47fSLisandro Dalcin 
1412f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1413f33bbcb6SJed Brown {
1414f33bbcb6SJed Brown   PetscFunctionBegin;
1415f33bbcb6SJed Brown   PetscFunctionReturn(0);
1416f33bbcb6SJed Brown }
1417f33bbcb6SJed Brown 
1418f33bbcb6SJed Brown /*MC
1419f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1420f33bbcb6SJed Brown 
1421f33bbcb6SJed Brown   Level: beginner
1422f33bbcb6SJed Brown 
14234eb428fdSBarry Smith   Notes:
1424c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14254eb428fdSBarry Smith 
14261566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14274eb428fdSBarry Smith 
1428db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1429f33bbcb6SJed Brown 
1430f33bbcb6SJed Brown M*/
14318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1432f33bbcb6SJed Brown {
1433f33bbcb6SJed Brown   PetscFunctionBegin;
14349566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14359566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts,1.0));
14369566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts,PETSC_FALSE));
14370c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1438f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1439f33bbcb6SJed Brown   PetscFunctionReturn(0);
1440f33bbcb6SJed Brown }
1441f33bbcb6SJed Brown 
14421566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14431566a47fSLisandro Dalcin {
14441566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14451566a47fSLisandro Dalcin 
14461566a47fSLisandro Dalcin   PetscFunctionBegin;
144708401ef6SPierre Jolivet   PetscCheck(th->Theta == 0.5,PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson");
14483c633725SBarry Smith   PetscCheck(th->endpoint,PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson");
14499566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14501566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14511566a47fSLisandro Dalcin }
14521566a47fSLisandro Dalcin 
1453f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1454f33bbcb6SJed Brown {
1455f33bbcb6SJed Brown   PetscFunctionBegin;
1456f33bbcb6SJed Brown   PetscFunctionReturn(0);
1457f33bbcb6SJed Brown }
1458f33bbcb6SJed Brown 
1459f33bbcb6SJed Brown /*MC
1460f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1461f33bbcb6SJed Brown 
1462f33bbcb6SJed Brown   Level: beginner
1463f33bbcb6SJed Brown 
1464f33bbcb6SJed Brown   Notes:
14657cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14667cf5af47SJed Brown 
14677cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1468f33bbcb6SJed Brown 
1469db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`
1470f33bbcb6SJed Brown 
1471f33bbcb6SJed Brown M*/
14728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1473f33bbcb6SJed Brown {
1474f33bbcb6SJed Brown   PetscFunctionBegin;
14759566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14769566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts,0.5));
14779566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts,PETSC_TRUE));
14780c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1479f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1480f33bbcb6SJed Brown   PetscFunctionReturn(0);
1481f33bbcb6SJed Brown }
1482