xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 6980b70547907d6716539c2ecc7abefc899c03de)
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;
171566a47fSLisandro Dalcin   PetscInt     order;
181566a47fSLisandro Dalcin   PetscBool    endpoint;
191566a47fSLisandro Dalcin   PetscBool    extrapolate;
20715f1b00SHong Zhang   TSStepStatus status;
21715f1b00SHong Zhang   Vec          VecCostIntegral0;         /* Backup for roll-backs due to events */
221566a47fSLisandro Dalcin 
23715f1b00SHong Zhang   /* context for sensitivity analysis */
24022c081aSHong Zhang   PetscInt     num_tlm;                  /* Total number of tangent linear equations */
25b5e0532dSHong Zhang   Vec          *VecsDeltaLam;            /* Increment of the adjoint sensitivity w.r.t IC at stage */
26b5e0532dSHong Zhang   Vec          *VecsDeltaMu;             /* Increment of the adjoint sensitivity w.r.t P at stage */
2713b1b0ffSHong Zhang   Vec          *VecsSensiTemp;           /* Vector to be multiplied with Jacobian transpose */
2813b1b0ffSHong Zhang   Mat          MatDeltaFwdSensip;        /* Increment of the forward sensitivity at stage */
29b5b4017aSHong Zhang   Vec          VecDeltaFwdSensipCol;     /* Working vector for holding one column of the sensitivity matrix */
3013b1b0ffSHong Zhang   Mat          MatFwdSensip0;            /* backup for roll-backs due to events */
317207e0fdSHong Zhang   Mat          MatIntegralSensipTemp;    /* Working vector for forward integral sensitivity */
327207e0fdSHong Zhang   Mat          MatIntegralSensip0;       /* backup for roll-backs due to events */
33b5b4017aSHong Zhang   Vec          *VecsDeltaLam2;           /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
34b5b4017aSHong Zhang   Vec          *VecsDeltaMu2;            /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
35b5b4017aSHong Zhang   Vec          *VecsSensi2Temp;          /* Working vectors that holds the residual for the second-order adjoint */
36b5b4017aSHong Zhang   Vec          *VecsAffine;              /* Working vectors to store residuals */
37715f1b00SHong Zhang   /* context for error estimation */
381566a47fSLisandro Dalcin   Vec          vec_sol_prev;
391566a47fSLisandro Dalcin   Vec          vec_lte_work;
40316643e7SJed Brown } TS_Theta;
41316643e7SJed Brown 
427445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
437445fe48SJed Brown {
447445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
457445fe48SJed Brown   PetscErrorCode ierr;
467445fe48SJed Brown 
477445fe48SJed Brown   PetscFunctionBegin;
487445fe48SJed Brown   if (X0) {
497445fe48SJed Brown     if (dm && dm != ts->dm) {
500d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
517445fe48SJed Brown     } else *X0 = ts->vec_sol;
527445fe48SJed Brown   }
537445fe48SJed Brown   if (Xdot) {
547445fe48SJed Brown     if (dm && dm != ts->dm) {
550d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
567445fe48SJed Brown     } else *Xdot = th->Xdot;
577445fe48SJed Brown   }
587445fe48SJed Brown   PetscFunctionReturn(0);
597445fe48SJed Brown }
607445fe48SJed Brown 
610d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
620d0b770aSPeter Brune {
630d0b770aSPeter Brune   PetscErrorCode ierr;
640d0b770aSPeter Brune 
650d0b770aSPeter Brune   PetscFunctionBegin;
660d0b770aSPeter Brune   if (X0) {
670d0b770aSPeter Brune     if (dm && dm != ts->dm) {
680d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
690d0b770aSPeter Brune     }
700d0b770aSPeter Brune   }
710d0b770aSPeter Brune   if (Xdot) {
720d0b770aSPeter Brune     if (dm && dm != ts->dm) {
730d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
740d0b770aSPeter Brune     }
750d0b770aSPeter Brune   }
760d0b770aSPeter Brune   PetscFunctionReturn(0);
770d0b770aSPeter Brune }
780d0b770aSPeter Brune 
797445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
807445fe48SJed Brown {
817445fe48SJed Brown   PetscFunctionBegin;
827445fe48SJed Brown   PetscFunctionReturn(0);
837445fe48SJed Brown }
847445fe48SJed Brown 
857445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
867445fe48SJed Brown {
877445fe48SJed Brown   TS             ts = (TS)ctx;
887445fe48SJed Brown   PetscErrorCode ierr;
897445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
907445fe48SJed Brown 
917445fe48SJed Brown   PetscFunctionBegin;
927445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
937445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
947445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
957445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
967445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
977445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
980d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
990d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
1007445fe48SJed Brown   PetscFunctionReturn(0);
1017445fe48SJed Brown }
1027445fe48SJed Brown 
103258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
104258e1594SPeter Brune {
105258e1594SPeter Brune   PetscFunctionBegin;
106258e1594SPeter Brune   PetscFunctionReturn(0);
107258e1594SPeter Brune }
108258e1594SPeter Brune 
109258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
110258e1594SPeter Brune {
111258e1594SPeter Brune   TS             ts = (TS)ctx;
112258e1594SPeter Brune   PetscErrorCode ierr;
113258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
114258e1594SPeter Brune 
115258e1594SPeter Brune   PetscFunctionBegin;
116258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
117258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
118258e1594SPeter Brune 
119258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune 
122258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
124258e1594SPeter Brune 
125258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
126258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
127258e1594SPeter Brune   PetscFunctionReturn(0);
128258e1594SPeter Brune }
129258e1594SPeter Brune 
130022c081aSHong Zhang static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
131022c081aSHong Zhang {
132022c081aSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
133cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
134022c081aSHong Zhang   PetscErrorCode ierr;
135022c081aSHong Zhang 
136022c081aSHong Zhang   PetscFunctionBegin;
137022c081aSHong Zhang   if (th->endpoint) {
138022c081aSHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
139022c081aSHong Zhang     if (th->Theta!=1.0) {
140cd4cee2dSHong Zhang       ierr = TSComputeRHSFunction(quadts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr);
141cd4cee2dSHong Zhang       ierr = VecAXPY(quadts->vec_sol,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr);
142022c081aSHong Zhang     }
143cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
144cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
145022c081aSHong Zhang   } else {
146cd4cee2dSHong Zhang     ierr = TSComputeRHSFunction(quadts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
147cd4cee2dSHong Zhang     ierr = VecAXPY(quadts->vec_sol,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
148022c081aSHong Zhang   }
149022c081aSHong Zhang   PetscFunctionReturn(0);
150022c081aSHong Zhang }
151022c081aSHong Zhang 
152b1cb36f3SHong Zhang static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
153b1cb36f3SHong Zhang {
154b1cb36f3SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
155cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
156b1cb36f3SHong Zhang   PetscErrorCode ierr;
157b1cb36f3SHong Zhang 
158b1cb36f3SHong Zhang   PetscFunctionBegin;
159b1cb36f3SHong Zhang   /* backup cost integral */
160cd4cee2dSHong Zhang   ierr = VecCopy(quadts->vec_sol,th->VecCostIntegral0);CHKERRQ(ierr);
161022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
162b1cb36f3SHong Zhang   PetscFunctionReturn(0);
163b1cb36f3SHong Zhang }
164b1cb36f3SHong Zhang 
165b1cb36f3SHong Zhang static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
166b1cb36f3SHong Zhang {
167b1cb36f3SHong Zhang   PetscErrorCode ierr;
168b1cb36f3SHong Zhang 
169b1cb36f3SHong Zhang   PetscFunctionBegin;
170022c081aSHong Zhang   ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr);
17124655328SShri   PetscFunctionReturn(0);
17224655328SShri }
17324655328SShri 
174470880abSPatrick Sanan static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x)
1751566a47fSLisandro Dalcin {
1761566a47fSLisandro Dalcin   PetscInt       nits,lits;
1771566a47fSLisandro Dalcin   PetscErrorCode ierr;
1781566a47fSLisandro Dalcin 
1791566a47fSLisandro Dalcin   PetscFunctionBegin;
1801566a47fSLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
1811566a47fSLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
1821566a47fSLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1831566a47fSLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
1841566a47fSLisandro Dalcin   PetscFunctionReturn(0);
1851566a47fSLisandro Dalcin }
1861566a47fSLisandro Dalcin 
187193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
188316643e7SJed Brown {
189316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1901566a47fSLisandro Dalcin   PetscInt       rejections = 0;
1914957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
1921566a47fSLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
193051f2191SLisandro Dalcin   PetscErrorCode ierr;
194316643e7SJed Brown 
195316643e7SJed Brown   PetscFunctionBegin;
1961566a47fSLisandro Dalcin   if (!ts->steprollback) {
1972ffb9264SLisandro Dalcin     if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
1983b1890cdSShri Abhyankar     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
1991566a47fSLisandro Dalcin   }
2001566a47fSLisandro Dalcin 
2011566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
2021566a47fSLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
203715f1b00SHong Zhang 
2041566a47fSLisandro Dalcin     PetscReal shift = 1/(th->Theta*ts->time_step);
2051566a47fSLisandro Dalcin     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step;
206316643e7SJed Brown 
207be5899b3SLisandro Dalcin     ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr);
208fecfb714SLisandro Dalcin     if (th->extrapolate && !ts->steprestart) {
209be5899b3SLisandro Dalcin       ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr);
210be5899b3SLisandro Dalcin     }
211eb284becSJed Brown     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
212eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
2131566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
2141566a47fSLisandro Dalcin       ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
2151566a47fSLisandro Dalcin       ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr);
2161566a47fSLisandro Dalcin     } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */
2171566a47fSLisandro Dalcin       ierr = VecZeroEntries(th->affine);CHKERRQ(ierr);
218eb284becSJed Brown     }
219be5899b3SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
220470880abSPatrick Sanan     ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr);
221be5899b3SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr);
2221566a47fSLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr);
223fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
224051f2191SLisandro Dalcin 
225051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2261566a47fSLisandro Dalcin     if (th->endpoint) {
2271566a47fSLisandro Dalcin       ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
2281566a47fSLisandro Dalcin     } else {
229cb3a5882SLisandro Dalcin       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr);
2301566a47fSLisandro Dalcin       ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
2311566a47fSLisandro Dalcin     }
2321566a47fSLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
2331566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2341566a47fSLisandro Dalcin     if (!accept) {
2351566a47fSLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
236be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
237be5899b3SLisandro Dalcin       goto reject_step;
238051f2191SLisandro Dalcin     }
2393b1890cdSShri Abhyankar 
240715f1b00SHong Zhang     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
241b1cb36f3SHong Zhang       th->ptime     = ts->ptime;
242b1cb36f3SHong Zhang       th->time_step = ts->time_step;
24317145e75SHong Zhang     }
2441566a47fSLisandro Dalcin 
2452b5a38e1SLisandro Dalcin     ts->ptime += ts->time_step;
246cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
247051f2191SLisandro Dalcin     break;
248051f2191SLisandro Dalcin 
249051f2191SLisandro Dalcin   reject_step:
250fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
2511566a47fSLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
252051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
2531566a47fSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
2543b1890cdSShri Abhyankar     }
2553b1890cdSShri Abhyankar   }
256316643e7SJed Brown   PetscFunctionReturn(0);
257316643e7SJed Brown }
258316643e7SJed Brown 
259*6980b705SHong Zhang /*
260*6980b705SHong Zhang   Use SNES to compute the Jacobian so that finite differencing could be used when TS Jacobian is not available.
261*6980b705SHong Zhang */
262*6980b705SHong Zhang static PetscErrorCode KSPTSFormOperator_Private(KSP ksp,Vec x,Mat J,Mat Jpre,TS ts)
263*6980b705SHong Zhang {
264*6980b705SHong Zhang   SNES           snes = ts->snes;
265*6980b705SHong Zhang   MatFDColoring  color;
266*6980b705SHong Zhang   PetscErrorCode ierr;
267*6980b705SHong Zhang 
268*6980b705SHong Zhang   PetscFunctionBegin;
269*6980b705SHong Zhang   /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
270*6980b705SHong Zhang   ierr = PetscObjectQuery((PetscObject)Jpre,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);
271*6980b705SHong Zhang   if (color) {
272*6980b705SHong Zhang     Vec f;
273*6980b705SHong Zhang     ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
274*6980b705SHong Zhang     ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
275*6980b705SHong Zhang   }
276*6980b705SHong Zhang   ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr);
277*6980b705SHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
278*6980b705SHong Zhang   PetscFunctionReturn(0);
279*6980b705SHong Zhang }
280*6980b705SHong Zhang 
281c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
282c9681e5dSHong Zhang {
283c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
284cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
285c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
286c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
287c9681e5dSHong Zhang   PetscInt       nadj;
2887207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
289c9681e5dSHong Zhang   KSP            ksp;
290c9681e5dSHong Zhang   PetscScalar    *xarr;
291077c4feaSHong Zhang   TSEquationType eqtype;
292077c4feaSHong Zhang   PetscBool      isexplicitode = PETSC_FALSE;
293c9681e5dSHong Zhang   PetscErrorCode ierr;
294c9681e5dSHong Zhang 
295c9681e5dSHong Zhang   PetscFunctionBegin;
296077c4feaSHong Zhang   ierr = TSGetEquationType(ts,&eqtype);CHKERRQ(ierr);
297077c4feaSHong Zhang   if (eqtype == TS_EQ_ODE_EXPLICIT) {
298077c4feaSHong Zhang     isexplicitode  = PETSC_TRUE;
299077c4feaSHong Zhang     VecsDeltaLam  = ts->vecs_sensi;
300077c4feaSHong Zhang     VecsDeltaLam2 = ts->vecs_sensi2;
301077c4feaSHong Zhang   }
302c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
303c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
304cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
3057207e0fdSHong Zhang   if (quadts) {
306cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
307cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
3087207e0fdSHong Zhang   }
309c9681e5dSHong Zhang 
310c9681e5dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
311c9681e5dSHong Zhang   th->stage_time = ts->ptime; /* time_step is negative*/
312c9681e5dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
313c9681e5dSHong Zhang   th->time_step  = -ts->time_step;
314c9681e5dSHong Zhang 
31533f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
3167207e0fdSHong Zhang   if (quadts) {
317cd4cee2dSHong Zhang     ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
3187207e0fdSHong Zhang   }
319cd4cee2dSHong Zhang 
320c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
321c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
322c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
323cd4cee2dSHong Zhang     if (quadJ) {
324cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
325cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
326cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
327cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
328cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
329c9681e5dSHong Zhang     }
330c9681e5dSHong Zhang   }
331c9681e5dSHong Zhang 
332c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
333*6980b705SHong Zhang   ierr = KSPTSFormOperator_Private(ksp,th->X,J,Jpre,ts);CHKERRQ(ierr);
334c9681e5dSHong Zhang 
335c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
336c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
337c9681e5dSHong Zhang     KSPConvergedReason kspreason;
338c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
339c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
340c9681e5dSHong Zhang     if (kspreason < 0) {
341c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
342c9681e5dSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
343c9681e5dSHong Zhang     }
344c9681e5dSHong Zhang   }
345c9681e5dSHong Zhang 
346c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
347c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
348c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
349c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
350c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
351b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
352c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
353b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
354c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
355c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
356*6980b705SHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],1./th->time_step);CHKERRQ(ierr);
35733f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
358c9681e5dSHong Zhang       if (ts->vecs_fup) {
35933f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
360c9681e5dSHong Zhang       }
361c9681e5dSHong Zhang     }
362c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
363c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
36405c9950eSHong Zhang       KSPConvergedReason kspreason;
365c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
36605c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
36705c9950eSHong Zhang       if (kspreason < 0) {
36805c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
36905c9950eSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
37005c9950eSHong Zhang       }
371c9681e5dSHong Zhang     }
372c9681e5dSHong Zhang   }
373c9681e5dSHong Zhang 
374c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
375077c4feaSHong Zhang   if (!isexplicitode) {
376*6980b705SHong Zhang     PetscReal shift = 0.0;
377cd4cee2dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */
378c9681e5dSHong Zhang     ierr  = MatScale(J,-1.);CHKERRQ(ierr);
379c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
38033f72c5dSHong Zhang       /* Add f_U \lambda_s to the original RHS */
381c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
382c9681e5dSHong Zhang       ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr);
383c9681e5dSHong Zhang       ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
384c9681e5dSHong Zhang       if (ts->vecs_sensi2) {
385c9681e5dSHong Zhang         ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
386c9681e5dSHong Zhang         ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr);
387c9681e5dSHong Zhang         ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
388c9681e5dSHong Zhang       }
389c9681e5dSHong Zhang     }
390077c4feaSHong Zhang   }
391c9681e5dSHong Zhang   if (ts->vecs_sensip) {
392*6980b705SHong Zhang     PetscReal shift = 1./th->time_step;;
393c9681e5dSHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
3947207e0fdSHong Zhang     if (quadts) {
395cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
3967207e0fdSHong Zhang     }
397c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
398c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
399b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
40033f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
401b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
402c9681e5dSHong Zhang     }
403cd4cee2dSHong Zhang 
404c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
405c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
406c9681e5dSHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
407cd4cee2dSHong Zhang       if (quadJp) {
408cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
409cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
410cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
411cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
412cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
41333f72c5dSHong Zhang       }
414c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
415c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
416c9681e5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
417c9681e5dSHong Zhang         if (ts->vecs_fpu) {
418c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
419c9681e5dSHong Zhang         }
420c9681e5dSHong Zhang         if (ts->vecs_fpp) {
421c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
422c9681e5dSHong Zhang         }
423c9681e5dSHong Zhang       }
424c9681e5dSHong Zhang     }
425c9681e5dSHong Zhang   }
426c9681e5dSHong Zhang 
427c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
428c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
429c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
430c9681e5dSHong Zhang   }
431c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
432c9681e5dSHong Zhang   PetscFunctionReturn(0);
433c9681e5dSHong Zhang }
434c9681e5dSHong Zhang 
43542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4362ca6e920SHong Zhang {
4372ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
438cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
439b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
440b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4412ca6e920SHong Zhang   PetscInt       nadj;
4427207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4432ca6e920SHong Zhang   KSP            ksp;
4442ca6e920SHong Zhang   PetscReal      shift;
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 */
488022c081aSHong Zhang   shift = 1./(th->Theta*th->time_step);
4893c54f92cSHong Zhang   if (th->endpoint) {
490cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
4913c54f92cSHong Zhang   } else {
492cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);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);
5211d14d03bSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],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 */
541022c081aSHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
542cd4cee2dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
54333f72c5dSHong Zhang     /* R_U at t_n */
5447207e0fdSHong Zhang     if (quadts) {
545cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
5467207e0fdSHong Zhang     }
547abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
548b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
549cd4cee2dSHong Zhang       if (quadJ) {
550cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
551cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
552cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr);
553cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
554cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
555b5b4017aSHong Zhang       }
5561d14d03bSHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
557b5b4017aSHong Zhang     }
5581d14d03bSHong Zhang 
5591d14d03bSHong Zhang     /* Second-order adjoint */
5601d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
561b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
562b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
563b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
564b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
565b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
566b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
56733f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
568b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
569b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
570b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
57133f72c5dSHong 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  */
572b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
57333f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
574b5b4017aSHong Zhang         if (ts->vecs_fup) {
57533f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
576b5b4017aSHong Zhang         }
5771d14d03bSHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr);
578b5b4017aSHong Zhang       }
579b5e0532dSHong Zhang     }
5803fd52205SHong Zhang 
5813fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
582b5b4017aSHong Zhang       /* U_{n+1} */
58333f72c5dSHong Zhang       shift = -1./(th->Theta*th->time_step);
58433f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5857207e0fdSHong Zhang       if (quadts) {
586cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
5877207e0fdSHong Zhang       }
588abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
589b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
59033f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5913fd52205SHong Zhang       }
59233f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
59333f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
59433f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
59533f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
596b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
597b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
59833f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
59933f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
60033f72c5dSHong Zhang 
601b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
602b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
603b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
60433f72c5dSHong 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)  */
605b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
60633f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
607b5b4017aSHong Zhang           if (ts->vecs_fpu) {
60833f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
609b5b4017aSHong Zhang           }
610b5b4017aSHong Zhang           if (ts->vecs_fpp) {
61133f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
612b5b4017aSHong Zhang           }
613b5b4017aSHong Zhang         }
614b5b4017aSHong Zhang       }
615b5b4017aSHong Zhang 
616b5b4017aSHong Zhang       /* U_s */
61733f72c5dSHong Zhang       shift = 1./((th->Theta-1.0)*th->time_step);
61833f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6197207e0fdSHong Zhang       if (quadts) {
620cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
6217207e0fdSHong Zhang       }
622abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
623b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
62433f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
62533f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
62633f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
62733f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
62833f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
629b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
630b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
63133f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
63233f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
633b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
634b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
635b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
63633f72c5dSHong 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) */
637b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
63833f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
639b5b4017aSHong Zhang             if (ts->vecs_fpu) {
64033f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
641b5e0532dSHong Zhang             }
642b5b4017aSHong Zhang             if (ts->vecs_fpp) {
64333f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
644b5b4017aSHong Zhang             }
64536eaed60SHong Zhang           }
64636eaed60SHong Zhang         }
6473fd52205SHong Zhang       }
648b5e0532dSHong Zhang     }
6493fd52205SHong Zhang   } else { /* one-stage case */
6503c54f92cSHong Zhang     shift = 0.0;
651cd4cee2dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
6527207e0fdSHong Zhang     if (quadts) {
653cd4cee2dSHong Zhang       ierr  = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
6547207e0fdSHong Zhang     }
655abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
656b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
657022c081aSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
658cd4cee2dSHong Zhang       if (quadJ) {
659cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
660cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
661cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr);
662cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
663cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
66436eaed60SHong Zhang       }
6652ca6e920SHong Zhang     }
6663fd52205SHong Zhang     if (ts->vecs_sensip) {
66733f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6687207e0fdSHong Zhang       if (quadts) {
669cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
6707207e0fdSHong Zhang       }
671abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
67233f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
67333f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
674cd4cee2dSHong Zhang         if (quadJp) {
675cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
676cd4cee2dSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
677cd4cee2dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
678cd4cee2dSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
679cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
68036eaed60SHong Zhang         }
68136eaed60SHong Zhang       }
6823fd52205SHong Zhang     }
6832ca6e920SHong Zhang   }
6842ca6e920SHong Zhang 
6852ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6862ca6e920SHong Zhang   PetscFunctionReturn(0);
6872ca6e920SHong Zhang }
6882ca6e920SHong Zhang 
689cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
690cd652676SJed Brown {
691cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
692be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
693cd652676SJed Brown   PetscErrorCode ierr;
694cd652676SJed Brown 
695cd652676SJed Brown   PetscFunctionBegin;
696a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
697be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
698be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
699cd652676SJed Brown   PetscFunctionReturn(0);
700cd652676SJed Brown }
701cd652676SJed Brown 
7021566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
7031566a47fSLisandro Dalcin {
7041566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7051566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
7061566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
7077453f775SEmil Constantinescu   PetscReal      wltea,wlter;
7081566a47fSLisandro Dalcin   PetscErrorCode ierr;
7091566a47fSLisandro Dalcin 
7101566a47fSLisandro Dalcin   PetscFunctionBegin;
7112ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
7121566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
713fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
7141566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
715fecfb714SLisandro Dalcin   {
716be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
717be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
7181566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
7191566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
7201566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
7211566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
7221566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
7237453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
7241566a47fSLisandro Dalcin   }
7251566a47fSLisandro Dalcin   if (order) *order = 2;
7261566a47fSLisandro Dalcin   PetscFunctionReturn(0);
7271566a47fSLisandro Dalcin }
7281566a47fSLisandro Dalcin 
7291566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7301566a47fSLisandro Dalcin {
7311566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7327207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7331566a47fSLisandro Dalcin   PetscErrorCode ierr;
7341566a47fSLisandro Dalcin 
7351566a47fSLisandro Dalcin   PetscFunctionBegin;
7361566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
737cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
738cd4cee2dSHong Zhang     ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr);
7391566a47fSLisandro Dalcin   }
740715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
74113b1b0ffSHong Zhang   if (ts->mat_sensip) {
74213b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7436f92202bSHong Zhang   }
7447207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7457207e0fdSHong Zhang     ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7466f92202bSHong Zhang   }
747715f1b00SHong Zhang   PetscFunctionReturn(0);
748715f1b00SHong Zhang }
749715f1b00SHong Zhang 
750715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
751715f1b00SHong Zhang {
752715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
753cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
75413b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
755b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7567207e0fdSHong Zhang   PetscInt       ntlm;
757715f1b00SHong Zhang   KSP            ksp;
7587207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
759715f1b00SHong Zhang   PetscReal      shift;
76013b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
761715f1b00SHong Zhang   PetscErrorCode ierr;
762715f1b00SHong Zhang 
763715f1b00SHong Zhang   PetscFunctionBegin;
76413b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
76513b1b0ffSHong Zhang 
7667207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7677207e0fdSHong Zhang     ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7686f92202bSHong Zhang   }
769715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
770cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
7717207e0fdSHong Zhang   if (quadts) {
772cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
773cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
7747207e0fdSHong Zhang   }
775715f1b00SHong Zhang 
776715f1b00SHong Zhang   /* Build RHS */
777715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
778715f1b00SHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
779cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
78013b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
78113b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
782715f1b00SHong Zhang 
783022c081aSHong Zhang     /* Add the f_p forcing terms */
7840e11c21fSHong Zhang     if (ts->Jacp) {
78533f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
78633f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
78733f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
78833f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7890e11c21fSHong Zhang     }
790715f1b00SHong Zhang   } else { /* 1-stage method */
791715f1b00SHong Zhang     shift = 0.0;
792cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
79313b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
79413b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
79513b1b0ffSHong Zhang 
796022c081aSHong Zhang     /* Add the f_p forcing terms */
7970e11c21fSHong Zhang     if (ts->Jacp) {
79833f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
79933f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
800715f1b00SHong Zhang     }
8010e11c21fSHong Zhang   }
802715f1b00SHong Zhang 
803715f1b00SHong Zhang   /* Build LHS */
804715f1b00SHong Zhang   shift = 1/(th->Theta*th->time_step);
805715f1b00SHong Zhang   if (th->endpoint) {
806cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
807715f1b00SHong Zhang   } else {
808cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
809715f1b00SHong Zhang   }
810cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
811715f1b00SHong Zhang 
812715f1b00SHong Zhang   /*
813715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
814c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
815715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
816715f1b00SHong Zhang   */
817715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
8187207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
819cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
820cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
8217207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8227207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8237207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
824715f1b00SHong Zhang     }
825715f1b00SHong Zhang   }
826715f1b00SHong Zhang 
827715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
828022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
829b5b4017aSHong Zhang     KSPConvergedReason kspreason;
83013b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
831b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
832715f1b00SHong Zhang     if (th->endpoint) {
83313b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
834b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
835b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
836b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
83713b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
838715f1b00SHong Zhang     } else {
839b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
840715f1b00SHong Zhang     }
841b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
842b5b4017aSHong Zhang     if (kspreason < 0) {
843b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
844b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
845b5b4017aSHong Zhang     }
846b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
84713b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
848715f1b00SHong Zhang   }
849715f1b00SHong Zhang 
850b5b4017aSHong Zhang 
85113b1b0ffSHong Zhang   /*
85213b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
853c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
85413b1b0ffSHong Zhang   */
8557207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
85613b1b0ffSHong Zhang     if (!th->endpoint) {
8577207e0fdSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */
858cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
859cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
8607207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8617207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8627207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
86313b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
86413b1b0ffSHong Zhang     } else {
865cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
866cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
8677207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8687207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8697207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
87013b1b0ffSHong Zhang     }
87113b1b0ffSHong Zhang   } else {
87213b1b0ffSHong Zhang     if (!th->endpoint) {
87313b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
874715f1b00SHong Zhang     }
875715f1b00SHong Zhang   }
8761566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8771566a47fSLisandro Dalcin }
8781566a47fSLisandro Dalcin 
8796affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
8806affb6f8SHong Zhang {
8816affb6f8SHong Zhang   TS_Theta *th = (TS_Theta*)ts->data;
8826affb6f8SHong Zhang 
8836affb6f8SHong Zhang   PetscFunctionBegin;
8846affb6f8SHong Zhang   if (ns) *ns = 1;
8856affb6f8SHong Zhang   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
8866affb6f8SHong Zhang   PetscFunctionReturn(0);
8876affb6f8SHong Zhang }
8886affb6f8SHong Zhang 
889316643e7SJed Brown /*------------------------------------------------------------*/
890277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
891316643e7SJed Brown {
892316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
893316643e7SJed Brown   PetscErrorCode ierr;
894316643e7SJed Brown 
895316643e7SJed Brown   PetscFunctionBegin;
8966bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
8976bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
8983b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
899eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
9001566a47fSLisandro Dalcin 
9011566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
9021566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
9031566a47fSLisandro Dalcin 
904b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
905277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
906277b19d0SLisandro Dalcin }
907277b19d0SLisandro Dalcin 
908ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
909ece46509SHong Zhang {
910ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
911ece46509SHong Zhang   PetscErrorCode ierr;
912ece46509SHong Zhang 
913ece46509SHong Zhang   PetscFunctionBegin;
914ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
915ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
916ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
917ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
918ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
919ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
920ece46509SHong Zhang   PetscFunctionReturn(0);
921ece46509SHong Zhang }
922ece46509SHong Zhang 
923277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
924277b19d0SLisandro Dalcin {
925277b19d0SLisandro Dalcin   PetscErrorCode ierr;
926277b19d0SLisandro Dalcin 
927277b19d0SLisandro Dalcin   PetscFunctionBegin;
928277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
929b3a6b972SJed Brown   if (ts->dm) {
930b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
931b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
932b3a6b972SJed Brown   }
933277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
934bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
935bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
936bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
937bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
938316643e7SJed Brown   PetscFunctionReturn(0);
939316643e7SJed Brown }
940316643e7SJed Brown 
941316643e7SJed Brown /*
942316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9432b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
944316643e7SJed Brown */
9450f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
946316643e7SJed Brown {
947316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
948316643e7SJed Brown   PetscErrorCode ierr;
9497445fe48SJed Brown   Vec            X0,Xdot;
9507445fe48SJed Brown   DM             dm,dmsave;
9511566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
952316643e7SJed Brown 
953316643e7SJed Brown   PetscFunctionBegin;
9547445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9555a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9567445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
957b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
9587445fe48SJed Brown 
9597445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9607445fe48SJed Brown   dmsave = ts->dm;
9617445fe48SJed Brown   ts->dm = dm;
9627445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
9637445fe48SJed Brown   ts->dm = dmsave;
9640d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
965316643e7SJed Brown   PetscFunctionReturn(0);
966316643e7SJed Brown }
967316643e7SJed Brown 
968d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
969316643e7SJed Brown {
970316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
971316643e7SJed Brown   PetscErrorCode ierr;
9727445fe48SJed Brown   Vec            Xdot;
9737445fe48SJed Brown   DM             dm,dmsave;
9741566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
975316643e7SJed Brown 
976316643e7SJed Brown   PetscFunctionBegin;
9777445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
978be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9790298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
9807445fe48SJed Brown 
9817445fe48SJed Brown   dmsave = ts->dm;
9827445fe48SJed Brown   ts->dm = dm;
983d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
9847445fe48SJed Brown   ts->dm = dmsave;
9850298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
986316643e7SJed Brown   PetscFunctionReturn(0);
987316643e7SJed Brown }
988316643e7SJed Brown 
989715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
990715f1b00SHong Zhang {
991715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9927207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
993715f1b00SHong Zhang   PetscErrorCode ierr;
994715f1b00SHong Zhang 
995715f1b00SHong Zhang   PetscFunctionBegin;
996022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
99713b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
99813b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
9997207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10007207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10017207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
1002715f1b00SHong Zhang   }
10036f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
100413b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
100513b1b0ffSHong Zhang 
1006b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
1007715f1b00SHong Zhang   PetscFunctionReturn(0);
1008715f1b00SHong Zhang }
1009715f1b00SHong Zhang 
10107adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
10117adebddeSHong Zhang {
10127adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
10137207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
10147adebddeSHong Zhang   PetscErrorCode ierr;
10157adebddeSHong Zhang 
10167adebddeSHong Zhang   PetscFunctionBegin;
10177207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
10187207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
10197207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
10207adebddeSHong Zhang   }
10217adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
10227adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
10237adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
10247adebddeSHong Zhang   PetscFunctionReturn(0);
10257adebddeSHong Zhang }
10267adebddeSHong Zhang 
1027316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
1028316643e7SJed Brown {
1029316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1030cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10312ffb9264SLisandro Dalcin   PetscBool      match;
1032316643e7SJed Brown   PetscErrorCode ierr;
1033316643e7SJed Brown 
1034316643e7SJed Brown   PetscFunctionBegin;
1035cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1036cd4cee2dSHong Zhang     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1037b1cb36f3SHong Zhang   }
103839d13666SHong Zhang   if (!th->X) {
1039316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
104039d13666SHong Zhang   }
104139d13666SHong Zhang   if (!th->Xdot) {
1042316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
104339d13666SHong Zhang   }
104439d13666SHong Zhang   if (!th->X0) {
10453b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
104639d13666SHong Zhang   }
10471566a47fSLisandro Dalcin   if (th->endpoint) {
10481566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
10497445fe48SJed Brown   }
10503b1890cdSShri Abhyankar 
10511566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
10521566a47fSLisandro Dalcin 
10531566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
10541566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10551566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10561566a47fSLisandro Dalcin 
10571566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
10581566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
10592ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
10602ffb9264SLisandro Dalcin   if (!match) {
10611566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
10621566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
10633b1890cdSShri Abhyankar   }
10641566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1065b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1066b4dd3bf9SBarry Smith }
10670c7376e5SHong Zhang 
1068b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1069b4dd3bf9SBarry Smith 
1070b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1071b4dd3bf9SBarry Smith {
1072b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1073b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1074b4dd3bf9SBarry Smith 
1075b4dd3bf9SBarry Smith   PetscFunctionBegin;
1076b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1077b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
10782ca6e920SHong Zhang   if (ts->vecs_sensip) {
1079b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
10802ca6e920SHong Zhang   }
1081b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1082b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1083b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
108467633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
108567633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
108667633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1087b5b4017aSHong Zhang   }
1088c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1089c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
109067633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
109167633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
109267633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1093b5b4017aSHong Zhang   }
1094316643e7SJed Brown   PetscFunctionReturn(0);
1095316643e7SJed Brown }
1096316643e7SJed Brown 
10974416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1098316643e7SJed Brown {
1099316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1100316643e7SJed Brown   PetscErrorCode ierr;
1101316643e7SJed Brown 
1102316643e7SJed Brown   PetscFunctionBegin;
1103e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1104316643e7SJed Brown   {
11050298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
11060298fd71SBarry 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);
110703842d09SLisandro 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);
1108316643e7SJed Brown   }
1109316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1110316643e7SJed Brown   PetscFunctionReturn(0);
1111316643e7SJed Brown }
1112316643e7SJed Brown 
1113316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1114316643e7SJed Brown {
1115316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1116ace3abfcSBarry Smith   PetscBool      iascii;
1117316643e7SJed Brown   PetscErrorCode ierr;
1118316643e7SJed Brown 
1119316643e7SJed Brown   PetscFunctionBegin;
1120251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1121316643e7SJed Brown   if (iascii) {
11227c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1123ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1124316643e7SJed Brown   }
1125316643e7SJed Brown   PetscFunctionReturn(0);
1126316643e7SJed Brown }
1127316643e7SJed Brown 
1128560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
11290de4c49aSJed Brown {
11300de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11310de4c49aSJed Brown 
11320de4c49aSJed Brown   PetscFunctionBegin;
11330de4c49aSJed Brown   *theta = th->Theta;
11340de4c49aSJed Brown   PetscFunctionReturn(0);
11350de4c49aSJed Brown }
11360de4c49aSJed Brown 
1137560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11380de4c49aSJed Brown {
11390de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11400de4c49aSJed Brown 
11410de4c49aSJed Brown   PetscFunctionBegin;
11427c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11430de4c49aSJed Brown   th->Theta = theta;
11441566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11450de4c49aSJed Brown   PetscFunctionReturn(0);
11460de4c49aSJed Brown }
1147eb284becSJed Brown 
1148560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
114926f2ff8fSLisandro Dalcin {
115026f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
115126f2ff8fSLisandro Dalcin 
115226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
115326f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
115426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
115526f2ff8fSLisandro Dalcin }
115626f2ff8fSLisandro Dalcin 
1157560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1158eb284becSJed Brown {
1159eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1160eb284becSJed Brown 
1161eb284becSJed Brown   PetscFunctionBegin;
1162eb284becSJed Brown   th->endpoint = flg;
1163eb284becSJed Brown   PetscFunctionReturn(0);
1164eb284becSJed Brown }
11650de4c49aSJed Brown 
1166f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1167f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1168f9c1d6abSBarry Smith {
1169f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1170f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11713fd8ae06SJed Brown   const PetscReal one = 1.0;
1172f9c1d6abSBarry Smith 
1173f9c1d6abSBarry Smith   PetscFunctionBegin;
11743fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1175f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1176f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1177f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1178f9c1d6abSBarry Smith }
1179f9c1d6abSBarry Smith #endif
1180f9c1d6abSBarry Smith 
118142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
118242682096SHong Zhang {
118342682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
118442682096SHong Zhang 
118542682096SHong Zhang   PetscFunctionBegin;
11861566a47fSLisandro Dalcin   if (ns) *ns = 1;
11871566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
118842682096SHong Zhang   PetscFunctionReturn(0);
118942682096SHong Zhang }
1190f9c1d6abSBarry Smith 
1191316643e7SJed Brown /* ------------------------------------------------------------ */
1192316643e7SJed Brown /*MC
119396f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1194316643e7SJed Brown 
1195316643e7SJed Brown    Level: beginner
1196316643e7SJed Brown 
11974eb428fdSBarry Smith    Options Database:
1198c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1199c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
120003842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
12014eb428fdSBarry Smith 
1202eb284becSJed Brown    Notes:
12030c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
12040c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
12054eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
12064eb428fdSBarry Smith 
1207eb284becSJed Brown    This method can be applied to DAE.
1208eb284becSJed Brown 
1209eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
1210eb284becSJed Brown 
1211eb284becSJed Brown .vb
1212eb284becSJed Brown   Theta | Theta
1213eb284becSJed Brown   -------------
1214eb284becSJed Brown         |  1
1215eb284becSJed Brown .ve
1216eb284becSJed Brown 
1217eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1218eb284becSJed Brown 
1219eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1220eb284becSJed Brown 
1221eb284becSJed Brown .vb
1222eb284becSJed Brown   0 | 0         0
1223eb284becSJed Brown   1 | 1-Theta   Theta
1224eb284becSJed Brown   -------------------
1225eb284becSJed Brown     | 1-Theta   Theta
1226eb284becSJed Brown .ve
1227eb284becSJed Brown 
1228eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1229eb284becSJed Brown 
1230eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1231eb284becSJed Brown 
1232eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1233eb284becSJed Brown 
12344eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1235eb284becSJed Brown 
1236eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1237316643e7SJed Brown 
1238316643e7SJed Brown M*/
12398cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1240316643e7SJed Brown {
1241316643e7SJed Brown   TS_Theta       *th;
1242316643e7SJed Brown   PetscErrorCode ierr;
1243316643e7SJed Brown 
1244316643e7SJed Brown   PetscFunctionBegin;
1245277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1246ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1247316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1248316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1249316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
125042f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1251ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1252316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1253cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12541566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
125524655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1256316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12570f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12580f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1259f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1260f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1261f9c1d6abSBarry Smith #endif
126242682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
126342f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1264b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1265b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12662ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12676affb6f8SHong Zhang 
1268715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12697adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1270715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12716affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1272316643e7SJed Brown 
1273efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1274efd4aadfSBarry Smith 
1275b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1276316643e7SJed Brown   ts->data = (void*)th;
1277316643e7SJed Brown 
1278715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1279715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1280715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1281b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1282715f1b00SHong Zhang 
12836f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1284316643e7SJed Brown   th->Theta       = 0.5;
12851566a47fSLisandro Dalcin   th->order       = 2;
1286bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1287bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1288bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1289bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1290316643e7SJed Brown   PetscFunctionReturn(0);
1291316643e7SJed Brown }
12920de4c49aSJed Brown 
12930de4c49aSJed Brown /*@
12940de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12950de4c49aSJed Brown 
12960de4c49aSJed Brown   Not Collective
12970de4c49aSJed Brown 
12980de4c49aSJed Brown   Input Parameter:
12990de4c49aSJed Brown .  ts - timestepping context
13000de4c49aSJed Brown 
13010de4c49aSJed Brown   Output Parameter:
13020de4c49aSJed Brown .  theta - stage abscissa
13030de4c49aSJed Brown 
13040de4c49aSJed Brown   Note:
13050de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
13060de4c49aSJed Brown 
13070de4c49aSJed Brown   Level: Advanced
13080de4c49aSJed Brown 
13090de4c49aSJed Brown .seealso: TSThetaSetTheta()
13100de4c49aSJed Brown @*/
13117087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
13120de4c49aSJed Brown {
13134ac538c5SBarry Smith   PetscErrorCode ierr;
13140de4c49aSJed Brown 
13150de4c49aSJed Brown   PetscFunctionBegin;
1316afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13170de4c49aSJed Brown   PetscValidPointer(theta,2);
13184ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
13190de4c49aSJed Brown   PetscFunctionReturn(0);
13200de4c49aSJed Brown }
13210de4c49aSJed Brown 
13220de4c49aSJed Brown /*@
13230de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
13240de4c49aSJed Brown 
13250de4c49aSJed Brown   Not Collective
13260de4c49aSJed Brown 
13270de4c49aSJed Brown   Input Parameter:
13280de4c49aSJed Brown +  ts - timestepping context
13290de4c49aSJed Brown -  theta - stage abscissa
13300de4c49aSJed Brown 
13310de4c49aSJed Brown   Options Database:
13320de4c49aSJed Brown .  -ts_theta_theta <theta>
13330de4c49aSJed Brown 
13340de4c49aSJed Brown   Level: Intermediate
13350de4c49aSJed Brown 
13360de4c49aSJed Brown .seealso: TSThetaGetTheta()
13370de4c49aSJed Brown @*/
13387087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13390de4c49aSJed Brown {
13404ac538c5SBarry Smith   PetscErrorCode ierr;
13410de4c49aSJed Brown 
13420de4c49aSJed Brown   PetscFunctionBegin;
1343afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13444ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
13450de4c49aSJed Brown   PetscFunctionReturn(0);
13460de4c49aSJed Brown }
1347f33bbcb6SJed Brown 
134826f2ff8fSLisandro Dalcin /*@
134926f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
135026f2ff8fSLisandro Dalcin 
135126f2ff8fSLisandro Dalcin   Not Collective
135226f2ff8fSLisandro Dalcin 
135326f2ff8fSLisandro Dalcin   Input Parameter:
135426f2ff8fSLisandro Dalcin .  ts - timestepping context
135526f2ff8fSLisandro Dalcin 
135626f2ff8fSLisandro Dalcin   Output Parameter:
135726f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
135826f2ff8fSLisandro Dalcin 
135926f2ff8fSLisandro Dalcin   Level: Advanced
136026f2ff8fSLisandro Dalcin 
136126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
136226f2ff8fSLisandro Dalcin @*/
136326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
136426f2ff8fSLisandro Dalcin {
136526f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
136626f2ff8fSLisandro Dalcin 
136726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
136826f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
136926f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1370163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
137126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
137226f2ff8fSLisandro Dalcin }
137326f2ff8fSLisandro Dalcin 
1374eb284becSJed Brown /*@
1375eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1376eb284becSJed Brown 
1377eb284becSJed Brown   Not Collective
1378eb284becSJed Brown 
1379eb284becSJed Brown   Input Parameter:
1380eb284becSJed Brown +  ts - timestepping context
1381eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1382eb284becSJed Brown 
1383eb284becSJed Brown   Options Database:
1384eb284becSJed Brown .  -ts_theta_endpoint <flg>
1385eb284becSJed Brown 
1386eb284becSJed Brown   Level: Intermediate
1387eb284becSJed Brown 
1388eb284becSJed Brown .seealso: TSTHETA, TSCN
1389eb284becSJed Brown @*/
1390eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1391eb284becSJed Brown {
1392eb284becSJed Brown   PetscErrorCode ierr;
1393eb284becSJed Brown 
1394eb284becSJed Brown   PetscFunctionBegin;
1395eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1396eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1397eb284becSJed Brown   PetscFunctionReturn(0);
1398eb284becSJed Brown }
1399eb284becSJed Brown 
1400f33bbcb6SJed Brown /*
1401f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1402f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1403f33bbcb6SJed Brown  */
1404f33bbcb6SJed Brown 
14051566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
14061566a47fSLisandro Dalcin {
14071566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14081566a47fSLisandro Dalcin   PetscErrorCode ierr;
14091566a47fSLisandro Dalcin 
14101566a47fSLisandro Dalcin   PetscFunctionBegin;
14111566a47fSLisandro 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");
14121566a47fSLisandro 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");
14131566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14141566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14151566a47fSLisandro Dalcin }
14161566a47fSLisandro Dalcin 
1417f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1418f33bbcb6SJed Brown {
1419f33bbcb6SJed Brown   PetscFunctionBegin;
1420f33bbcb6SJed Brown   PetscFunctionReturn(0);
1421f33bbcb6SJed Brown }
1422f33bbcb6SJed Brown 
1423f33bbcb6SJed Brown /*MC
1424f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1425f33bbcb6SJed Brown 
1426f33bbcb6SJed Brown   Level: beginner
1427f33bbcb6SJed Brown 
14284eb428fdSBarry Smith   Notes:
1429c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14304eb428fdSBarry Smith 
14311566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14324eb428fdSBarry Smith 
1433f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1434f33bbcb6SJed Brown 
1435f33bbcb6SJed Brown M*/
14368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1437f33bbcb6SJed Brown {
1438f33bbcb6SJed Brown   PetscErrorCode ierr;
1439f33bbcb6SJed Brown 
1440f33bbcb6SJed Brown   PetscFunctionBegin;
1441f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1442f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
14431566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
14440c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1445f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1446f33bbcb6SJed Brown   PetscFunctionReturn(0);
1447f33bbcb6SJed Brown }
1448f33bbcb6SJed Brown 
14491566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14501566a47fSLisandro Dalcin {
14511566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14521566a47fSLisandro Dalcin   PetscErrorCode ierr;
14531566a47fSLisandro Dalcin 
14541566a47fSLisandro Dalcin   PetscFunctionBegin;
14551566a47fSLisandro 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");
14561566a47fSLisandro 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");
14571566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14581566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14591566a47fSLisandro Dalcin }
14601566a47fSLisandro Dalcin 
1461f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1462f33bbcb6SJed Brown {
1463f33bbcb6SJed Brown   PetscFunctionBegin;
1464f33bbcb6SJed Brown   PetscFunctionReturn(0);
1465f33bbcb6SJed Brown }
1466f33bbcb6SJed Brown 
1467f33bbcb6SJed Brown /*MC
1468f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1469f33bbcb6SJed Brown 
1470f33bbcb6SJed Brown   Level: beginner
1471f33bbcb6SJed Brown 
1472f33bbcb6SJed Brown   Notes:
14737cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14747cf5af47SJed Brown 
14757cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1476f33bbcb6SJed Brown 
1477f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1478f33bbcb6SJed Brown 
1479f33bbcb6SJed Brown M*/
14808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1481f33bbcb6SJed Brown {
1482f33bbcb6SJed Brown   PetscErrorCode ierr;
1483f33bbcb6SJed Brown 
1484f33bbcb6SJed Brown   PetscFunctionBegin;
1485f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1486f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1487eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
14880c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1489f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1490f33bbcb6SJed Brown   PetscFunctionReturn(0);
1491f33bbcb6SJed Brown }
1492