xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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;
2479566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D 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) */
2867207e0fdSHong Zhang   if (quadts) {
2879566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
2887207e0fdSHong Zhang   }
289cd4cee2dSHong Zhang 
290c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
2919566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]));
2929566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj],1./adjoint_time_step)); /* lambda_{n+1}/h */
293cd4cee2dSHong Zhang     if (quadJ) {
2949566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
2959566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
2969566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col));
2979566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
2989566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
299c9681e5dSHong Zhang     }
300c9681e5dSHong Zhang   }
301c9681e5dSHong Zhang 
302c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
303c87ba875SIvan Yashchuk   th->shift = 1./adjoint_time_step;
3049566063dSJacob Faibussowitsch   PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
3059566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
306c9681e5dSHong Zhang 
307c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
308c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
309c9681e5dSHong Zhang     KSPConvergedReason kspreason;
3109566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]));
3119566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
312c9681e5dSHong Zhang     if (kspreason < 0) {
313c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
3149566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj));
315c9681e5dSHong Zhang     }
316c9681e5dSHong Zhang   }
317c9681e5dSHong Zhang 
318c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
319c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
3209566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
3219566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
322c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
3239566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
324c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
3259566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
326c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
3279566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
3289566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj],1./adjoint_time_step));
3299566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]));
330c9681e5dSHong Zhang       if (ts->vecs_fup) {
3319566063dSJacob Faibussowitsch         PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]));
332c9681e5dSHong Zhang       }
333c9681e5dSHong Zhang     }
334c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
335c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
33605c9950eSHong Zhang       KSPConvergedReason kspreason;
3379566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]));
3389566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp,&kspreason));
33905c9950eSHong Zhang       if (kspreason < 0) {
34005c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
3419566063dSJacob Faibussowitsch         PetscCall(PetscInfo(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj));
34205c9950eSHong Zhang       }
343c9681e5dSHong Zhang     }
344c9681e5dSHong Zhang   }
345c9681e5dSHong Zhang 
346c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
347077c4feaSHong Zhang   if (!isexplicitode) {
3483e05ceb1SHong Zhang     th->shift = 0.0;
3499566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
3509566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
351c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
35233f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
3539566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj],-1.));
3549566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]));
3559566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensiTemp[nadj],-adjoint_time_step));
3569566063dSJacob Faibussowitsch       PetscCall(VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]));
357c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
3589566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]));
3599566063dSJacob Faibussowitsch         PetscCall(VecScale(VecsSensi2Temp[nadj],-adjoint_time_step));
3609566063dSJacob Faibussowitsch         PetscCall(VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]));
361c9681e5dSHong Zhang       }
362c9681e5dSHong Zhang     }
363077c4feaSHong Zhang   }
364c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3659566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->Xdot,-1./adjoint_time_step,1.0/adjoint_time_step,0,th->X0,ts->vec_sol));
3669566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,1./adjoint_time_step,ts->Jacp,PETSC_FALSE)); /* get -f_p */
3677207e0fdSHong Zhang     if (quadts) {
3689566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
3697207e0fdSHong Zhang     }
370c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
371c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
3729566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
37333f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
3749566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
375c9681e5dSHong Zhang     }
376cd4cee2dSHong Zhang 
377c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
3789566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
3799566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]));
380cd4cee2dSHong Zhang       if (quadJp) {
3819566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
3829566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
3839566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col));
3849566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdp_col));
3859566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
38633f72c5dSHong Zhang       }
387c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
3889566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
3899566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,VecsDeltaMu2[nadj]));
390c9681e5dSHong Zhang         if (ts->vecs_fpu) {
3919566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpu[nadj]));
392c9681e5dSHong Zhang         }
393c9681e5dSHong Zhang         if (ts->vecs_fpp) {
3949566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpp[nadj]));
395c9681e5dSHong Zhang         }
396c9681e5dSHong Zhang       }
397c9681e5dSHong Zhang     }
398c9681e5dSHong Zhang   }
399c9681e5dSHong Zhang 
400c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
4019566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4029566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
403c9681e5dSHong Zhang   }
404c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
405c9681e5dSHong Zhang   PetscFunctionReturn(0);
406c9681e5dSHong Zhang }
407c9681e5dSHong Zhang 
40842f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4092ca6e920SHong Zhang {
4102ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
411cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
412b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
413b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4142ca6e920SHong Zhang   PetscInt       nadj;
4157207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4162ca6e920SHong Zhang   KSP            ksp;
417b5b4017aSHong Zhang   PetscScalar    *xarr;
4181a352fa8SHong Zhang   PetscReal      adjoint_time_step;
4191a352fa8SHong 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) */
4202ca6e920SHong Zhang 
4212ca6e920SHong Zhang   PetscFunctionBegin;
422c9681e5dSHong Zhang   if (th->Theta == 1.) {
4239566063dSJacob Faibussowitsch     PetscCall(TSAdjointStepBEuler_Private(ts));
424c9681e5dSHong Zhang     PetscFunctionReturn(0);
425c9681e5dSHong Zhang   }
4262ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
4279566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes,&ksp));
4289566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,&Jpre,NULL,NULL));
4297207e0fdSHong Zhang   if (quadts) {
4309566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL));
4319566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL));
4327207e0fdSHong Zhang   }
433b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4341a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step);
4351a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4361a352fa8SHong Zhang   adjoint_time_step = -ts->time_step;  /* always positive since time_step is negative */
4372ca6e920SHong Zhang 
43882ad9101SHong Zhang   if (!th->endpoint) {
4395072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4409566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(th->X0,1.0/(1.0-th->Theta),th->Theta/(th->Theta-1.0),0,th->X,ts->vec_sol));
44182ad9101SHong Zhang   }
44282ad9101SHong Zhang 
443b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
44433f72c5dSHong Zhang   /* Cost function has an integral term */
4457207e0fdSHong Zhang   if (quadts) {
44605755b9cSHong Zhang     if (th->endpoint) {
4479566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
44805755b9cSHong Zhang     } else {
4499566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
45005755b9cSHong Zhang     }
4517207e0fdSHong Zhang   }
45233f72c5dSHong Zhang 
453abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4549566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]));
4559566063dSJacob Faibussowitsch     PetscCall(VecScale(VecsSensiTemp[nadj],1./(th->Theta*adjoint_time_step)));
456cd4cee2dSHong Zhang     if (quadJ) {
4579566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
4589566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
4599566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col));
4609566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_drdu_col));
4619566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
46236eaed60SHong Zhang     }
4632ca6e920SHong Zhang   }
4643c54f92cSHong Zhang 
465b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4661a352fa8SHong Zhang   th->shift = 1./(th->Theta*adjoint_time_step);
4673c54f92cSHong Zhang   if (th->endpoint) {
4689566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre));
4693c54f92cSHong Zhang   } else {
4709566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X,J,Jpre));
4713c54f92cSHong Zhang   }
4729566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
4732ca6e920SHong Zhang 
474b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
475abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4761d14d03bSHong Zhang     KSPConvergedReason kspreason;
4779566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]));
4789566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
4791d14d03bSHong Zhang     if (kspreason < 0) {
4801d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
4819566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj));
4821d14d03bSHong Zhang     }
4832ca6e920SHong Zhang   }
4843c54f92cSHong Zhang 
4851d14d03bSHong Zhang   /* Second-order adjoint */
486b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
4873c633725SBarry Smith     PetscCheck(th->endpoint,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
488b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
4899566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
4909566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
491b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
4929566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
4939566063dSJacob Faibussowitsch     PetscCall(VecResetArray(ts->vec_sensip_col));
4949566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
495b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
4969566063dSJacob Faibussowitsch     PetscCall(TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
497b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
4989566063dSJacob Faibussowitsch       PetscCall(VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]));
4999566063dSJacob Faibussowitsch       PetscCall(VecScale(VecsSensi2Temp[nadj],th->shift));
5009566063dSJacob Faibussowitsch       PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]));
501b5b4017aSHong Zhang       if (ts->vecs_fup) {
5029566063dSJacob Faibussowitsch         PetscCall(VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]));
503b5b4017aSHong Zhang       }
504b5b4017aSHong Zhang     }
505b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
506b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
5071d14d03bSHong Zhang       KSPConvergedReason kspreason;
5089566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]));
5099566063dSJacob Faibussowitsch       PetscCall(KSPGetConvergedReason(ksp,&kspreason));
5101d14d03bSHong Zhang       if (kspreason < 0) {
5111d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
5129566063dSJacob Faibussowitsch         PetscCall(PetscInfo(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj));
5131d14d03bSHong Zhang       }
514b5b4017aSHong Zhang     }
515b5b4017aSHong Zhang   }
516b5b4017aSHong Zhang 
51736eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5181d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5191a352fa8SHong Zhang     th->shift      = 1./((th->Theta-1.)*adjoint_time_step);
5201a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
5219566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X0,J,Jpre));
5229566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
52333f72c5dSHong Zhang     /* R_U at t_n */
5247207e0fdSHong Zhang     if (quadts) {
5259566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,adjoint_ptime,th->X0,quadJ,NULL));
5267207e0fdSHong Zhang     }
527abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
5289566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]));
529cd4cee2dSHong Zhang       if (quadJ) {
5309566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
5319566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
5329566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col));
5339566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
5349566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
535b5b4017aSHong Zhang       }
5369566063dSJacob Faibussowitsch       PetscCall(VecScale(ts->vecs_sensi[nadj],1./th->shift));
537b5b4017aSHong Zhang     }
5381d14d03bSHong Zhang 
5391d14d03bSHong Zhang     /* Second-order adjoint */
5401d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
541b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
5429566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(th->MatFwdSensip0,0,&xarr));
5439566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
544b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5459566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu));
5469566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
5479566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0,&xarr));
548b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5499566063dSJacob Faibussowitsch       PetscCall(TSComputeIHessianProductFunctionUP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup));
550b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
55133f72c5dSHong 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  */
5529566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]));
5539566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]));
554b5b4017aSHong Zhang         if (ts->vecs_fup) {
5559566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]));
556b5b4017aSHong Zhang         }
5579566063dSJacob Faibussowitsch         PetscCall(VecScale(ts->vecs_sensi2[nadj],1./th->shift));
558b5b4017aSHong Zhang       }
559b5e0532dSHong Zhang     }
5603fd52205SHong Zhang 
561c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
562c586f2e8SHong Zhang 
5633fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
564b5b4017aSHong Zhang       /* U_{n+1} */
56582ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
5669566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol));
5679566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,-1./(th->Theta*adjoint_time_step),ts->Jacp,PETSC_FALSE));
5687207e0fdSHong Zhang       if (quadts) {
5699566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
5707207e0fdSHong Zhang       }
571abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
5729566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
5739566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu[nadj]));
5743b711c3fSHong Zhang         if (quadJp) {
5759566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
5769566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
5779566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*th->Theta,ts->vec_drdp_col));
5789566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
5799566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
5803b711c3fSHong Zhang         }
5813fd52205SHong Zhang       }
58233f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
58333f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
5849566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(ts->mat_sensip,0,&xarr));
5859566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
586b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
5879566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
5889566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_sensip_col));
5899566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
59033f72c5dSHong Zhang 
591b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
5929566063dSJacob Faibussowitsch         PetscCall(TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
593b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
59433f72c5dSHong 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)  */
5959566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
5969566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu2[nadj]));
597b5b4017aSHong Zhang           if (ts->vecs_fpu) {
5989566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpu[nadj]));
599b5b4017aSHong Zhang           }
600b5b4017aSHong Zhang           if (ts->vecs_fpp) {
6019566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpp[nadj]));
602b5b4017aSHong Zhang           }
603b5b4017aSHong Zhang         }
604b5b4017aSHong Zhang       }
605b5b4017aSHong Zhang 
606b5b4017aSHong Zhang       /* U_s */
6079566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
6089566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,adjoint_ptime,th->X0,th->Xdot,1./((th->Theta-1.0)*adjoint_time_step),ts->Jacp,PETSC_FALSE));
6097207e0fdSHong Zhang       if (quadts) {
6109566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSJacobianP(quadts,adjoint_ptime,th->X0,quadJp));
6117207e0fdSHong Zhang       }
612abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
6139566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
6149566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu[nadj]));
6153b711c3fSHong Zhang         if (quadJp) {
6169566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
6179566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
6189566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*(1.0-th->Theta),ts->vec_drdp_col));
6199566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6209566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
6213b711c3fSHong Zhang         }
62233f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
62333f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
6249566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(th->MatFwdSensip0,0,&xarr));
6259566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
626b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6279566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu));
6289566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_sensip_col));
6299566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0,&xarr));
630b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6319566063dSJacob Faibussowitsch           PetscCall(TSComputeIHessianProductFunctionPP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp));
632b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
63333f72c5dSHong 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) */
6349566063dSJacob Faibussowitsch             PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]));
6359566063dSJacob Faibussowitsch             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]));
636b5b4017aSHong Zhang             if (ts->vecs_fpu) {
6379566063dSJacob Faibussowitsch               PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]));
638b5e0532dSHong Zhang             }
639b5b4017aSHong Zhang             if (ts->vecs_fpp) {
6409566063dSJacob Faibussowitsch               PetscCall(VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]));
641b5b4017aSHong Zhang             }
64236eaed60SHong Zhang           }
64336eaed60SHong Zhang         }
6443fd52205SHong Zhang       }
645b5e0532dSHong Zhang     }
6463fd52205SHong Zhang   } else { /* one-stage case */
6473e05ceb1SHong Zhang     th->shift = 0.0;
6489566063dSJacob Faibussowitsch     PetscCall(TSComputeSNESJacobian(ts,th->X,J,Jpre)); /* get -f_y */
6499566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp,J,Jpre));
6507207e0fdSHong Zhang     if (quadts) {
6519566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
6527207e0fdSHong Zhang     }
653abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
6549566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]));
6559566063dSJacob Faibussowitsch       PetscCall(VecAXPY(ts->vecs_sensi[nadj],-adjoint_time_step,VecsSensiTemp[nadj]));
656cd4cee2dSHong Zhang       if (quadJ) {
6579566063dSJacob Faibussowitsch         PetscCall(MatDenseGetColumn(quadJ,nadj,&xarr));
6589566063dSJacob Faibussowitsch         PetscCall(VecPlaceArray(ts->vec_drdu_col,xarr));
6599566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensi[nadj],adjoint_time_step,ts->vec_drdu_col));
6609566063dSJacob Faibussowitsch         PetscCall(VecResetArray(ts->vec_drdu_col));
6619566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreColumn(quadJ,&xarr));
66236eaed60SHong Zhang       }
6632ca6e920SHong Zhang     }
6643fd52205SHong Zhang     if (ts->vecs_sensip) {
66582ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
6669566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X));
6679566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
6687207e0fdSHong Zhang       if (quadts) {
6699566063dSJacob Faibussowitsch         PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp));
6707207e0fdSHong Zhang       }
671abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
6729566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]));
6739566063dSJacob Faibussowitsch         PetscCall(VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]));
674cd4cee2dSHong Zhang         if (quadJp) {
6759566063dSJacob Faibussowitsch           PetscCall(MatDenseGetColumn(quadJp,nadj,&xarr));
6769566063dSJacob Faibussowitsch           PetscCall(VecPlaceArray(ts->vec_drdp_col,xarr));
6779566063dSJacob Faibussowitsch           PetscCall(VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col));
6789566063dSJacob Faibussowitsch           PetscCall(VecResetArray(ts->vec_drdp_col));
6799566063dSJacob Faibussowitsch           PetscCall(MatDenseRestoreColumn(quadJp,&xarr));
68036eaed60SHong Zhang         }
68136eaed60SHong Zhang       }
6823fd52205SHong Zhang     }
6832ca6e920SHong Zhang   }
6842ca6e920SHong Zhang 
6852ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6862ca6e920SHong Zhang   PetscFunctionReturn(0);
6872ca6e920SHong Zhang }
6882ca6e920SHong Zhang 
689cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
690cd652676SJed Brown {
691cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
692be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
693cd652676SJed Brown 
694cd652676SJed Brown   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,th->X));
696be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
6979566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(X,dt,th->Xdot,th->X));
698cd652676SJed Brown   PetscFunctionReturn(0);
699cd652676SJed Brown }
700cd652676SJed Brown 
7011566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
7021566a47fSLisandro Dalcin {
7031566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7041566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
7051566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
7067453f775SEmil Constantinescu   PetscReal      wltea,wlter;
7071566a47fSLisandro Dalcin 
7081566a47fSLisandro Dalcin   PetscFunctionBegin;
7092ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
7101566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
711fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
7121566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
713fecfb714SLisandro Dalcin   {
714be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
715be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
7161566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
7171566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
7181566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
7199566063dSJacob Faibussowitsch     PetscCall(VecCopy(X,Y));
7209566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y,3,scal,vecs));
7219566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter));
7221566a47fSLisandro Dalcin   }
7231566a47fSLisandro Dalcin   if (order) *order = 2;
7241566a47fSLisandro Dalcin   PetscFunctionReturn(0);
7251566a47fSLisandro Dalcin }
7261566a47fSLisandro Dalcin 
7271566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7281566a47fSLisandro Dalcin {
7291566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7307207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7311566a47fSLisandro Dalcin 
7321566a47fSLisandro Dalcin   PetscFunctionBegin;
7339566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0,ts->vec_sol));
734cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
7359566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->VecCostIntegral0,quadts->vec_sol));
7361566a47fSLisandro Dalcin   }
737715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
73813b1b0ffSHong Zhang   if (ts->mat_sensip) {
7399566063dSJacob Faibussowitsch     PetscCall(MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN));
7406f92202bSHong Zhang   }
7417207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7429566063dSJacob Faibussowitsch     PetscCall(MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN));
7436f92202bSHong Zhang   }
744715f1b00SHong Zhang   PetscFunctionReturn(0);
745715f1b00SHong Zhang }
746715f1b00SHong Zhang 
747715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
748715f1b00SHong Zhang {
749715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
750cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
75113b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
752b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7537207e0fdSHong Zhang   PetscInt       ntlm;
754715f1b00SHong Zhang   KSP            ksp;
7557207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
75613b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
7571a352fa8SHong Zhang   PetscReal      previous_shift;
758715f1b00SHong Zhang 
759715f1b00SHong Zhang   PetscFunctionBegin;
7601a352fa8SHong Zhang   previous_shift = th->shift;
7619566063dSJacob Faibussowitsch   PetscCall(MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN));
76213b1b0ffSHong Zhang 
7637207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7649566063dSJacob Faibussowitsch     PetscCall(MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN));
7656f92202bSHong Zhang   }
7669566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(ts->snes,&ksp));
7679566063dSJacob Faibussowitsch   PetscCall(TSGetIJacobian(ts,&J,&Jpre,NULL,NULL));
7687207e0fdSHong Zhang   if (quadts) {
7699566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL));
7709566063dSJacob Faibussowitsch     PetscCall(TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL));
7717207e0fdSHong Zhang   }
772715f1b00SHong Zhang 
773715f1b00SHong Zhang   /* Build RHS */
774715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7751a352fa8SHong Zhang     th->shift = 1./((th->Theta-1.)*th->time_step0);
7769566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->ptime0,th->X0,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,(th->Theta-1.)/th->Theta));
779715f1b00SHong Zhang 
780022c081aSHong Zhang     /* Add the f_p forcing terms */
7810e11c21fSHong Zhang     if (ts->Jacp) {
7829566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(th->Xdot));
7839566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->ptime0,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
7849566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN));
78582ad9101SHong Zhang       th->shift = previous_shift;
7869566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol));
7879566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
7889566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN));
7890e11c21fSHong Zhang     }
790715f1b00SHong Zhang   } else { /* 1-stage method */
7913e05ceb1SHong Zhang     th->shift = 0.0;
7929566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
7939566063dSJacob Faibussowitsch     PetscCall(MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip));
7949566063dSJacob Faibussowitsch     PetscCall(MatScale(MatDeltaFwdSensip,-1.));
79513b1b0ffSHong Zhang 
796022c081aSHong Zhang     /* Add the f_p forcing terms */
7970e11c21fSHong Zhang     if (ts->Jacp) {
79882ad9101SHong Zhang       th->shift = previous_shift;
7999566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X));
8009566063dSJacob Faibussowitsch       PetscCall(TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE));
8019566063dSJacob Faibussowitsch       PetscCall(MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN));
802715f1b00SHong Zhang     }
8030e11c21fSHong Zhang   }
804715f1b00SHong Zhang 
805715f1b00SHong Zhang   /* Build LHS */
8061a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
807715f1b00SHong Zhang   if (th->endpoint) {
8089566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
809715f1b00SHong Zhang   } else {
8109566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE));
811715f1b00SHong Zhang   }
8129566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp,J,Jpre));
813715f1b00SHong Zhang 
814715f1b00SHong Zhang   /*
815715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
816c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
817715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
818715f1b00SHong Zhang   */
819715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8207207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8219566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->ptime0,th->X0,quadJ,NULL));
8229566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->ptime0,th->X0,quadJp));
8239566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8249566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8259566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
826715f1b00SHong Zhang     }
827715f1b00SHong Zhang   }
828715f1b00SHong Zhang 
829715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
830022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
831b5b4017aSHong Zhang     KSPConvergedReason kspreason;
8329566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr));
8339566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol,barr));
834715f1b00SHong Zhang     if (th->endpoint) {
8359566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr));
8369566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ts->vec_sensip_col,xarr));
8379566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col));
8389566063dSJacob Faibussowitsch       PetscCall(VecResetArray(ts->vec_sensip_col));
8399566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumn(ts->mat_sensip,&xarr));
840715f1b00SHong Zhang     } else {
8419566063dSJacob Faibussowitsch       PetscCall(KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol));
842715f1b00SHong Zhang     }
8439566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(ksp,&kspreason));
844b5b4017aSHong Zhang     if (kspreason < 0) {
845b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
8469566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm));
847b5b4017aSHong Zhang     }
8489566063dSJacob Faibussowitsch     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
8499566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip,&barr));
850715f1b00SHong Zhang   }
851715f1b00SHong Zhang 
85213b1b0ffSHong Zhang   /*
85313b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
854c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
85513b1b0ffSHong Zhang   */
8567207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
85713b1b0ffSHong Zhang     if (!th->endpoint) {
8589566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN)); /* stage sensitivity */
8599566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL));
8609566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp));
8619566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8629566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8639566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
8649566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN));
86513b1b0ffSHong Zhang     } else {
8669566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL));
8679566063dSJacob Faibussowitsch       PetscCall(TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp));
8689566063dSJacob Faibussowitsch       PetscCall(MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp));
8699566063dSJacob Faibussowitsch       PetscCall(MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN));
8709566063dSJacob Faibussowitsch       PetscCall(MatAXPY(quadts->mat_sensip,th->time_step0*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN));
87113b1b0ffSHong Zhang     }
87213b1b0ffSHong Zhang   } else {
87313b1b0ffSHong Zhang     if (!th->endpoint) {
8749566063dSJacob Faibussowitsch       PetscCall(MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN));
875715f1b00SHong Zhang     }
876715f1b00SHong Zhang   }
8771566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8781566a47fSLisandro Dalcin }
8791566a47fSLisandro Dalcin 
880cc4f23bcSHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat *stagesensip[])
8816affb6f8SHong Zhang {
8826affb6f8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
8836affb6f8SHong Zhang 
8846affb6f8SHong Zhang   PetscFunctionBegin;
8857409eceaSHong Zhang   if (ns) {
8867409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
8877409eceaSHong Zhang     else *ns = 2; /* endpoint form */
8887409eceaSHong Zhang   }
8895072d23cSHong Zhang   if (stagesensip) {
890cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
8917409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
892cc4f23bcSHong Zhang     } else {
893cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
894cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
895cc4f23bcSHong Zhang     }
896cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
8975072d23cSHong Zhang   }
8986affb6f8SHong Zhang   PetscFunctionReturn(0);
8996affb6f8SHong Zhang }
9006affb6f8SHong Zhang 
901316643e7SJed Brown /*------------------------------------------------------------*/
902277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
903316643e7SJed Brown {
904316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
905316643e7SJed Brown 
906316643e7SJed Brown   PetscFunctionBegin;
9079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X));
9089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xdot));
9099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
9109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->affine));
9111566a47fSLisandro Dalcin 
9129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
9139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
9141566a47fSLisandro Dalcin 
9159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecCostIntegral0));
916277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
917277b19d0SLisandro Dalcin }
918277b19d0SLisandro Dalcin 
919ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
920ece46509SHong Zhang {
921ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
922ece46509SHong Zhang 
923ece46509SHong Zhang   PetscFunctionBegin;
9249566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaLam));
9259566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaMu));
9269566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2));
9279566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2));
9289566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsSensiTemp));
9299566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp));
930ece46509SHong Zhang   PetscFunctionReturn(0);
931ece46509SHong Zhang }
932ece46509SHong Zhang 
933277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
934277b19d0SLisandro Dalcin {
935277b19d0SLisandro Dalcin   PetscFunctionBegin;
9369566063dSJacob Faibussowitsch   PetscCall(TSReset_Theta(ts));
937b3a6b972SJed Brown   if (ts->dm) {
9389566063dSJacob Faibussowitsch     PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts));
9399566063dSJacob Faibussowitsch     PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts));
940b3a6b972SJed Brown   }
9419566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
9429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL));
9439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL));
9449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL));
9459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL));
946316643e7SJed Brown   PetscFunctionReturn(0);
947316643e7SJed Brown }
948316643e7SJed Brown 
949316643e7SJed Brown /*
950316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9512b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9520fd10508SMatthew Knepley 
9530fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9540fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9550fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
956316643e7SJed Brown */
9570f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
958316643e7SJed Brown {
959316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
9607445fe48SJed Brown   Vec            X0,Xdot;
9617445fe48SJed Brown   DM             dm,dmsave;
9623e05ceb1SHong Zhang   PetscReal      shift = th->shift;
963316643e7SJed Brown 
964316643e7SJed Brown   PetscFunctionBegin;
9659566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
9665a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9679566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot));
968494ce9fcSHong Zhang   if (x != X0) {
9699566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x));
970494ce9fcSHong Zhang   } else {
9719566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Xdot));
972494ce9fcSHong Zhang   }
9737445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9747445fe48SJed Brown   dmsave = ts->dm;
9757445fe48SJed Brown   ts->dm = dm;
9769566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE));
9777445fe48SJed Brown   ts->dm = dmsave;
9789566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot));
979316643e7SJed Brown   PetscFunctionReturn(0);
980316643e7SJed Brown }
981316643e7SJed Brown 
982d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
983316643e7SJed Brown {
984316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
9857445fe48SJed Brown   Vec            Xdot;
9867445fe48SJed Brown   DM             dm,dmsave;
9873e05ceb1SHong Zhang   PetscReal      shift = th->shift;
988316643e7SJed Brown 
989316643e7SJed Brown   PetscFunctionBegin;
9909566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
991be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9929566063dSJacob Faibussowitsch   PetscCall(TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot));
9937445fe48SJed Brown 
9947445fe48SJed Brown   dmsave = ts->dm;
9957445fe48SJed Brown   ts->dm = dm;
9969566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE));
9977445fe48SJed Brown   ts->dm = dmsave;
9989566063dSJacob Faibussowitsch   PetscCall(TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot));
999316643e7SJed Brown   PetscFunctionReturn(0);
1000316643e7SJed Brown }
1001316643e7SJed Brown 
1002715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
1003715f1b00SHong Zhang {
1004715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10057207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
1006715f1b00SHong Zhang 
1007715f1b00SHong Zhang   PetscFunctionBegin;
1008022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
100913b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
10109566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip));
10117207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10129566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp));
10139566063dSJacob Faibussowitsch     PetscCall(MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0));
1014715f1b00SHong Zhang   }
10156f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
10169566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0));
101713b1b0ffSHong Zhang 
10189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol));
1019715f1b00SHong Zhang   PetscFunctionReturn(0);
1020715f1b00SHong Zhang }
1021715f1b00SHong Zhang 
10227adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
10237adebddeSHong Zhang {
10247adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10257207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10267adebddeSHong Zhang 
10277adebddeSHong Zhang   PetscFunctionBegin;
10287207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10299566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
10309566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&th->MatIntegralSensip0));
10317adebddeSHong Zhang   }
10329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
10339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
10349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&th->MatFwdSensip0));
10357adebddeSHong Zhang   PetscFunctionReturn(0);
10367adebddeSHong Zhang }
10377adebddeSHong Zhang 
1038316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1039316643e7SJed Brown {
1040316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1041cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10422ffb9264SLisandro Dalcin   PetscBool      match;
1043316643e7SJed Brown 
1044316643e7SJed Brown   PetscFunctionBegin;
1045cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
10469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0));
1047b1cb36f3SHong Zhang   }
104839d13666SHong Zhang   if (!th->X) {
10499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->X));
105039d13666SHong Zhang   }
105139d13666SHong Zhang   if (!th->Xdot) {
10529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->Xdot));
105339d13666SHong Zhang   }
105439d13666SHong Zhang   if (!th->X0) {
10559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->X0));
105639d13666SHong Zhang   }
10571566a47fSLisandro Dalcin   if (th->endpoint) {
10589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->affine));
10597445fe48SJed Brown   }
10603b1890cdSShri Abhyankar 
10611566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
106220d49056SMatthew G. Knepley   th->shift = 1/(th->Theta*ts->time_step);
10631566a47fSLisandro Dalcin 
10649566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&ts->dm));
10659566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts));
10669566063dSJacob Faibussowitsch   PetscCall(DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts));
10671566a47fSLisandro Dalcin 
10689566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&ts->adapt));
10699566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
10709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match));
10712ffb9264SLisandro Dalcin   if (!match) {
10729566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->vec_sol_prev));
10739566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&th->vec_lte_work));
10743b1890cdSShri Abhyankar   }
10759566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&ts->snes));
1076cc4f23bcSHong Zhang 
1077cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1078b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1079b4dd3bf9SBarry Smith }
10800c7376e5SHong Zhang 
1081b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1082b4dd3bf9SBarry Smith 
1083b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1084b4dd3bf9SBarry Smith {
1085b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1086b4dd3bf9SBarry Smith 
1087b4dd3bf9SBarry Smith   PetscFunctionBegin;
10889566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam));
10899566063dSJacob Faibussowitsch   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp));
10902ca6e920SHong Zhang   if (ts->vecs_sensip) {
10919566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu));
10922ca6e920SHong Zhang   }
1093b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
10949566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2));
10959566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp));
109667633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
109767633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
109867633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1099b5b4017aSHong Zhang   }
1100c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
11019566063dSJacob Faibussowitsch     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2));
110267633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
110367633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
110467633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1105b5b4017aSHong Zhang   }
1106316643e7SJed Brown   PetscFunctionReturn(0);
1107316643e7SJed Brown }
1108316643e7SJed Brown 
11094416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1110316643e7SJed Brown {
1111316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1112316643e7SJed Brown 
1113316643e7SJed Brown   PetscFunctionBegin;
11149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options"));
1115316643e7SJed Brown   {
11169566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL));
11179566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL));
11189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL));
1119316643e7SJed Brown   }
11209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
1121316643e7SJed Brown   PetscFunctionReturn(0);
1122316643e7SJed Brown }
1123316643e7SJed Brown 
1124316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1125316643e7SJed Brown {
1126316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1127ace3abfcSBarry Smith   PetscBool      iascii;
1128316643e7SJed Brown 
1129316643e7SJed Brown   PetscFunctionBegin;
11309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1131316643e7SJed Brown   if (iascii) {
11329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta));
11339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no"));
1134316643e7SJed Brown   }
1135316643e7SJed Brown   PetscFunctionReturn(0);
1136316643e7SJed Brown }
1137316643e7SJed Brown 
1138560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11390de4c49aSJed Brown {
11400de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11410de4c49aSJed Brown 
11420de4c49aSJed Brown   PetscFunctionBegin;
11430de4c49aSJed Brown   *theta = th->Theta;
11440de4c49aSJed Brown   PetscFunctionReturn(0);
11450de4c49aSJed Brown }
11460de4c49aSJed Brown 
1147560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11480de4c49aSJed Brown {
11490de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11500de4c49aSJed Brown 
11510de4c49aSJed Brown   PetscFunctionBegin;
11522c71b3e2SJacob Faibussowitsch   PetscCheckFalse(theta <= 0 || 1 < theta,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11530de4c49aSJed Brown   th->Theta = theta;
11541566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11550de4c49aSJed Brown   PetscFunctionReturn(0);
11560de4c49aSJed Brown }
1157eb284becSJed Brown 
1158560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
115926f2ff8fSLisandro Dalcin {
116026f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
116126f2ff8fSLisandro Dalcin 
116226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
116326f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
116426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
116526f2ff8fSLisandro Dalcin }
116626f2ff8fSLisandro Dalcin 
1167560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1168eb284becSJed Brown {
1169eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1170eb284becSJed Brown 
1171eb284becSJed Brown   PetscFunctionBegin;
1172eb284becSJed Brown   th->endpoint = flg;
1173eb284becSJed Brown   PetscFunctionReturn(0);
1174eb284becSJed Brown }
11750de4c49aSJed Brown 
1176f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1177f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1178f9c1d6abSBarry Smith {
1179f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1180f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11813fd8ae06SJed Brown   const PetscReal one = 1.0;
1182f9c1d6abSBarry Smith 
1183f9c1d6abSBarry Smith   PetscFunctionBegin;
11843fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1185f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1186f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1187f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1188f9c1d6abSBarry Smith }
1189f9c1d6abSBarry Smith #endif
1190f9c1d6abSBarry Smith 
1191cc4f23bcSHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec *Y[])
119242682096SHong Zhang {
119342682096SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
119442682096SHong Zhang 
119542682096SHong Zhang   PetscFunctionBegin;
11967409eceaSHong Zhang   if (ns) {
11977409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
11987409eceaSHong Zhang     else *ns = 2; /* endpoint form */
11997409eceaSHong Zhang   }
12005072d23cSHong Zhang   if (Y) {
1201cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
12027409eceaSHong Zhang       th->Stages[0] = th->X;
1203cc4f23bcSHong Zhang     } else {
1204cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1205cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1206cc4f23bcSHong Zhang     }
1207cc4f23bcSHong Zhang     *Y = th->Stages;
12085072d23cSHong Zhang   }
120942682096SHong Zhang   PetscFunctionReturn(0);
121042682096SHong Zhang }
1211f9c1d6abSBarry Smith 
1212316643e7SJed Brown /* ------------------------------------------------------------ */
1213316643e7SJed Brown /*MC
121496f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1215316643e7SJed Brown 
1216316643e7SJed Brown    Level: beginner
1217316643e7SJed Brown 
12184eb428fdSBarry Smith    Options Database:
1219c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1220c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
122103842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
12224eb428fdSBarry Smith 
1223eb284becSJed Brown    Notes:
12240c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
12250c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
12264eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
12274eb428fdSBarry Smith 
12287409eceaSHong 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.
1229eb284becSJed Brown 
12307409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1231eb284becSJed Brown 
1232eb284becSJed Brown .vb
1233eb284becSJed Brown   Theta | Theta
1234eb284becSJed Brown   -------------
1235eb284becSJed Brown         |  1
1236eb284becSJed Brown .ve
1237eb284becSJed Brown 
1238eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1239eb284becSJed Brown 
1240eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1241eb284becSJed Brown 
1242eb284becSJed Brown .vb
1243eb284becSJed Brown   0 | 0         0
1244eb284becSJed Brown   1 | 1-Theta   Theta
1245eb284becSJed Brown   -------------------
1246eb284becSJed Brown     | 1-Theta   Theta
1247eb284becSJed Brown .ve
1248eb284becSJed Brown 
1249eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1250eb284becSJed Brown 
1251eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1252eb284becSJed Brown 
1253eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1254eb284becSJed Brown 
12554eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1256eb284becSJed Brown 
1257eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1258316643e7SJed Brown 
1259316643e7SJed Brown M*/
12608cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1261316643e7SJed Brown {
1262316643e7SJed Brown   TS_Theta       *th;
1263316643e7SJed Brown 
1264316643e7SJed Brown   PetscFunctionBegin;
1265277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1266ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1267316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1268316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1269316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
127042f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1271ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1272316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1273cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12741566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
127524655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1276316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12770f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12780f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1279f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1280f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1281f9c1d6abSBarry Smith #endif
128242682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
128342f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1284b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1285b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12862ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12876affb6f8SHong Zhang 
1288715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12897adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1290715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12916affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1292316643e7SJed Brown 
1293efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1294efd4aadfSBarry Smith 
12959566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&th));
1296316643e7SJed Brown   ts->data = (void*)th;
1297316643e7SJed Brown 
1298715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1299715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1300715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1301b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1302715f1b00SHong Zhang 
13036f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1304316643e7SJed Brown   th->Theta       = 0.5;
13051566a47fSLisandro Dalcin   th->order       = 2;
13069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta));
13079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta));
13089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta));
13099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta));
1310316643e7SJed Brown   PetscFunctionReturn(0);
1311316643e7SJed Brown }
13120de4c49aSJed Brown 
13130de4c49aSJed Brown /*@
13140de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
13150de4c49aSJed Brown 
13160de4c49aSJed Brown   Not Collective
13170de4c49aSJed Brown 
13180de4c49aSJed Brown   Input Parameter:
13190de4c49aSJed Brown .  ts - timestepping context
13200de4c49aSJed Brown 
13210de4c49aSJed Brown   Output Parameter:
13220de4c49aSJed Brown .  theta - stage abscissa
13230de4c49aSJed Brown 
13240de4c49aSJed Brown   Note:
13250de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
13260de4c49aSJed Brown 
13270de4c49aSJed Brown   Level: Advanced
13280de4c49aSJed Brown 
13290de4c49aSJed Brown .seealso: TSThetaSetTheta()
13300de4c49aSJed Brown @*/
13317087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13320de4c49aSJed Brown {
13330de4c49aSJed Brown   PetscFunctionBegin;
1334afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1335dadcf809SJacob Faibussowitsch   PetscValidRealPointer(theta,2);
1336*cac4c232SBarry Smith   PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));
13370de4c49aSJed Brown   PetscFunctionReturn(0);
13380de4c49aSJed Brown }
13390de4c49aSJed Brown 
13400de4c49aSJed Brown /*@
13410de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13420de4c49aSJed Brown 
13430de4c49aSJed Brown   Not Collective
13440de4c49aSJed Brown 
1345d8d19677SJose E. Roman   Input Parameters:
13460de4c49aSJed Brown +  ts - timestepping context
13470de4c49aSJed Brown -  theta - stage abscissa
13480de4c49aSJed Brown 
13490de4c49aSJed Brown   Options Database:
135067b8a455SSatish Balay .  -ts_theta_theta <theta> - set theta
13510de4c49aSJed Brown 
13520de4c49aSJed Brown   Level: Intermediate
13530de4c49aSJed Brown 
13540de4c49aSJed Brown .seealso: TSThetaGetTheta()
13550de4c49aSJed Brown @*/
13567087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13570de4c49aSJed Brown {
13580de4c49aSJed Brown   PetscFunctionBegin;
1359afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1360*cac4c232SBarry Smith   PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));
13610de4c49aSJed Brown   PetscFunctionReturn(0);
13620de4c49aSJed Brown }
1363f33bbcb6SJed Brown 
136426f2ff8fSLisandro Dalcin /*@
136526f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
136626f2ff8fSLisandro Dalcin 
136726f2ff8fSLisandro Dalcin   Not Collective
136826f2ff8fSLisandro Dalcin 
136926f2ff8fSLisandro Dalcin   Input Parameter:
137026f2ff8fSLisandro Dalcin .  ts - timestepping context
137126f2ff8fSLisandro Dalcin 
137226f2ff8fSLisandro Dalcin   Output Parameter:
137326f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
137426f2ff8fSLisandro Dalcin 
137526f2ff8fSLisandro Dalcin   Level: Advanced
137626f2ff8fSLisandro Dalcin 
137726f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
137826f2ff8fSLisandro Dalcin @*/
137926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
138026f2ff8fSLisandro Dalcin {
138126f2ff8fSLisandro Dalcin   PetscFunctionBegin;
138226f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1383dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(endpoint,2);
1384*cac4c232SBarry Smith   PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));
138526f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
138626f2ff8fSLisandro Dalcin }
138726f2ff8fSLisandro Dalcin 
1388eb284becSJed Brown /*@
1389eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1390eb284becSJed Brown 
1391eb284becSJed Brown   Not Collective
1392eb284becSJed Brown 
1393d8d19677SJose E. Roman   Input Parameters:
1394eb284becSJed Brown +  ts - timestepping context
1395eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1396eb284becSJed Brown 
1397eb284becSJed Brown   Options Database:
139867b8a455SSatish Balay .  -ts_theta_endpoint <flg> - use the endpoint variant
1399eb284becSJed Brown 
1400eb284becSJed Brown   Level: Intermediate
1401eb284becSJed Brown 
1402eb284becSJed Brown .seealso: TSTHETA, TSCN
1403eb284becSJed Brown @*/
1404eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1405eb284becSJed Brown {
1406eb284becSJed Brown   PetscFunctionBegin;
1407eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1408*cac4c232SBarry Smith   PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));
1409eb284becSJed Brown   PetscFunctionReturn(0);
1410eb284becSJed Brown }
1411eb284becSJed Brown 
1412f33bbcb6SJed Brown /*
1413f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1414f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1415f33bbcb6SJed Brown  */
1416f33bbcb6SJed Brown 
14171566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
14181566a47fSLisandro Dalcin {
14191566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14201566a47fSLisandro Dalcin 
14211566a47fSLisandro Dalcin   PetscFunctionBegin;
14222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(th->Theta != 1.0,PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler");
142328b400f6SJacob 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");
14249566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14251566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14261566a47fSLisandro Dalcin }
14271566a47fSLisandro Dalcin 
1428f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1429f33bbcb6SJed Brown {
1430f33bbcb6SJed Brown   PetscFunctionBegin;
1431f33bbcb6SJed Brown   PetscFunctionReturn(0);
1432f33bbcb6SJed Brown }
1433f33bbcb6SJed Brown 
1434f33bbcb6SJed Brown /*MC
1435f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1436f33bbcb6SJed Brown 
1437f33bbcb6SJed Brown   Level: beginner
1438f33bbcb6SJed Brown 
14394eb428fdSBarry Smith   Notes:
1440c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14414eb428fdSBarry Smith 
14421566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14434eb428fdSBarry Smith 
1444f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1445f33bbcb6SJed Brown 
1446f33bbcb6SJed Brown M*/
14478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1448f33bbcb6SJed Brown {
1449f33bbcb6SJed Brown   PetscFunctionBegin;
14509566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14519566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts,1.0));
14529566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts,PETSC_FALSE));
14530c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1454f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1455f33bbcb6SJed Brown   PetscFunctionReturn(0);
1456f33bbcb6SJed Brown }
1457f33bbcb6SJed Brown 
14581566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14591566a47fSLisandro Dalcin {
14601566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14611566a47fSLisandro Dalcin 
14621566a47fSLisandro Dalcin   PetscFunctionBegin;
14632c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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");
14643c633725SBarry 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");
14659566063dSJacob Faibussowitsch   PetscCall(TSSetUp_Theta(ts));
14661566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14671566a47fSLisandro Dalcin }
14681566a47fSLisandro Dalcin 
1469f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1470f33bbcb6SJed Brown {
1471f33bbcb6SJed Brown   PetscFunctionBegin;
1472f33bbcb6SJed Brown   PetscFunctionReturn(0);
1473f33bbcb6SJed Brown }
1474f33bbcb6SJed Brown 
1475f33bbcb6SJed Brown /*MC
1476f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1477f33bbcb6SJed Brown 
1478f33bbcb6SJed Brown   Level: beginner
1479f33bbcb6SJed Brown 
1480f33bbcb6SJed Brown   Notes:
14817cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14827cf5af47SJed Brown 
14837cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1484f33bbcb6SJed Brown 
1485f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1486f33bbcb6SJed Brown 
1487f33bbcb6SJed Brown M*/
14888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1489f33bbcb6SJed Brown {
1490f33bbcb6SJed Brown   PetscFunctionBegin;
14919566063dSJacob Faibussowitsch   PetscCall(TSCreate_Theta(ts));
14929566063dSJacob Faibussowitsch   PetscCall(TSThetaSetTheta(ts,0.5));
14939566063dSJacob Faibussowitsch   PetscCall(TSThetaSetEndpoint(ts,PETSC_TRUE));
14940c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1495f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1496f33bbcb6SJed Brown   PetscFunctionReturn(0);
1497f33bbcb6SJed Brown }
1498