xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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   PetscErrorCode ierr;
497445fe48SJed Brown 
507445fe48SJed Brown   PetscFunctionBegin;
517445fe48SJed Brown   if (X0) {
527445fe48SJed Brown     if (dm && dm != ts->dm) {
530d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
547445fe48SJed Brown     } else *X0 = ts->vec_sol;
557445fe48SJed Brown   }
567445fe48SJed Brown   if (Xdot) {
577445fe48SJed Brown     if (dm && dm != ts->dm) {
580d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
597445fe48SJed Brown     } else *Xdot = th->Xdot;
607445fe48SJed Brown   }
617445fe48SJed Brown   PetscFunctionReturn(0);
627445fe48SJed Brown }
637445fe48SJed Brown 
640d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
650d0b770aSPeter Brune {
660d0b770aSPeter Brune   PetscErrorCode ierr;
670d0b770aSPeter Brune 
680d0b770aSPeter Brune   PetscFunctionBegin;
690d0b770aSPeter Brune   if (X0) {
700d0b770aSPeter Brune     if (dm && dm != ts->dm) {
710d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
720d0b770aSPeter Brune     }
730d0b770aSPeter Brune   }
740d0b770aSPeter Brune   if (Xdot) {
750d0b770aSPeter Brune     if (dm && dm != ts->dm) {
760d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
770d0b770aSPeter Brune     }
780d0b770aSPeter Brune   }
790d0b770aSPeter Brune   PetscFunctionReturn(0);
800d0b770aSPeter Brune }
810d0b770aSPeter Brune 
827445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
837445fe48SJed Brown {
847445fe48SJed Brown   PetscFunctionBegin;
857445fe48SJed Brown   PetscFunctionReturn(0);
867445fe48SJed Brown }
877445fe48SJed Brown 
887445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
897445fe48SJed Brown {
907445fe48SJed Brown   TS             ts = (TS)ctx;
917445fe48SJed Brown   PetscErrorCode ierr;
927445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
937445fe48SJed Brown 
947445fe48SJed Brown   PetscFunctionBegin;
957445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
977445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
987445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
997445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
1007445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
1010d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
1020d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1037445fe48SJed Brown   PetscFunctionReturn(0);
1047445fe48SJed Brown }
1057445fe48SJed Brown 
106258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
107258e1594SPeter Brune {
108258e1594SPeter Brune   PetscFunctionBegin;
109258e1594SPeter Brune   PetscFunctionReturn(0);
110258e1594SPeter Brune }
111258e1594SPeter Brune 
112258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
113258e1594SPeter Brune {
114258e1594SPeter Brune   TS             ts = (TS)ctx;
115258e1594SPeter Brune   PetscErrorCode ierr;
116258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
117258e1594SPeter Brune 
118258e1594SPeter Brune   PetscFunctionBegin;
119258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
120258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
121258e1594SPeter Brune 
122258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune 
125258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127258e1594SPeter Brune 
128258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
129258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
130258e1594SPeter Brune   PetscFunctionReturn(0);
131258e1594SPeter Brune }
132258e1594SPeter Brune 
133022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
134022c081aSHong Zhang {
135022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
136cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
137022c081aSHong Zhang   PetscErrorCode ierr;
138022c081aSHong Zhang 
139022c081aSHong Zhang   PetscFunctionBegin;
140022c081aSHong Zhang   if (th->endpoint) {
141022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
142022c081aSHong Zhang     if (th->Theta!=1.0) {
1431a352fa8SHong Zhang       ierr = TSComputeRHSFunction(quadts,th->ptime0,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
1441a352fa8SHong Zhang       ierr = VecAXPY(quadts->vec_sol,th->time_step0*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
145022c081aSHong Zhang     }
146cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
1471a352fa8SHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step0*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
148022c081aSHong Zhang   } else {
149cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
1501a352fa8SHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step0,ts->vec_costintegrand);CHKERRQ(ierr);
151022c081aSHong Zhang   }
152022c081aSHong Zhang   PetscFunctionReturn(0);
153022c081aSHong Zhang }
154022c081aSHong Zhang 
155b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
156b1cb36f3SHong Zhang {
157b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
158cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
159b1cb36f3SHong Zhang   PetscErrorCode ierr;
160b1cb36f3SHong Zhang 
161b1cb36f3SHong Zhang   PetscFunctionBegin;
162b1cb36f3SHong Zhang   /* backup cost integral */
163cd4cee2dSHong Zhang   ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr);
164022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
165b1cb36f3SHong Zhang   PetscFunctionReturn(0);
166b1cb36f3SHong Zhang }
167b1cb36f3SHong Zhang 
168b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
169b1cb36f3SHong Zhang {
1701a352fa8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
171b1cb36f3SHong Zhang   PetscErrorCode ierr;
172b1cb36f3SHong Zhang 
173b1cb36f3SHong Zhang   PetscFunctionBegin;
1741a352fa8SHong Zhang   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
1751a352fa8SHong Zhang   th->ptime0     = ts->ptime + ts->time_step;
1761a352fa8SHong Zhang   th->time_step0 = -ts->time_step;
177022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
17824655328SShri   PetscFunctionReturn(0);
17924655328SShri }
18024655328SShri 
181470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1821566a47fSLisandro Dalcin {
1831566a47fSLisandro Dalcin   PetscInt       nits,lits;
1841566a47fSLisandro Dalcin   PetscErrorCode ierr;
1851566a47fSLisandro Dalcin 
1861566a47fSLisandro Dalcin   PetscFunctionBegin;
1871566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1881566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1891566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1901566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1911566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1921566a47fSLisandro Dalcin }
1931566a47fSLisandro Dalcin 
194193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
195316643e7SJed Brown {
196316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1971566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1984957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1991566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
200051f2191SLisandro Dalcin   PetscErrorCode ierr;
201316643e7SJed Brown 
202316643e7SJed Brown   PetscFunctionBegin;
2031566a47fSLisandro Dalcin   if (!ts->steprollback) {
2042ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
2053b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
2061566a47fSLisandro Dalcin   }
2071566a47fSLisandro Dalcin 
2081566a47fSLisandro Dalcin   th->status     = TS_STEP_INCOMPLETE;
20999260352SHong Zhang   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
2103e05ceb1SHong Zhang     th->shift      = 1/(th->Theta*ts->time_step);
2111566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
212be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
213fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
2143e05ceb1SHong Zhang       ierr = VecAXPY(th->X,1/th->shift,th->Xdot);CHKERRQ(ierr);
215be5899b3SLisandro Dalcin     }
216eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
217eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2181566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2191566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2201566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2211566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2221566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
223eb284becSJed Brown     }
224be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
225470880abSPatrick Sanan     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
226be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2271566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
228fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
229051f2191SLisandro Dalcin 
230051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2311566a47fSLisandro Dalcin     if (th->endpoint) {
2321566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2331566a47fSLisandro Dalcin     } else {
2343e05ceb1SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr);
2351566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2361566a47fSLisandro Dalcin     }
2371566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2381566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2391566a47fSLisandro Dalcin     if (!accept) {
2401566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
241be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
242be5899b3SLisandro Dalcin       goto reject_step;
243051f2191SLisandro Dalcin     }
2443b1890cdSShri Abhyankar 
245715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
2461a352fa8SHong Zhang       th->ptime0     = ts->ptime;
2471a352fa8SHong Zhang       th->time_step0 = ts->time_step;
24817145e75SHong Zhang     }
2492b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
250cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
251051f2191SLisandro Dalcin     break;
252051f2191SLisandro Dalcin 
253051f2191SLisandro Dalcin   reject_step:
254fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2551566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
256051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2571566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2583b1890cdSShri Abhyankar     }
2593b1890cdSShri Abhyankar   }
260316643e7SJed Brown   PetscFunctionReturn(0);
261316643e7SJed Brown }
262316643e7SJed Brown 
263c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
264c9681e5dSHong Zhang {
265c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
266cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
267c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
268c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
269c9681e5dSHong Zhang   PetscInt       nadj;
2707207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
271c9681e5dSHong Zhang   KSP            ksp;
272c9681e5dSHong Zhang   PetscScalar    *xarr;
273077c4feaSHong Zhang   TSEquationType eqtype;
274077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
2751a352fa8SHong Zhang   PetscReal      adjoint_time_step;
276c9681e5dSHong Zhang   PetscErrorCode ierr;
277c9681e5dSHong Zhang 
278c9681e5dSHong Zhang   PetscFunctionBegin;
279077c4feaSHong Zhang   ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr);
280077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
281077c4feaSHong Zhang     isexplicitode  = PETSC_TRUE;
282077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
283077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
284077c4feaSHong Zhang   }
285c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
286c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
287cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
2887207e0fdSHong Zhang   if (quadts) {
289cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
290cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
2917207e0fdSHong Zhang   }
292c9681e5dSHong Zhang 
2931a352fa8SHong Zhang   th->stage_time    = ts->ptime;
2941a352fa8SHong Zhang   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */
295c9681e5dSHong Zhang 
29633f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2977207e0fdSHong Zhang   if (quadts) {
29882ad9101SHong Zhang     ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
2997207e0fdSHong Zhang   }
300cd4cee2dSHong Zhang 
301c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
302c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
3031a352fa8SHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./adjoint_time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
304cd4cee2dSHong Zhang     if (quadJ) {
305cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
306cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
307cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
308cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
309cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
310c9681e5dSHong Zhang     }
311c9681e5dSHong Zhang   }
312c9681e5dSHong Zhang 
313c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
314c87ba875SIvan Yashchuk   th->shift = 1./adjoint_time_step;
31582ad9101SHong Zhang   ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr);
31604f7482eSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
317c9681e5dSHong Zhang 
318c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
319c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
320c9681e5dSHong Zhang     KSPConvergedReason kspreason;
321c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
322c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
323c9681e5dSHong Zhang     if (kspreason < 0) {
324c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
32529b3c8a1SHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
326c9681e5dSHong Zhang     }
327c9681e5dSHong Zhang   }
328c9681e5dSHong Zhang 
329c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
330c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
331c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
332c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
333c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
33482ad9101SHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
335c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
33682ad9101SHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
337c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
338c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
3391a352fa8SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],1./adjoint_time_step);CHKERRQ(ierr);
34033f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
341c9681e5dSHong Zhang       if (ts->vecs_fup) {
34233f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
343c9681e5dSHong Zhang       }
344c9681e5dSHong Zhang     }
345c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
346c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
34705c9950eSHong Zhang       KSPConvergedReason kspreason;
348c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
34905c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
35005c9950eSHong Zhang       if (kspreason < 0) {
35105c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
35229b3c8a1SHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
35305c9950eSHong Zhang       }
354c9681e5dSHong Zhang     }
355c9681e5dSHong Zhang   }
356c9681e5dSHong Zhang 
357c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
358077c4feaSHong Zhang   if (!isexplicitode) {
3593e05ceb1SHong Zhang     th->shift = 0.0;
36082ad9101SHong Zhang     ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr);
36104f7482eSHong Zhang     ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
362c9681e5dSHong Zhang     ierr = MatScale(J,-1.);CHKERRQ(ierr);
363c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
36433f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
365c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
3661a352fa8SHong Zhang       ierr = VecScale(VecsSensiTemp[nadj],adjoint_time_step);CHKERRQ(ierr);
367c9681e5dSHong Zhang       ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
368c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
369c9681e5dSHong Zhang         ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
3701a352fa8SHong Zhang         ierr = VecScale(VecsSensi2Temp[nadj],adjoint_time_step);CHKERRQ(ierr);
371c9681e5dSHong Zhang         ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
372c9681e5dSHong Zhang       }
373c9681e5dSHong Zhang     }
374077c4feaSHong Zhang   }
375c9681e5dSHong Zhang   if (ts->vecs_sensip) {
37682ad9101SHong Zhang     ierr = VecAXPBYPCZ(th->Xdot,-1./adjoint_time_step,1.0/adjoint_time_step,0,th->X0,ts->vec_sol);CHKERRQ(ierr);
37782ad9101SHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,1./adjoint_time_step,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
3787207e0fdSHong Zhang     if (quadts) {
37982ad9101SHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
3807207e0fdSHong Zhang     }
381c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
382c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
38382ad9101SHong Zhang       ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
38433f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
38582ad9101SHong Zhang       ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
386c9681e5dSHong Zhang     }
387cd4cee2dSHong Zhang 
388c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
389c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
3901a352fa8SHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
391cd4cee2dSHong Zhang       if (quadJp) {
392cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
393cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
3941a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr);
395cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
396cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
39733f72c5dSHong Zhang       }
398c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
399c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
4001a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
401c9681e5dSHong Zhang         if (ts->vecs_fpu) {
4021a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
403c9681e5dSHong Zhang         }
404c9681e5dSHong Zhang         if (ts->vecs_fpp) {
4051a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
406c9681e5dSHong Zhang         }
407c9681e5dSHong Zhang       }
408c9681e5dSHong Zhang     }
409c9681e5dSHong Zhang   }
410c9681e5dSHong Zhang 
411c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
412c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
413c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
414c9681e5dSHong Zhang   }
415c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
416c9681e5dSHong Zhang   PetscFunctionReturn(0);
417c9681e5dSHong Zhang }
418c9681e5dSHong Zhang 
41942f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4202ca6e920SHong Zhang {
4212ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
422cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
423b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
424b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4252ca6e920SHong Zhang   PetscInt       nadj;
4267207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4272ca6e920SHong Zhang   KSP            ksp;
428b5b4017aSHong Zhang   PetscScalar    *xarr;
4291a352fa8SHong Zhang   PetscReal      adjoint_time_step;
4301a352fa8SHong 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) */
431b5b4017aSHong Zhang   PetscErrorCode ierr;
4322ca6e920SHong Zhang 
4332ca6e920SHong Zhang   PetscFunctionBegin;
434c9681e5dSHong Zhang   if (th->Theta == 1.) {
435c9681e5dSHong Zhang     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
436c9681e5dSHong Zhang     PetscFunctionReturn(0);
437c9681e5dSHong Zhang   }
4382ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
439302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
440cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
4417207e0fdSHong Zhang   if (quadts) {
442cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
443cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
4447207e0fdSHong Zhang   }
445b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
4461a352fa8SHong Zhang   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step);
4471a352fa8SHong Zhang   adjoint_ptime     = ts->ptime + ts->time_step;
4481a352fa8SHong Zhang   adjoint_time_step = -ts->time_step;  /* always positive since time_step is negative */
4492ca6e920SHong Zhang 
45082ad9101SHong Zhang   if (!th->endpoint) {
4515072d23cSHong Zhang     /* recover th->X0 using vec_sol and the stage value th->X */
4525072d23cSHong Zhang     ierr = VecAXPBYPCZ(th->X0,1.0/(1.0-th->Theta),th->Theta/(th->Theta-1.0),0,th->X,ts->vec_sol);CHKERRQ(ierr);
45382ad9101SHong Zhang   }
45482ad9101SHong Zhang 
455b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
45633f72c5dSHong Zhang   /* Cost function has an integral term */
4577207e0fdSHong Zhang   if (quadts) {
45805755b9cSHong Zhang     if (th->endpoint) {
459cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
46005755b9cSHong Zhang     } else {
461cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
46205755b9cSHong Zhang     }
4637207e0fdSHong Zhang   }
46433f72c5dSHong Zhang 
465abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
466b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
4671a352fa8SHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*adjoint_time_step));CHKERRQ(ierr);
468cd4cee2dSHong Zhang     if (quadJ) {
469cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
470cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
471cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
472cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
473cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
47436eaed60SHong Zhang     }
4752ca6e920SHong Zhang   }
4763c54f92cSHong Zhang 
477b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4781a352fa8SHong Zhang   th->shift = 1./(th->Theta*adjoint_time_step);
4793c54f92cSHong Zhang   if (th->endpoint) {
48004f7482eSHong Zhang     ierr = TSComputeSNESJacobian(ts,ts->vec_sol,J,Jpre);CHKERRQ(ierr);
4813c54f92cSHong Zhang   } else {
48204f7482eSHong Zhang     ierr = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr);
4833c54f92cSHong Zhang   }
484cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
4852ca6e920SHong Zhang 
486b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
487abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4881d14d03bSHong Zhang     KSPConvergedReason kspreason;
489b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
4901d14d03bSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
4911d14d03bSHong Zhang     if (kspreason < 0) {
4921d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
49329b3c8a1SHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
4941d14d03bSHong Zhang     }
4952ca6e920SHong Zhang   }
4963c54f92cSHong Zhang 
4971d14d03bSHong Zhang   /* Second-order adjoint */
498b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
499b5b4017aSHong Zhang     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
500b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
501b5b4017aSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
502b5b4017aSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
503b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
504b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
505b5b4017aSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
506b5b4017aSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
507b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
508b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
509b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
510b5b4017aSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
5113e05ceb1SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],th->shift);CHKERRQ(ierr);
51233f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
513b5b4017aSHong Zhang       if (ts->vecs_fup) {
51433f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
515b5b4017aSHong Zhang       }
516b5b4017aSHong Zhang     }
517b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
518b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
5191d14d03bSHong Zhang       KSPConvergedReason kspreason;
5201d14d03bSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
5211d14d03bSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
5221d14d03bSHong Zhang       if (kspreason < 0) {
5231d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
52429b3c8a1SHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
5251d14d03bSHong Zhang       }
526b5b4017aSHong Zhang     }
527b5b4017aSHong Zhang   }
528b5b4017aSHong Zhang 
52936eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5301d14d03bSHong Zhang   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5311a352fa8SHong Zhang     th->shift      = 1./((th->Theta-1.)*adjoint_time_step);
5321a352fa8SHong Zhang     th->stage_time = adjoint_ptime;
53304f7482eSHong Zhang     ierr           = TSComputeSNESJacobian(ts,th->X0,J,Jpre);CHKERRQ(ierr);
53404f7482eSHong Zhang     ierr           = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
53533f72c5dSHong Zhang     /* R_U at t_n */
5367207e0fdSHong Zhang     if (quadts) {
5371a352fa8SHong Zhang       ierr = TSComputeRHSJacobian(quadts,adjoint_ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
5387207e0fdSHong Zhang     }
539abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
540b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
541cd4cee2dSHong Zhang       if (quadJ) {
542cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
543cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
544cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr);
545cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
546cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
547b5b4017aSHong Zhang       }
5483e05ceb1SHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./th->shift);CHKERRQ(ierr);
549b5b4017aSHong Zhang     }
5501d14d03bSHong Zhang 
5511d14d03bSHong Zhang     /* Second-order adjoint */
5521d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
553b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
554b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
555b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
556b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
5571a352fa8SHong Zhang       ierr = TSComputeIHessianProductFunctionUU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
558b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
55933f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
560b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
5611a352fa8SHong Zhang       ierr = TSComputeIHessianProductFunctionUP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
562b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
56333f72c5dSHong 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  */
564b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
56533f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
566b5b4017aSHong Zhang         if (ts->vecs_fup) {
56733f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
568b5b4017aSHong Zhang         }
5693e05ceb1SHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./th->shift);CHKERRQ(ierr);
570b5b4017aSHong Zhang       }
571b5e0532dSHong Zhang     }
5723fd52205SHong Zhang 
573c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
574c586f2e8SHong Zhang 
5753fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
576b5b4017aSHong Zhang       /* U_{n+1} */
57782ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
57882ad9101SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol);CHKERRQ(ierr);
5791a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,-1./(th->Theta*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5807207e0fdSHong Zhang       if (quadts) {
581cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
5827207e0fdSHong Zhang       }
583abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
584b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
5851a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5863b711c3fSHong Zhang         if (quadJp) {
5873b711c3fSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
5883b711c3fSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
5893b711c3fSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*th->Theta,ts->vec_drdp_col);CHKERRQ(ierr);
5903b711c3fSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
5913b711c3fSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
5923b711c3fSHong Zhang         }
5933fd52205SHong Zhang       }
59433f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
59533f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
59633f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
59733f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
598b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
599b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
60033f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
60133f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
60233f72c5dSHong Zhang 
603b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
604b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
605b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
60633f72c5dSHong 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)  */
607b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
6081a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
609b5b4017aSHong Zhang           if (ts->vecs_fpu) {
6101a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
611b5b4017aSHong Zhang           }
612b5b4017aSHong Zhang           if (ts->vecs_fpp) {
6131a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
614b5b4017aSHong Zhang           }
615b5b4017aSHong Zhang         }
616b5b4017aSHong Zhang       }
617b5b4017aSHong Zhang 
618b5b4017aSHong Zhang       /* U_s */
61982ad9101SHong Zhang       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
6201a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,adjoint_ptime,th->X0,th->Xdot,1./((th->Theta-1.0)*adjoint_time_step),ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6217207e0fdSHong Zhang       if (quadts) {
6221a352fa8SHong Zhang         ierr = TSComputeRHSJacobianP(quadts,adjoint_ptime,th->X0,quadJp);CHKERRQ(ierr);
6237207e0fdSHong Zhang       }
624abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
625b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
6261a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
6273b711c3fSHong Zhang         if (quadJp) {
6283b711c3fSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
6293b711c3fSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
6303b711c3fSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step*(1.0-th->Theta),ts->vec_drdp_col);CHKERRQ(ierr);
6313b711c3fSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
6323b711c3fSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
6333b711c3fSHong Zhang         }
63433f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
63533f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
63633f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
63733f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
638b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
6391a352fa8SHong Zhang           ierr = TSComputeIHessianProductFunctionPU(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
64033f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
64133f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
642b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
6431a352fa8SHong Zhang           ierr = TSComputeIHessianProductFunctionPP(ts,adjoint_ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
644b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
64533f72c5dSHong 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) */
646b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
6471a352fa8SHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
648b5b4017aSHong Zhang             if (ts->vecs_fpu) {
6491a352fa8SHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
650b5e0532dSHong Zhang             }
651b5b4017aSHong Zhang             if (ts->vecs_fpp) {
6521a352fa8SHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-adjoint_time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
653b5b4017aSHong Zhang             }
65436eaed60SHong Zhang           }
65536eaed60SHong Zhang         }
6563fd52205SHong Zhang       }
657b5e0532dSHong Zhang     }
6583fd52205SHong Zhang   } else { /* one-stage case */
6593e05ceb1SHong Zhang     th->shift = 0.0;
66004f7482eSHong Zhang     ierr      = TSComputeSNESJacobian(ts,th->X,J,Jpre);CHKERRQ(ierr); /* get -f_y */
66104f7482eSHong Zhang     ierr      = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
6627207e0fdSHong Zhang     if (quadts) {
663cd4cee2dSHong Zhang       ierr  = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
6647207e0fdSHong Zhang     }
665abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
666b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
6671a352fa8SHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-adjoint_time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
668cd4cee2dSHong Zhang       if (quadJ) {
669cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
670cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
6711a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],adjoint_time_step,ts->vec_drdu_col);CHKERRQ(ierr);
672cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
673cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
67436eaed60SHong Zhang       }
6752ca6e920SHong Zhang     }
6763fd52205SHong Zhang     if (ts->vecs_sensip) {
67782ad9101SHong Zhang       th->shift = 1.0/(adjoint_time_step*th->Theta);
67882ad9101SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr);
6793e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6807207e0fdSHong Zhang       if (quadts) {
681cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
6827207e0fdSHong Zhang       }
683abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
68433f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
6851a352fa8SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-adjoint_time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
686cd4cee2dSHong Zhang         if (quadJp) {
687cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
688cd4cee2dSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
6891a352fa8SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],adjoint_time_step,ts->vec_drdp_col);CHKERRQ(ierr);
690cd4cee2dSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
691cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
69236eaed60SHong Zhang         }
69336eaed60SHong Zhang       }
6943fd52205SHong Zhang     }
6952ca6e920SHong Zhang   }
6962ca6e920SHong Zhang 
6972ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6982ca6e920SHong Zhang   PetscFunctionReturn(0);
6992ca6e920SHong Zhang }
7002ca6e920SHong Zhang 
701cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
702cd652676SJed Brown {
703cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
704be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
705cd652676SJed Brown   PetscErrorCode ierr;
706cd652676SJed Brown 
707cd652676SJed Brown   PetscFunctionBegin;
708a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
709be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
710be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
711cd652676SJed Brown   PetscFunctionReturn(0);
712cd652676SJed Brown }
713cd652676SJed Brown 
7141566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
7151566a47fSLisandro Dalcin {
7161566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7171566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
7181566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
7197453f775SEmil Constantinescu   PetscReal      wltea,wlter;
7201566a47fSLisandro Dalcin   PetscErrorCode ierr;
7211566a47fSLisandro Dalcin 
7221566a47fSLisandro Dalcin   PetscFunctionBegin;
7232ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
7241566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
725fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
7261566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
727fecfb714SLisandro Dalcin   {
728be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
729be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
7301566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
7311566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
7321566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
7331566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
7341566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
7357453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
7361566a47fSLisandro Dalcin   }
7371566a47fSLisandro Dalcin   if (order) *order = 2;
7381566a47fSLisandro Dalcin   PetscFunctionReturn(0);
7391566a47fSLisandro Dalcin }
7401566a47fSLisandro Dalcin 
7411566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7421566a47fSLisandro Dalcin {
7431566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7447207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7451566a47fSLisandro Dalcin   PetscErrorCode ierr;
7461566a47fSLisandro Dalcin 
7471566a47fSLisandro Dalcin   PetscFunctionBegin;
7481566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
749cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
750cd4cee2dSHong Zhang     ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr);
7511566a47fSLisandro Dalcin   }
752715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
75313b1b0ffSHong Zhang   if (ts->mat_sensip) {
75413b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7556f92202bSHong Zhang   }
7567207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7577207e0fdSHong Zhang     ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7586f92202bSHong Zhang   }
759715f1b00SHong Zhang   PetscFunctionReturn(0);
760715f1b00SHong Zhang }
761715f1b00SHong Zhang 
762715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
763715f1b00SHong Zhang {
764715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
765cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
76613b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
767b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7687207e0fdSHong Zhang   PetscInt       ntlm;
769715f1b00SHong Zhang   KSP            ksp;
7707207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
77113b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
7721a352fa8SHong Zhang   PetscReal      previous_shift;
773715f1b00SHong Zhang   PetscErrorCode ierr;
774715f1b00SHong Zhang 
775715f1b00SHong Zhang   PetscFunctionBegin;
7761a352fa8SHong Zhang   previous_shift = th->shift;
77713b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
77813b1b0ffSHong Zhang 
7797207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7807207e0fdSHong Zhang     ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7816f92202bSHong Zhang   }
782715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
783cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
7847207e0fdSHong Zhang   if (quadts) {
785cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
786cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
7877207e0fdSHong Zhang   }
788715f1b00SHong Zhang 
789715f1b00SHong Zhang   /* Build RHS */
790715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7911a352fa8SHong Zhang     th->shift = 1./((th->Theta-1.)*th->time_step0);
7921a352fa8SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime0,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
79313b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
79413b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
795715f1b00SHong Zhang 
796022c081aSHong Zhang     /* Add the f_p forcing terms */
7970e11c21fSHong Zhang     if (ts->Jacp) {
79882ad9101SHong Zhang       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
7991a352fa8SHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime0,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
80033f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
80182ad9101SHong Zhang       th->shift = previous_shift;
80282ad9101SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,ts->vec_sol);CHKERRQ(ierr);
8033e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
80433f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
8050e11c21fSHong Zhang     }
806715f1b00SHong Zhang   } else { /* 1-stage method */
8073e05ceb1SHong Zhang     th->shift = 0.0;
8083e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
80913b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
81013b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
81113b1b0ffSHong Zhang 
812022c081aSHong Zhang     /* Add the f_p forcing terms */
8130e11c21fSHong Zhang     if (ts->Jacp) {
81482ad9101SHong Zhang       th->shift = previous_shift;
81582ad9101SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr);
8163e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
81733f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
818715f1b00SHong Zhang     }
8190e11c21fSHong Zhang   }
820715f1b00SHong Zhang 
821715f1b00SHong Zhang   /* Build LHS */
8221a352fa8SHong Zhang   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
823715f1b00SHong Zhang   if (th->endpoint) {
8243e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
825715f1b00SHong Zhang   } else {
8263e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
827715f1b00SHong Zhang   }
828cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
829715f1b00SHong Zhang 
830715f1b00SHong Zhang   /*
831715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
832c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
833715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
834715f1b00SHong Zhang   */
835715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8367207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
8371a352fa8SHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime0,th->X0,quadJ,NULL);CHKERRQ(ierr);
8381a352fa8SHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->ptime0,th->X0,quadJp);CHKERRQ(ierr);
8397207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8407207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8411a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
842715f1b00SHong Zhang     }
843715f1b00SHong Zhang   }
844715f1b00SHong Zhang 
845715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
846022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
847b5b4017aSHong Zhang     KSPConvergedReason kspreason;
84813b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
849b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
850715f1b00SHong Zhang     if (th->endpoint) {
85113b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
852b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
853b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
854b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
85513b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
856715f1b00SHong Zhang     } else {
857b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
858715f1b00SHong Zhang     }
859b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
860b5b4017aSHong Zhang     if (kspreason < 0) {
861b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
862b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
863b5b4017aSHong Zhang     }
864b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
86513b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
866715f1b00SHong Zhang   }
867715f1b00SHong Zhang 
86813b1b0ffSHong Zhang   /*
86913b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
870c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
87113b1b0ffSHong Zhang   */
8727207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
87313b1b0ffSHong Zhang     if (!th->endpoint) {
8747207e0fdSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */
875cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
876cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
8777207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8787207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8791a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
88013b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
88113b1b0ffSHong Zhang     } else {
882cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
883cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
8847207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8857207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8861a352fa8SHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step0*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
88713b1b0ffSHong Zhang     }
88813b1b0ffSHong Zhang   } else {
88913b1b0ffSHong Zhang     if (!th->endpoint) {
89013b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
891715f1b00SHong Zhang     }
892715f1b00SHong Zhang   }
8931566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8941566a47fSLisandro Dalcin }
8951566a47fSLisandro Dalcin 
896cc4f23bcSHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat *stagesensip[])
8976affb6f8SHong Zhang {
8986affb6f8SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
8996affb6f8SHong Zhang 
9006affb6f8SHong Zhang   PetscFunctionBegin;
9017409eceaSHong Zhang   if (ns) {
9027409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
9037409eceaSHong Zhang     else *ns = 2; /* endpoint form */
9047409eceaSHong Zhang   }
9055072d23cSHong Zhang   if (stagesensip) {
906cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
9077409eceaSHong Zhang       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
908cc4f23bcSHong Zhang     } else {
909cc4f23bcSHong Zhang       th->MatFwdStages[0] = th->MatFwdSensip0;
910cc4f23bcSHong Zhang       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
911cc4f23bcSHong Zhang     }
912cc4f23bcSHong Zhang     *stagesensip = th->MatFwdStages;
9135072d23cSHong Zhang   }
9146affb6f8SHong Zhang   PetscFunctionReturn(0);
9156affb6f8SHong Zhang }
9166affb6f8SHong Zhang 
917316643e7SJed Brown /*------------------------------------------------------------*/
918277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
919316643e7SJed Brown {
920316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
921316643e7SJed Brown   PetscErrorCode ierr;
922316643e7SJed Brown 
923316643e7SJed Brown   PetscFunctionBegin;
9246bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
9256bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
9263b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
927eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
9281566a47fSLisandro Dalcin 
9291566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
9301566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
9311566a47fSLisandro Dalcin 
932b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
933277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
934277b19d0SLisandro Dalcin }
935277b19d0SLisandro Dalcin 
936ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
937ece46509SHong Zhang {
938ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
939ece46509SHong Zhang   PetscErrorCode ierr;
940ece46509SHong Zhang 
941ece46509SHong Zhang   PetscFunctionBegin;
942ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
943ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
944ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
945ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
946ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
947ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
948ece46509SHong Zhang   PetscFunctionReturn(0);
949ece46509SHong Zhang }
950ece46509SHong Zhang 
951277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
952277b19d0SLisandro Dalcin {
953277b19d0SLisandro Dalcin   PetscErrorCode ierr;
954277b19d0SLisandro Dalcin 
955277b19d0SLisandro Dalcin   PetscFunctionBegin;
956277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
957b3a6b972SJed Brown   if (ts->dm) {
958b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
959b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
960b3a6b972SJed Brown   }
961277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
962bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
963bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
964bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
965bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
966316643e7SJed Brown   PetscFunctionReturn(0);
967316643e7SJed Brown }
968316643e7SJed Brown 
969316643e7SJed Brown /*
970316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9712b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9720fd10508SMatthew Knepley 
9730fd10508SMatthew Knepley   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
9740fd10508SMatthew Knepley   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
9750fd10508SMatthew Knepley   U = (U_{n+1} + U0)/2
976316643e7SJed Brown */
9770f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
978316643e7SJed Brown {
979316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
980316643e7SJed Brown   PetscErrorCode ierr;
9817445fe48SJed Brown   Vec            X0,Xdot;
9827445fe48SJed Brown   DM             dm,dmsave;
9833e05ceb1SHong Zhang   PetscReal      shift = th->shift;
984316643e7SJed Brown 
985316643e7SJed Brown   PetscFunctionBegin;
9867445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9875a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9887445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
989494ce9fcSHong Zhang   if (x != X0) {
990b296d7d5SJed Brown     ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
991494ce9fcSHong Zhang   } else {
992494ce9fcSHong Zhang     ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
993494ce9fcSHong Zhang   }
9947445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9957445fe48SJed Brown   dmsave = ts->dm;
9967445fe48SJed Brown   ts->dm = dm;
9977445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
9987445fe48SJed Brown   ts->dm = dmsave;
9990d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
1000316643e7SJed Brown   PetscFunctionReturn(0);
1001316643e7SJed Brown }
1002316643e7SJed Brown 
1003d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
1004316643e7SJed Brown {
1005316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1006316643e7SJed Brown   PetscErrorCode ierr;
10077445fe48SJed Brown   Vec            Xdot;
10087445fe48SJed Brown   DM             dm,dmsave;
10093e05ceb1SHong Zhang   PetscReal      shift = th->shift;
1010316643e7SJed Brown 
1011316643e7SJed Brown   PetscFunctionBegin;
10127445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1013be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
10140298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
10157445fe48SJed Brown 
10167445fe48SJed Brown   dmsave = ts->dm;
10177445fe48SJed Brown   ts->dm = dm;
1018d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
10197445fe48SJed Brown   ts->dm = dmsave;
10200298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
1021316643e7SJed Brown   PetscFunctionReturn(0);
1022316643e7SJed Brown }
1023316643e7SJed Brown 
1024715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
1025715f1b00SHong Zhang {
1026715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10277207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
1028715f1b00SHong Zhang   PetscErrorCode ierr;
1029715f1b00SHong Zhang 
1030715f1b00SHong Zhang   PetscFunctionBegin;
1031022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
103213b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
103313b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10347207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10357207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10367207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
1037715f1b00SHong Zhang   }
10386f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
103913b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
104013b1b0ffSHong Zhang 
1041b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
1042715f1b00SHong Zhang   PetscFunctionReturn(0);
1043715f1b00SHong Zhang }
1044715f1b00SHong Zhang 
10457adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
10467adebddeSHong Zhang {
10477adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10487207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10497adebddeSHong Zhang   PetscErrorCode ierr;
10507adebddeSHong Zhang 
10517adebddeSHong Zhang   PetscFunctionBegin;
10527207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10537207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10547207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
10557adebddeSHong Zhang   }
10567adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
10577adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10587adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
10597adebddeSHong Zhang   PetscFunctionReturn(0);
10607adebddeSHong Zhang }
10617adebddeSHong Zhang 
1062316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1063316643e7SJed Brown {
1064316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1065cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10662ffb9264SLisandro Dalcin   PetscBool      match;
1067316643e7SJed Brown   PetscErrorCode ierr;
1068316643e7SJed Brown 
1069316643e7SJed Brown   PetscFunctionBegin;
1070cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1071cd4cee2dSHong Zhang     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1072b1cb36f3SHong Zhang   }
107339d13666SHong Zhang   if (!th->X) {
1074316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
107539d13666SHong Zhang   }
107639d13666SHong Zhang   if (!th->Xdot) {
1077316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
107839d13666SHong Zhang   }
107939d13666SHong Zhang   if (!th->X0) {
10803b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
108139d13666SHong Zhang   }
10821566a47fSLisandro Dalcin   if (th->endpoint) {
10831566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
10847445fe48SJed Brown   }
10853b1890cdSShri Abhyankar 
10861566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
108720d49056SMatthew G. Knepley   th->shift = 1/(th->Theta*ts->time_step);
10881566a47fSLisandro Dalcin 
10891566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
10901566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10911566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10921566a47fSLisandro Dalcin 
10931566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
10941566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
10952ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
10962ffb9264SLisandro Dalcin   if (!match) {
10971566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
10981566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
10993b1890cdSShri Abhyankar   }
11001566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1101cc4f23bcSHong Zhang 
1102cc4f23bcSHong Zhang   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1103b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1104b4dd3bf9SBarry Smith }
11050c7376e5SHong Zhang 
1106b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1107b4dd3bf9SBarry Smith 
1108b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1109b4dd3bf9SBarry Smith {
1110b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1111b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1112b4dd3bf9SBarry Smith 
1113b4dd3bf9SBarry Smith   PetscFunctionBegin;
1114b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1115b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
11162ca6e920SHong Zhang   if (ts->vecs_sensip) {
1117b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
11182ca6e920SHong Zhang   }
1119b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1120b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1121b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
112267633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
112367633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
112467633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1125b5b4017aSHong Zhang   }
1126c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1127c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
112867633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
112967633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
113067633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1131b5b4017aSHong Zhang   }
1132316643e7SJed Brown   PetscFunctionReturn(0);
1133316643e7SJed Brown }
1134316643e7SJed Brown 
11354416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1136316643e7SJed Brown {
1137316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1138316643e7SJed Brown   PetscErrorCode ierr;
1139316643e7SJed Brown 
1140316643e7SJed Brown   PetscFunctionBegin;
1141e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1142316643e7SJed Brown   {
11430298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
11440298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr);
114503842d09SLisandro Dalcin     ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
1146316643e7SJed Brown   }
1147316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1148316643e7SJed Brown   PetscFunctionReturn(0);
1149316643e7SJed Brown }
1150316643e7SJed Brown 
1151316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1152316643e7SJed Brown {
1153316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1154ace3abfcSBarry Smith   PetscBool      iascii;
1155316643e7SJed Brown   PetscErrorCode ierr;
1156316643e7SJed Brown 
1157316643e7SJed Brown   PetscFunctionBegin;
1158251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1159316643e7SJed Brown   if (iascii) {
11607c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1161ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1162316643e7SJed Brown   }
1163316643e7SJed Brown   PetscFunctionReturn(0);
1164316643e7SJed Brown }
1165316643e7SJed Brown 
1166560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11670de4c49aSJed Brown {
11680de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11690de4c49aSJed Brown 
11700de4c49aSJed Brown   PetscFunctionBegin;
11710de4c49aSJed Brown   *theta = th->Theta;
11720de4c49aSJed Brown   PetscFunctionReturn(0);
11730de4c49aSJed Brown }
11740de4c49aSJed Brown 
1175560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11760de4c49aSJed Brown {
11770de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11780de4c49aSJed Brown 
11790de4c49aSJed Brown   PetscFunctionBegin;
11807c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11810de4c49aSJed Brown   th->Theta = theta;
11821566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11830de4c49aSJed Brown   PetscFunctionReturn(0);
11840de4c49aSJed Brown }
1185eb284becSJed Brown 
1186560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
118726f2ff8fSLisandro Dalcin {
118826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
118926f2ff8fSLisandro Dalcin 
119026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
119126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
119226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
119326f2ff8fSLisandro Dalcin }
119426f2ff8fSLisandro Dalcin 
1195560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1196eb284becSJed Brown {
1197eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1198eb284becSJed Brown 
1199eb284becSJed Brown   PetscFunctionBegin;
1200eb284becSJed Brown   th->endpoint = flg;
1201eb284becSJed Brown   PetscFunctionReturn(0);
1202eb284becSJed Brown }
12030de4c49aSJed Brown 
1204f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1205f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1206f9c1d6abSBarry Smith {
1207f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1208f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
12093fd8ae06SJed Brown   const PetscReal one = 1.0;
1210f9c1d6abSBarry Smith 
1211f9c1d6abSBarry Smith   PetscFunctionBegin;
12123fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1213f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1214f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1215f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1216f9c1d6abSBarry Smith }
1217f9c1d6abSBarry Smith #endif
1218f9c1d6abSBarry Smith 
1219cc4f23bcSHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec *Y[])
122042682096SHong Zhang {
122142682096SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
122242682096SHong Zhang 
122342682096SHong Zhang   PetscFunctionBegin;
12247409eceaSHong Zhang   if (ns) {
12257409eceaSHong Zhang     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
12267409eceaSHong Zhang     else *ns = 2; /* endpoint form */
12277409eceaSHong Zhang   }
12285072d23cSHong Zhang   if (Y) {
1229cc4f23bcSHong Zhang     if (!th->endpoint && th->Theta != 1.0) {
12307409eceaSHong Zhang       th->Stages[0] = th->X;
1231cc4f23bcSHong Zhang     } else {
1232cc4f23bcSHong Zhang       th->Stages[0] = th->X0;
1233cc4f23bcSHong Zhang       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1234cc4f23bcSHong Zhang     }
1235cc4f23bcSHong Zhang     *Y = th->Stages;
12365072d23cSHong Zhang   }
123742682096SHong Zhang   PetscFunctionReturn(0);
123842682096SHong Zhang }
1239f9c1d6abSBarry Smith 
1240316643e7SJed Brown /* ------------------------------------------------------------ */
1241316643e7SJed Brown /*MC
124296f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1243316643e7SJed Brown 
1244316643e7SJed Brown    Level: beginner
1245316643e7SJed Brown 
12464eb428fdSBarry Smith    Options Database:
1247c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1248c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
124903842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
12504eb428fdSBarry Smith 
1251eb284becSJed Brown    Notes:
12520c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
12530c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
12544eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
12554eb428fdSBarry Smith 
12567409eceaSHong 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.
1257eb284becSJed Brown 
12587409eceaSHong Zhang    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.
1259eb284becSJed Brown 
1260eb284becSJed Brown .vb
1261eb284becSJed Brown   Theta | Theta
1262eb284becSJed Brown   -------------
1263eb284becSJed Brown         |  1
1264eb284becSJed Brown .ve
1265eb284becSJed Brown 
1266eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1267eb284becSJed Brown 
1268eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1269eb284becSJed Brown 
1270eb284becSJed Brown .vb
1271eb284becSJed Brown   0 | 0         0
1272eb284becSJed Brown   1 | 1-Theta   Theta
1273eb284becSJed Brown   -------------------
1274eb284becSJed Brown     | 1-Theta   Theta
1275eb284becSJed Brown .ve
1276eb284becSJed Brown 
1277eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1278eb284becSJed Brown 
1279eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1280eb284becSJed Brown 
1281eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1282eb284becSJed Brown 
12834eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1284eb284becSJed Brown 
1285eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1286316643e7SJed Brown 
1287316643e7SJed Brown M*/
12888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1289316643e7SJed Brown {
1290316643e7SJed Brown   TS_Theta       *th;
1291316643e7SJed Brown   PetscErrorCode ierr;
1292316643e7SJed Brown 
1293316643e7SJed Brown   PetscFunctionBegin;
1294277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1295ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1296316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1297316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1298316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
129942f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1300ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1301316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1302cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
13031566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
130424655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1305316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
13060f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
13070f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1308f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1309f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1310f9c1d6abSBarry Smith #endif
131142682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
131242f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1313b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1314b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
13152ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
13166affb6f8SHong Zhang 
1317715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
13187adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1319715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
13206affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1321316643e7SJed Brown 
1322efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1323efd4aadfSBarry Smith 
1324b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1325316643e7SJed Brown   ts->data = (void*)th;
1326316643e7SJed Brown 
1327715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1328715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1329715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1330b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1331715f1b00SHong Zhang 
13326f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1333316643e7SJed Brown   th->Theta       = 0.5;
13341566a47fSLisandro Dalcin   th->order       = 2;
1335bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1336bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1337bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1338bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1339316643e7SJed Brown   PetscFunctionReturn(0);
1340316643e7SJed Brown }
13410de4c49aSJed Brown 
13420de4c49aSJed Brown /*@
13430de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
13440de4c49aSJed Brown 
13450de4c49aSJed Brown   Not Collective
13460de4c49aSJed Brown 
13470de4c49aSJed Brown   Input Parameter:
13480de4c49aSJed Brown .  ts - timestepping context
13490de4c49aSJed Brown 
13500de4c49aSJed Brown   Output Parameter:
13510de4c49aSJed Brown .  theta - stage abscissa
13520de4c49aSJed Brown 
13530de4c49aSJed Brown   Note:
13540de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
13550de4c49aSJed Brown 
13560de4c49aSJed Brown   Level: Advanced
13570de4c49aSJed Brown 
13580de4c49aSJed Brown .seealso: TSThetaSetTheta()
13590de4c49aSJed Brown @*/
13607087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13610de4c49aSJed Brown {
13624ac538c5SBarry Smith   PetscErrorCode ierr;
13630de4c49aSJed Brown 
13640de4c49aSJed Brown   PetscFunctionBegin;
1365afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13660de4c49aSJed Brown   PetscValidPointer(theta,2);
13674ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
13680de4c49aSJed Brown   PetscFunctionReturn(0);
13690de4c49aSJed Brown }
13700de4c49aSJed Brown 
13710de4c49aSJed Brown /*@
13720de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13730de4c49aSJed Brown 
13740de4c49aSJed Brown   Not Collective
13750de4c49aSJed Brown 
1376*d8d19677SJose E. Roman   Input Parameters:
13770de4c49aSJed Brown +  ts - timestepping context
13780de4c49aSJed Brown -  theta - stage abscissa
13790de4c49aSJed Brown 
13800de4c49aSJed Brown   Options Database:
13810de4c49aSJed Brown .  -ts_theta_theta <theta>
13820de4c49aSJed Brown 
13830de4c49aSJed Brown   Level: Intermediate
13840de4c49aSJed Brown 
13850de4c49aSJed Brown .seealso: TSThetaGetTheta()
13860de4c49aSJed Brown @*/
13877087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13880de4c49aSJed Brown {
13894ac538c5SBarry Smith   PetscErrorCode ierr;
13900de4c49aSJed Brown 
13910de4c49aSJed Brown   PetscFunctionBegin;
1392afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13934ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
13940de4c49aSJed Brown   PetscFunctionReturn(0);
13950de4c49aSJed Brown }
1396f33bbcb6SJed Brown 
139726f2ff8fSLisandro Dalcin /*@
139826f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
139926f2ff8fSLisandro Dalcin 
140026f2ff8fSLisandro Dalcin   Not Collective
140126f2ff8fSLisandro Dalcin 
140226f2ff8fSLisandro Dalcin   Input Parameter:
140326f2ff8fSLisandro Dalcin .  ts - timestepping context
140426f2ff8fSLisandro Dalcin 
140526f2ff8fSLisandro Dalcin   Output Parameter:
140626f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
140726f2ff8fSLisandro Dalcin 
140826f2ff8fSLisandro Dalcin   Level: Advanced
140926f2ff8fSLisandro Dalcin 
141026f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
141126f2ff8fSLisandro Dalcin @*/
141226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
141326f2ff8fSLisandro Dalcin {
141426f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
141526f2ff8fSLisandro Dalcin 
141626f2ff8fSLisandro Dalcin   PetscFunctionBegin;
141726f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
141826f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1419163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
142026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
142126f2ff8fSLisandro Dalcin }
142226f2ff8fSLisandro Dalcin 
1423eb284becSJed Brown /*@
1424eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1425eb284becSJed Brown 
1426eb284becSJed Brown   Not Collective
1427eb284becSJed Brown 
1428*d8d19677SJose E. Roman   Input Parameters:
1429eb284becSJed Brown +  ts - timestepping context
1430eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1431eb284becSJed Brown 
1432eb284becSJed Brown   Options Database:
1433eb284becSJed Brown .  -ts_theta_endpoint <flg>
1434eb284becSJed Brown 
1435eb284becSJed Brown   Level: Intermediate
1436eb284becSJed Brown 
1437eb284becSJed Brown .seealso: TSTHETA, TSCN
1438eb284becSJed Brown @*/
1439eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1440eb284becSJed Brown {
1441eb284becSJed Brown   PetscErrorCode ierr;
1442eb284becSJed Brown 
1443eb284becSJed Brown   PetscFunctionBegin;
1444eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1445eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1446eb284becSJed Brown   PetscFunctionReturn(0);
1447eb284becSJed Brown }
1448eb284becSJed Brown 
1449f33bbcb6SJed Brown /*
1450f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1451f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1452f33bbcb6SJed Brown  */
1453f33bbcb6SJed Brown 
14541566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
14551566a47fSLisandro Dalcin {
14561566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14571566a47fSLisandro Dalcin   PetscErrorCode ierr;
14581566a47fSLisandro Dalcin 
14591566a47fSLisandro Dalcin   PetscFunctionBegin;
14601566a47fSLisandro Dalcin   if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n");
14611566a47fSLisandro Dalcin   if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n");
14621566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14631566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14641566a47fSLisandro Dalcin }
14651566a47fSLisandro Dalcin 
1466f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1467f33bbcb6SJed Brown {
1468f33bbcb6SJed Brown   PetscFunctionBegin;
1469f33bbcb6SJed Brown   PetscFunctionReturn(0);
1470f33bbcb6SJed Brown }
1471f33bbcb6SJed Brown 
1472f33bbcb6SJed Brown /*MC
1473f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1474f33bbcb6SJed Brown 
1475f33bbcb6SJed Brown   Level: beginner
1476f33bbcb6SJed Brown 
14774eb428fdSBarry Smith   Notes:
1478c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14794eb428fdSBarry Smith 
14801566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14814eb428fdSBarry Smith 
1482f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1483f33bbcb6SJed Brown 
1484f33bbcb6SJed Brown M*/
14858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1486f33bbcb6SJed Brown {
1487f33bbcb6SJed Brown   PetscErrorCode ierr;
1488f33bbcb6SJed Brown 
1489f33bbcb6SJed Brown   PetscFunctionBegin;
1490f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1491f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
14921566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
14930c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1494f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1495f33bbcb6SJed Brown   PetscFunctionReturn(0);
1496f33bbcb6SJed Brown }
1497f33bbcb6SJed Brown 
14981566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14991566a47fSLisandro Dalcin {
15001566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
15011566a47fSLisandro Dalcin   PetscErrorCode ierr;
15021566a47fSLisandro Dalcin 
15031566a47fSLisandro Dalcin   PetscFunctionBegin;
15041566a47fSLisandro Dalcin   if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n");
15051566a47fSLisandro Dalcin   if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n");
15061566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
15071566a47fSLisandro Dalcin   PetscFunctionReturn(0);
15081566a47fSLisandro Dalcin }
15091566a47fSLisandro Dalcin 
1510f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1511f33bbcb6SJed Brown {
1512f33bbcb6SJed Brown   PetscFunctionBegin;
1513f33bbcb6SJed Brown   PetscFunctionReturn(0);
1514f33bbcb6SJed Brown }
1515f33bbcb6SJed Brown 
1516f33bbcb6SJed Brown /*MC
1517f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1518f33bbcb6SJed Brown 
1519f33bbcb6SJed Brown   Level: beginner
1520f33bbcb6SJed Brown 
1521f33bbcb6SJed Brown   Notes:
15227cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
15237cf5af47SJed Brown 
15247cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1525f33bbcb6SJed Brown 
1526f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1527f33bbcb6SJed Brown 
1528f33bbcb6SJed Brown M*/
15298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1530f33bbcb6SJed Brown {
1531f33bbcb6SJed Brown   PetscErrorCode ierr;
1532f33bbcb6SJed Brown 
1533f33bbcb6SJed Brown   PetscFunctionBegin;
1534f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1535f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1536eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
15370c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1538f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1539f33bbcb6SJed Brown   PetscFunctionReturn(0);
1540f33bbcb6SJed Brown }
1541