xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 2ca6e9204ea7b65f3b887030ac693c318e313a26)
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>
7316643e7SJed Brown 
8316643e7SJed Brown typedef struct {
9316643e7SJed Brown   Vec          X,Xdot;                   /* Storage for one stage */
103b1890cdSShri Abhyankar   Vec          X0;                       /* work vector to store X0 */
11eb284becSJed Brown   Vec          affine;                   /* Affine vector needed for residual at beginning of step */
12*2ca6e920SHong Zhang   Vec          *VecDeltaLam;             /* Increment of the adjoint sensitivity w.r.t IC at stage*/
13*2ca6e920SHong Zhang   Vec          *VecDeltaMu;              /* Increment of the adjoint sensitivity w.r.t P at stage*/
14*2ca6e920SHong Zhang   Vec          *VecSensiTemp;            /* Vector to be timed with Jacobian transpose*/
15ace3abfcSBarry Smith   PetscBool    extrapolate;
16eb284becSJed Brown   PetscBool    endpoint;
17316643e7SJed Brown   PetscReal    Theta;
18316643e7SJed Brown   PetscReal    stage_time;
193b1890cdSShri Abhyankar   TSStepStatus status;
203b1890cdSShri Abhyankar   char         *name;
213b1890cdSShri Abhyankar   PetscInt     order;
223b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
233b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
24316643e7SJed Brown } TS_Theta;
25316643e7SJed Brown 
26316643e7SJed Brown #undef __FUNCT__
277445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
287445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
297445fe48SJed Brown {
307445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
317445fe48SJed Brown   PetscErrorCode ierr;
327445fe48SJed Brown 
337445fe48SJed Brown   PetscFunctionBegin;
347445fe48SJed Brown   if (X0) {
357445fe48SJed Brown     if (dm && dm != ts->dm) {
360d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
377445fe48SJed Brown     } else *X0 = ts->vec_sol;
387445fe48SJed Brown   }
397445fe48SJed Brown   if (Xdot) {
407445fe48SJed Brown     if (dm && dm != ts->dm) {
410d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
427445fe48SJed Brown     } else *Xdot = th->Xdot;
437445fe48SJed Brown   }
447445fe48SJed Brown   PetscFunctionReturn(0);
457445fe48SJed Brown }
467445fe48SJed Brown 
470d0b770aSPeter Brune 
480d0b770aSPeter Brune #undef __FUNCT__
490d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
500d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
510d0b770aSPeter Brune {
520d0b770aSPeter Brune   PetscErrorCode ierr;
530d0b770aSPeter Brune 
540d0b770aSPeter Brune   PetscFunctionBegin;
550d0b770aSPeter Brune   if (X0) {
560d0b770aSPeter Brune     if (dm && dm != ts->dm) {
570d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
580d0b770aSPeter Brune     }
590d0b770aSPeter Brune   }
600d0b770aSPeter Brune   if (Xdot) {
610d0b770aSPeter Brune     if (dm && dm != ts->dm) {
620d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
630d0b770aSPeter Brune     }
640d0b770aSPeter Brune   }
650d0b770aSPeter Brune   PetscFunctionReturn(0);
660d0b770aSPeter Brune }
670d0b770aSPeter Brune 
687445fe48SJed Brown #undef __FUNCT__
697445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
707445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
717445fe48SJed Brown {
727445fe48SJed Brown 
737445fe48SJed Brown   PetscFunctionBegin;
747445fe48SJed Brown   PetscFunctionReturn(0);
757445fe48SJed Brown }
767445fe48SJed Brown 
777445fe48SJed Brown #undef __FUNCT__
787445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
797445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
807445fe48SJed Brown {
817445fe48SJed Brown   TS             ts = (TS)ctx;
827445fe48SJed Brown   PetscErrorCode ierr;
837445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
847445fe48SJed Brown 
857445fe48SJed Brown   PetscFunctionBegin;
867445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
877445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
897445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
907445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
917445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
920d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
930d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
947445fe48SJed Brown   PetscFunctionReturn(0);
957445fe48SJed Brown }
967445fe48SJed Brown 
977445fe48SJed Brown #undef __FUNCT__
98258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
99258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
100258e1594SPeter Brune {
101258e1594SPeter Brune 
102258e1594SPeter Brune   PetscFunctionBegin;
103258e1594SPeter Brune   PetscFunctionReturn(0);
104258e1594SPeter Brune }
105258e1594SPeter Brune 
106258e1594SPeter Brune #undef __FUNCT__
107258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
108258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
109258e1594SPeter Brune {
110258e1594SPeter Brune   TS             ts = (TS)ctx;
111258e1594SPeter Brune   PetscErrorCode ierr;
112258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
113258e1594SPeter Brune 
114258e1594SPeter Brune   PetscFunctionBegin;
115258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
116258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
117258e1594SPeter Brune 
118258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune 
121258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
123258e1594SPeter Brune 
124258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
125258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
126258e1594SPeter Brune   PetscFunctionReturn(0);
127258e1594SPeter Brune }
128258e1594SPeter Brune 
1293b1890cdSShri Abhyankar #undef __FUNCT__
1303b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1313b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1323b1890cdSShri Abhyankar {
1333b1890cdSShri Abhyankar   PetscErrorCode ierr;
1343b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1353b1890cdSShri Abhyankar 
1363b1890cdSShri Abhyankar   PetscFunctionBegin;
137ce94432eSBarry 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");
1383b1890cdSShri Abhyankar   if (order == th->order) {
1393b1890cdSShri Abhyankar     if (th->endpoint) {
1403b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1413b1890cdSShri Abhyankar     } else {
1423b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1433b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1443b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1453b1890cdSShri Abhyankar     }
1463b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1473b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1483b1890cdSShri Abhyankar   }
1493b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1503b1890cdSShri Abhyankar }
151258e1594SPeter Brune 
152258e1594SPeter Brune #undef __FUNCT__
15324655328SShri #define __FUNCT__ "TSRollBack_Theta"
15424655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15524655328SShri {
15624655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
15724655328SShri   PetscErrorCode ierr;
15824655328SShri 
15924655328SShri   PetscFunctionBegin;
16024655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
16124655328SShri   th->status    = TS_STEP_INCOMPLETE;
16224655328SShri   PetscFunctionReturn(0);
16324655328SShri }
16424655328SShri 
16524655328SShri #undef __FUNCT__
166316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
167193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
168316643e7SJed Brown {
169316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1703b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
1713b1890cdSShri Abhyankar   PetscReal      next_time_step;
1723b1890cdSShri Abhyankar   TSAdapt        adapt;
1734957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
174051f2191SLisandro Dalcin   PetscErrorCode ierr;
175316643e7SJed Brown 
176316643e7SJed Brown   PetscFunctionBegin;
1773b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1783b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
179051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
180b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
181eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
182b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
183b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
184316643e7SJed Brown 
185eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
186eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
187eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
188eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
189eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
190eb284becSJed Brown     }
191ace68cafSJed Brown     if (th->extrapolate) {
192b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
193ace68cafSJed Brown     } else {
1942b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
195ace68cafSJed Brown     }
196eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
197316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
198316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1995ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
200051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
201552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2024957b756SLisandro Dalcin     ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr);
203051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
204051f2191SLisandro Dalcin 
2050298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
206051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2073b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
208552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2093b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2100298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2113b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
212051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
213051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
214051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
215051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
216051f2191SLisandro Dalcin       goto reject_step;
217051f2191SLisandro Dalcin     }
2183b1890cdSShri Abhyankar 
2193b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2202b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
221cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
222316643e7SJed Brown     ts->steps++;
2233b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
224051f2191SLisandro Dalcin     break;
225051f2191SLisandro Dalcin 
226051f2191SLisandro Dalcin reject_step:
227051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
228051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
229051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2303b1890cdSShri Abhyankar     }
231051f2191SLisandro Dalcin     continue;
2323b1890cdSShri Abhyankar   }
233316643e7SJed Brown   PetscFunctionReturn(0);
234316643e7SJed Brown }
235316643e7SJed Brown 
236cd652676SJed Brown #undef __FUNCT__
237*2ca6e920SHong Zhang #define __FUNCT__ "TSStepAdj_Theta"
238*2ca6e920SHong Zhang static PetscErrorCode TSStepAdj_Theta(TS ts)
239*2ca6e920SHong Zhang {
240*2ca6e920SHong Zhang   TS_Theta            *th = (TS_Theta*)ts->data;
241*2ca6e920SHong Zhang   Vec                 *VecDeltaLam = th->VecDeltaLam,*VecDeltaMu = th->VecDeltaMu,*VecSensiTemp = th->VecSensiTemp;
242*2ca6e920SHong Zhang   PetscInt            nadj;
243*2ca6e920SHong Zhang   PetscErrorCode      ierr;
244*2ca6e920SHong Zhang   Mat                 J,Jp;
245*2ca6e920SHong Zhang   KSP                 ksp;
246*2ca6e920SHong Zhang   PetscReal           shift;
247*2ca6e920SHong Zhang 
248*2ca6e920SHong Zhang   PetscFunctionBegin;
249*2ca6e920SHong Zhang 
250*2ca6e920SHong Zhang   th->status = TS_STEP_INCOMPLETE;
251*2ca6e920SHong Zhang   ierr = SNESGetKSP(ts->snes,&ksp);
252*2ca6e920SHong Zhang   ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
253*2ca6e920SHong Zhang   th->stage_time = ts->ptime + (th->endpoint ? 0. : (1-th->Theta)*ts->time_step); /* time_step is negative*/
254*2ca6e920SHong Zhang 
255*2ca6e920SHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
256*2ca6e920SHong Zhang   /* Build RHS */
257*2ca6e920SHong Zhang   shift = -1./((th->Theta-1)*ts->time_step);
258*2ca6e920SHong Zhang   ierr  = TSComputeIJacobian(ts,ts->ptime,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
259*2ca6e920SHong Zhang 
260*2ca6e920SHong Zhang   if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
261*2ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
262*2ca6e920SHong Zhang       ierr  = MatMultTranspose(J,ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
263*2ca6e920SHong Zhang     }
264*2ca6e920SHong Zhang   }else {
265*2ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
266*2ca6e920SHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
267*2ca6e920SHong Zhang       ierr = VecScale(VecSensiTemp[nadj],-1./(th->Theta*ts->time_step));CHKERRQ(ierr);
268*2ca6e920SHong Zhang     }
269*2ca6e920SHong Zhang   }
270*2ca6e920SHong Zhang   /* Build LHS */
271*2ca6e920SHong Zhang   shift = -1./(th->Theta*ts->time_step);
272*2ca6e920SHong Zhang   ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr);
273*2ca6e920SHong Zhang   ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr);
274*2ca6e920SHong Zhang 
275*2ca6e920SHong Zhang   /* Solve LHS X = RHS */
276*2ca6e920SHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
277*2ca6e920SHong Zhang     ierr = KSPSolveTranspose(ksp,VecSensiTemp[nadj],VecDeltaLam[nadj]);CHKERRQ(ierr);
278*2ca6e920SHong Zhang   }
279*2ca6e920SHong Zhang   if(th->endpoint) {
280*2ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
281*2ca6e920SHong Zhang       ierr = VecCopy(VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
282*2ca6e920SHong Zhang     }
283*2ca6e920SHong Zhang   }else {
284*2ca6e920SHong Zhang     shift = -1./(th->Theta*ts->time_step);
285*2ca6e920SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
286*2ca6e920SHong Zhang       ierr = VecAXPBYPCZ(VecSensiTemp[nadj],-shift,shift,0,VecDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr);
287*2ca6e920SHong Zhang       ierr = VecAXPY(ts->vecs_sensi[nadj],-ts->time_step,VecSensiTemp[nadj]);CHKERRQ(ierr);
288*2ca6e920SHong Zhang     }
289*2ca6e920SHong Zhang   }
290*2ca6e920SHong Zhang 
291*2ca6e920SHong Zhang   ts->ptime += ts->time_step;
292*2ca6e920SHong Zhang   ts->steps++;
293*2ca6e920SHong Zhang   th->status = TS_STEP_COMPLETE;
294*2ca6e920SHong Zhang 
295*2ca6e920SHong Zhang   PetscFunctionReturn(0);
296*2ca6e920SHong Zhang }
297*2ca6e920SHong Zhang 
298*2ca6e920SHong Zhang #undef __FUNCT__
299cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
300cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
301cd652676SJed Brown {
302cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
3035a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
304cd652676SJed Brown   PetscErrorCode ierr;
305cd652676SJed Brown 
306cd652676SJed Brown   PetscFunctionBegin;
307a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
3085a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
3095a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
310cd652676SJed Brown   PetscFunctionReturn(0);
311cd652676SJed Brown }
312cd652676SJed Brown 
313316643e7SJed Brown /*------------------------------------------------------------*/
314316643e7SJed Brown #undef __FUNCT__
315277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
316277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
317316643e7SJed Brown {
318316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
319316643e7SJed Brown   PetscErrorCode ierr;
320316643e7SJed Brown 
321316643e7SJed Brown   PetscFunctionBegin;
3226bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
3236bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
3243b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
325eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
326*2ca6e920SHong Zhang   if(ts->reverse_mode) {
327*2ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
328*2ca6e920SHong Zhang     if(th->VecDeltaMu) {
329*2ca6e920SHong Zhang       ierr = VecDestroyVecs(ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
330*2ca6e920SHong Zhang     }
331*2ca6e920SHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
332*2ca6e920SHong Zhang   }
333277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
334277b19d0SLisandro Dalcin }
335277b19d0SLisandro Dalcin 
336277b19d0SLisandro Dalcin #undef __FUNCT__
337277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
338277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
339277b19d0SLisandro Dalcin {
340277b19d0SLisandro Dalcin   PetscErrorCode ierr;
341277b19d0SLisandro Dalcin 
342277b19d0SLisandro Dalcin   PetscFunctionBegin;
343277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
344277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
345bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
346bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
347bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
348bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
349316643e7SJed Brown   PetscFunctionReturn(0);
350316643e7SJed Brown }
351316643e7SJed Brown 
352316643e7SJed Brown /*
353316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
3542b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
355316643e7SJed Brown */
356316643e7SJed Brown #undef __FUNCT__
3570f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
3580f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
359316643e7SJed Brown {
360316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
361316643e7SJed Brown   PetscErrorCode ierr;
3627445fe48SJed Brown   Vec            X0,Xdot;
3637445fe48SJed Brown   DM             dm,dmsave;
364b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
365316643e7SJed Brown 
366316643e7SJed Brown   PetscFunctionBegin;
3677445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3685a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
3697445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
370b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
3717445fe48SJed Brown 
3727445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
3737445fe48SJed Brown   dmsave = ts->dm;
3747445fe48SJed Brown   ts->dm = dm;
3757445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
3767445fe48SJed Brown   ts->dm = dmsave;
3770d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
378316643e7SJed Brown   PetscFunctionReturn(0);
379316643e7SJed Brown }
380316643e7SJed Brown 
381316643e7SJed Brown #undef __FUNCT__
3820f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
383d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
384316643e7SJed Brown {
385316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
386316643e7SJed Brown   PetscErrorCode ierr;
3877445fe48SJed Brown   Vec            Xdot;
3887445fe48SJed Brown   DM             dm,dmsave;
389b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
390316643e7SJed Brown 
391316643e7SJed Brown   PetscFunctionBegin;
3927445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3937445fe48SJed Brown 
3940f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
3950298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
3967445fe48SJed Brown 
3977445fe48SJed Brown   dmsave = ts->dm;
3987445fe48SJed Brown   ts->dm = dm;
399d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
4007445fe48SJed Brown   ts->dm = dmsave;
4010298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
402316643e7SJed Brown   PetscFunctionReturn(0);
403316643e7SJed Brown }
404316643e7SJed Brown 
405316643e7SJed Brown #undef __FUNCT__
406316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
407316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
408316643e7SJed Brown {
409316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
410316643e7SJed Brown   PetscErrorCode ierr;
4117445fe48SJed Brown   SNES           snes;
412ef749922SLisandro Dalcin   TSAdapt        adapt;
4137445fe48SJed Brown   DM             dm;
414316643e7SJed Brown 
415316643e7SJed Brown   PetscFunctionBegin;
416316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
417316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
4183b1890cdSShri Abhyankar   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
4197445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
4207445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
4217445fe48SJed Brown   if (dm) {
4227445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
423258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
4247445fe48SJed Brown   }
4253b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
4263b1890cdSShri Abhyankar   else th->order = 1;
4273b1890cdSShri Abhyankar 
428552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
429ef749922SLisandro Dalcin   if (!th->adapt) {
4303b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
4313b1890cdSShri Abhyankar   }
432*2ca6e920SHong Zhang   if (ts->reverse_mode) {
433*2ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecDeltaLam);CHKERRQ(ierr);
434*2ca6e920SHong Zhang     if(ts->vecs_sensip) {
435*2ca6e920SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numberadjs,&th->VecDeltaMu);CHKERRQ(ierr);
436*2ca6e920SHong Zhang     }
437*2ca6e920SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&th->VecSensiTemp);CHKERRQ(ierr);
438*2ca6e920SHong Zhang   }
439316643e7SJed Brown   PetscFunctionReturn(0);
440316643e7SJed Brown }
441316643e7SJed Brown /*------------------------------------------------------------*/
442316643e7SJed Brown 
443316643e7SJed Brown #undef __FUNCT__
444316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
4458c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
446316643e7SJed Brown {
447316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
448316643e7SJed Brown   PetscErrorCode ierr;
449316643e7SJed Brown 
450316643e7SJed Brown   PetscFunctionBegin;
451e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
452316643e7SJed Brown   {
4530298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
4540298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
4550298fd71SBarry 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);
4560298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
457d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
458316643e7SJed Brown   }
459316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
460316643e7SJed Brown   PetscFunctionReturn(0);
461316643e7SJed Brown }
462316643e7SJed Brown 
463316643e7SJed Brown #undef __FUNCT__
464316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
465316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
466316643e7SJed Brown {
467316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
468ace3abfcSBarry Smith   PetscBool      iascii;
469316643e7SJed Brown   PetscErrorCode ierr;
470316643e7SJed Brown 
471316643e7SJed Brown   PetscFunctionBegin;
472251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
473316643e7SJed Brown   if (iascii) {
4747c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
475ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
476316643e7SJed Brown   }
477ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
478316643e7SJed Brown   PetscFunctionReturn(0);
479316643e7SJed Brown }
480316643e7SJed Brown 
4810de4c49aSJed Brown #undef __FUNCT__
4820de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
4837087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
4840de4c49aSJed Brown {
4850de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4860de4c49aSJed Brown 
4870de4c49aSJed Brown   PetscFunctionBegin;
4880de4c49aSJed Brown   *theta = th->Theta;
4890de4c49aSJed Brown   PetscFunctionReturn(0);
4900de4c49aSJed Brown }
4910de4c49aSJed Brown 
4920de4c49aSJed Brown #undef __FUNCT__
4930de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
4947087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
4950de4c49aSJed Brown {
4960de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4970de4c49aSJed Brown 
4980de4c49aSJed Brown   PetscFunctionBegin;
4997c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
5000de4c49aSJed Brown   th->Theta = theta;
5010de4c49aSJed Brown   PetscFunctionReturn(0);
5020de4c49aSJed Brown }
503eb284becSJed Brown 
504eb284becSJed Brown #undef __FUNCT__
50578e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
50626f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
50726f2ff8fSLisandro Dalcin {
50826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
50926f2ff8fSLisandro Dalcin 
51026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
51126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
51226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
51326f2ff8fSLisandro Dalcin }
51426f2ff8fSLisandro Dalcin 
51526f2ff8fSLisandro Dalcin #undef __FUNCT__
51626f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
517eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
518eb284becSJed Brown {
519eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
520eb284becSJed Brown 
521eb284becSJed Brown   PetscFunctionBegin;
522eb284becSJed Brown   th->endpoint = flg;
523eb284becSJed Brown   PetscFunctionReturn(0);
524eb284becSJed Brown }
5250de4c49aSJed Brown 
526f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
527f9c1d6abSBarry Smith #undef __FUNCT__
528f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
529f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
530f9c1d6abSBarry Smith {
531f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
532f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
5333fd8ae06SJed Brown   const PetscReal one = 1.0;
534f9c1d6abSBarry Smith 
535f9c1d6abSBarry Smith   PetscFunctionBegin;
5363fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
537f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
538f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
539f9c1d6abSBarry Smith   PetscFunctionReturn(0);
540f9c1d6abSBarry Smith }
541f9c1d6abSBarry Smith #endif
542f9c1d6abSBarry Smith 
54342682096SHong Zhang #undef __FUNCT__
54442682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
54542682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
54642682096SHong Zhang {
54742682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
54842682096SHong Zhang 
54942682096SHong Zhang   PetscFunctionBegin;
550*2ca6e920SHong Zhang   *ns = 1;
551*2ca6e920SHong Zhang   if(Y) {
552*2ca6e920SHong Zhang     if(th->endpoint) { /* return the first (explicit) stage X0 for checkpointing */
553*2ca6e920SHong Zhang       *Y  = &(th->X0);
554*2ca6e920SHong Zhang     }else { /* return the stage value*/
555*2ca6e920SHong Zhang       *Y  = &(th->X);
556*2ca6e920SHong Zhang     }
557*2ca6e920SHong Zhang   }
55842682096SHong Zhang   PetscFunctionReturn(0);
55942682096SHong Zhang }
560f9c1d6abSBarry Smith 
561316643e7SJed Brown /* ------------------------------------------------------------ */
562316643e7SJed Brown /*MC
56396f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
564316643e7SJed Brown 
565316643e7SJed Brown    Level: beginner
566316643e7SJed Brown 
5674eb428fdSBarry Smith    Options Database:
5680c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
5694eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
5700c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
5714eb428fdSBarry Smith 
572eb284becSJed Brown    Notes:
5730c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
5740c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
5754eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
5764eb428fdSBarry Smith 
5774eb428fdSBarry Smith 
5784eb428fdSBarry Smith 
579eb284becSJed Brown    This method can be applied to DAE.
580eb284becSJed Brown 
581eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
582eb284becSJed Brown 
583eb284becSJed Brown .vb
584eb284becSJed Brown   Theta | Theta
585eb284becSJed Brown   -------------
586eb284becSJed Brown         |  1
587eb284becSJed Brown .ve
588eb284becSJed Brown 
589eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
590eb284becSJed Brown 
591eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
592eb284becSJed Brown 
593eb284becSJed Brown .vb
594eb284becSJed Brown   0 | 0         0
595eb284becSJed Brown   1 | 1-Theta   Theta
596eb284becSJed Brown   -------------------
597eb284becSJed Brown     | 1-Theta   Theta
598eb284becSJed Brown .ve
599eb284becSJed Brown 
600eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
601eb284becSJed Brown 
602eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
603eb284becSJed Brown 
604eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
605eb284becSJed Brown 
6064eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
607eb284becSJed Brown 
608eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
609316643e7SJed Brown 
610316643e7SJed Brown M*/
611316643e7SJed Brown #undef __FUNCT__
612316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
6138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
614316643e7SJed Brown {
615316643e7SJed Brown   TS_Theta       *th;
616316643e7SJed Brown   PetscErrorCode ierr;
617316643e7SJed Brown 
618316643e7SJed Brown   PetscFunctionBegin;
619277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
620316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
621316643e7SJed Brown   ts->ops->view           = TSView_Theta;
622316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
623316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
624cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
6253b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
62624655328SShri   ts->ops->rollback       = TSRollBack_Theta;
627316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
6280f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
6290f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
630f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
631f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
632f9c1d6abSBarry Smith #endif
63342682096SHong Zhang   ts->ops->getstages      = TSGetStages_Theta;
634*2ca6e920SHong Zhang   ts->ops->stepadj        = TSStepAdj_Theta;
635316643e7SJed Brown 
636b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
637316643e7SJed Brown   ts->data = (void*)th;
638316643e7SJed Brown 
6396f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
640316643e7SJed Brown   th->Theta       = 0.5;
6413b1890cdSShri Abhyankar   th->ccfl        = 1.0;
6423b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
643bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
644bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
645bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
646bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
647316643e7SJed Brown   PetscFunctionReturn(0);
648316643e7SJed Brown }
6490de4c49aSJed Brown 
6500de4c49aSJed Brown #undef __FUNCT__
6510de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
6520de4c49aSJed Brown /*@
6530de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
6540de4c49aSJed Brown 
6550de4c49aSJed Brown   Not Collective
6560de4c49aSJed Brown 
6570de4c49aSJed Brown   Input Parameter:
6580de4c49aSJed Brown .  ts - timestepping context
6590de4c49aSJed Brown 
6600de4c49aSJed Brown   Output Parameter:
6610de4c49aSJed Brown .  theta - stage abscissa
6620de4c49aSJed Brown 
6630de4c49aSJed Brown   Note:
6640de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
6650de4c49aSJed Brown 
6660de4c49aSJed Brown   Level: Advanced
6670de4c49aSJed Brown 
6680de4c49aSJed Brown .seealso: TSThetaSetTheta()
6690de4c49aSJed Brown @*/
6707087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
6710de4c49aSJed Brown {
6724ac538c5SBarry Smith   PetscErrorCode ierr;
6730de4c49aSJed Brown 
6740de4c49aSJed Brown   PetscFunctionBegin;
675afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6760de4c49aSJed Brown   PetscValidPointer(theta,2);
6774ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
6780de4c49aSJed Brown   PetscFunctionReturn(0);
6790de4c49aSJed Brown }
6800de4c49aSJed Brown 
6810de4c49aSJed Brown #undef __FUNCT__
6820de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
6830de4c49aSJed Brown /*@
6840de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
6850de4c49aSJed Brown 
6860de4c49aSJed Brown   Not Collective
6870de4c49aSJed Brown 
6880de4c49aSJed Brown   Input Parameter:
6890de4c49aSJed Brown +  ts - timestepping context
6900de4c49aSJed Brown -  theta - stage abscissa
6910de4c49aSJed Brown 
6920de4c49aSJed Brown   Options Database:
6930de4c49aSJed Brown .  -ts_theta_theta <theta>
6940de4c49aSJed Brown 
6950de4c49aSJed Brown   Level: Intermediate
6960de4c49aSJed Brown 
6970de4c49aSJed Brown .seealso: TSThetaGetTheta()
6980de4c49aSJed Brown @*/
6997087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
7000de4c49aSJed Brown {
7014ac538c5SBarry Smith   PetscErrorCode ierr;
7020de4c49aSJed Brown 
7030de4c49aSJed Brown   PetscFunctionBegin;
704afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
7054ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
7060de4c49aSJed Brown   PetscFunctionReturn(0);
7070de4c49aSJed Brown }
708f33bbcb6SJed Brown 
709eb284becSJed Brown #undef __FUNCT__
71026f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
71126f2ff8fSLisandro Dalcin /*@
71226f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
71326f2ff8fSLisandro Dalcin 
71426f2ff8fSLisandro Dalcin   Not Collective
71526f2ff8fSLisandro Dalcin 
71626f2ff8fSLisandro Dalcin   Input Parameter:
71726f2ff8fSLisandro Dalcin .  ts - timestepping context
71826f2ff8fSLisandro Dalcin 
71926f2ff8fSLisandro Dalcin   Output Parameter:
72026f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
72126f2ff8fSLisandro Dalcin 
72226f2ff8fSLisandro Dalcin   Level: Advanced
72326f2ff8fSLisandro Dalcin 
72426f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
72526f2ff8fSLisandro Dalcin @*/
72626f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
72726f2ff8fSLisandro Dalcin {
72826f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
72926f2ff8fSLisandro Dalcin 
73026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
73126f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
73226f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
73326f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
73426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
73526f2ff8fSLisandro Dalcin }
73626f2ff8fSLisandro Dalcin 
73726f2ff8fSLisandro Dalcin #undef __FUNCT__
738eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
739eb284becSJed Brown /*@
740eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
741eb284becSJed Brown 
742eb284becSJed Brown   Not Collective
743eb284becSJed Brown 
744eb284becSJed Brown   Input Parameter:
745eb284becSJed Brown +  ts - timestepping context
746eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
747eb284becSJed Brown 
748eb284becSJed Brown   Options Database:
749eb284becSJed Brown .  -ts_theta_endpoint <flg>
750eb284becSJed Brown 
751eb284becSJed Brown   Level: Intermediate
752eb284becSJed Brown 
753eb284becSJed Brown .seealso: TSTHETA, TSCN
754eb284becSJed Brown @*/
755eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
756eb284becSJed Brown {
757eb284becSJed Brown   PetscErrorCode ierr;
758eb284becSJed Brown 
759eb284becSJed Brown   PetscFunctionBegin;
760eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
761eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
762eb284becSJed Brown   PetscFunctionReturn(0);
763eb284becSJed Brown }
764eb284becSJed Brown 
765f33bbcb6SJed Brown /*
766f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
767f33bbcb6SJed Brown  * The creation functions for these specializations are below.
768f33bbcb6SJed Brown  */
769f33bbcb6SJed Brown 
770f33bbcb6SJed Brown #undef __FUNCT__
771f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
772f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
773f33bbcb6SJed Brown {
774d52bd9f3SBarry Smith   PetscErrorCode ierr;
775d52bd9f3SBarry Smith 
776f33bbcb6SJed Brown   PetscFunctionBegin;
777d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
778f33bbcb6SJed Brown   PetscFunctionReturn(0);
779f33bbcb6SJed Brown }
780f33bbcb6SJed Brown 
781f33bbcb6SJed Brown /*MC
782f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
783f33bbcb6SJed Brown 
784f33bbcb6SJed Brown   Level: beginner
785f33bbcb6SJed Brown 
7864eb428fdSBarry Smith   Notes:
787c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
7884eb428fdSBarry Smith 
7894eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
7904eb428fdSBarry Smith 
791f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
792f33bbcb6SJed Brown 
793f33bbcb6SJed Brown M*/
794f33bbcb6SJed Brown #undef __FUNCT__
795f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
7968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
797f33bbcb6SJed Brown {
798f33bbcb6SJed Brown   PetscErrorCode ierr;
799f33bbcb6SJed Brown 
800f33bbcb6SJed Brown   PetscFunctionBegin;
801f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
802f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
803f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
804f33bbcb6SJed Brown   PetscFunctionReturn(0);
805f33bbcb6SJed Brown }
806f33bbcb6SJed Brown 
807f33bbcb6SJed Brown #undef __FUNCT__
808f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
809f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
810f33bbcb6SJed Brown {
811d52bd9f3SBarry Smith   PetscErrorCode ierr;
812d52bd9f3SBarry Smith 
813f33bbcb6SJed Brown   PetscFunctionBegin;
814d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
815f33bbcb6SJed Brown   PetscFunctionReturn(0);
816f33bbcb6SJed Brown }
817f33bbcb6SJed Brown 
818f33bbcb6SJed Brown /*MC
819f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
820f33bbcb6SJed Brown 
821f33bbcb6SJed Brown   Level: beginner
822f33bbcb6SJed Brown 
823f33bbcb6SJed Brown   Notes:
8247cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
8257cf5af47SJed Brown 
8267cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
827f33bbcb6SJed Brown 
828f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
829f33bbcb6SJed Brown 
830f33bbcb6SJed Brown M*/
831f33bbcb6SJed Brown #undef __FUNCT__
832f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
8338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
834f33bbcb6SJed Brown {
835f33bbcb6SJed Brown   PetscErrorCode ierr;
836f33bbcb6SJed Brown 
837f33bbcb6SJed Brown   PetscFunctionBegin;
838f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
839f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
840eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
841f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
842f33bbcb6SJed Brown   PetscFunctionReturn(0);
843f33bbcb6SJed Brown }
844