xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 3c54f92c5fb7121d6e9d66350d795c36b24d7f9f)
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>
7*3c54f92cSHong 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;
242*3c54f92cSHong 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);
254a4cab896SHong 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 */
2592ca6e920SHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
2602ca6e920SHong Zhang     ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
2612ca6e920SHong Zhang     ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
2622ca6e920SHong Zhang   }
263*3c54f92cSHong Zhang 
2642ca6e920SHong Zhang   /* Build LHS */
2652ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
266*3c54f92cSHong Zhang   if (th->endpoint) {
267*3c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
268*3c54f92cSHong Zhang   }else {
269*3c54f92cSHong Zhang     ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
270*3c54f92cSHong Zhang   }
2712ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
2722ca6e920SHong Zhang 
2732ca6e920SHong Zhang   /* Solve LHS X = RHS */
2742ca6e920SHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
2752ca6e920SHong Zhang     ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr);
2762ca6e920SHong Zhang   }
277*3c54f92cSHong Zhang 
278*3c54f92cSHong Zhang   /* Update sensitivities */
279*3c54f92cSHong Zhang   if(th->endpoint && th->Theta!=1.0) {
280*3c54f92cSHong Zhang     shift = -1./((th->Theta-1.0)*ts->time_step);
281*3c54f92cSHong Zhang     ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
2822ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
283*3c54f92cSHong Zhang       ierr = MatMultTranspose(J,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
284*3c54f92cSHong Zhang       ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr);
2852ca6e920SHong Zhang     }
2862ca6e920SHong Zhang   }else {
287*3c54f92cSHong Zhang     shift = 0.0;
288*3c54f92cSHong Zhang     if(th->endpoint) {
289*3c54f92cSHong Zhang       ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
290*3c54f92cSHong Zhang     }else {
291*3c54f92cSHong Zhang       ierr  = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
292*3c54f92cSHong Zhang     }
2932ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
294*3c54f92cSHong Zhang       ierr = MatMultTranspose(J,VecDeltaLam[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
295*3c54f92cSHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr);
2962ca6e920SHong Zhang     }
2972ca6e920SHong Zhang   }
2982ca6e920SHong Zhang 
2992ca6e920SHong Zhang   ts->ptime += ts->time_step;
3002ca6e920SHong Zhang   ts->steps++;
3012ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
3022ca6e920SHong Zhang   PetscFunctionReturn(0);
3032ca6e920SHong Zhang }
3042ca6e920SHong Zhang 
3052ca6e920SHong Zhang #undef __FUNCT__
306cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
307cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
308cd652676SJed Brown {
309cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
3105a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
311cd652676SJed Brown   PetscErrorCode ierr;
312cd652676SJed Brown 
313cd652676SJed Brown   PetscFunctionBegin;
314a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
3155a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
3165a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
317cd652676SJed Brown   PetscFunctionReturn(0);
318cd652676SJed Brown }
319cd652676SJed Brown 
320316643e7SJed Brown /*------------------------------------------------------------*/
321316643e7SJed Brown #undef __FUNCT__
322277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
323277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
324316643e7SJed Brown {
325316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
326316643e7SJed Brown   PetscErrorCode ierr;
327316643e7SJed Brown 
328316643e7SJed Brown   PetscFunctionBegin;
3296bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
3306bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
3313b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
332eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
3332ca6e920SHong Zhang   if(ts->reverse_mode) {
3342ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
3352ca6e920SHong Zhang     if(th->VecDeltaMu) {
3362ca6e920SHong Zhang       ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
3372ca6e920SHong Zhang     }
3382ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
3392ca6e920SHong Zhang   }
340277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
341277b19d0SLisandro Dalcin }
342277b19d0SLisandro Dalcin 
343277b19d0SLisandro Dalcin #undef __FUNCT__
344277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
345277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
346277b19d0SLisandro Dalcin {
347277b19d0SLisandro Dalcin   PetscErrorCode ierr;
348277b19d0SLisandro Dalcin 
349277b19d0SLisandro Dalcin   PetscFunctionBegin;
350277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
351277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
352bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
353bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
354bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
355bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
356316643e7SJed Brown   PetscFunctionReturn(0);
357316643e7SJed Brown }
358316643e7SJed Brown 
359316643e7SJed Brown /*
360316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
3612b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
362316643e7SJed Brown */
363316643e7SJed Brown #undef __FUNCT__
3640f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
3650f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
366316643e7SJed Brown {
367316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
368316643e7SJed Brown   PetscErrorCode ierr;
3697445fe48SJed Brown   Vec            X0,Xdot;
3707445fe48SJed Brown   DM             dm,dmsave;
371b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
372316643e7SJed Brown 
373316643e7SJed Brown   PetscFunctionBegin;
3747445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3755a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
3767445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
377b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
3787445fe48SJed Brown 
3797445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
3807445fe48SJed Brown   dmsave = ts->dm;
3817445fe48SJed Brown   ts->dm = dm;
3827445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
3837445fe48SJed Brown   ts->dm = dmsave;
3840d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
385316643e7SJed Brown   PetscFunctionReturn(0);
386316643e7SJed Brown }
387316643e7SJed Brown 
388316643e7SJed Brown #undef __FUNCT__
3890f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
390d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
391316643e7SJed Brown {
392316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
393316643e7SJed Brown   PetscErrorCode ierr;
3947445fe48SJed Brown   Vec            Xdot;
3957445fe48SJed Brown   DM             dm,dmsave;
396b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
397316643e7SJed Brown 
398316643e7SJed Brown   PetscFunctionBegin;
3997445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
4007445fe48SJed Brown 
4010f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
4020298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
4037445fe48SJed Brown 
4047445fe48SJed Brown   dmsave = ts->dm;
4057445fe48SJed Brown   ts->dm = dm;
406d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
4077445fe48SJed Brown   ts->dm = dmsave;
4080298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
409316643e7SJed Brown   PetscFunctionReturn(0);
410316643e7SJed Brown }
411316643e7SJed Brown 
412316643e7SJed Brown #undef __FUNCT__
413316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
414316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
415316643e7SJed Brown {
416316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
417316643e7SJed Brown   PetscErrorCode ierr;
4187445fe48SJed Brown   SNES           snes;
419ef749922SLisandro Dalcin   TSAdapt        adapt;
4207445fe48SJed Brown   DM             dm;
421316643e7SJed Brown 
422316643e7SJed Brown   PetscFunctionBegin;
423316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
424316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
4253b1890cdSShri Abhyankar   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
4267445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4277445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
4287445fe48SJed Brown   if (dm) {
4297445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
430258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
4317445fe48SJed Brown   }
4323b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
4333b1890cdSShri Abhyankar   else th->order = 1;
4343b1890cdSShri Abhyankar 
435552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
436ef749922SLisandro Dalcin   if (!th->adapt) {
4373b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
4383b1890cdSShri Abhyankar   }
4392ca6e920SHong Zhang   if (ts->reverse_mode) {
4402ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
4412ca6e920SHong Zhang     if(ts->vecs_sensip) {
4422ca6e920SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
4432ca6e920SHong Zhang     }
4442ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
4452ca6e920SHong Zhang   }
446316643e7SJed Brown   PetscFunctionReturn(0);
447316643e7SJed Brown }
448316643e7SJed Brown /*------------------------------------------------------------*/
449316643e7SJed Brown 
450316643e7SJed Brown #undef __FUNCT__
451316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
4528c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
453316643e7SJed Brown {
454316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
455316643e7SJed Brown   PetscErrorCode ierr;
456316643e7SJed Brown 
457316643e7SJed Brown   PetscFunctionBegin;
458e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
459316643e7SJed Brown   {
4600298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
4610298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
4620298fd71SBarry 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);
4630298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
464d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
465316643e7SJed Brown   }
466316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
467316643e7SJed Brown   PetscFunctionReturn(0);
468316643e7SJed Brown }
469316643e7SJed Brown 
470316643e7SJed Brown #undef __FUNCT__
471316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
472316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
473316643e7SJed Brown {
474316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
475ace3abfcSBarry Smith   PetscBool      iascii;
476316643e7SJed Brown   PetscErrorCode ierr;
477316643e7SJed Brown 
478316643e7SJed Brown   PetscFunctionBegin;
479251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
480316643e7SJed Brown   if (iascii) {
4817c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
482ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
483316643e7SJed Brown   }
484ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
485316643e7SJed Brown   PetscFunctionReturn(0);
486316643e7SJed Brown }
487316643e7SJed Brown 
4880de4c49aSJed Brown #undef __FUNCT__
4890de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
4907087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
4910de4c49aSJed Brown {
4920de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4930de4c49aSJed Brown 
4940de4c49aSJed Brown   PetscFunctionBegin;
4950de4c49aSJed Brown   *theta = th->Theta;
4960de4c49aSJed Brown   PetscFunctionReturn(0);
4970de4c49aSJed Brown }
4980de4c49aSJed Brown 
4990de4c49aSJed Brown #undef __FUNCT__
5000de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
5017087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
5020de4c49aSJed Brown {
5030de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
5040de4c49aSJed Brown 
5050de4c49aSJed Brown   PetscFunctionBegin;
5067c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
5070de4c49aSJed Brown   th->Theta = theta;
5080de4c49aSJed Brown   PetscFunctionReturn(0);
5090de4c49aSJed Brown }
510eb284becSJed Brown 
511eb284becSJed Brown #undef __FUNCT__
51278e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
51326f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
51426f2ff8fSLisandro Dalcin {
51526f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
51626f2ff8fSLisandro Dalcin 
51726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
51826f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
51926f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
52026f2ff8fSLisandro Dalcin }
52126f2ff8fSLisandro Dalcin 
52226f2ff8fSLisandro Dalcin #undef __FUNCT__
52326f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
524eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
525eb284becSJed Brown {
526eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
527eb284becSJed Brown 
528eb284becSJed Brown   PetscFunctionBegin;
529eb284becSJed Brown   th->endpoint = flg;
530eb284becSJed Brown   PetscFunctionReturn(0);
531eb284becSJed Brown }
5320de4c49aSJed Brown 
533f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
534f9c1d6abSBarry Smith #undef __FUNCT__
535f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
536f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
537f9c1d6abSBarry Smith {
538f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
539f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
5403fd8ae06SJed Brown   const PetscReal one = 1.0;
541f9c1d6abSBarry Smith 
542f9c1d6abSBarry Smith   PetscFunctionBegin;
5433fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
544f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
545f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
546f9c1d6abSBarry Smith   PetscFunctionReturn(0);
547f9c1d6abSBarry Smith }
548f9c1d6abSBarry Smith #endif
549f9c1d6abSBarry Smith 
55042682096SHong Zhang #undef __FUNCT__
55142682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
55242682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
55342682096SHong Zhang {
55442682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
55542682096SHong Zhang 
55642682096SHong Zhang   PetscFunctionBegin;
5572ca6e920SHong Zhang   *ns = 1;
5582ca6e920SHong Zhang   if(Y) {
5592ca6e920SHong Zhang     *Y  = &(th->X);
5602ca6e920SHong Zhang   }
56142682096SHong Zhang   PetscFunctionReturn(0);
56242682096SHong Zhang }
563f9c1d6abSBarry Smith 
564316643e7SJed Brown /* ------------------------------------------------------------ */
565316643e7SJed Brown /*MC
56696f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
567316643e7SJed Brown 
568316643e7SJed Brown    Level: beginner
569316643e7SJed Brown 
5704eb428fdSBarry Smith    Options Database:
5710c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
5724eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
5730c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
5744eb428fdSBarry Smith 
575eb284becSJed Brown    Notes:
5760c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
5770c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
5784eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
5794eb428fdSBarry Smith 
5804eb428fdSBarry Smith 
5814eb428fdSBarry Smith 
582eb284becSJed Brown    This method can be applied to DAE.
583eb284becSJed Brown 
584eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
585eb284becSJed Brown 
586eb284becSJed Brown .vb
587eb284becSJed Brown   Theta | Theta
588eb284becSJed Brown   -------------
589eb284becSJed Brown         |  1
590eb284becSJed Brown .ve
591eb284becSJed Brown 
592eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
593eb284becSJed Brown 
594eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
595eb284becSJed Brown 
596eb284becSJed Brown .vb
597eb284becSJed Brown   0 | 0         0
598eb284becSJed Brown   1 | 1-Theta   Theta
599eb284becSJed Brown   -------------------
600eb284becSJed Brown     | 1-Theta   Theta
601eb284becSJed Brown .ve
602eb284becSJed Brown 
603eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
604eb284becSJed Brown 
605eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
606eb284becSJed Brown 
607eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
608eb284becSJed Brown 
6094eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
610eb284becSJed Brown 
611eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
612316643e7SJed Brown 
613316643e7SJed Brown M*/
614316643e7SJed Brown #undef __FUNCT__
615316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
6168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
617316643e7SJed Brown {
618316643e7SJed Brown   TS_Theta       *th;
619316643e7SJed Brown   PetscErrorCode ierr;
620316643e7SJed Brown 
621316643e7SJed Brown   PetscFunctionBegin;
622277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
623316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
624316643e7SJed Brown   ts->ops->view           = TSView_Theta;
625316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
626316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
627cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
6283b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
62924655328SShri   ts->ops->rollback       = TSRollBack_Theta;
630316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
6310f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
6320f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
633f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
634f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
635f9c1d6abSBarry Smith #endif
63642682096SHong Zhang   ts->ops->getstages      = TSGetStages_Theta;
6372ca6e920SHong Zhang   ts->ops->stepadj        = TSStepAdj_Theta;
638316643e7SJed Brown 
639b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
640316643e7SJed Brown   ts->data = (void*)th;
641316643e7SJed Brown 
6426f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
643316643e7SJed Brown   th->Theta       = 0.5;
6443b1890cdSShri Abhyankar   th->ccfl        = 1.0;
6453b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
646bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
647bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
648bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
649bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
650316643e7SJed Brown   PetscFunctionReturn(0);
651316643e7SJed Brown }
6520de4c49aSJed Brown 
6530de4c49aSJed Brown #undef __FUNCT__
6540de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
6550de4c49aSJed Brown /*@
6560de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
6570de4c49aSJed Brown 
6580de4c49aSJed Brown   Not Collective
6590de4c49aSJed Brown 
6600de4c49aSJed Brown   Input Parameter:
6610de4c49aSJed Brown .  ts - timestepping context
6620de4c49aSJed Brown 
6630de4c49aSJed Brown   Output Parameter:
6640de4c49aSJed Brown .  theta - stage abscissa
6650de4c49aSJed Brown 
6660de4c49aSJed Brown   Note:
6670de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
6680de4c49aSJed Brown 
6690de4c49aSJed Brown   Level: Advanced
6700de4c49aSJed Brown 
6710de4c49aSJed Brown .seealso: TSThetaSetTheta()
6720de4c49aSJed Brown @*/
6737087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
6740de4c49aSJed Brown {
6754ac538c5SBarry Smith   PetscErrorCode ierr;
6760de4c49aSJed Brown 
6770de4c49aSJed Brown   PetscFunctionBegin;
678afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6790de4c49aSJed Brown   PetscValidPointer(theta,2);
6804ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
6810de4c49aSJed Brown   PetscFunctionReturn(0);
6820de4c49aSJed Brown }
6830de4c49aSJed Brown 
6840de4c49aSJed Brown #undef __FUNCT__
6850de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
6860de4c49aSJed Brown /*@
6870de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
6880de4c49aSJed Brown 
6890de4c49aSJed Brown   Not Collective
6900de4c49aSJed Brown 
6910de4c49aSJed Brown   Input Parameter:
6920de4c49aSJed Brown +  ts - timestepping context
6930de4c49aSJed Brown -  theta - stage abscissa
6940de4c49aSJed Brown 
6950de4c49aSJed Brown   Options Database:
6960de4c49aSJed Brown .  -ts_theta_theta <theta>
6970de4c49aSJed Brown 
6980de4c49aSJed Brown   Level: Intermediate
6990de4c49aSJed Brown 
7000de4c49aSJed Brown .seealso: TSThetaGetTheta()
7010de4c49aSJed Brown @*/
7027087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
7030de4c49aSJed Brown {
7044ac538c5SBarry Smith   PetscErrorCode ierr;
7050de4c49aSJed Brown 
7060de4c49aSJed Brown   PetscFunctionBegin;
707afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7084ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
7090de4c49aSJed Brown   PetscFunctionReturn(0);
7100de4c49aSJed Brown }
711f33bbcb6SJed Brown 
712eb284becSJed Brown #undef __FUNCT__
71326f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
71426f2ff8fSLisandro Dalcin /*@
71526f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
71626f2ff8fSLisandro Dalcin 
71726f2ff8fSLisandro Dalcin   Not Collective
71826f2ff8fSLisandro Dalcin 
71926f2ff8fSLisandro Dalcin   Input Parameter:
72026f2ff8fSLisandro Dalcin .  ts - timestepping context
72126f2ff8fSLisandro Dalcin 
72226f2ff8fSLisandro Dalcin   Output Parameter:
72326f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
72426f2ff8fSLisandro Dalcin 
72526f2ff8fSLisandro Dalcin   Level: Advanced
72626f2ff8fSLisandro Dalcin 
72726f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
72826f2ff8fSLisandro Dalcin @*/
72926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
73026f2ff8fSLisandro Dalcin {
73126f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
73226f2ff8fSLisandro Dalcin 
73326f2ff8fSLisandro Dalcin   PetscFunctionBegin;
73426f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
73526f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
73626f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
73726f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
73826f2ff8fSLisandro Dalcin }
73926f2ff8fSLisandro Dalcin 
74026f2ff8fSLisandro Dalcin #undef __FUNCT__
741eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
742eb284becSJed Brown /*@
743eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
744eb284becSJed Brown 
745eb284becSJed Brown   Not Collective
746eb284becSJed Brown 
747eb284becSJed Brown   Input Parameter:
748eb284becSJed Brown +  ts - timestepping context
749eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
750eb284becSJed Brown 
751eb284becSJed Brown   Options Database:
752eb284becSJed Brown .  -ts_theta_endpoint <flg>
753eb284becSJed Brown 
754eb284becSJed Brown   Level: Intermediate
755eb284becSJed Brown 
756eb284becSJed Brown .seealso: TSTHETA, TSCN
757eb284becSJed Brown @*/
758eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
759eb284becSJed Brown {
760eb284becSJed Brown   PetscErrorCode ierr;
761eb284becSJed Brown 
762eb284becSJed Brown   PetscFunctionBegin;
763eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
764eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
765eb284becSJed Brown   PetscFunctionReturn(0);
766eb284becSJed Brown }
767eb284becSJed Brown 
768f33bbcb6SJed Brown /*
769f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
770f33bbcb6SJed Brown  * The creation functions for these specializations are below.
771f33bbcb6SJed Brown  */
772f33bbcb6SJed Brown 
773f33bbcb6SJed Brown #undef __FUNCT__
774f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
775f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
776f33bbcb6SJed Brown {
777d52bd9f3SBarry Smith   PetscErrorCode ierr;
778d52bd9f3SBarry Smith 
779f33bbcb6SJed Brown   PetscFunctionBegin;
780d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
781f33bbcb6SJed Brown   PetscFunctionReturn(0);
782f33bbcb6SJed Brown }
783f33bbcb6SJed Brown 
784f33bbcb6SJed Brown /*MC
785f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
786f33bbcb6SJed Brown 
787f33bbcb6SJed Brown   Level: beginner
788f33bbcb6SJed Brown 
7894eb428fdSBarry Smith   Notes:
790c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
7914eb428fdSBarry Smith 
7924eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
7934eb428fdSBarry Smith 
794f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
795f33bbcb6SJed Brown 
796f33bbcb6SJed Brown M*/
797f33bbcb6SJed Brown #undef __FUNCT__
798f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
7998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
800f33bbcb6SJed Brown {
801f33bbcb6SJed Brown   PetscErrorCode ierr;
802f33bbcb6SJed Brown 
803f33bbcb6SJed Brown   PetscFunctionBegin;
804f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
805f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
806f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
807f33bbcb6SJed Brown   PetscFunctionReturn(0);
808f33bbcb6SJed Brown }
809f33bbcb6SJed Brown 
810f33bbcb6SJed Brown #undef __FUNCT__
811f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
812f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
813f33bbcb6SJed Brown {
814d52bd9f3SBarry Smith   PetscErrorCode ierr;
815d52bd9f3SBarry Smith 
816f33bbcb6SJed Brown   PetscFunctionBegin;
817d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
818f33bbcb6SJed Brown   PetscFunctionReturn(0);
819f33bbcb6SJed Brown }
820f33bbcb6SJed Brown 
821f33bbcb6SJed Brown /*MC
822f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
823f33bbcb6SJed Brown 
824f33bbcb6SJed Brown   Level: beginner
825f33bbcb6SJed Brown 
826f33bbcb6SJed Brown   Notes:
8277cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
8287cf5af47SJed Brown 
8297cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
830f33bbcb6SJed Brown 
831f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
832f33bbcb6SJed Brown 
833f33bbcb6SJed Brown M*/
834f33bbcb6SJed Brown #undef __FUNCT__
835f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
8368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
837f33bbcb6SJed Brown {
838f33bbcb6SJed Brown   PetscErrorCode ierr;
839f33bbcb6SJed Brown 
840f33bbcb6SJed Brown   PetscFunctionBegin;
841f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
842f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
843eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
844f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
845f33bbcb6SJed Brown   PetscFunctionReturn(0);
846f33bbcb6SJed Brown }
847