xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision d4aa0a5823bf4d2cb392c6170f9569e0e6833bf1)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4b45d2f2cSJed Brown #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 {
10316643e7SJed Brown   Vec          X,Xdot;                   /* Storage for one stage */
113b1890cdSShri Abhyankar   Vec          X0;                       /* work vector to store X0 */
12eb284becSJed Brown   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
132ca6e920SHong Zhang   Vec          *VecDeltaLam;             /* Increment of the adjoint sensitivity w.r.t IC at stage*/
142ca6e920SHong Zhang   Vec          *VecDeltaMu;              /* Increment of the adjoint sensitivity w.r.t P at stage*/
152ca6e920SHong Zhang   Vec          *VecSensiTemp;            /* Vector to be timed with Jacobian transpose*/
16ace3abfcSBarry Smith   PetscBool    extrapolate;
17eb284becSJed Brown   PetscBool    endpoint;
18316643e7SJed Brown   PetscReal    Theta;
19316643e7SJed Brown   PetscReal    stage_time;
203b1890cdSShri Abhyankar   TSStepStatus status;
213b1890cdSShri Abhyankar   char         *name;
223b1890cdSShri Abhyankar   PetscInt     order;
233b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
243b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
25316643e7SJed Brown } TS_Theta;
26316643e7SJed Brown 
27316643e7SJed Brown #undef __FUNCT__
287445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
297445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
307445fe48SJed Brown {
317445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
327445fe48SJed Brown   PetscErrorCode ierr;
337445fe48SJed Brown 
347445fe48SJed Brown   PetscFunctionBegin;
357445fe48SJed Brown   if (X0) {
367445fe48SJed Brown     if (dm && dm != ts->dm) {
370d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
387445fe48SJed Brown     } else *X0 = ts->vec_sol;
397445fe48SJed Brown   }
407445fe48SJed Brown   if (Xdot) {
417445fe48SJed Brown     if (dm && dm != ts->dm) {
420d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
437445fe48SJed Brown     } else *Xdot = th->Xdot;
447445fe48SJed Brown   }
457445fe48SJed Brown   PetscFunctionReturn(0);
467445fe48SJed Brown }
477445fe48SJed Brown 
480d0b770aSPeter Brune 
490d0b770aSPeter Brune #undef __FUNCT__
500d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
510d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
520d0b770aSPeter Brune {
530d0b770aSPeter Brune   PetscErrorCode ierr;
540d0b770aSPeter Brune 
550d0b770aSPeter Brune   PetscFunctionBegin;
560d0b770aSPeter Brune   if (X0) {
570d0b770aSPeter Brune     if (dm && dm != ts->dm) {
580d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
590d0b770aSPeter Brune     }
600d0b770aSPeter Brune   }
610d0b770aSPeter Brune   if (Xdot) {
620d0b770aSPeter Brune     if (dm && dm != ts->dm) {
630d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
640d0b770aSPeter Brune     }
650d0b770aSPeter Brune   }
660d0b770aSPeter Brune   PetscFunctionReturn(0);
670d0b770aSPeter Brune }
680d0b770aSPeter Brune 
697445fe48SJed Brown #undef __FUNCT__
707445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
717445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
727445fe48SJed Brown {
737445fe48SJed Brown 
747445fe48SJed Brown   PetscFunctionBegin;
757445fe48SJed Brown   PetscFunctionReturn(0);
767445fe48SJed Brown }
777445fe48SJed Brown 
787445fe48SJed Brown #undef __FUNCT__
797445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
807445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
817445fe48SJed Brown {
827445fe48SJed Brown   TS             ts = (TS)ctx;
837445fe48SJed Brown   PetscErrorCode ierr;
847445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
857445fe48SJed Brown 
867445fe48SJed Brown   PetscFunctionBegin;
877445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
897445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
907445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
917445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
927445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
930d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
940d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
957445fe48SJed Brown   PetscFunctionReturn(0);
967445fe48SJed Brown }
977445fe48SJed Brown 
987445fe48SJed Brown #undef __FUNCT__
99258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
100258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
101258e1594SPeter Brune {
102258e1594SPeter Brune 
103258e1594SPeter Brune   PetscFunctionBegin;
104258e1594SPeter Brune   PetscFunctionReturn(0);
105258e1594SPeter Brune }
106258e1594SPeter Brune 
107258e1594SPeter Brune #undef __FUNCT__
108258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
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 
1303b1890cdSShri Abhyankar #undef __FUNCT__
1313b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1323b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1333b1890cdSShri Abhyankar {
1343b1890cdSShri Abhyankar   PetscErrorCode ierr;
1353b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1363b1890cdSShri Abhyankar 
1373b1890cdSShri Abhyankar   PetscFunctionBegin;
138ce94432eSBarry Smith   if (order == 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"No time-step adaptivity implemented for 1st order theta method; Run with -ts_adapt_type none");
1393b1890cdSShri Abhyankar   if (order == th->order) {
1403b1890cdSShri Abhyankar     if (th->endpoint) {
1413b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1423b1890cdSShri Abhyankar     } else {
1433b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1443b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1453b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1463b1890cdSShri Abhyankar     }
1473b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1483b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1493b1890cdSShri Abhyankar   }
1503b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1513b1890cdSShri Abhyankar }
152258e1594SPeter Brune 
153258e1594SPeter Brune #undef __FUNCT__
15424655328SShri #define __FUNCT__ "TSRollBack_Theta"
15524655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15624655328SShri {
15724655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
15824655328SShri   PetscErrorCode ierr;
15924655328SShri 
16024655328SShri   PetscFunctionBegin;
16124655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
16224655328SShri   th->status    = TS_STEP_INCOMPLETE;
16324655328SShri   PetscFunctionReturn(0);
16424655328SShri }
16524655328SShri 
16624655328SShri #undef __FUNCT__
167316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
168193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
169316643e7SJed Brown {
170316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1713b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
1723b1890cdSShri Abhyankar   PetscReal      next_time_step;
1733b1890cdSShri Abhyankar   TSAdapt        adapt;
1744957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
175051f2191SLisandro Dalcin   PetscErrorCode ierr;
176316643e7SJed Brown 
177316643e7SJed Brown   PetscFunctionBegin;
1783b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1793b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
180051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
181b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
182eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
183b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
184b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
185316643e7SJed Brown 
186eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
187eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
188eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
189eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
190eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
191eb284becSJed Brown     }
192ace68cafSJed Brown     if (th->extrapolate) {
193b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
194ace68cafSJed Brown     } else {
1952b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
196ace68cafSJed Brown     }
197eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
198316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
199316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
2005ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
201051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
202552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2034957b756SLisandro Dalcin     ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr);
204051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
205051f2191SLisandro Dalcin 
2060298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
207051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2083b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
209552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2103b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2110298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2123b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
213051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
214051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
215051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
216051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
217051f2191SLisandro Dalcin       goto reject_step;
218051f2191SLisandro Dalcin     }
2193b1890cdSShri Abhyankar 
2203b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2212b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
222cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
223316643e7SJed Brown     ts->steps++;
2243b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
225051f2191SLisandro Dalcin     break;
226051f2191SLisandro Dalcin 
227051f2191SLisandro Dalcin reject_step:
228051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
229051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
230051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2313b1890cdSShri Abhyankar     }
232051f2191SLisandro Dalcin     continue;
2333b1890cdSShri Abhyankar   }
234316643e7SJed Brown   PetscFunctionReturn(0);
235316643e7SJed Brown }
236316643e7SJed Brown 
237cd652676SJed Brown #undef __FUNCT__
2382ca6e920SHong Zhang #define __FUNCT__ "TSStepAdj_Theta"
2392ca6e920SHong Zhang static PetscErrorCode TSStepAdj_Theta(TS ts)
2402ca6e920SHong Zhang {
2412ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
2423c54f92cSHong Zhang   Vec                 *VecDeltaLam = th->VecDeltaLam,*VecDeltaMu = th->VecDeltaMu,*VecSensiTemp = th->VecSensiTemp;
2432ca6e920SHong Zhang   PetscInt            nadj;
2442ca6e920SHong Zhang   PetscErrorCode      ierr;
2452ca6e920SHong Zhang   Mat                 J,Jp;
2462ca6e920SHong Zhang   KSP                 ksp;
2472ca6e920SHong Zhang   PetscReal           shift;
2482ca6e920SHong Zhang 
2492ca6e920SHong Zhang   PetscFunctionBegin;
2502ca6e920SHong Zhang 
2512ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
2522ca6e920SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);
2532ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
2543fd52205SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? ts->time_step : (1.-th->Theta)*ts->time_step); /* time_step is negative*/
2552ca6e920SHong Zhang 
2562ca6e920SHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
2572ca6e920SHong Zhang 
258a4cab896SHong Zhang   /* Build RHS */
25905755b9cSHong Zhang   if (ts->vec_costquad) { /* Cost function has an integral (quadrature) term */
26005755b9cSHong Zhang     if (th->endpoint) {
261*d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdy);CHKERRQ(ierr);
26205755b9cSHong Zhang     }else {
263*d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
26405755b9cSHong Zhang     }
26505755b9cSHong Zhang   }
2662ca6e920SHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
2672ca6e920SHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
2682ca6e920SHong Zhang     ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
26936eaed60SHong Zhang     if (ts->vec_costquad) {
27036eaed60SHong Zhang       ierr = VecAXPY(VecSensiTemp[nadj],1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
27136eaed60SHong Zhang     }
2722ca6e920SHong Zhang   }
2733c54f92cSHong Zhang 
2742ca6e920SHong Zhang   /* Build LHS */
2752ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
2763c54f92cSHong Zhang   if (th->endpoint) {
2773c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2783c54f92cSHong Zhang   }else {
2793c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2803c54f92cSHong Zhang   }
2812ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
2822ca6e920SHong Zhang 
2832ca6e920SHong Zhang   /* Solve LHS X = RHS */
2842ca6e920SHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
2852ca6e920SHong Zhang     ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr);
2862ca6e920SHong Zhang   }
2873c54f92cSHong Zhang 
28836eaed60SHong Zhang   /* Update sensitivities, and evaluate integrals if there is any */
2893fd52205SHong Zhang   if(th->endpoint && th->Theta!=1.) { /* two-stage case */
2903fd52205SHong Zhang     shift = -1./((th->Theta-1.)*ts->time_step);
2913c54f92cSHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
29236eaed60SHong Zhang     if (ts->vec_costquad) {
293*d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
29436eaed60SHong Zhang       /* Evolve ts->vec_costquad to compute integrals */
295*d4aa0a58SBarry Smith       ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr);
29636eaed60SHong Zhang       ierr = VecAXPY(ts->vec_costquad,-ts->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr);
297*d4aa0a58SBarry Smith       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
29836eaed60SHong Zhang       ierr = VecAXPY(ts->vec_costquad,ts->time_step*(th->Theta-1.),ts->vec_costintegrand);CHKERRQ(ierr);
29936eaed60SHong Zhang     }
3002ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
3013c54f92cSHong Zhang       ierr = MatMultTranspose(J,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
30236eaed60SHong Zhang       if (ts->vec_costquad) {
30336eaed60SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-1.,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
30436eaed60SHong Zhang       }
3053c54f92cSHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
3062ca6e920SHong Zhang     }
3073fd52205SHong Zhang 
3083fd52205SHong Zhang     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
309*d4aa0a58SBarry Smith       ierr = TSAdjointComputeRHSJacobianP(ts,ts->ptime,ts->vec_sol,ts->Jacp);CHKERRQ(ierr);
3103fd52205SHong Zhang       for (nadj=0; nadj<ts->numberadjs; nadj++) {
3113fd52205SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
3123fd52205SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,VecDeltaMu[nadj]);CHKERRQ(ierr);
3133fd52205SHong Zhang       }
314*d4aa0a58SBarry Smith       ierr = TSAdjointComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
3153fd52205SHong Zhang       for (nadj=0; nadj<ts->numberadjs; nadj++) {
3163fd52205SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
3173fd52205SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),VecDeltaMu[nadj]);CHKERRQ(ierr);
3183fd52205SHong Zhang       }
31936eaed60SHong Zhang       if (ts->vec_costquad) {
320*d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,ts->ptime,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr);
32136eaed60SHong Zhang         for (nadj=0; nadj<ts->numberadjs; nadj++) {
32236eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
32336eaed60SHong Zhang         }
324*d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
32536eaed60SHong Zhang         for (nadj=0; nadj<ts->numberadjs; nadj++) {
32636eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr);
32736eaed60SHong Zhang         }
32836eaed60SHong Zhang       }
3293fd52205SHong Zhang     }
3303fd52205SHong Zhang   }else { /* one-stage case */
3313c54f92cSHong Zhang     shift = 0.0;
332a2ae3dd2SHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */
33336eaed60SHong Zhang     if (ts->vec_costquad) {
334*d4aa0a58SBarry Smith       ierr = TSAdjointComputeDRDYFunction(ts,th->stage_time,th->X,ts->vecs_drdy);CHKERRQ(ierr);
33536eaed60SHong Zhang       /* Evolve ts->vec_costquad to compute integrals */
336*d4aa0a58SBarry Smith       ierr = TSAdjointComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr);
33736eaed60SHong Zhang       ierr = VecAXPY(ts->vec_costquad,-ts->time_step,ts->vec_costintegrand);CHKERRQ(ierr);
33836eaed60SHong Zhang     }
3393fd52205SHong Zhang      /* When th->endpoint is true and th->Theta==1 (beuler method), the Jacobian is supposed to be evaluated at ts->ptime like this:
3403c54f92cSHong Zhang     if(th->endpoint) {
3413c54f92cSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
3423c54f92cSHong Zhang     }
3433fd52205SHong Zhang     but ts->ptime and ts->vec_sol have the same values as th->stage_time and th->X in this case. So the code is simplified here.
3443fd52205SHong Zhang     */
3452ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
3463c54f92cSHong Zhang       ierr = MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
3473c54f92cSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr);
34836eaed60SHong Zhang       if (ts->vec_costquad) {
34936eaed60SHong Zhang         ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,ts->vecs_drdy[nadj]);CHKERRQ(ierr);
35036eaed60SHong Zhang       }
3512ca6e920SHong Zhang     }
3523fd52205SHong Zhang     if (ts->vecs_sensip) {
353*d4aa0a58SBarry Smith       ierr = TSAdjointComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr);
3543fd52205SHong Zhang       for (nadj=0; nadj<ts->numberadjs; nadj++) {
3553fd52205SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecDeltaLam[nadj],VecDeltaMu[nadj]);CHKERRQ(ierr);
3563fd52205SHong Zhang         ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,VecDeltaMu[nadj]);CHKERRQ(ierr);
3573fd52205SHong Zhang       }
35836eaed60SHong Zhang       if (ts->vec_costquad) {
359*d4aa0a58SBarry Smith         ierr = TSAdjointComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr);
36036eaed60SHong Zhang         for (nadj=0; nadj<ts->numberadjs; nadj++) {
36136eaed60SHong Zhang           ierr = VecAXPY(ts->vecs_sensip[nadj],-ts->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr);
36236eaed60SHong Zhang         }
36336eaed60SHong Zhang       }
3643fd52205SHong Zhang     }
3652ca6e920SHong Zhang   }
3662ca6e920SHong Zhang 
3672ca6e920SHong Zhang   ts->ptime += ts->time_step;
3682ca6e920SHong Zhang   ts->steps++;
3692ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
3702ca6e920SHong Zhang   PetscFunctionReturn(0);
3712ca6e920SHong Zhang }
3722ca6e920SHong Zhang 
3732ca6e920SHong Zhang #undef __FUNCT__
374cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
375cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
376cd652676SJed Brown {
377cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
3785a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
379cd652676SJed Brown   PetscErrorCode ierr;
380cd652676SJed Brown 
381cd652676SJed Brown   PetscFunctionBegin;
382a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
3835a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
3845a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
385cd652676SJed Brown   PetscFunctionReturn(0);
386cd652676SJed Brown }
387cd652676SJed Brown 
388316643e7SJed Brown /*------------------------------------------------------------*/
389316643e7SJed Brown #undef __FUNCT__
390277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
391277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
392316643e7SJed Brown {
393316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
394316643e7SJed Brown   PetscErrorCode ierr;
395316643e7SJed Brown 
396316643e7SJed Brown   PetscFunctionBegin;
3976bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
3986bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
3993b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
400eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
4012ca6e920SHong Zhang   if(ts->reverse_mode) {
4022ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
4032ca6e920SHong Zhang     if(th->VecDeltaMu) {
4042ca6e920SHong Zhang       ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
4052ca6e920SHong Zhang     }
4062ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
4072ca6e920SHong Zhang   }
408277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
409277b19d0SLisandro Dalcin }
410277b19d0SLisandro Dalcin 
411277b19d0SLisandro Dalcin #undef __FUNCT__
412277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
413277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
414277b19d0SLisandro Dalcin {
415277b19d0SLisandro Dalcin   PetscErrorCode ierr;
416277b19d0SLisandro Dalcin 
417277b19d0SLisandro Dalcin   PetscFunctionBegin;
418277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
419277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
420bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
421bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
422bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
423bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
424316643e7SJed Brown   PetscFunctionReturn(0);
425316643e7SJed Brown }
426316643e7SJed Brown 
427316643e7SJed Brown /*
428316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
4292b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
430316643e7SJed Brown */
431316643e7SJed Brown #undef __FUNCT__
4320f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
4330f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
434316643e7SJed Brown {
435316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
436316643e7SJed Brown   PetscErrorCode ierr;
4377445fe48SJed Brown   Vec            X0,Xdot;
4387445fe48SJed Brown   DM             dm,dmsave;
439b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
440316643e7SJed Brown 
441316643e7SJed Brown   PetscFunctionBegin;
4427445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
4435a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
4447445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
445b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
4467445fe48SJed Brown 
4477445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
4487445fe48SJed Brown   dmsave = ts->dm;
4497445fe48SJed Brown   ts->dm = dm;
4507445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
4517445fe48SJed Brown   ts->dm = dmsave;
4520d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
453316643e7SJed Brown   PetscFunctionReturn(0);
454316643e7SJed Brown }
455316643e7SJed Brown 
456316643e7SJed Brown #undef __FUNCT__
4570f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
458d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
459316643e7SJed Brown {
460316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
461316643e7SJed Brown   PetscErrorCode ierr;
4627445fe48SJed Brown   Vec            Xdot;
4637445fe48SJed Brown   DM             dm,dmsave;
464b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
465316643e7SJed Brown 
466316643e7SJed Brown   PetscFunctionBegin;
4677445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
4687445fe48SJed Brown 
4690f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
4700298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
4717445fe48SJed Brown 
4727445fe48SJed Brown   dmsave = ts->dm;
4737445fe48SJed Brown   ts->dm = dm;
474d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
4757445fe48SJed Brown   ts->dm = dmsave;
4760298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
477316643e7SJed Brown   PetscFunctionReturn(0);
478316643e7SJed Brown }
479316643e7SJed Brown 
480316643e7SJed Brown #undef __FUNCT__
481316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
482316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
483316643e7SJed Brown {
484316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
485316643e7SJed Brown   PetscErrorCode ierr;
4867445fe48SJed Brown   SNES           snes;
487ef749922SLisandro Dalcin   TSAdapt        adapt;
4887445fe48SJed Brown   DM             dm;
489316643e7SJed Brown 
490316643e7SJed Brown   PetscFunctionBegin;
49139d13666SHong Zhang   if (!th->X) {
492316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
49339d13666SHong Zhang   }
49439d13666SHong Zhang   if (!th->Xdot) {
495316643e7SJed Brown     ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
49639d13666SHong Zhang   }
49739d13666SHong Zhang   if (!th->X0) {
4983b1890cdSShri Abhyankar     ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
49939d13666SHong Zhang   }
5007445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
5017445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
5027445fe48SJed Brown   if (dm) {
5037445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
504258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
5057445fe48SJed Brown   }
5063b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
5073b1890cdSShri Abhyankar   else th->order = 1;
5083b1890cdSShri Abhyankar 
509552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
510ef749922SLisandro Dalcin   if (!th->adapt) {
5113b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
5123b1890cdSShri Abhyankar   }
5132ca6e920SHong Zhang   if (ts->reverse_mode) {
5142ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
5152ca6e920SHong Zhang     if(ts->vecs_sensip) {
5162ca6e920SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
5172ca6e920SHong Zhang     }
5182ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
5192ca6e920SHong Zhang   }
520316643e7SJed Brown   PetscFunctionReturn(0);
521316643e7SJed Brown }
522316643e7SJed Brown /*------------------------------------------------------------*/
523316643e7SJed Brown 
524316643e7SJed Brown #undef __FUNCT__
525316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
5268c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
527316643e7SJed Brown {
528316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
529316643e7SJed Brown   PetscErrorCode ierr;
530316643e7SJed Brown 
531316643e7SJed Brown   PetscFunctionBegin;
532e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
533316643e7SJed Brown   {
5340298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
5350298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
5360298fd71SBarry 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);
5370298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
538d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
539316643e7SJed Brown   }
540316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
541316643e7SJed Brown   PetscFunctionReturn(0);
542316643e7SJed Brown }
543316643e7SJed Brown 
544316643e7SJed Brown #undef __FUNCT__
545316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
546316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
547316643e7SJed Brown {
548316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
549ace3abfcSBarry Smith   PetscBool      iascii;
550316643e7SJed Brown   PetscErrorCode ierr;
551316643e7SJed Brown 
552316643e7SJed Brown   PetscFunctionBegin;
553251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
554316643e7SJed Brown   if (iascii) {
5557c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
556ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
557316643e7SJed Brown   }
558ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
559316643e7SJed Brown   PetscFunctionReturn(0);
560316643e7SJed Brown }
561316643e7SJed Brown 
5620de4c49aSJed Brown #undef __FUNCT__
5630de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
5647087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
5650de4c49aSJed Brown {
5660de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
5670de4c49aSJed Brown 
5680de4c49aSJed Brown   PetscFunctionBegin;
5690de4c49aSJed Brown   *theta = th->Theta;
5700de4c49aSJed Brown   PetscFunctionReturn(0);
5710de4c49aSJed Brown }
5720de4c49aSJed Brown 
5730de4c49aSJed Brown #undef __FUNCT__
5740de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
5757087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
5760de4c49aSJed Brown {
5770de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
5780de4c49aSJed Brown 
5790de4c49aSJed Brown   PetscFunctionBegin;
5807c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
5810de4c49aSJed Brown   th->Theta = theta;
5820de4c49aSJed Brown   PetscFunctionReturn(0);
5830de4c49aSJed Brown }
584eb284becSJed Brown 
585eb284becSJed Brown #undef __FUNCT__
58678e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
58726f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
58826f2ff8fSLisandro Dalcin {
58926f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
59026f2ff8fSLisandro Dalcin 
59126f2ff8fSLisandro Dalcin   PetscFunctionBegin;
59226f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
59326f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
59426f2ff8fSLisandro Dalcin }
59526f2ff8fSLisandro Dalcin 
59626f2ff8fSLisandro Dalcin #undef __FUNCT__
59726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
598eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
599eb284becSJed Brown {
600eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
601eb284becSJed Brown 
602eb284becSJed Brown   PetscFunctionBegin;
603eb284becSJed Brown   th->endpoint = flg;
604eb284becSJed Brown   PetscFunctionReturn(0);
605eb284becSJed Brown }
6060de4c49aSJed Brown 
607f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
608f9c1d6abSBarry Smith #undef __FUNCT__
609f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
610f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
611f9c1d6abSBarry Smith {
612f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
613f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
6143fd8ae06SJed Brown   const PetscReal one = 1.0;
615f9c1d6abSBarry Smith 
616f9c1d6abSBarry Smith   PetscFunctionBegin;
6173fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
618f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
619f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
620f9c1d6abSBarry Smith   PetscFunctionReturn(0);
621f9c1d6abSBarry Smith }
622f9c1d6abSBarry Smith #endif
623f9c1d6abSBarry Smith 
62442682096SHong Zhang #undef __FUNCT__
62542682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
62642682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
62742682096SHong Zhang {
62842682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
62942682096SHong Zhang 
63042682096SHong Zhang   PetscFunctionBegin;
6312ca6e920SHong Zhang   *ns = 1;
6322ca6e920SHong Zhang   if(Y) {
6332ca6e920SHong Zhang     *Y  = &(th->X);
6342ca6e920SHong Zhang   }
63542682096SHong Zhang   PetscFunctionReturn(0);
63642682096SHong Zhang }
637f9c1d6abSBarry Smith 
638316643e7SJed Brown /* ------------------------------------------------------------ */
639316643e7SJed Brown /*MC
64096f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
641316643e7SJed Brown 
642316643e7SJed Brown    Level: beginner
643316643e7SJed Brown 
6444eb428fdSBarry Smith    Options Database:
6450c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
6464eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
6470c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
6484eb428fdSBarry Smith 
649eb284becSJed Brown    Notes:
6500c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
6510c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
6524eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
6534eb428fdSBarry Smith 
6544eb428fdSBarry Smith 
6554eb428fdSBarry Smith 
656eb284becSJed Brown    This method can be applied to DAE.
657eb284becSJed Brown 
658eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
659eb284becSJed Brown 
660eb284becSJed Brown .vb
661eb284becSJed Brown   Theta | Theta
662eb284becSJed Brown   -------------
663eb284becSJed Brown         |  1
664eb284becSJed Brown .ve
665eb284becSJed Brown 
666eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
667eb284becSJed Brown 
668eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
669eb284becSJed Brown 
670eb284becSJed Brown .vb
671eb284becSJed Brown   0 | 0         0
672eb284becSJed Brown   1 | 1-Theta   Theta
673eb284becSJed Brown   -------------------
674eb284becSJed Brown     | 1-Theta   Theta
675eb284becSJed Brown .ve
676eb284becSJed Brown 
677eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
678eb284becSJed Brown 
679eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
680eb284becSJed Brown 
681eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
682eb284becSJed Brown 
6834eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
684eb284becSJed Brown 
685eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
686316643e7SJed Brown 
687316643e7SJed Brown M*/
688316643e7SJed Brown #undef __FUNCT__
689316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
6908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
691316643e7SJed Brown {
692316643e7SJed Brown   TS_Theta       *th;
693316643e7SJed Brown   PetscErrorCode ierr;
694316643e7SJed Brown 
695316643e7SJed Brown   PetscFunctionBegin;
696277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
697316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
698316643e7SJed Brown   ts->ops->view           = TSView_Theta;
699316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
700316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
701cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
7023b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
70324655328SShri   ts->ops->rollback       = TSRollBack_Theta;
704316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
7050f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
7060f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
707f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
708f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
709f9c1d6abSBarry Smith #endif
71042682096SHong Zhang   ts->ops->getstages      = TSGetStages_Theta;
7112ca6e920SHong Zhang   ts->ops->stepadj        = TSStepAdj_Theta;
712316643e7SJed Brown 
713b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
714316643e7SJed Brown   ts->data = (void*)th;
715316643e7SJed Brown 
7166f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
717316643e7SJed Brown   th->Theta       = 0.5;
7183b1890cdSShri Abhyankar   th->ccfl        = 1.0;
7193b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
720bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
721bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
722bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
723bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
724316643e7SJed Brown   PetscFunctionReturn(0);
725316643e7SJed Brown }
7260de4c49aSJed Brown 
7270de4c49aSJed Brown #undef __FUNCT__
7280de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
7290de4c49aSJed Brown /*@
7300de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
7310de4c49aSJed Brown 
7320de4c49aSJed Brown   Not Collective
7330de4c49aSJed Brown 
7340de4c49aSJed Brown   Input Parameter:
7350de4c49aSJed Brown .  ts - timestepping context
7360de4c49aSJed Brown 
7370de4c49aSJed Brown   Output Parameter:
7380de4c49aSJed Brown .  theta - stage abscissa
7390de4c49aSJed Brown 
7400de4c49aSJed Brown   Note:
7410de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
7420de4c49aSJed Brown 
7430de4c49aSJed Brown   Level: Advanced
7440de4c49aSJed Brown 
7450de4c49aSJed Brown .seealso: TSThetaSetTheta()
7460de4c49aSJed Brown @*/
7477087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
7480de4c49aSJed Brown {
7494ac538c5SBarry Smith   PetscErrorCode ierr;
7500de4c49aSJed Brown 
7510de4c49aSJed Brown   PetscFunctionBegin;
752afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7530de4c49aSJed Brown   PetscValidPointer(theta,2);
7544ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
7550de4c49aSJed Brown   PetscFunctionReturn(0);
7560de4c49aSJed Brown }
7570de4c49aSJed Brown 
7580de4c49aSJed Brown #undef __FUNCT__
7590de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
7600de4c49aSJed Brown /*@
7610de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
7620de4c49aSJed Brown 
7630de4c49aSJed Brown   Not Collective
7640de4c49aSJed Brown 
7650de4c49aSJed Brown   Input Parameter:
7660de4c49aSJed Brown +  ts - timestepping context
7670de4c49aSJed Brown -  theta - stage abscissa
7680de4c49aSJed Brown 
7690de4c49aSJed Brown   Options Database:
7700de4c49aSJed Brown .  -ts_theta_theta <theta>
7710de4c49aSJed Brown 
7720de4c49aSJed Brown   Level: Intermediate
7730de4c49aSJed Brown 
7740de4c49aSJed Brown .seealso: TSThetaGetTheta()
7750de4c49aSJed Brown @*/
7767087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
7770de4c49aSJed Brown {
7784ac538c5SBarry Smith   PetscErrorCode ierr;
7790de4c49aSJed Brown 
7800de4c49aSJed Brown   PetscFunctionBegin;
781afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7824ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
7830de4c49aSJed Brown   PetscFunctionReturn(0);
7840de4c49aSJed Brown }
785f33bbcb6SJed Brown 
786eb284becSJed Brown #undef __FUNCT__
78726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
78826f2ff8fSLisandro Dalcin /*@
78926f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
79026f2ff8fSLisandro Dalcin 
79126f2ff8fSLisandro Dalcin   Not Collective
79226f2ff8fSLisandro Dalcin 
79326f2ff8fSLisandro Dalcin   Input Parameter:
79426f2ff8fSLisandro Dalcin .  ts - timestepping context
79526f2ff8fSLisandro Dalcin 
79626f2ff8fSLisandro Dalcin   Output Parameter:
79726f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
79826f2ff8fSLisandro Dalcin 
79926f2ff8fSLisandro Dalcin   Level: Advanced
80026f2ff8fSLisandro Dalcin 
80126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
80226f2ff8fSLisandro Dalcin @*/
80326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
80426f2ff8fSLisandro Dalcin {
80526f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
80626f2ff8fSLisandro Dalcin 
80726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
80826f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
80926f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
81026f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
81126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
81226f2ff8fSLisandro Dalcin }
81326f2ff8fSLisandro Dalcin 
81426f2ff8fSLisandro Dalcin #undef __FUNCT__
815eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
816eb284becSJed Brown /*@
817eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
818eb284becSJed Brown 
819eb284becSJed Brown   Not Collective
820eb284becSJed Brown 
821eb284becSJed Brown   Input Parameter:
822eb284becSJed Brown +  ts - timestepping context
823eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
824eb284becSJed Brown 
825eb284becSJed Brown   Options Database:
826eb284becSJed Brown .  -ts_theta_endpoint <flg>
827eb284becSJed Brown 
828eb284becSJed Brown   Level: Intermediate
829eb284becSJed Brown 
830eb284becSJed Brown .seealso: TSTHETA, TSCN
831eb284becSJed Brown @*/
832eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
833eb284becSJed Brown {
834eb284becSJed Brown   PetscErrorCode ierr;
835eb284becSJed Brown 
836eb284becSJed Brown   PetscFunctionBegin;
837eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
838eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
839eb284becSJed Brown   PetscFunctionReturn(0);
840eb284becSJed Brown }
841eb284becSJed Brown 
842f33bbcb6SJed Brown /*
843f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
844f33bbcb6SJed Brown  * The creation functions for these specializations are below.
845f33bbcb6SJed Brown  */
846f33bbcb6SJed Brown 
847f33bbcb6SJed Brown #undef __FUNCT__
848f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
849f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
850f33bbcb6SJed Brown {
851d52bd9f3SBarry Smith   PetscErrorCode ierr;
852d52bd9f3SBarry Smith 
853f33bbcb6SJed Brown   PetscFunctionBegin;
854d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
855f33bbcb6SJed Brown   PetscFunctionReturn(0);
856f33bbcb6SJed Brown }
857f33bbcb6SJed Brown 
858f33bbcb6SJed Brown /*MC
859f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
860f33bbcb6SJed Brown 
861f33bbcb6SJed Brown   Level: beginner
862f33bbcb6SJed Brown 
8634eb428fdSBarry Smith   Notes:
864c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
8654eb428fdSBarry Smith 
8664eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
8674eb428fdSBarry Smith 
868f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
869f33bbcb6SJed Brown 
870f33bbcb6SJed Brown M*/
871f33bbcb6SJed Brown #undef __FUNCT__
872f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
8738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
874f33bbcb6SJed Brown {
875f33bbcb6SJed Brown   PetscErrorCode ierr;
876f33bbcb6SJed Brown 
877f33bbcb6SJed Brown   PetscFunctionBegin;
878f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
879f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
880f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
881f33bbcb6SJed Brown   PetscFunctionReturn(0);
882f33bbcb6SJed Brown }
883f33bbcb6SJed Brown 
884f33bbcb6SJed Brown #undef __FUNCT__
885f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
886f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
887f33bbcb6SJed Brown {
888d52bd9f3SBarry Smith   PetscErrorCode ierr;
889d52bd9f3SBarry Smith 
890f33bbcb6SJed Brown   PetscFunctionBegin;
891d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
892f33bbcb6SJed Brown   PetscFunctionReturn(0);
893f33bbcb6SJed Brown }
894f33bbcb6SJed Brown 
895f33bbcb6SJed Brown /*MC
896f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
897f33bbcb6SJed Brown 
898f33bbcb6SJed Brown   Level: beginner
899f33bbcb6SJed Brown 
900f33bbcb6SJed Brown   Notes:
9017cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
9027cf5af47SJed Brown 
9037cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
904f33bbcb6SJed Brown 
905f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
906f33bbcb6SJed Brown 
907f33bbcb6SJed Brown M*/
908f33bbcb6SJed Brown #undef __FUNCT__
909f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
9108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
911f33bbcb6SJed Brown {
912f33bbcb6SJed Brown   PetscErrorCode ierr;
913f33bbcb6SJed Brown 
914f33bbcb6SJed Brown   PetscFunctionBegin;
915f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
916f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
917eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
918f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
919f33bbcb6SJed Brown   PetscFunctionReturn(0);
920f33bbcb6SJed Brown }
921