xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision b1e111eb0d8a8c8f79792012adb4e9519947d865)
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 
259c9681e5dSHong Zhang static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
260c9681e5dSHong Zhang {
261c9681e5dSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
262cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
263c9681e5dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
264c9681e5dSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
265c9681e5dSHong Zhang   PetscInt       nadj;
2667207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
267c9681e5dSHong Zhang   KSP            ksp;
268c9681e5dSHong Zhang   PetscReal      shift;
269c9681e5dSHong Zhang   PetscScalar    *xarr;
270c9681e5dSHong Zhang   PetscErrorCode ierr;
271c9681e5dSHong Zhang 
272c9681e5dSHong Zhang   PetscFunctionBegin;
273c9681e5dSHong Zhang   th->status = TS_STEP_INCOMPLETE;
274c9681e5dSHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
275cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
2767207e0fdSHong Zhang   if (quadts) {
277cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
278cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
2797207e0fdSHong Zhang   }
280c9681e5dSHong Zhang 
281c9681e5dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
282c9681e5dSHong Zhang   th->stage_time = ts->ptime; /* time_step is negative*/
283c9681e5dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
284c9681e5dSHong Zhang   th->time_step  = -ts->time_step;
285c9681e5dSHong Zhang 
28633f72c5dSHong Zhang   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
2877207e0fdSHong Zhang   if (quadts) {
288cd4cee2dSHong Zhang     ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
2897207e0fdSHong Zhang   }
290cd4cee2dSHong Zhang 
291c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
292c9681e5dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
293c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./th->time_step);CHKERRQ(ierr); /* lambda_{n+1}/h */
294cd4cee2dSHong Zhang     if (quadJ) {
295cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
296cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
297cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
298cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
299cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
300c9681e5dSHong Zhang     }
301c9681e5dSHong Zhang   }
302c9681e5dSHong Zhang 
303c9681e5dSHong Zhang   /* Build LHS for first-order adjoint */
304c9681e5dSHong Zhang   shift = 1./th->time_step;
305cd4cee2dSHong Zhang   ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
306cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
307c9681e5dSHong Zhang 
308c9681e5dSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
309c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
310c9681e5dSHong Zhang     KSPConvergedReason kspreason;
311c9681e5dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
312c9681e5dSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
313c9681e5dSHong Zhang     if (kspreason < 0) {
314c9681e5dSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
315c9681e5dSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
316c9681e5dSHong Zhang     }
317c9681e5dSHong Zhang   }
318c9681e5dSHong Zhang 
319c9681e5dSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
320c9681e5dSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
321c9681e5dSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
322c9681e5dSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
323c9681e5dSHong Zhang     /* lambda_s^T F_UU w_1 */
324*b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
325c9681e5dSHong Zhang     /* lambda_s^T F_UP w_2 */
326*b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
327c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
328c9681e5dSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
329c9681e5dSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
33033f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
331c9681e5dSHong Zhang       if (ts->vecs_fup) {
33233f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
333c9681e5dSHong Zhang       }
334c9681e5dSHong Zhang     }
335c9681e5dSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
336c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
33705c9950eSHong Zhang       KSPConvergedReason kspreason;
338c9681e5dSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
33905c9950eSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
34005c9950eSHong Zhang       if (kspreason < 0) {
34105c9950eSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
34205c9950eSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
34305c9950eSHong Zhang       }
344c9681e5dSHong Zhang     }
345c9681e5dSHong Zhang   }
346c9681e5dSHong Zhang 
347c9681e5dSHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
348c9681e5dSHong Zhang   shift = 0.0;
349cd4cee2dSHong Zhang   ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_U */
350c9681e5dSHong Zhang   ierr  = MatScale(J,-1.);CHKERRQ(ierr);
351c9681e5dSHong Zhang   for (nadj=0; nadj<ts->numcost; nadj++) {
35233f72c5dSHong Zhang     /* Add f_U \lambda_s to the original RHS */
353c9681e5dSHong Zhang     ierr = MatMultTransposeAdd(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
354c9681e5dSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],th->time_step);CHKERRQ(ierr);
355c9681e5dSHong Zhang     ierr = VecCopy(VecsSensiTemp[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
356c9681e5dSHong Zhang     if (ts->vecs_sensi2) {
357c9681e5dSHong Zhang       ierr = MatMultTransposeAdd(J,VecsDeltaLam2[nadj],VecsSensi2Temp[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
358c9681e5dSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],th->time_step);CHKERRQ(ierr);
359c9681e5dSHong Zhang       ierr = VecCopy(VecsSensi2Temp[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
360c9681e5dSHong Zhang     }
361c9681e5dSHong Zhang   }
362c9681e5dSHong Zhang   if (ts->vecs_sensip) {
363c9681e5dSHong Zhang     ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_p */
3647207e0fdSHong Zhang     if (quadts) {
365cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
3667207e0fdSHong Zhang     }
367c9681e5dSHong Zhang     if (ts->vecs_sensi2p) {
368c9681e5dSHong Zhang       /* lambda_s^T F_PU w_1 */
369*b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
37033f72c5dSHong Zhang       /* lambda_s^T F_PP w_2 */
371*b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,th->X,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
372c9681e5dSHong Zhang     }
373cd4cee2dSHong Zhang 
374c9681e5dSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
375c9681e5dSHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
376c9681e5dSHong Zhang       ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
377cd4cee2dSHong Zhang       if (quadJp) {
378cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
379cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
380cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
381cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
382cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
38333f72c5dSHong Zhang       }
384c9681e5dSHong Zhang       if (ts->vecs_sensi2p) {
385c9681e5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
386c9681e5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
387c9681e5dSHong Zhang         if (ts->vecs_fpu) {
388c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
389c9681e5dSHong Zhang         }
390c9681e5dSHong Zhang         if (ts->vecs_fpp) {
391c9681e5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
392c9681e5dSHong Zhang         }
393c9681e5dSHong Zhang       }
394c9681e5dSHong Zhang     }
395c9681e5dSHong Zhang   }
396c9681e5dSHong Zhang 
397c9681e5dSHong Zhang   if (ts->vecs_sensi2) {
398c9681e5dSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
399c9681e5dSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
400c9681e5dSHong Zhang   }
401c9681e5dSHong Zhang   th->status = TS_STEP_COMPLETE;
402c9681e5dSHong Zhang   PetscFunctionReturn(0);
403c9681e5dSHong Zhang }
404c9681e5dSHong Zhang 
40542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_Theta(TS ts)
4062ca6e920SHong Zhang {
4072ca6e920SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
408cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
409b5e0532dSHong Zhang   Vec            *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp;
410b5b4017aSHong Zhang   Vec            *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp;
4112ca6e920SHong Zhang   PetscInt       nadj;
4127207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
4132ca6e920SHong Zhang   KSP            ksp;
4142ca6e920SHong Zhang   PetscReal      shift;
415b5b4017aSHong Zhang   PetscScalar    *xarr;
416b5b4017aSHong Zhang   PetscErrorCode ierr;
4172ca6e920SHong Zhang 
4182ca6e920SHong Zhang   PetscFunctionBegin;
419c9681e5dSHong Zhang   if (th->Theta == 1.) {
420c9681e5dSHong Zhang     ierr = TSAdjointStepBEuler_Private(ts);CHKERRQ(ierr);
421c9681e5dSHong Zhang     PetscFunctionReturn(0);
422c9681e5dSHong Zhang   }
4232ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
424302440fdSBarry Smith   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
425cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
4267207e0fdSHong Zhang   if (quadts) {
427cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
428cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
4297207e0fdSHong Zhang   }
430b5e0532dSHong Zhang   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
431022c081aSHong Zhang   th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/
432b5e0532dSHong Zhang   th->ptime      = ts->ptime + ts->time_step;
433022c081aSHong Zhang   th->time_step  = -ts->time_step;
4342ca6e920SHong Zhang 
435b5b4017aSHong Zhang   /* Build RHS for first-order adjoint */
43633f72c5dSHong Zhang   /* Cost function has an integral term */
4377207e0fdSHong Zhang   if (quadts) {
43805755b9cSHong Zhang     if (th->endpoint) {
439cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
44005755b9cSHong Zhang     } else {
441cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
44205755b9cSHong Zhang     }
4437207e0fdSHong Zhang   }
44433f72c5dSHong Zhang 
445abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
446b5e0532dSHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
447022c081aSHong Zhang     ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr);
448cd4cee2dSHong Zhang     if (quadJ) {
449cd4cee2dSHong Zhang       ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
450cd4cee2dSHong Zhang       ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
451cd4cee2dSHong Zhang       ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vec_drdu_col);CHKERRQ(ierr);
452cd4cee2dSHong Zhang       ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
453cd4cee2dSHong Zhang       ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
45436eaed60SHong Zhang     }
4552ca6e920SHong Zhang   }
4563c54f92cSHong Zhang 
457b5b4017aSHong Zhang   /* Build LHS for first-order adjoint */
458022c081aSHong Zhang   shift = 1./(th->Theta*th->time_step);
4593c54f92cSHong Zhang   if (th->endpoint) {
460cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
4613c54f92cSHong Zhang   } else {
462cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
4633c54f92cSHong Zhang   }
464cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
4652ca6e920SHong Zhang 
466b5b4017aSHong Zhang   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
467abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
4681d14d03bSHong Zhang     KSPConvergedReason kspreason;
469b5e0532dSHong Zhang     ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr);
4701d14d03bSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
4711d14d03bSHong Zhang     if (kspreason < 0) {
4721d14d03bSHong Zhang       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
4731d14d03bSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
4741d14d03bSHong Zhang     }
4752ca6e920SHong Zhang   }
4763c54f92cSHong Zhang 
4771d14d03bSHong Zhang   /* Second-order adjoint */
478b5b4017aSHong Zhang   if (ts->vecs_sensi2) { /* U_{n+1} */
479b5b4017aSHong Zhang     if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta");
480b5b4017aSHong Zhang     /* Get w1 at t_{n+1} from TLM matrix */
481b5b4017aSHong Zhang     ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
482b5b4017aSHong Zhang     ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
483b5b4017aSHong Zhang     /* lambda_s^T F_UU w_1 */
484*b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
485b5b4017aSHong Zhang     ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
486b5b4017aSHong Zhang     ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
487b5b4017aSHong Zhang     /* lambda_s^T F_UP w_2 */
488*b1e111ebSHong Zhang     ierr = TSComputeIHessianProductFunctionUP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
489b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */
490b5b4017aSHong Zhang       ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr);
4911d14d03bSHong Zhang       ierr = VecScale(VecsSensi2Temp[nadj],shift);CHKERRQ(ierr);
49233f72c5dSHong Zhang       ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
493b5b4017aSHong Zhang       if (ts->vecs_fup) {
49433f72c5dSHong Zhang         ierr = VecAXPY(VecsSensi2Temp[nadj],-1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
495b5b4017aSHong Zhang       }
496b5b4017aSHong Zhang     }
497b5b4017aSHong Zhang     /* Solve stage equation LHS X = RHS for second-order adjoint */
498b5b4017aSHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
4991d14d03bSHong Zhang       KSPConvergedReason kspreason;
5001d14d03bSHong Zhang       ierr = KSPSolveTranspose(ksp,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr);
5011d14d03bSHong Zhang       ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
5021d14d03bSHong Zhang       if (kspreason < 0) {
5031d14d03bSHong Zhang         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
5041d14d03bSHong Zhang         ierr = PetscInfo2(ts,"Step=%D, %Dth cost function, transposed linear solve fails, stopping adjoint solve\n",ts->steps,nadj);CHKERRQ(ierr);
5051d14d03bSHong Zhang       }
506b5b4017aSHong Zhang     }
507b5b4017aSHong Zhang   }
508b5b4017aSHong Zhang 
50936eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
5101d14d03bSHong Zhang   if(th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
511022c081aSHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
512cd4cee2dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
51333f72c5dSHong Zhang     /* R_U at t_n */
5147207e0fdSHong Zhang     if (quadts) {
515cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
5167207e0fdSHong Zhang     }
517abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
518b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
519cd4cee2dSHong Zhang       if (quadJ) {
520cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
521cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
522cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vec_drdu_col);CHKERRQ(ierr);
523cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
524cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
525b5b4017aSHong Zhang       }
5261d14d03bSHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
527b5b4017aSHong Zhang     }
5281d14d03bSHong Zhang 
5291d14d03bSHong Zhang     /* Second-order adjoint */
5301d14d03bSHong Zhang     if (ts->vecs_sensi2) { /* U_n */
531b5b4017aSHong Zhang       /* Get w1 at t_n from TLM matrix */
532b5b4017aSHong Zhang       ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
533b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
534b5b4017aSHong Zhang       /* lambda_s^T F_UU w_1 */
535*b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr);
536b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
53733f72c5dSHong Zhang       ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
538b5b4017aSHong Zhang       /* lambda_s^T F_UU w_2 */
539*b1e111ebSHong Zhang       ierr = TSComputeIHessianProductFunctionUP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr);
540b5b4017aSHong Zhang       for (nadj=0; nadj<ts->numcost; nadj++) {
54133f72c5dSHong 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  */
542b5b4017aSHong Zhang         ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr);
54333f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr);
544b5b4017aSHong Zhang         if (ts->vecs_fup) {
54533f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr);
546b5b4017aSHong Zhang         }
5471d14d03bSHong Zhang         ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr);
548b5b4017aSHong Zhang       }
549b5e0532dSHong Zhang     }
5503fd52205SHong Zhang 
5513fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
552b5b4017aSHong Zhang       /* U_{n+1} */
55333f72c5dSHong Zhang       shift = -1./(th->Theta*th->time_step);
55433f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5557207e0fdSHong Zhang       if (quadts) {
556cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
5577207e0fdSHong Zhang       }
558abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
559b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
56033f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr);
5613fd52205SHong Zhang       }
56233f72c5dSHong Zhang       if (ts->vecs_sensi2p) { /* second-order */
56333f72c5dSHong Zhang         /* Get w1 at t_{n+1} from TLM matrix */
56433f72c5dSHong Zhang         ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr);
56533f72c5dSHong Zhang         ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
566b5b4017aSHong Zhang         /* lambda_s^T F_PU w_1 */
567*b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPU(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
56833f72c5dSHong Zhang         ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
56933f72c5dSHong Zhang         ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
57033f72c5dSHong Zhang 
571b5b4017aSHong Zhang         /* lambda_s^T F_PP w_2 */
572*b1e111ebSHong Zhang         ierr = TSComputeIHessianProductFunctionPP(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
573b5b4017aSHong Zhang         for (nadj=0; nadj<ts->numcost; nadj++) {
57433f72c5dSHong 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)  */
575b5b4017aSHong Zhang           ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
57633f72c5dSHong Zhang           ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr);
577b5b4017aSHong Zhang           if (ts->vecs_fpu) {
57833f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr);
579b5b4017aSHong Zhang           }
580b5b4017aSHong Zhang           if (ts->vecs_fpp) {
58133f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr);
582b5b4017aSHong Zhang           }
583b5b4017aSHong Zhang         }
584b5b4017aSHong Zhang       }
585b5b4017aSHong Zhang 
586b5b4017aSHong Zhang       /* U_s */
58733f72c5dSHong Zhang       shift = 1./((th->Theta-1.0)*th->time_step);
58833f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
5897207e0fdSHong Zhang       if (quadts) {
590cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
5917207e0fdSHong Zhang       }
592abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
593b5e0532dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
59433f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr);
59533f72c5dSHong Zhang         if (ts->vecs_sensi2p) { /* second-order */
59633f72c5dSHong Zhang           /* Get w1 at t_n from TLM matrix */
59733f72c5dSHong Zhang           ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr);
59833f72c5dSHong Zhang           ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
599b5b4017aSHong Zhang           /* lambda_s^T F_PU w_1 */
600*b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPU(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr);
60133f72c5dSHong Zhang           ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
60233f72c5dSHong Zhang           ierr = MatDenseRestoreColumn(th->MatFwdSensip0,&xarr);CHKERRQ(ierr);
603b5b4017aSHong Zhang           /* lambda_s^T F_PP w_2 */
604*b1e111ebSHong Zhang           ierr = TSComputeIHessianProductFunctionPP(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr);
605b5b4017aSHong Zhang           for (nadj=0; nadj<ts->numcost; nadj++) {
60633f72c5dSHong 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) */
607b5b4017aSHong Zhang             ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr);
60833f72c5dSHong Zhang             ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr);
609b5b4017aSHong Zhang             if (ts->vecs_fpu) {
61033f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr);
611b5e0532dSHong Zhang             }
612b5b4017aSHong Zhang             if (ts->vecs_fpp) {
61333f72c5dSHong Zhang               ierr = VecAXPY(ts->vecs_sensi2p[nadj],-th->time_step*(1.0-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr);
614b5b4017aSHong Zhang             }
61536eaed60SHong Zhang           }
61636eaed60SHong Zhang         }
6173fd52205SHong Zhang       }
618b5e0532dSHong Zhang     }
6193fd52205SHong Zhang   } else { /* one-stage case */
6203c54f92cSHong Zhang     shift = 0.0;
621cd4cee2dSHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
6227207e0fdSHong Zhang     if (quadts) {
623cd4cee2dSHong Zhang       ierr  = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
6247207e0fdSHong Zhang     }
625abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
626b5e0532dSHong Zhang       ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr);
627022c081aSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr);
628cd4cee2dSHong Zhang       if (quadJ) {
629cd4cee2dSHong Zhang         ierr = MatDenseGetColumn(quadJ,nadj,&xarr);CHKERRQ(ierr);
630cd4cee2dSHong Zhang         ierr = VecPlaceArray(ts->vec_drdu_col,xarr);CHKERRQ(ierr);
631cd4cee2dSHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vec_drdu_col);CHKERRQ(ierr);
632cd4cee2dSHong Zhang         ierr = VecResetArray(ts->vec_drdu_col);CHKERRQ(ierr);
633cd4cee2dSHong Zhang         ierr = MatDenseRestoreColumn(quadJ,&xarr);CHKERRQ(ierr);
63436eaed60SHong Zhang       }
6352ca6e920SHong Zhang     }
6363fd52205SHong Zhang     if (ts->vecs_sensip) {
63733f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
6387207e0fdSHong Zhang       if (quadts) {
639cd4cee2dSHong Zhang         ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
6407207e0fdSHong Zhang       }
641abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
64233f72c5dSHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr);
64333f72c5dSHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr);
644cd4cee2dSHong Zhang         if (quadJp) {
645cd4cee2dSHong Zhang           ierr = MatDenseGetColumn(quadJp,nadj,&xarr);CHKERRQ(ierr);
646cd4cee2dSHong Zhang           ierr = VecPlaceArray(ts->vec_drdp_col,xarr);CHKERRQ(ierr);
647cd4cee2dSHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vec_drdp_col);CHKERRQ(ierr);
648cd4cee2dSHong Zhang           ierr = VecResetArray(ts->vec_drdp_col);CHKERRQ(ierr);
649cd4cee2dSHong Zhang           ierr = MatDenseRestoreColumn(quadJp,&xarr);CHKERRQ(ierr);
65036eaed60SHong Zhang         }
65136eaed60SHong Zhang       }
6523fd52205SHong Zhang     }
6532ca6e920SHong Zhang   }
6542ca6e920SHong Zhang 
6552ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
6562ca6e920SHong Zhang   PetscFunctionReturn(0);
6572ca6e920SHong Zhang }
6582ca6e920SHong Zhang 
659cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
660cd652676SJed Brown {
661cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
662be5899b3SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
663cd652676SJed Brown   PetscErrorCode ierr;
664cd652676SJed Brown 
665cd652676SJed Brown   PetscFunctionBegin;
666a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
667be5899b3SLisandro Dalcin   if (th->endpoint) dt *= th->Theta;
668be5899b3SLisandro Dalcin   ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr);
669cd652676SJed Brown   PetscFunctionReturn(0);
670cd652676SJed Brown }
671cd652676SJed Brown 
6721566a47fSLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
6731566a47fSLisandro Dalcin {
6741566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
6751566a47fSLisandro Dalcin   Vec            X = ts->vec_sol;      /* X = solution */
6761566a47fSLisandro Dalcin   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
6777453f775SEmil Constantinescu   PetscReal      wltea,wlter;
6781566a47fSLisandro Dalcin   PetscErrorCode ierr;
6791566a47fSLisandro Dalcin 
6801566a47fSLisandro Dalcin   PetscFunctionBegin;
6812ffb9264SLisandro Dalcin   if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);}
6821566a47fSLisandro Dalcin   /* Cannot compute LTE in first step or in restart after event */
683fecfb714SLisandro Dalcin   if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);}
6841566a47fSLisandro Dalcin   /* Compute LTE using backward differences with non-constant time step */
685fecfb714SLisandro Dalcin   {
686be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
687be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
6881566a47fSLisandro Dalcin     PetscScalar scal[3]; Vec vecs[3];
6891566a47fSLisandro Dalcin     scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
6901566a47fSLisandro Dalcin     vecs[0] = X;    vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
6911566a47fSLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
6921566a47fSLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
6937453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
6941566a47fSLisandro Dalcin   }
6951566a47fSLisandro Dalcin   if (order) *order = 2;
6961566a47fSLisandro Dalcin   PetscFunctionReturn(0);
6971566a47fSLisandro Dalcin }
6981566a47fSLisandro Dalcin 
6991566a47fSLisandro Dalcin static PetscErrorCode TSRollBack_Theta(TS ts)
7001566a47fSLisandro Dalcin {
7011566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
7027207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
7031566a47fSLisandro Dalcin   PetscErrorCode ierr;
7041566a47fSLisandro Dalcin 
7051566a47fSLisandro Dalcin   PetscFunctionBegin;
7061566a47fSLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
707cd4cee2dSHong Zhang   if (quadts && ts->costintegralfwd) {
708cd4cee2dSHong Zhang     ierr = VecCopy(th->VecCostIntegral0,quadts->vec_sol);CHKERRQ(ierr);
7091566a47fSLisandro Dalcin   }
710715f1b00SHong Zhang   th->status = TS_STEP_INCOMPLETE;
71113b1b0ffSHong Zhang   if (ts->mat_sensip) {
71213b1b0ffSHong Zhang     ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7136f92202bSHong Zhang   }
7147207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7157207e0fdSHong Zhang     ierr = MatCopy(th->MatIntegralSensip0,quadts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7166f92202bSHong Zhang   }
717715f1b00SHong Zhang   PetscFunctionReturn(0);
718715f1b00SHong Zhang }
719715f1b00SHong Zhang 
720715f1b00SHong Zhang static PetscErrorCode TSForwardStep_Theta(TS ts)
721715f1b00SHong Zhang {
722715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
723cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
72413b1b0ffSHong Zhang   Mat            MatDeltaFwdSensip = th->MatDeltaFwdSensip;
725b5b4017aSHong Zhang   Vec            VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
7267207e0fdSHong Zhang   PetscInt       ntlm;
727715f1b00SHong Zhang   KSP            ksp;
7287207e0fdSHong Zhang   Mat            J,Jpre,quadJ = NULL,quadJp = NULL;
729715f1b00SHong Zhang   PetscReal      shift;
73013b1b0ffSHong Zhang   PetscScalar    *barr,*xarr;
731715f1b00SHong Zhang   PetscErrorCode ierr;
732715f1b00SHong Zhang 
733715f1b00SHong Zhang   PetscFunctionBegin;
73413b1b0ffSHong Zhang   ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
73513b1b0ffSHong Zhang 
7367207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
7377207e0fdSHong Zhang     ierr = MatCopy(quadts->mat_sensip,th->MatIntegralSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7386f92202bSHong Zhang   }
739715f1b00SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr);
740cd4cee2dSHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jpre,NULL,NULL);CHKERRQ(ierr);
7417207e0fdSHong Zhang   if (quadts) {
742cd4cee2dSHong Zhang     ierr = TSGetRHSJacobian(quadts,&quadJ,NULL,NULL,NULL);CHKERRQ(ierr);
743cd4cee2dSHong Zhang     ierr = TSGetRHSJacobianP(quadts,&quadJp,NULL,NULL);CHKERRQ(ierr);
7447207e0fdSHong Zhang   }
745715f1b00SHong Zhang 
746715f1b00SHong Zhang   /* Build RHS */
747715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method*/
748715f1b00SHong Zhang     shift = 1./((th->Theta-1.)*th->time_step);
749cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
75013b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
75113b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
752715f1b00SHong Zhang 
753022c081aSHong Zhang     /* Add the f_p forcing terms */
7540e11c21fSHong Zhang     if (ts->Jacp) {
75533f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->ptime,th->X0,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
75633f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
75733f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
75833f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
7590e11c21fSHong Zhang     }
760715f1b00SHong Zhang   } else { /* 1-stage method */
761715f1b00SHong Zhang     shift = 0.0;
762cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
76313b1b0ffSHong Zhang     ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr);
76413b1b0ffSHong Zhang     ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr);
76513b1b0ffSHong Zhang 
766022c081aSHong Zhang     /* Add the f_p forcing terms */
7670e11c21fSHong Zhang     if (ts->Jacp) {
76833f72c5dSHong Zhang       ierr = TSComputeIJacobianP(ts,th->stage_time,th->X,th->Xdot,shift,ts->Jacp,PETSC_FALSE);CHKERRQ(ierr);
76933f72c5dSHong Zhang       ierr = MatAXPY(MatDeltaFwdSensip,-1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr);
770715f1b00SHong Zhang     }
7710e11c21fSHong Zhang   }
772715f1b00SHong Zhang 
773715f1b00SHong Zhang   /* Build LHS */
774715f1b00SHong Zhang   shift = 1/(th->Theta*th->time_step);
775715f1b00SHong Zhang   if (th->endpoint) {
776cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
777715f1b00SHong Zhang   } else {
778cd4cee2dSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jpre,PETSC_FALSE);CHKERRQ(ierr);
779715f1b00SHong Zhang   }
780cd4cee2dSHong Zhang   ierr = KSPSetOperators(ksp,J,Jpre);CHKERRQ(ierr);
781715f1b00SHong Zhang 
782715f1b00SHong Zhang   /*
783715f1b00SHong Zhang     Evaluate the first stage of integral gradients with the 2-stage method:
784c9ad9de0SHong Zhang     drdu|t_n*S(t_n) + drdp|t_n
785715f1b00SHong Zhang     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
786715f1b00SHong Zhang   */
787715f1b00SHong Zhang   if (th->endpoint) { /* 2-stage method only */
7887207e0fdSHong Zhang     if (quadts && quadts->mat_sensip) {
789cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->ptime,th->X0,quadJ,NULL);CHKERRQ(ierr);
790cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->ptime,th->X0,quadJp);CHKERRQ(ierr);
7917207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
7927207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7937207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*(1.-th->Theta),th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
794715f1b00SHong Zhang     }
795715f1b00SHong Zhang   }
796715f1b00SHong Zhang 
797715f1b00SHong Zhang   /* Solve the tangent linear equation for forward sensitivities to parameters */
798022c081aSHong Zhang   for (ntlm=0; ntlm<th->num_tlm; ntlm++) {
799b5b4017aSHong Zhang     KSPConvergedReason kspreason;
80013b1b0ffSHong Zhang     ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr);
801b5b4017aSHong Zhang     ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr);
802715f1b00SHong Zhang     if (th->endpoint) {
80313b1b0ffSHong Zhang       ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr);
804b5b4017aSHong Zhang       ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr);
805b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr);
806b5b4017aSHong Zhang       ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr);
80713b1b0ffSHong Zhang       ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr);
808715f1b00SHong Zhang     } else {
809b5b4017aSHong Zhang       ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr);
810715f1b00SHong Zhang     }
811b5b4017aSHong Zhang     ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr);
812b5b4017aSHong Zhang     if (kspreason < 0) {
813b5b4017aSHong Zhang       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
814b5b4017aSHong Zhang       ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr);
815b5b4017aSHong Zhang     }
816b5b4017aSHong Zhang     ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr);
81713b1b0ffSHong Zhang     ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr);
818715f1b00SHong Zhang   }
819715f1b00SHong Zhang 
820b5b4017aSHong Zhang 
82113b1b0ffSHong Zhang   /*
82213b1b0ffSHong Zhang     Evaluate the second stage of integral gradients with the 2-stage method:
823c9ad9de0SHong Zhang     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
82413b1b0ffSHong Zhang   */
8257207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
82613b1b0ffSHong Zhang     if (!th->endpoint) {
8277207e0fdSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); /* stage sensitivity */
828cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,th->X,quadJ,NULL);CHKERRQ(ierr);
829cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,th->X,quadJp);CHKERRQ(ierr);
8307207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8317207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8327207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
83313b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
83413b1b0ffSHong Zhang     } else {
835cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobian(quadts,th->stage_time,ts->vec_sol,quadJ,NULL);CHKERRQ(ierr);
836cd4cee2dSHong Zhang       ierr = TSComputeRHSJacobianP(quadts,th->stage_time,ts->vec_sol,quadJp);CHKERRQ(ierr);
8377207e0fdSHong Zhang       ierr = MatTransposeMatMult(ts->mat_sensip,quadJ,MAT_REUSE_MATRIX,PETSC_DEFAULT,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
8387207e0fdSHong Zhang       ierr = MatAXPY(th->MatIntegralSensipTemp,1,quadJp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
8397207e0fdSHong Zhang       ierr = MatAXPY(quadts->mat_sensip,th->time_step*th->Theta,th->MatIntegralSensipTemp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
84013b1b0ffSHong Zhang     }
84113b1b0ffSHong Zhang   } else {
84213b1b0ffSHong Zhang     if (!th->endpoint) {
84313b1b0ffSHong Zhang       ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
844715f1b00SHong Zhang     }
845715f1b00SHong Zhang   }
8461566a47fSLisandro Dalcin   PetscFunctionReturn(0);
8471566a47fSLisandro Dalcin }
8481566a47fSLisandro Dalcin 
8496affb6f8SHong Zhang static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip)
8506affb6f8SHong Zhang {
8516affb6f8SHong Zhang   TS_Theta *th = (TS_Theta*)ts->data;
8526affb6f8SHong Zhang 
8536affb6f8SHong Zhang   PetscFunctionBegin;
8546affb6f8SHong Zhang   if (ns) *ns = 1;
8556affb6f8SHong Zhang   if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip);
8566affb6f8SHong Zhang   PetscFunctionReturn(0);
8576affb6f8SHong Zhang }
8586affb6f8SHong Zhang 
859316643e7SJed Brown /*------------------------------------------------------------*/
860277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
861316643e7SJed Brown {
862316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
863316643e7SJed Brown   PetscErrorCode ierr;
864316643e7SJed Brown 
865316643e7SJed Brown   PetscFunctionBegin;
8666bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
8676bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
8683b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
869eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
8701566a47fSLisandro Dalcin 
8711566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
8721566a47fSLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
8731566a47fSLisandro Dalcin 
874b1cb36f3SHong Zhang   ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr);
875277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
876277b19d0SLisandro Dalcin }
877277b19d0SLisandro Dalcin 
878ece46509SHong Zhang static PetscErrorCode TSAdjointReset_Theta(TS ts)
879ece46509SHong Zhang {
880ece46509SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
881ece46509SHong Zhang   PetscErrorCode ierr;
882ece46509SHong Zhang 
883ece46509SHong Zhang   PetscFunctionBegin;
884ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
885ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
886ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
887ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
888ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
889ece46509SHong Zhang   ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
890ece46509SHong Zhang   PetscFunctionReturn(0);
891ece46509SHong Zhang }
892ece46509SHong Zhang 
893277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
894277b19d0SLisandro Dalcin {
895277b19d0SLisandro Dalcin   PetscErrorCode ierr;
896277b19d0SLisandro Dalcin 
897277b19d0SLisandro Dalcin   PetscFunctionBegin;
898277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
899b3a6b972SJed Brown   if (ts->dm) {
900b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
901b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
902b3a6b972SJed Brown   }
903277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
904bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
905bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
906bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
907bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
908316643e7SJed Brown   PetscFunctionReturn(0);
909316643e7SJed Brown }
910316643e7SJed Brown 
911316643e7SJed Brown /*
912316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9132b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
914316643e7SJed Brown */
9150f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
916316643e7SJed Brown {
917316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
918316643e7SJed Brown   PetscErrorCode ierr;
9197445fe48SJed Brown   Vec            X0,Xdot;
9207445fe48SJed Brown   DM             dm,dmsave;
9211566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
922316643e7SJed Brown 
923316643e7SJed Brown   PetscFunctionBegin;
9247445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9255a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
9267445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
927b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
9287445fe48SJed Brown 
9297445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
9307445fe48SJed Brown   dmsave = ts->dm;
9317445fe48SJed Brown   ts->dm = dm;
9327445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
9337445fe48SJed Brown   ts->dm = dmsave;
9340d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
935316643e7SJed Brown   PetscFunctionReturn(0);
936316643e7SJed Brown }
937316643e7SJed Brown 
938d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
939316643e7SJed Brown {
940316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
941316643e7SJed Brown   PetscErrorCode ierr;
9427445fe48SJed Brown   Vec            Xdot;
9437445fe48SJed Brown   DM             dm,dmsave;
9441566a47fSLisandro Dalcin   PetscReal      shift = 1/(th->Theta*ts->time_step);
945316643e7SJed Brown 
946316643e7SJed Brown   PetscFunctionBegin;
9477445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
948be5899b3SLisandro Dalcin   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
9490298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
9507445fe48SJed Brown 
9517445fe48SJed Brown   dmsave = ts->dm;
9527445fe48SJed Brown   ts->dm = dm;
953d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
9547445fe48SJed Brown   ts->dm = dmsave;
9550298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
956316643e7SJed Brown   PetscFunctionReturn(0);
957316643e7SJed Brown }
958316643e7SJed Brown 
959715f1b00SHong Zhang static PetscErrorCode TSForwardSetUp_Theta(TS ts)
960715f1b00SHong Zhang {
961715f1b00SHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9627207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
963715f1b00SHong Zhang   PetscErrorCode ierr;
964715f1b00SHong Zhang 
965715f1b00SHong Zhang   PetscFunctionBegin;
966022c081aSHong Zhang   /* combine sensitivities to parameters and sensitivities to initial values into one array */
96713b1b0ffSHong Zhang   th->num_tlm = ts->num_parameters;
96813b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr);
9697207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9707207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensipTemp);CHKERRQ(ierr);
9717207e0fdSHong Zhang     ierr = MatDuplicate(quadts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatIntegralSensip0);CHKERRQ(ierr);
972715f1b00SHong Zhang   }
9736f92202bSHong Zhang   /* backup sensitivity results for roll-backs */
97413b1b0ffSHong Zhang   ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr);
97513b1b0ffSHong Zhang 
976b5b4017aSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
977715f1b00SHong Zhang   PetscFunctionReturn(0);
978715f1b00SHong Zhang }
979715f1b00SHong Zhang 
9807adebddeSHong Zhang static PetscErrorCode TSForwardReset_Theta(TS ts)
9817adebddeSHong Zhang {
9827adebddeSHong Zhang   TS_Theta       *th = (TS_Theta*)ts->data;
9837207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
9847adebddeSHong Zhang   PetscErrorCode ierr;
9857adebddeSHong Zhang 
9867adebddeSHong Zhang   PetscFunctionBegin;
9877207e0fdSHong Zhang   if (quadts && quadts->mat_sensip) {
9887207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensipTemp);CHKERRQ(ierr);
9897207e0fdSHong Zhang     ierr = MatDestroy(&th->MatIntegralSensip0);CHKERRQ(ierr);
9907adebddeSHong Zhang   }
9917adebddeSHong Zhang   ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr);
9927adebddeSHong Zhang   ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr);
9937adebddeSHong Zhang   ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr);
9947adebddeSHong Zhang   PetscFunctionReturn(0);
9957adebddeSHong Zhang }
9967adebddeSHong Zhang 
997316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
998316643e7SJed Brown {
999316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1000cd4cee2dSHong Zhang   TS             quadts = ts->quadraturets;
10012ffb9264SLisandro Dalcin   PetscBool      match;
1002316643e7SJed Brown   PetscErrorCode ierr;
1003316643e7SJed Brown 
1004316643e7SJed Brown   PetscFunctionBegin;
1005cd4cee2dSHong Zhang   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1006cd4cee2dSHong Zhang     ierr = VecDuplicate(quadts->vec_sol,&th->VecCostIntegral0);CHKERRQ(ierr);
1007b1cb36f3SHong Zhang   }
100839d13666SHong Zhang   if (!th->X) {
1009316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
101039d13666SHong Zhang   }
101139d13666SHong Zhang   if (!th->Xdot) {
1012316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
101339d13666SHong Zhang   }
101439d13666SHong Zhang   if (!th->X0) {
10153b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
101639d13666SHong Zhang   }
10171566a47fSLisandro Dalcin   if (th->endpoint) {
10181566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);
10197445fe48SJed Brown   }
10203b1890cdSShri Abhyankar 
10211566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
10221566a47fSLisandro Dalcin 
10231566a47fSLisandro Dalcin   ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr);
10241566a47fSLisandro Dalcin   ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10251566a47fSLisandro Dalcin   ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
10261566a47fSLisandro Dalcin 
10271566a47fSLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
10281566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
10292ffb9264SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr);
10302ffb9264SLisandro Dalcin   if (!match) {
10311566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
10321566a47fSLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
10333b1890cdSShri Abhyankar   }
10341566a47fSLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
1035b4dd3bf9SBarry Smith   PetscFunctionReturn(0);
1036b4dd3bf9SBarry Smith }
10370c7376e5SHong Zhang 
1038b4dd3bf9SBarry Smith /*------------------------------------------------------------*/
1039b4dd3bf9SBarry Smith 
1040b4dd3bf9SBarry Smith static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1041b4dd3bf9SBarry Smith {
1042b4dd3bf9SBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
1043b4dd3bf9SBarry Smith   PetscErrorCode ierr;
1044b4dd3bf9SBarry Smith 
1045b4dd3bf9SBarry Smith   PetscFunctionBegin;
1046b5e0532dSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr);
1047b5b4017aSHong Zhang   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr);
10482ca6e920SHong Zhang   if (ts->vecs_sensip) {
1049b5e0532dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr);
10502ca6e920SHong Zhang   }
1051b5b4017aSHong Zhang   if (ts->vecs_sensi2) {
1052b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr);
1053b5b4017aSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr);
105467633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
105567633408SHong Zhang     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
105667633408SHong Zhang     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1057b5b4017aSHong Zhang   }
1058c9681e5dSHong Zhang   if (ts->vecs_sensi2p) {
1059c9681e5dSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi2p[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr);
106067633408SHong Zhang     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
106167633408SHong Zhang     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
106267633408SHong Zhang     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1063b5b4017aSHong Zhang   }
1064316643e7SJed Brown   PetscFunctionReturn(0);
1065316643e7SJed Brown }
1066316643e7SJed Brown 
10674416b707SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts)
1068316643e7SJed Brown {
1069316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1070316643e7SJed Brown   PetscErrorCode ierr;
1071316643e7SJed Brown 
1072316643e7SJed Brown   PetscFunctionBegin;
1073e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
1074316643e7SJed Brown   {
10750298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
10760298fd71SBarry 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);
107703842d09SLisandro 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);
1078316643e7SJed Brown   }
1079316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1080316643e7SJed Brown   PetscFunctionReturn(0);
1081316643e7SJed Brown }
1082316643e7SJed Brown 
1083316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
1084316643e7SJed Brown {
1085316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1086ace3abfcSBarry Smith   PetscBool      iascii;
1087316643e7SJed Brown   PetscErrorCode ierr;
1088316643e7SJed Brown 
1089316643e7SJed Brown   PetscFunctionBegin;
1090251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1091316643e7SJed Brown   if (iascii) {
10927c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
1093ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
1094316643e7SJed Brown   }
1095316643e7SJed Brown   PetscFunctionReturn(0);
1096316643e7SJed Brown }
1097316643e7SJed Brown 
1098560360afSLisandro Dalcin static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
10990de4c49aSJed Brown {
11000de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11010de4c49aSJed Brown 
11020de4c49aSJed Brown   PetscFunctionBegin;
11030de4c49aSJed Brown   *theta = th->Theta;
11040de4c49aSJed Brown   PetscFunctionReturn(0);
11050de4c49aSJed Brown }
11060de4c49aSJed Brown 
1107560360afSLisandro Dalcin static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta)
11080de4c49aSJed Brown {
11090de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
11100de4c49aSJed Brown 
11110de4c49aSJed Brown   PetscFunctionBegin;
11127c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
11130de4c49aSJed Brown   th->Theta = theta;
11141566a47fSLisandro Dalcin   th->order = (th->Theta == 0.5) ? 2 : 1;
11150de4c49aSJed Brown   PetscFunctionReturn(0);
11160de4c49aSJed Brown }
1117eb284becSJed Brown 
1118560360afSLisandro Dalcin static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
111926f2ff8fSLisandro Dalcin {
112026f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
112126f2ff8fSLisandro Dalcin 
112226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
112326f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
112426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
112526f2ff8fSLisandro Dalcin }
112626f2ff8fSLisandro Dalcin 
1127560360afSLisandro Dalcin static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
1128eb284becSJed Brown {
1129eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
1130eb284becSJed Brown 
1131eb284becSJed Brown   PetscFunctionBegin;
1132eb284becSJed Brown   th->endpoint = flg;
1133eb284becSJed Brown   PetscFunctionReturn(0);
1134eb284becSJed Brown }
11350de4c49aSJed Brown 
1136f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1137f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
1138f9c1d6abSBarry Smith {
1139f9c1d6abSBarry Smith   PetscComplex   z   = xr + xi*PETSC_i,f;
1140f9c1d6abSBarry Smith   TS_Theta       *th = (TS_Theta*)ts->data;
11413fd8ae06SJed Brown   const PetscReal one = 1.0;
1142f9c1d6abSBarry Smith 
1143f9c1d6abSBarry Smith   PetscFunctionBegin;
11443fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
1145f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
1146f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
1147f9c1d6abSBarry Smith   PetscFunctionReturn(0);
1148f9c1d6abSBarry Smith }
1149f9c1d6abSBarry Smith #endif
1150f9c1d6abSBarry Smith 
115142682096SHong Zhang static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
115242682096SHong Zhang {
115342682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
115442682096SHong Zhang 
115542682096SHong Zhang   PetscFunctionBegin;
11561566a47fSLisandro Dalcin   if (ns) *ns = 1;
11571566a47fSLisandro Dalcin   if (Y)  *Y  = th->endpoint ? &(th->X0) : &(th->X);
115842682096SHong Zhang   PetscFunctionReturn(0);
115942682096SHong Zhang }
1160f9c1d6abSBarry Smith 
1161316643e7SJed Brown /* ------------------------------------------------------------ */
1162316643e7SJed Brown /*MC
116396f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
1164316643e7SJed Brown 
1165316643e7SJed Brown    Level: beginner
1166316643e7SJed Brown 
11674eb428fdSBarry Smith    Options Database:
1168c82ed962SBarry Smith +  -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
1169c82ed962SBarry Smith .  -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
117003842d09SLisandro Dalcin -  -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable)
11714eb428fdSBarry Smith 
1172eb284becSJed Brown    Notes:
11730c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
11740c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
11754eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
11764eb428fdSBarry Smith 
1177eb284becSJed Brown    This method can be applied to DAE.
1178eb284becSJed Brown 
1179eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
1180eb284becSJed Brown 
1181eb284becSJed Brown .vb
1182eb284becSJed Brown   Theta | Theta
1183eb284becSJed Brown   -------------
1184eb284becSJed Brown         |  1
1185eb284becSJed Brown .ve
1186eb284becSJed Brown 
1187eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
1188eb284becSJed Brown 
1189eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
1190eb284becSJed Brown 
1191eb284becSJed Brown .vb
1192eb284becSJed Brown   0 | 0         0
1193eb284becSJed Brown   1 | 1-Theta   Theta
1194eb284becSJed Brown   -------------------
1195eb284becSJed Brown     | 1-Theta   Theta
1196eb284becSJed Brown .ve
1197eb284becSJed Brown 
1198eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
1199eb284becSJed Brown 
1200eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
1201eb284becSJed Brown 
1202eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
1203eb284becSJed Brown 
12044eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
1205eb284becSJed Brown 
1206eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
1207316643e7SJed Brown 
1208316643e7SJed Brown M*/
12098cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1210316643e7SJed Brown {
1211316643e7SJed Brown   TS_Theta       *th;
1212316643e7SJed Brown   PetscErrorCode ierr;
1213316643e7SJed Brown 
1214316643e7SJed Brown   PetscFunctionBegin;
1215277b19d0SLisandro Dalcin   ts->ops->reset           = TSReset_Theta;
1216ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1217316643e7SJed Brown   ts->ops->destroy         = TSDestroy_Theta;
1218316643e7SJed Brown   ts->ops->view            = TSView_Theta;
1219316643e7SJed Brown   ts->ops->setup           = TSSetUp_Theta;
122042f2b339SBarry Smith   ts->ops->adjointsetup    = TSAdjointSetUp_Theta;
1221ece46509SHong Zhang   ts->ops->adjointreset    = TSAdjointReset_Theta;
1222316643e7SJed Brown   ts->ops->step            = TSStep_Theta;
1223cd652676SJed Brown   ts->ops->interpolate     = TSInterpolate_Theta;
12241566a47fSLisandro Dalcin   ts->ops->evaluatewlte    = TSEvaluateWLTE_Theta;
122524655328SShri   ts->ops->rollback        = TSRollBack_Theta;
1226316643e7SJed Brown   ts->ops->setfromoptions  = TSSetFromOptions_Theta;
12270f5c6efeSJed Brown   ts->ops->snesfunction    = SNESTSFormFunction_Theta;
12280f5c6efeSJed Brown   ts->ops->snesjacobian    = SNESTSFormJacobian_Theta;
1229f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
1230f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
1231f9c1d6abSBarry Smith #endif
123242682096SHong Zhang   ts->ops->getstages       = TSGetStages_Theta;
123342f2b339SBarry Smith   ts->ops->adjointstep     = TSAdjointStep_Theta;
1234b1cb36f3SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1235b1cb36f3SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
12362ffb9264SLisandro Dalcin   ts->default_adapt_type   = TSADAPTNONE;
12376affb6f8SHong Zhang 
1238715f1b00SHong Zhang   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
12397adebddeSHong Zhang   ts->ops->forwardreset     = TSForwardReset_Theta;
1240715f1b00SHong Zhang   ts->ops->forwardstep      = TSForwardStep_Theta;
12416affb6f8SHong Zhang   ts->ops->forwardgetstages = TSForwardGetStages_Theta;
1242316643e7SJed Brown 
1243efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1244efd4aadfSBarry Smith 
1245b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1246316643e7SJed Brown   ts->data = (void*)th;
1247316643e7SJed Brown 
1248715f1b00SHong Zhang   th->VecsDeltaLam    = NULL;
1249715f1b00SHong Zhang   th->VecsDeltaMu     = NULL;
1250715f1b00SHong Zhang   th->VecsSensiTemp   = NULL;
1251b5b4017aSHong Zhang   th->VecsSensi2Temp  = NULL;
1252715f1b00SHong Zhang 
12536f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
1254316643e7SJed Brown   th->Theta       = 0.5;
12551566a47fSLisandro Dalcin   th->order       = 2;
1256bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
1257bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
1258bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
1259bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
1260316643e7SJed Brown   PetscFunctionReturn(0);
1261316643e7SJed Brown }
12620de4c49aSJed Brown 
12630de4c49aSJed Brown /*@
12640de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
12650de4c49aSJed Brown 
12660de4c49aSJed Brown   Not Collective
12670de4c49aSJed Brown 
12680de4c49aSJed Brown   Input Parameter:
12690de4c49aSJed Brown .  ts - timestepping context
12700de4c49aSJed Brown 
12710de4c49aSJed Brown   Output Parameter:
12720de4c49aSJed Brown .  theta - stage abscissa
12730de4c49aSJed Brown 
12740de4c49aSJed Brown   Note:
12750de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
12760de4c49aSJed Brown 
12770de4c49aSJed Brown   Level: Advanced
12780de4c49aSJed Brown 
12790de4c49aSJed Brown .seealso: TSThetaSetTheta()
12800de4c49aSJed Brown @*/
12817087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
12820de4c49aSJed Brown {
12834ac538c5SBarry Smith   PetscErrorCode ierr;
12840de4c49aSJed Brown 
12850de4c49aSJed Brown   PetscFunctionBegin;
1286afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12870de4c49aSJed Brown   PetscValidPointer(theta,2);
12884ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
12890de4c49aSJed Brown   PetscFunctionReturn(0);
12900de4c49aSJed Brown }
12910de4c49aSJed Brown 
12920de4c49aSJed Brown /*@
12930de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
12940de4c49aSJed Brown 
12950de4c49aSJed Brown   Not Collective
12960de4c49aSJed Brown 
12970de4c49aSJed Brown   Input Parameter:
12980de4c49aSJed Brown +  ts - timestepping context
12990de4c49aSJed Brown -  theta - stage abscissa
13000de4c49aSJed Brown 
13010de4c49aSJed Brown   Options Database:
13020de4c49aSJed Brown .  -ts_theta_theta <theta>
13030de4c49aSJed Brown 
13040de4c49aSJed Brown   Level: Intermediate
13050de4c49aSJed Brown 
13060de4c49aSJed Brown .seealso: TSThetaGetTheta()
13070de4c49aSJed Brown @*/
13087087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
13090de4c49aSJed Brown {
13104ac538c5SBarry Smith   PetscErrorCode ierr;
13110de4c49aSJed Brown 
13120de4c49aSJed Brown   PetscFunctionBegin;
1313afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13144ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
13150de4c49aSJed Brown   PetscFunctionReturn(0);
13160de4c49aSJed Brown }
1317f33bbcb6SJed Brown 
131826f2ff8fSLisandro Dalcin /*@
131926f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
132026f2ff8fSLisandro Dalcin 
132126f2ff8fSLisandro Dalcin   Not Collective
132226f2ff8fSLisandro Dalcin 
132326f2ff8fSLisandro Dalcin   Input Parameter:
132426f2ff8fSLisandro Dalcin .  ts - timestepping context
132526f2ff8fSLisandro Dalcin 
132626f2ff8fSLisandro Dalcin   Output Parameter:
132726f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
132826f2ff8fSLisandro Dalcin 
132926f2ff8fSLisandro Dalcin   Level: Advanced
133026f2ff8fSLisandro Dalcin 
133126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
133226f2ff8fSLisandro Dalcin @*/
133326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
133426f2ff8fSLisandro Dalcin {
133526f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
133626f2ff8fSLisandro Dalcin 
133726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
133826f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
133926f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
1340163d334eSBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
134126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
134226f2ff8fSLisandro Dalcin }
134326f2ff8fSLisandro Dalcin 
1344eb284becSJed Brown /*@
1345eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
1346eb284becSJed Brown 
1347eb284becSJed Brown   Not Collective
1348eb284becSJed Brown 
1349eb284becSJed Brown   Input Parameter:
1350eb284becSJed Brown +  ts - timestepping context
1351eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
1352eb284becSJed Brown 
1353eb284becSJed Brown   Options Database:
1354eb284becSJed Brown .  -ts_theta_endpoint <flg>
1355eb284becSJed Brown 
1356eb284becSJed Brown   Level: Intermediate
1357eb284becSJed Brown 
1358eb284becSJed Brown .seealso: TSTHETA, TSCN
1359eb284becSJed Brown @*/
1360eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
1361eb284becSJed Brown {
1362eb284becSJed Brown   PetscErrorCode ierr;
1363eb284becSJed Brown 
1364eb284becSJed Brown   PetscFunctionBegin;
1365eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1366eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
1367eb284becSJed Brown   PetscFunctionReturn(0);
1368eb284becSJed Brown }
1369eb284becSJed Brown 
1370f33bbcb6SJed Brown /*
1371f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1372f33bbcb6SJed Brown  * The creation functions for these specializations are below.
1373f33bbcb6SJed Brown  */
1374f33bbcb6SJed Brown 
13751566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_BEuler(TS ts)
13761566a47fSLisandro Dalcin {
13771566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
13781566a47fSLisandro Dalcin   PetscErrorCode ierr;
13791566a47fSLisandro Dalcin 
13801566a47fSLisandro Dalcin   PetscFunctionBegin;
13811566a47fSLisandro 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");
13821566a47fSLisandro 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");
13831566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
13841566a47fSLisandro Dalcin   PetscFunctionReturn(0);
13851566a47fSLisandro Dalcin }
13861566a47fSLisandro Dalcin 
1387f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
1388f33bbcb6SJed Brown {
1389f33bbcb6SJed Brown   PetscFunctionBegin;
1390f33bbcb6SJed Brown   PetscFunctionReturn(0);
1391f33bbcb6SJed Brown }
1392f33bbcb6SJed Brown 
1393f33bbcb6SJed Brown /*MC
1394f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
1395f33bbcb6SJed Brown 
1396f33bbcb6SJed Brown   Level: beginner
1397f33bbcb6SJed Brown 
13984eb428fdSBarry Smith   Notes:
1399c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
14004eb428fdSBarry Smith 
14011566a47fSLisandro Dalcin $  -ts_type theta -ts_theta_theta 1.0
14024eb428fdSBarry Smith 
1403f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
1404f33bbcb6SJed Brown 
1405f33bbcb6SJed Brown M*/
14068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1407f33bbcb6SJed Brown {
1408f33bbcb6SJed Brown   PetscErrorCode ierr;
1409f33bbcb6SJed Brown 
1410f33bbcb6SJed Brown   PetscFunctionBegin;
1411f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1412f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
14131566a47fSLisandro Dalcin   ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr);
14140c7376e5SHong Zhang   ts->ops->setup = TSSetUp_BEuler;
1415f33bbcb6SJed Brown   ts->ops->view  = TSView_BEuler;
1416f33bbcb6SJed Brown   PetscFunctionReturn(0);
1417f33bbcb6SJed Brown }
1418f33bbcb6SJed Brown 
14191566a47fSLisandro Dalcin static PetscErrorCode TSSetUp_CN(TS ts)
14201566a47fSLisandro Dalcin {
14211566a47fSLisandro Dalcin   TS_Theta       *th = (TS_Theta*)ts->data;
14221566a47fSLisandro Dalcin   PetscErrorCode ierr;
14231566a47fSLisandro Dalcin 
14241566a47fSLisandro Dalcin   PetscFunctionBegin;
14251566a47fSLisandro 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");
14261566a47fSLisandro 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");
14271566a47fSLisandro Dalcin   ierr = TSSetUp_Theta(ts);CHKERRQ(ierr);
14281566a47fSLisandro Dalcin   PetscFunctionReturn(0);
14291566a47fSLisandro Dalcin }
14301566a47fSLisandro Dalcin 
1431f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
1432f33bbcb6SJed Brown {
1433f33bbcb6SJed Brown   PetscFunctionBegin;
1434f33bbcb6SJed Brown   PetscFunctionReturn(0);
1435f33bbcb6SJed Brown }
1436f33bbcb6SJed Brown 
1437f33bbcb6SJed Brown /*MC
1438f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
1439f33bbcb6SJed Brown 
1440f33bbcb6SJed Brown   Level: beginner
1441f33bbcb6SJed Brown 
1442f33bbcb6SJed Brown   Notes:
14437cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
14447cf5af47SJed Brown 
14457cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
1446f33bbcb6SJed Brown 
1447f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
1448f33bbcb6SJed Brown 
1449f33bbcb6SJed Brown M*/
14508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1451f33bbcb6SJed Brown {
1452f33bbcb6SJed Brown   PetscErrorCode ierr;
1453f33bbcb6SJed Brown 
1454f33bbcb6SJed Brown   PetscFunctionBegin;
1455f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
1456f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
1457eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
14580c7376e5SHong Zhang   ts->ops->setup = TSSetUp_CN;
1459f33bbcb6SJed Brown   ts->ops->view  = TSView_CN;
1460f33bbcb6SJed Brown   PetscFunctionReturn(0);
1461f33bbcb6SJed Brown }
1462