xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 494ce9fcc6b758f1f9163de9f43b431662f5beef)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
5a055c8caSBarry Smith #include <petscsnes.h>
61e25c274SJed Brown #include <petscdm.h>
73c54f92cSHong Zhang #include <petscmat.h>
8316643e7SJed Brown 
9316643e7SJed Brown typedef struct {
10715f1b00SHong Zhang   /* context for time stepping */
111566a47fSLisandro Dalcin   PetscReal    stage_time;
121566a47fSLisandro Dalcin   Vec          X0,X,Xdot;                /* Storage for stages and time derivative */
131566a47fSLisandro Dalcin   Vec          affine;                   /* Affine vector needed for residual at beginning of step in endpoint formulation */
141566a47fSLisandro Dalcin   PetscReal    Theta;
15715f1b00SHong Zhang   PetscReal    ptime;
16715f1b00SHong Zhang   PetscReal    time_step;
173e05ceb1SHong Zhang   PetscReal    shift;
181566a47fSLisandro Dalcin   PetscInt     order;
191566a47fSLisandro Dalcin   PetscBool    endpoint;
201566a47fSLisandro Dalcin   PetscBool    extrapolate;
21715f1b00SHong Zhang   TSStepStatus status;
22715f1b00SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
231566a47fSLisandro Dalcin 
24715f1b00SHong Zhang   /* context for sensitivity analysis */
25022c081aSHong Zhang   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
26b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
27b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
2813b1b0ffSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
2913b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
30b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
3113b1b0ffSHong Zhang   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
327207e0fdSHong Zhang   Mat          MatIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
337207e0fdSHong Zhang   Mat          MatIntegralSensip0;       /* backup for roll-backs due to events */
34b5b4017aSHong Zhang   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
35b5b4017aSHong Zhang   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
36b5b4017aSHong Zhang   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
37b5b4017aSHong Zhang   Vec          *VecsAffine;              /* Working vectors to store residuals */
38715f1b00SHong Zhang   /* context for error estimation */
391566a47fSLisandro Dalcin   Vec          vec_sol_prev;
401566a47fSLisandro Dalcin   Vec          vec_lte_work;
41316643e7SJed Brown } TS_Theta;
42316643e7SJed Brown 
437445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
447445fe48SJed Brown {
457445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
467445fe48SJed Brown   PetscErrorCode ierr;
477445fe48SJed Brown 
487445fe48SJed Brown   PetscFunctionBegin;
497445fe48SJed Brown   if (X0) {
507445fe48SJed Brown     if (dm && dm != ts->dm) {
510d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
527445fe48SJed Brown     } else *X0 = ts->vec_sol;
537445fe48SJed Brown   }
547445fe48SJed Brown   if (Xdot) {
557445fe48SJed Brown     if (dm && dm != ts->dm) {
560d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
577445fe48SJed Brown     } else *Xdot = th->Xdot;
587445fe48SJed Brown   }
597445fe48SJed Brown   PetscFunctionReturn(0);
607445fe48SJed Brown }
617445fe48SJed Brown 
620d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
630d0b770aSPeter Brune {
640d0b770aSPeter Brune   PetscErrorCode ierr;
650d0b770aSPeter Brune 
660d0b770aSPeter Brune   PetscFunctionBegin;
670d0b770aSPeter Brune   if (X0) {
680d0b770aSPeter Brune     if (dm && dm != ts->dm) {
690d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
700d0b770aSPeter Brune     }
710d0b770aSPeter Brune   }
720d0b770aSPeter Brune   if (Xdot) {
730d0b770aSPeter Brune     if (dm && dm != ts->dm) {
740d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
750d0b770aSPeter Brune     }
760d0b770aSPeter Brune   }
770d0b770aSPeter Brune   PetscFunctionReturn(0);
780d0b770aSPeter Brune }
790d0b770aSPeter Brune 
807445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
817445fe48SJed Brown {
827445fe48SJed Brown   PetscFunctionBegin;
837445fe48SJed Brown   PetscFunctionReturn(0);
847445fe48SJed Brown }
857445fe48SJed Brown 
867445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
877445fe48SJed Brown {
887445fe48SJed Brown   TS             ts = (TS)ctx;
897445fe48SJed Brown   PetscErrorCode ierr;
907445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
917445fe48SJed Brown 
927445fe48SJed Brown   PetscFunctionBegin;
937445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
977445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
987445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
990d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
1000d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1017445fe48SJed Brown   PetscFunctionReturn(0);
1027445fe48SJed Brown }
1037445fe48SJed Brown 
104258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
105258e1594SPeter Brune {
106258e1594SPeter Brune   PetscFunctionBegin;
107258e1594SPeter Brune   PetscFunctionReturn(0);
108258e1594SPeter Brune }
109258e1594SPeter Brune 
110258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
111258e1594SPeter Brune {
112258e1594SPeter Brune   TS             ts = (TS)ctx;
113258e1594SPeter Brune   PetscErrorCode ierr;
114258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
115258e1594SPeter Brune 
116258e1594SPeter Brune   PetscFunctionBegin;
117258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
118258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
119258e1594SPeter Brune 
120258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122258e1594SPeter Brune 
123258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
125258e1594SPeter Brune 
126258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
127258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
128258e1594SPeter Brune   PetscFunctionReturn(0);
129258e1594SPeter Brune }
130258e1594SPeter Brune 
131022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
132022c081aSHong Zhang {
133022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
134cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
135022c081aSHong Zhang   PetscErrorCode ierr;
136022c081aSHong Zhang 
137022c081aSHong Zhang   PetscFunctionBegin;
138022c081aSHong Zhang   if (th->endpoint) {
139022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
140022c081aSHong Zhang     if (th->Theta!=1.0) {
141cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
142cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
143022c081aSHong Zhang     }
144cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
145cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
146022c081aSHong Zhang   } else {
147cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
148cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
149022c081aSHong Zhang   }
150022c081aSHong Zhang   PetscFunctionReturn(0);
151022c081aSHong Zhang }
152022c081aSHong Zhang 
153b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
154b1cb36f3SHong Zhang {
155b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
156cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
157b1cb36f3SHong Zhang   PetscErrorCode ierr;
158b1cb36f3SHong Zhang 
159b1cb36f3SHong Zhang   PetscFunctionBegin;
160b1cb36f3SHong Zhang   /* backup cost integral */
161cd4cee2dSHong Zhang   ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr);
162022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
163b1cb36f3SHong Zhang   PetscFunctionReturn(0);
164b1cb36f3SHong Zhang }
165b1cb36f3SHong Zhang 
166b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
167b1cb36f3SHong Zhang {
168b1cb36f3SHong Zhang   PetscErrorCode ierr;
169b1cb36f3SHong Zhang 
170b1cb36f3SHong Zhang   PetscFunctionBegin;
171022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
17224655328SShri   PetscFunctionReturn(0);
17324655328SShri }
17424655328SShri 
175470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1761566a47fSLisandro Dalcin {
1771566a47fSLisandro Dalcin   PetscInt       nits,lits;
1781566a47fSLisandro Dalcin   PetscErrorCode ierr;
1791566a47fSLisandro Dalcin 
1801566a47fSLisandro Dalcin   PetscFunctionBegin;
1811566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1821566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1831566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1841566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1851566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1861566a47fSLisandro Dalcin }
1871566a47fSLisandro Dalcin 
188193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
189316643e7SJed Brown {
190316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1911566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1924957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1931566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
194051f2191SLisandro Dalcin   PetscErrorCode ierr;
195316643e7SJed Brown 
196316643e7SJed Brown   PetscFunctionBegin;
1971566a47fSLisandro Dalcin   if (!ts->steprollback) {
1982ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
1993b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
2001566a47fSLisandro Dalcin   }
2011566a47fSLisandro Dalcin 
2021566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
2031566a47fSLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
204715f1b00SHong Zhang 
2053e05ceb1SHong Zhang     th->shift = 1/(th->Theta*ts->time_step);
2061566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
207316643e7SJed Brown 
208be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
209fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
2103e05ceb1SHong Zhang       ierr = VecAXPY(th->X,1/th->shift,th->Xdot);CHKERRQ(ierr);
211be5899b3SLisandro Dalcin     }
212eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
213eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2141566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2151566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2161566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2171566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2181566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
219eb284becSJed Brown     }
220be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
221470880abSPatrick Sanan     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
222be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2231566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
224fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
225051f2191SLisandro Dalcin 
226051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2271566a47fSLisandro Dalcin     if (th->endpoint) {
2281566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2291566a47fSLisandro Dalcin     } else {
2303e05ceb1SHong Zhang       ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->X0,th->X);CHKERRQ(ierr);
2311566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2321566a47fSLisandro Dalcin     }
2331566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2341566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2351566a47fSLisandro Dalcin     if (!accept) {
2361566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
237be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
238be5899b3SLisandro Dalcin       goto reject_step;
239051f2191SLisandro Dalcin     }
2403b1890cdSShri Abhyankar 
241715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
242b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
243b1cb36f3SHong Zhang       th->time_step = ts->time_step;
24417145e75SHong Zhang     }
2451566a47fSLisandro Dalcin 
2462b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
247cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
248051f2191SLisandro Dalcin     break;
249051f2191SLisandro Dalcin 
250051f2191SLisandro Dalcin   reject_step:
251fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2521566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
253051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2541566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2553b1890cdSShri Abhyankar     }
2563b1890cdSShri Abhyankar   }
257316643e7SJed Brown   PetscFunctionReturn(0);
258316643e7SJed Brown }
259316643e7SJed Brown 
2606980b705SHong Zhang /*
2616980b705SHong Zhang   Use SNES to compute the Jacobian so that finite differencing could be used when TS Jacobian is not available.
2626980b705SHong Zhang */
2636980b705SHong Zhang static PetscErrorCode KSPTSFormOperator_Private(KSP ksp,Vec x,Mat J,Mat Jpre,TS ts)
2646980b705SHong Zhang {
2656980b705SHong Zhang   SNES           snes = ts->snes;
2666980b705SHong Zhang   MatFDColoring  color;
2676980b705SHong Zhang   PetscErrorCode ierr;
2686980b705SHong Zhang 
2696980b705SHong Zhang   PetscFunctionBegin;
2706980b705SHong Zhang   /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
2716980b705SHong Zhang   ierr = PetscObjectQuery((PetscObject)Jpre,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);
2726980b705SHong Zhang   if (color) {
2736980b705SHong Zhang     Vec f;
2746980b705SHong Zhang     ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
2756980b705SHong Zhang     ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
2766980b705SHong Zhang   }
2776980b705SHong Zhang   ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr);
2786980b705SHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
2796980b705SHong Zhang   PetscFunctionReturn(0);
2806980b705SHong Zhang }
2816980b705SHong Zhang 
282c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
283c9681e5dSHong Zhang {
284c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
285cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
286c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
287c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
288c9681e5dSHong Zhang   PetscInt       nadj;
2897207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
290c9681e5dSHong Zhang   KSP            ksp;
291c9681e5dSHong Zhang   PetscScalar    *xarr;
292077c4feaSHong Zhang   TSEquationType eqtype;
293077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
294c9681e5dSHong Zhang   PetscErrorCode ierr;
295c9681e5dSHong Zhang 
296c9681e5dSHong Zhang   PetscFunctionBegin;
297077c4feaSHong Zhang   ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr);
298077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
299077c4feaSHong Zhang     isexplicitode  = PETSC_TRUE;
300077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
301077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
302077c4feaSHong Zhang   }
303c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
304c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
305cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
3067207e0fdSHong Zhang   if (quadts) {
307cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
308cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
3097207e0fdSHong Zhang   }
310c9681e5dSHong Zhang 
311c9681e5dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
312c9681e5dSHong Zhang   th->stage_time = ts->ptime; /* time_step is negative*/
313c9681e5dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
314c9681e5dSHong Zhang   th->time_step  = -ts->time_step;
315c9681e5dSHong Zhang 
31633f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
3177207e0fdSHong Zhang   if (quadts) {
318cd4cee2dSHong Zhang     ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
3197207e0fdSHong Zhang   }
320cd4cee2dSHong Zhang 
321c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
322c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
323c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
324cd4cee2dSHong Zhang     if (quadJ) {
325cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
326cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
327cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
328cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
329cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
330c9681e5dSHong Zhang     }
331c9681e5dSHong Zhang   }
332c9681e5dSHong Zhang 
333c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
3346980b705SHong Zhang   ierr = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr);
335c9681e5dSHong Zhang 
336c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
337c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
338c9681e5dSHong Zhang     KSPConvergedReason kspreason;
339c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
340c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
341c9681e5dSHong Zhang     if (kspreason < 0) {
342c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
343c9681e5dSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
344c9681e5dSHong Zhang     }
345c9681e5dSHong Zhang   }
346c9681e5dSHong Zhang 
347c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
348c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
349c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
350c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
351c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
352b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
353c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
354b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
355c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
356c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
3576980b705SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],1./th->time_step);CHKERRQ(ierr);
35833f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
359c9681e5dSHong Zhang       if (ts->vecs_fup) {
36033f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
361c9681e5dSHong Zhang       }
362c9681e5dSHong Zhang     }
363c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
364c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
36505c9950eSHong Zhang       KSPConvergedReason kspreason;
366c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
36705c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
36805c9950eSHong Zhang       if (kspreason < 0) {
36905c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
37005c9950eSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
37105c9950eSHong Zhang       }
372c9681e5dSHong Zhang     }
373c9681e5dSHong Zhang   }
374c9681e5dSHong Zhang 
375c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
376077c4feaSHong Zhang   if (!isexplicitode) {
3773e05ceb1SHong Zhang     th->shift = 0.0;
3783e05ceb1SHong Zhang     ierr  = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr);
379c9681e5dSHong Zhang     ierr  = MatScale(J,-1.);CHKERRQ(ierr);
380c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
38133f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
382c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
383c9681e5dSHong Zhang       ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr);
384c9681e5dSHong Zhang       ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
385c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
386c9681e5dSHong Zhang         ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
387c9681e5dSHong Zhang         ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr);
388c9681e5dSHong Zhang         ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
389c9681e5dSHong Zhang       }
390c9681e5dSHong Zhang     }
391077c4feaSHong Zhang   }
392c9681e5dSHong Zhang   if (ts->vecs_sensip) {
3933e05ceb1SHong Zhang     th->shift = 1./th->time_step;;
3943e05ceb1SHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
3957207e0fdSHong Zhang     if (quadts) {
396cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
3977207e0fdSHong Zhang     }
398c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
399c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
400b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
40133f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
402b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
403c9681e5dSHong Zhang     }
404cd4cee2dSHong Zhang 
405c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
406c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
407c9681e5dSHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
408cd4cee2dSHong Zhang       if (quadJp) {
409cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
410cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
411cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
412cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
413cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
41433f72c5dSHong Zhang       }
415c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
416c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
417c9681e5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
418c9681e5dSHong Zhang         if (ts->vecs_fpu) {
419c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
420c9681e5dSHong Zhang         }
421c9681e5dSHong Zhang         if (ts->vecs_fpp) {
422c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
423c9681e5dSHong Zhang         }
424c9681e5dSHong Zhang       }
425c9681e5dSHong Zhang     }
426c9681e5dSHong Zhang   }
427c9681e5dSHong Zhang 
428c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
429c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
430c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
431c9681e5dSHong Zhang   }
432c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
433c9681e5dSHong Zhang   PetscFunctionReturn(0);
434c9681e5dSHong Zhang }
435c9681e5dSHong Zhang 
43642f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4372ca6e920SHong Zhang {
4382ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
439cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
440b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
441b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4422ca6e920SHong Zhang   PetscInt       nadj;
4437207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4442ca6e920SHong Zhang   KSP            ksp;
445b5b4017aSHong Zhang   PetscScalar    *xarr;
446b5b4017aSHong Zhang   PetscErrorCode ierr;
4472ca6e920SHong Zhang 
4482ca6e920SHong Zhang   PetscFunctionBegin;
449c9681e5dSHong Zhang   if (th->Theta == 1.) {
450c9681e5dSHong Zhang     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
451c9681e5dSHong Zhang     PetscFunctionReturn(0);
452c9681e5dSHong Zhang   }
4532ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
454302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
455cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
4567207e0fdSHong Zhang   if (quadts) {
457cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
458cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
4597207e0fdSHong Zhang   }
460b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
461022c081aSHong Zhang   th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/
462b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
463022c081aSHong Zhang   th->time_step  = -ts->time_step;
4642ca6e920SHong Zhang 
465b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
46633f72c5dSHong Zhang   /* Cost function has an integral term */
4677207e0fdSHong Zhang   if (quadts) {
46805755b9cSHong Zhang     if (th->endpoint) {
469cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
47005755b9cSHong Zhang     } else {
471cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
47205755b9cSHong Zhang     }
4737207e0fdSHong Zhang   }
47433f72c5dSHong Zhang 
475abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
476b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
477022c081aSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr);
478cd4cee2dSHong Zhang     if (quadJ) {
479cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
480cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
481cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
482cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
483cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
48436eaed60SHong Zhang     }
4852ca6e920SHong Zhang   }
4863c54f92cSHong Zhang 
487b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
4883e05ceb1SHong Zhang   th->shift = 1./(th->Theta*th->time_step);
4893c54f92cSHong Zhang   if (th->endpoint) {
490c586f2e8SHong Zhang     ierr = KSPTSFormOperator_Private(ksp,ts->vec_sol,J,Jpre,ts);CHKERRQ(ierr);
4913c54f92cSHong Zhang   } else {
492c586f2e8SHong Zhang     ierr = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr);
4933c54f92cSHong Zhang   }
494cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
4952ca6e920SHong Zhang 
496b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
497abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4981d14d03bSHong Zhang     KSPConvergedReason kspreason;
499b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
5001d14d03bSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
5011d14d03bSHong Zhang     if (kspreason < 0) {
5021d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
5031d14d03bSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
5041d14d03bSHong Zhang     }
5052ca6e920SHong Zhang   }
5063c54f92cSHong Zhang 
5071d14d03bSHong Zhang   /* Second-order adjoint */
508b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
509b5b4017aSHong Zhang     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
510b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
511b5b4017aSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
512b5b4017aSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
513b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
514b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
515b5b4017aSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
516b5b4017aSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
517b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
518b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
519b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
520b5b4017aSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
5213e05ceb1SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],th->shift);CHKERRQ(ierr);
52233f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
523b5b4017aSHong Zhang       if (ts->vecs_fup) {
52433f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
525b5b4017aSHong Zhang       }
526b5b4017aSHong Zhang     }
527b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
528b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
5291d14d03bSHong Zhang       KSPConvergedReason kspreason;
5301d14d03bSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
5311d14d03bSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
5321d14d03bSHong Zhang       if (kspreason < 0) {
5331d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
5341d14d03bSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
5351d14d03bSHong Zhang       }
536b5b4017aSHong Zhang     }
537b5b4017aSHong Zhang   }
538b5b4017aSHong Zhang 
53936eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5401d14d03bSHong Zhang   if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
5413e05ceb1SHong Zhang     th->shift      = 1./((th->Theta-1.)*th->time_step);
542c586f2e8SHong Zhang     th->stage_time = th->ptime;
543c586f2e8SHong Zhang     ierr           = KSPTSFormOperator_Private(ksp,th->X0,J,Jpre,ts);CHKERRQ(ierr);
54433f72c5dSHong Zhang     /* R_U at t_n */
5457207e0fdSHong Zhang     if (quadts) {
546cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
5477207e0fdSHong Zhang     }
548abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
549b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
550cd4cee2dSHong Zhang       if (quadJ) {
551cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
552cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
553cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr);
554cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
555cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
556b5b4017aSHong Zhang       }
5573e05ceb1SHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./th->shift);CHKERRQ(ierr);
558b5b4017aSHong Zhang     }
5591d14d03bSHong Zhang 
5601d14d03bSHong Zhang     /* Second-order adjoint */
5611d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
562b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
563b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
564b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
565b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
566b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
567b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
56833f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
569b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
570b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
571b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
57233f72c5dSHong 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  */
573b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
57433f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
575b5b4017aSHong Zhang         if (ts->vecs_fup) {
57633f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
577b5b4017aSHong Zhang         }
5783e05ceb1SHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./th->shift);CHKERRQ(ierr);
579b5b4017aSHong Zhang       }
580b5e0532dSHong Zhang     }
5813fd52205SHong Zhang 
582c586f2e8SHong Zhang     th->stage_time = ts->ptime; /* recover the old value */
583c586f2e8SHong Zhang 
5843fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
585b5b4017aSHong Zhang       /* U_{n+1} */
5863e05ceb1SHong Zhang       th->shift = -1./(th->Theta*th->time_step);
5873e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5887207e0fdSHong Zhang       if (quadts) {
589cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
5907207e0fdSHong Zhang       }
591abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
592b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
59333f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5943fd52205SHong Zhang       }
59533f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
59633f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
59733f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
59833f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
599b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
600b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
60133f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
60233f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
60333f72c5dSHong Zhang 
604b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
605b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
606b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
60733f72c5dSHong 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)  */
608b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
60933f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
610b5b4017aSHong Zhang           if (ts->vecs_fpu) {
61133f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
612b5b4017aSHong Zhang           }
613b5b4017aSHong Zhang           if (ts->vecs_fpp) {
61433f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
615b5b4017aSHong Zhang           }
616b5b4017aSHong Zhang         }
617b5b4017aSHong Zhang       }
618b5b4017aSHong Zhang 
619b5b4017aSHong Zhang       /* U_s */
6203e05ceb1SHong Zhang       th->shift = 1./((th->Theta-1.0)*th->time_step);
6213e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6227207e0fdSHong Zhang       if (quadts) {
623cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
6247207e0fdSHong Zhang       }
625abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
626b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
62733f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
62833f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
62933f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
63033f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
63133f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
632b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
633b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
63433f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
63533f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
636b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
637b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
638b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
63933f72c5dSHong 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) */
640b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
64133f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
642b5b4017aSHong Zhang             if (ts->vecs_fpu) {
64333f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
644b5e0532dSHong Zhang             }
645b5b4017aSHong Zhang             if (ts->vecs_fpp) {
64633f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
647b5b4017aSHong Zhang             }
64836eaed60SHong Zhang           }
64936eaed60SHong Zhang         }
6503fd52205SHong Zhang       }
651b5e0532dSHong Zhang     }
6523fd52205SHong Zhang   } else { /* one-stage case */
6533e05ceb1SHong Zhang     th->shift = 0.0;
654c586f2e8SHong Zhang     ierr      = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr); /* get -f_y*/
6557207e0fdSHong Zhang     if (quadts) {
656cd4cee2dSHong Zhang       ierr  = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
6577207e0fdSHong Zhang     }
658abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
659b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
660022c081aSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
661cd4cee2dSHong Zhang       if (quadJ) {
662cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
663cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
664cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr);
665cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
666cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
66736eaed60SHong Zhang       }
6682ca6e920SHong Zhang     }
6693fd52205SHong Zhang     if (ts->vecs_sensip) {
6703e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6717207e0fdSHong Zhang       if (quadts) {
672cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
6737207e0fdSHong Zhang       }
674abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
67533f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
67633f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
677cd4cee2dSHong Zhang         if (quadJp) {
678cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
679cd4cee2dSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
680cd4cee2dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
681cd4cee2dSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
682cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
68336eaed60SHong Zhang         }
68436eaed60SHong Zhang       }
6853fd52205SHong Zhang     }
6862ca6e920SHong Zhang   }
6872ca6e920SHong Zhang 
6882ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6892ca6e920SHong Zhang   PetscFunctionReturn(0);
6902ca6e920SHong Zhang }
6912ca6e920SHong Zhang 
692cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
693cd652676SJed Brown {
694cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
695be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
696cd652676SJed Brown   PetscErrorCode ierr;
697cd652676SJed Brown 
698cd652676SJed Brown   PetscFunctionBegin;
699a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
700be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
701be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
702cd652676SJed Brown   PetscFunctionReturn(0);
703cd652676SJed Brown }
704cd652676SJed Brown 
7051566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
7061566a47fSLisandro Dalcin {
7071566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7081566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
7091566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
7107453f775SEmil Constantinescu   PetscReal      wltea,wlter;
7111566a47fSLisandro Dalcin   PetscErrorCode ierr;
7121566a47fSLisandro Dalcin 
7131566a47fSLisandro Dalcin   PetscFunctionBegin;
7142ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
7151566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
716fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
7171566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
718fecfb714SLisandro Dalcin   {
719be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
720be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
7211566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
7221566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
7231566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
7241566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
7251566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
7267453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
7271566a47fSLisandro Dalcin   }
7281566a47fSLisandro Dalcin   if (order) *order = 2;
7291566a47fSLisandro Dalcin   PetscFunctionReturn(0);
7301566a47fSLisandro Dalcin }
7311566a47fSLisandro Dalcin 
7321566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7331566a47fSLisandro Dalcin {
7341566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7357207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7361566a47fSLisandro Dalcin   PetscErrorCode ierr;
7371566a47fSLisandro Dalcin 
7381566a47fSLisandro Dalcin   PetscFunctionBegin;
7391566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
740cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
741cd4cee2dSHong Zhang     ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr);
7421566a47fSLisandro Dalcin   }
743715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
74413b1b0ffSHong Zhang   if (ts->mat_sensip) {
74513b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7466f92202bSHong Zhang   }
7477207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7487207e0fdSHong Zhang     ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7496f92202bSHong Zhang   }
750715f1b00SHong Zhang   PetscFunctionReturn(0);
751715f1b00SHong Zhang }
752715f1b00SHong Zhang 
753715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
754715f1b00SHong Zhang {
755715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
756cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
75713b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
758b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7597207e0fdSHong Zhang   PetscInt       ntlm;
760715f1b00SHong Zhang   KSP            ksp;
7617207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
76213b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
763715f1b00SHong Zhang   PetscErrorCode ierr;
764715f1b00SHong Zhang 
765715f1b00SHong Zhang   PetscFunctionBegin;
76613b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
76713b1b0ffSHong Zhang 
7687207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7697207e0fdSHong Zhang     ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7706f92202bSHong Zhang   }
771715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
772cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
7737207e0fdSHong Zhang   if (quadts) {
774cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
775cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
7767207e0fdSHong Zhang   }
777715f1b00SHong Zhang 
778715f1b00SHong Zhang   /* Build RHS */
779715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
7803e05ceb1SHong Zhang     th->shift = 1./((th->Theta-1.)*th->time_step);
7813e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
78213b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
78313b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
784715f1b00SHong Zhang 
785022c081aSHong Zhang     /* Add the f_p forcing terms */
7860e11c21fSHong Zhang     if (ts->Jacp) {
7873e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
78833f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7893e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
79033f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7910e11c21fSHong Zhang     }
792715f1b00SHong Zhang   } else { /* 1-stage method */
7933e05ceb1SHong Zhang     th->shift = 0.0;
7943e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
79513b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
79613b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
79713b1b0ffSHong Zhang 
798022c081aSHong Zhang     /* Add the f_p forcing terms */
7990e11c21fSHong Zhang     if (ts->Jacp) {
8003e05ceb1SHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,th->shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
80133f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
802715f1b00SHong Zhang     }
8030e11c21fSHong Zhang   }
804715f1b00SHong Zhang 
805715f1b00SHong Zhang   /* Build LHS */
8063e05ceb1SHong Zhang   th->shift = 1/(th->Theta*th->time_step);
807715f1b00SHong Zhang   if (th->endpoint) {
8083e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
809715f1b00SHong Zhang   } else {
8103e05ceb1SHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,th->shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
811715f1b00SHong Zhang   }
812cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
813715f1b00SHong Zhang 
814715f1b00SHong Zhang   /*
815715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
816c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
817715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
818715f1b00SHong Zhang   */
819715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8207207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
821cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
822cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
8237207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8247207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8257207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
826715f1b00SHong Zhang     }
827715f1b00SHong Zhang   }
828715f1b00SHong Zhang 
829715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
830022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
831b5b4017aSHong Zhang     KSPConvergedReason kspreason;
83213b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
833b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
834715f1b00SHong Zhang     if (th->endpoint) {
83513b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
836b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
837b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
838b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
83913b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
840715f1b00SHong Zhang     } else {
841b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
842715f1b00SHong Zhang     }
843b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
844b5b4017aSHong Zhang     if (kspreason < 0) {
845b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
846b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
847b5b4017aSHong Zhang     }
848b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
84913b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
850715f1b00SHong Zhang   }
851715f1b00SHong Zhang 
852b5b4017aSHong Zhang 
85313b1b0ffSHong Zhang   /*
85413b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
855c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
85613b1b0ffSHong Zhang   */
8577207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
85813b1b0ffSHong Zhang     if (!th->endpoint) {
8597207e0fdSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */
860cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
861cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
8627207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8637207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8647207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
86513b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
86613b1b0ffSHong Zhang     } else {
867cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
868cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
8697207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8707207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8717207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
87213b1b0ffSHong Zhang     }
87313b1b0ffSHong Zhang   } else {
87413b1b0ffSHong Zhang     if (!th->endpoint) {
87513b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
876715f1b00SHong Zhang     }
877715f1b00SHong Zhang   }
8781566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8791566a47fSLisandro Dalcin }
8801566a47fSLisandro Dalcin 
8816affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
8826affb6f8SHong Zhang {
8836affb6f8SHong Zhang   TS_Theta *th = (TS_Theta*)ts->data;
8846affb6f8SHong Zhang 
8856affb6f8SHong Zhang   PetscFunctionBegin;
8866affb6f8SHong Zhang   if (ns) *ns = 1;
8876affb6f8SHong Zhang   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
8886affb6f8SHong Zhang   PetscFunctionReturn(0);
8896affb6f8SHong Zhang }
8906affb6f8SHong Zhang 
891316643e7SJed Brown /*------------------------------------------------------------*/
892277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
893316643e7SJed Brown {
894316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
895316643e7SJed Brown   PetscErrorCode ierr;
896316643e7SJed Brown 
897316643e7SJed Brown   PetscFunctionBegin;
8986bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
8996bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
9003b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
901eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
9021566a47fSLisandro Dalcin 
9031566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
9041566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
9051566a47fSLisandro Dalcin 
906b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
907277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
908277b19d0SLisandro Dalcin }
909277b19d0SLisandro Dalcin 
910ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
911ece46509SHong Zhang {
912ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
913ece46509SHong Zhang   PetscErrorCode ierr;
914ece46509SHong Zhang 
915ece46509SHong Zhang   PetscFunctionBegin;
916ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
917ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
918ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
919ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
920ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
921ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
922ece46509SHong Zhang   PetscFunctionReturn(0);
923ece46509SHong Zhang }
924ece46509SHong Zhang 
925277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
926277b19d0SLisandro Dalcin {
927277b19d0SLisandro Dalcin   PetscErrorCode ierr;
928277b19d0SLisandro Dalcin 
929277b19d0SLisandro Dalcin   PetscFunctionBegin;
930277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
931b3a6b972SJed Brown   if (ts->dm) {
932b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
933b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
934b3a6b972SJed Brown   }
935277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
936bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
937bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
938bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
939bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
940316643e7SJed Brown   PetscFunctionReturn(0);
941316643e7SJed Brown }
942316643e7SJed Brown 
943316643e7SJed Brown /*
944316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9452b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
946316643e7SJed Brown */
9470f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
948316643e7SJed Brown {
949316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
950316643e7SJed Brown   PetscErrorCode ierr;
9517445fe48SJed Brown   Vec            X0,Xdot;
9527445fe48SJed Brown   DM             dm,dmsave;
9533e05ceb1SHong Zhang   PetscReal      shift = th->shift;
954316643e7SJed Brown 
955316643e7SJed Brown   PetscFunctionBegin;
9567445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9575a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9587445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
959*494ce9fcSHong Zhang   if (x != X0) {
960b296d7d5SJed Brown     ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
961*494ce9fcSHong Zhang   } else {
962*494ce9fcSHong Zhang     ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
963*494ce9fcSHong Zhang   }
9647445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9657445fe48SJed Brown   dmsave = ts->dm;
9667445fe48SJed Brown   ts->dm = dm;
9677445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
9687445fe48SJed Brown   ts->dm = dmsave;
9690d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
970316643e7SJed Brown   PetscFunctionReturn(0);
971316643e7SJed Brown }
972316643e7SJed Brown 
973d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
974316643e7SJed Brown {
975316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
976316643e7SJed Brown   PetscErrorCode ierr;
9777445fe48SJed Brown   Vec            Xdot;
9787445fe48SJed Brown   DM             dm,dmsave;
9793e05ceb1SHong Zhang   PetscReal      shift = th->shift;
980316643e7SJed Brown 
981316643e7SJed Brown   PetscFunctionBegin;
9827445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
983be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9840298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
9857445fe48SJed Brown 
9867445fe48SJed Brown   dmsave = ts->dm;
9877445fe48SJed Brown   ts->dm = dm;
988d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
9897445fe48SJed Brown   ts->dm = dmsave;
9900298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
991316643e7SJed Brown   PetscFunctionReturn(0);
992316643e7SJed Brown }
993316643e7SJed Brown 
994715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
995715f1b00SHong Zhang {
996715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9977207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
998715f1b00SHong Zhang   PetscErrorCode ierr;
999715f1b00SHong Zhang 
1000715f1b00SHong Zhang   PetscFunctionBegin;
1001022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
100213b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
100313b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10047207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10057207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10067207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
1007715f1b00SHong Zhang   }
10086f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
100913b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
101013b1b0ffSHong Zhang 
1011b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
1012715f1b00SHong Zhang   PetscFunctionReturn(0);
1013715f1b00SHong Zhang }
1014715f1b00SHong Zhang 
10157adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
10167adebddeSHong Zhang {
10177adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10187207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10197adebddeSHong Zhang   PetscErrorCode ierr;
10207adebddeSHong Zhang 
10217adebddeSHong Zhang   PetscFunctionBegin;
10227207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10237207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10247207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
10257adebddeSHong Zhang   }
10267adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
10277adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10287adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
10297adebddeSHong Zhang   PetscFunctionReturn(0);
10307adebddeSHong Zhang }
10317adebddeSHong Zhang 
1032316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1033316643e7SJed Brown {
1034316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1035cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10362ffb9264SLisandro Dalcin   PetscBool      match;
1037316643e7SJed Brown   PetscErrorCode ierr;
1038316643e7SJed Brown 
1039316643e7SJed Brown   PetscFunctionBegin;
1040cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1041cd4cee2dSHong Zhang     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1042b1cb36f3SHong Zhang   }
104339d13666SHong Zhang   if (!th->X) {
1044316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
104539d13666SHong Zhang   }
104639d13666SHong Zhang   if (!th->Xdot) {
1047316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
104839d13666SHong Zhang   }
104939d13666SHong Zhang   if (!th->X0) {
10503b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
105139d13666SHong Zhang   }
10521566a47fSLisandro Dalcin   if (th->endpoint) {
10531566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
10547445fe48SJed Brown   }
10553b1890cdSShri Abhyankar 
10561566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
10571566a47fSLisandro Dalcin 
10581566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
10591566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10601566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10611566a47fSLisandro Dalcin 
10621566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
10631566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
10642ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
10652ffb9264SLisandro Dalcin   if (!match) {
10661566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
10671566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
10683b1890cdSShri Abhyankar   }
10691566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1070b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1071b4dd3bf9SBarry Smith }
10720c7376e5SHong Zhang 
1073b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1074b4dd3bf9SBarry Smith 
1075b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1076b4dd3bf9SBarry Smith {
1077b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1078b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1079b4dd3bf9SBarry Smith 
1080b4dd3bf9SBarry Smith   PetscFunctionBegin;
1081b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1082b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
10832ca6e920SHong Zhang   if (ts->vecs_sensip) {
1084b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
10852ca6e920SHong Zhang   }
1086b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1087b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1088b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
108967633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
109067633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
109167633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1092b5b4017aSHong Zhang   }
1093c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1094c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
109567633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
109667633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
109767633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1098b5b4017aSHong Zhang   }
1099316643e7SJed Brown   PetscFunctionReturn(0);
1100316643e7SJed Brown }
1101316643e7SJed Brown 
11024416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1103316643e7SJed Brown {
1104316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1105316643e7SJed Brown   PetscErrorCode ierr;
1106316643e7SJed Brown 
1107316643e7SJed Brown   PetscFunctionBegin;
1108e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1109316643e7SJed Brown   {
11100298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
11110298fd71SBarry 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);
111203842d09SLisandro 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);
1113316643e7SJed Brown   }
1114316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1115316643e7SJed Brown   PetscFunctionReturn(0);
1116316643e7SJed Brown }
1117316643e7SJed Brown 
1118316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1119316643e7SJed Brown {
1120316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1121ace3abfcSBarry Smith   PetscBool      iascii;
1122316643e7SJed Brown   PetscErrorCode ierr;
1123316643e7SJed Brown 
1124316643e7SJed Brown   PetscFunctionBegin;
1125251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1126316643e7SJed Brown   if (iascii) {
11277c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1128ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1129316643e7SJed Brown   }
1130316643e7SJed Brown   PetscFunctionReturn(0);
1131316643e7SJed Brown }
1132316643e7SJed Brown 
1133560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11340de4c49aSJed Brown {
11350de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11360de4c49aSJed Brown 
11370de4c49aSJed Brown   PetscFunctionBegin;
11380de4c49aSJed Brown   *theta = th->Theta;
11390de4c49aSJed Brown   PetscFunctionReturn(0);
11400de4c49aSJed Brown }
11410de4c49aSJed Brown 
1142560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11430de4c49aSJed Brown {
11440de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11450de4c49aSJed Brown 
11460de4c49aSJed Brown   PetscFunctionBegin;
11477c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11480de4c49aSJed Brown   th->Theta = theta;
11491566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11500de4c49aSJed Brown   PetscFunctionReturn(0);
11510de4c49aSJed Brown }
1152eb284becSJed Brown 
1153560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
115426f2ff8fSLisandro Dalcin {
115526f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
115626f2ff8fSLisandro Dalcin 
115726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
115826f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
115926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
116026f2ff8fSLisandro Dalcin }
116126f2ff8fSLisandro Dalcin 
1162560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1163eb284becSJed Brown {
1164eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1165eb284becSJed Brown 
1166eb284becSJed Brown   PetscFunctionBegin;
1167eb284becSJed Brown   th->endpoint = flg;
1168eb284becSJed Brown   PetscFunctionReturn(0);
1169eb284becSJed Brown }
11700de4c49aSJed Brown 
1171f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1172f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1173f9c1d6abSBarry Smith {
1174f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1175f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11763fd8ae06SJed Brown   const PetscReal one = 1.0;
1177f9c1d6abSBarry Smith 
1178f9c1d6abSBarry Smith   PetscFunctionBegin;
11793fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1180f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1181f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1182f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1183f9c1d6abSBarry Smith }
1184f9c1d6abSBarry Smith #endif
1185f9c1d6abSBarry Smith 
118642682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
118742682096SHong Zhang {
118842682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
118942682096SHong Zhang 
119042682096SHong Zhang   PetscFunctionBegin;
11911566a47fSLisandro Dalcin   if (ns) *ns = 1;
11921566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
119342682096SHong Zhang   PetscFunctionReturn(0);
119442682096SHong Zhang }
1195f9c1d6abSBarry Smith 
1196316643e7SJed Brown /* ------------------------------------------------------------ */
1197316643e7SJed Brown /*MC
119896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1199316643e7SJed Brown 
1200316643e7SJed Brown    Level: beginner
1201316643e7SJed Brown 
12024eb428fdSBarry Smith    Options Database:
1203c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1204c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
120503842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
12064eb428fdSBarry Smith 
1207eb284becSJed Brown    Notes:
12080c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
12090c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
12104eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
12114eb428fdSBarry Smith 
1212eb284becSJed Brown    This method can be applied to DAE.
1213eb284becSJed Brown 
1214eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
1215eb284becSJed Brown 
1216eb284becSJed Brown .vb
1217eb284becSJed Brown   Theta | Theta
1218eb284becSJed Brown   -------------
1219eb284becSJed Brown         |  1
1220eb284becSJed Brown .ve
1221eb284becSJed Brown 
1222eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1223eb284becSJed Brown 
1224eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1225eb284becSJed Brown 
1226eb284becSJed Brown .vb
1227eb284becSJed Brown   0 | 0         0
1228eb284becSJed Brown   1 | 1-Theta   Theta
1229eb284becSJed Brown   -------------------
1230eb284becSJed Brown     | 1-Theta   Theta
1231eb284becSJed Brown .ve
1232eb284becSJed Brown 
1233eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1234eb284becSJed Brown 
1235eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1236eb284becSJed Brown 
1237eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1238eb284becSJed Brown 
12394eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1240eb284becSJed Brown 
1241eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1242316643e7SJed Brown 
1243316643e7SJed Brown M*/
12448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1245316643e7SJed Brown {
1246316643e7SJed Brown   TS_Theta       *th;
1247316643e7SJed Brown   PetscErrorCode ierr;
1248316643e7SJed Brown 
1249316643e7SJed Brown   PetscFunctionBegin;
1250277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1251ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1252316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1253316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1254316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
125542f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1256ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1257316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1258cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12591566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
126024655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1261316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12620f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12630f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1264f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1265f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1266f9c1d6abSBarry Smith #endif
126742682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
126842f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1269b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1270b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12712ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12726affb6f8SHong Zhang 
1273715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12747adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1275715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12766affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1277316643e7SJed Brown 
1278efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1279efd4aadfSBarry Smith 
1280b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1281316643e7SJed Brown   ts->data = (void*)th;
1282316643e7SJed Brown 
1283715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1284715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1285715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1286b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1287715f1b00SHong Zhang 
12886f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1289316643e7SJed Brown   th->Theta       = 0.5;
12901566a47fSLisandro Dalcin   th->order       = 2;
1291bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1292bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1293bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1294bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1295316643e7SJed Brown   PetscFunctionReturn(0);
1296316643e7SJed Brown }
12970de4c49aSJed Brown 
12980de4c49aSJed Brown /*@
12990de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
13000de4c49aSJed Brown 
13010de4c49aSJed Brown   Not Collective
13020de4c49aSJed Brown 
13030de4c49aSJed Brown   Input Parameter:
13040de4c49aSJed Brown .  ts - timestepping context
13050de4c49aSJed Brown 
13060de4c49aSJed Brown   Output Parameter:
13070de4c49aSJed Brown .  theta - stage abscissa
13080de4c49aSJed Brown 
13090de4c49aSJed Brown   Note:
13100de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
13110de4c49aSJed Brown 
13120de4c49aSJed Brown   Level: Advanced
13130de4c49aSJed Brown 
13140de4c49aSJed Brown .seealso: TSThetaSetTheta()
13150de4c49aSJed Brown @*/
13167087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13170de4c49aSJed Brown {
13184ac538c5SBarry Smith   PetscErrorCode ierr;
13190de4c49aSJed Brown 
13200de4c49aSJed Brown   PetscFunctionBegin;
1321afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13220de4c49aSJed Brown   PetscValidPointer(theta,2);
13234ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
13240de4c49aSJed Brown   PetscFunctionReturn(0);
13250de4c49aSJed Brown }
13260de4c49aSJed Brown 
13270de4c49aSJed Brown /*@
13280de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13290de4c49aSJed Brown 
13300de4c49aSJed Brown   Not Collective
13310de4c49aSJed Brown 
13320de4c49aSJed Brown   Input Parameter:
13330de4c49aSJed Brown +  ts - timestepping context
13340de4c49aSJed Brown -  theta - stage abscissa
13350de4c49aSJed Brown 
13360de4c49aSJed Brown   Options Database:
13370de4c49aSJed Brown .  -ts_theta_theta <theta>
13380de4c49aSJed Brown 
13390de4c49aSJed Brown   Level: Intermediate
13400de4c49aSJed Brown 
13410de4c49aSJed Brown .seealso: TSThetaGetTheta()
13420de4c49aSJed Brown @*/
13437087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13440de4c49aSJed Brown {
13454ac538c5SBarry Smith   PetscErrorCode ierr;
13460de4c49aSJed Brown 
13470de4c49aSJed Brown   PetscFunctionBegin;
1348afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13494ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
13500de4c49aSJed Brown   PetscFunctionReturn(0);
13510de4c49aSJed Brown }
1352f33bbcb6SJed Brown 
135326f2ff8fSLisandro Dalcin /*@
135426f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
135526f2ff8fSLisandro Dalcin 
135626f2ff8fSLisandro Dalcin   Not Collective
135726f2ff8fSLisandro Dalcin 
135826f2ff8fSLisandro Dalcin   Input Parameter:
135926f2ff8fSLisandro Dalcin .  ts - timestepping context
136026f2ff8fSLisandro Dalcin 
136126f2ff8fSLisandro Dalcin   Output Parameter:
136226f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
136326f2ff8fSLisandro Dalcin 
136426f2ff8fSLisandro Dalcin   Level: Advanced
136526f2ff8fSLisandro Dalcin 
136626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
136726f2ff8fSLisandro Dalcin @*/
136826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
136926f2ff8fSLisandro Dalcin {
137026f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
137126f2ff8fSLisandro Dalcin 
137226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
137326f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
137426f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1375163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
137626f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
137726f2ff8fSLisandro Dalcin }
137826f2ff8fSLisandro Dalcin 
1379eb284becSJed Brown /*@
1380eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1381eb284becSJed Brown 
1382eb284becSJed Brown   Not Collective
1383eb284becSJed Brown 
1384eb284becSJed Brown   Input Parameter:
1385eb284becSJed Brown +  ts - timestepping context
1386eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1387eb284becSJed Brown 
1388eb284becSJed Brown   Options Database:
1389eb284becSJed Brown .  -ts_theta_endpoint <flg>
1390eb284becSJed Brown 
1391eb284becSJed Brown   Level: Intermediate
1392eb284becSJed Brown 
1393eb284becSJed Brown .seealso: TSTHETA, TSCN
1394eb284becSJed Brown @*/
1395eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1396eb284becSJed Brown {
1397eb284becSJed Brown   PetscErrorCode ierr;
1398eb284becSJed Brown 
1399eb284becSJed Brown   PetscFunctionBegin;
1400eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1401eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1402eb284becSJed Brown   PetscFunctionReturn(0);
1403eb284becSJed Brown }
1404eb284becSJed Brown 
1405f33bbcb6SJed Brown /*
1406f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1407f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1408f33bbcb6SJed Brown  */
1409f33bbcb6SJed Brown 
14101566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
14111566a47fSLisandro Dalcin {
14121566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14131566a47fSLisandro Dalcin   PetscErrorCode ierr;
14141566a47fSLisandro Dalcin 
14151566a47fSLisandro Dalcin   PetscFunctionBegin;
14161566a47fSLisandro 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");
14171566a47fSLisandro 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");
14181566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14191566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14201566a47fSLisandro Dalcin }
14211566a47fSLisandro Dalcin 
1422f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1423f33bbcb6SJed Brown {
1424f33bbcb6SJed Brown   PetscFunctionBegin;
1425f33bbcb6SJed Brown   PetscFunctionReturn(0);
1426f33bbcb6SJed Brown }
1427f33bbcb6SJed Brown 
1428f33bbcb6SJed Brown /*MC
1429f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1430f33bbcb6SJed Brown 
1431f33bbcb6SJed Brown   Level: beginner
1432f33bbcb6SJed Brown 
14334eb428fdSBarry Smith   Notes:
1434c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14354eb428fdSBarry Smith 
14361566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14374eb428fdSBarry Smith 
1438f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1439f33bbcb6SJed Brown 
1440f33bbcb6SJed Brown M*/
14418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1442f33bbcb6SJed Brown {
1443f33bbcb6SJed Brown   PetscErrorCode ierr;
1444f33bbcb6SJed Brown 
1445f33bbcb6SJed Brown   PetscFunctionBegin;
1446f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1447f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
14481566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
14490c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1450f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1451f33bbcb6SJed Brown   PetscFunctionReturn(0);
1452f33bbcb6SJed Brown }
1453f33bbcb6SJed Brown 
14541566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(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 != 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");
14611566a47fSLisandro 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");
14621566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14631566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14641566a47fSLisandro Dalcin }
14651566a47fSLisandro Dalcin 
1466f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1467f33bbcb6SJed Brown {
1468f33bbcb6SJed Brown   PetscFunctionBegin;
1469f33bbcb6SJed Brown   PetscFunctionReturn(0);
1470f33bbcb6SJed Brown }
1471f33bbcb6SJed Brown 
1472f33bbcb6SJed Brown /*MC
1473f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1474f33bbcb6SJed Brown 
1475f33bbcb6SJed Brown   Level: beginner
1476f33bbcb6SJed Brown 
1477f33bbcb6SJed Brown   Notes:
14787cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14797cf5af47SJed Brown 
14807cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1481f33bbcb6SJed Brown 
1482f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1483f33bbcb6SJed Brown 
1484f33bbcb6SJed Brown M*/
14858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(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,0.5);CHKERRQ(ierr);
1492eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
14930c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1494f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1495f33bbcb6SJed Brown   PetscFunctionReturn(0);
1496f33bbcb6SJed Brown }
1497