xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 42682096d3b4a0b3863ee209ed6652af8e092da7)
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 */
12ace3abfcSBarry Smith   PetscBool    extrapolate;
13eb284becSJed Brown   PetscBool    endpoint;
14316643e7SJed Brown   PetscReal    Theta;
15316643e7SJed Brown   PetscReal    stage_time;
163b1890cdSShri Abhyankar   TSStepStatus status;
173b1890cdSShri Abhyankar   char         *name;
183b1890cdSShri Abhyankar   PetscInt     order;
193b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
203b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
21316643e7SJed Brown } TS_Theta;
22316643e7SJed Brown 
23316643e7SJed Brown #undef __FUNCT__
247445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
257445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
267445fe48SJed Brown {
277445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
287445fe48SJed Brown   PetscErrorCode ierr;
297445fe48SJed Brown 
307445fe48SJed Brown   PetscFunctionBegin;
317445fe48SJed Brown   if (X0) {
327445fe48SJed Brown     if (dm && dm != ts->dm) {
330d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
347445fe48SJed Brown     } else *X0 = ts->vec_sol;
357445fe48SJed Brown   }
367445fe48SJed Brown   if (Xdot) {
377445fe48SJed Brown     if (dm && dm != ts->dm) {
380d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
397445fe48SJed Brown     } else *Xdot = th->Xdot;
407445fe48SJed Brown   }
417445fe48SJed Brown   PetscFunctionReturn(0);
427445fe48SJed Brown }
437445fe48SJed Brown 
440d0b770aSPeter Brune 
450d0b770aSPeter Brune #undef __FUNCT__
460d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
470d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
480d0b770aSPeter Brune {
490d0b770aSPeter Brune   PetscErrorCode ierr;
500d0b770aSPeter Brune 
510d0b770aSPeter Brune   PetscFunctionBegin;
520d0b770aSPeter Brune   if (X0) {
530d0b770aSPeter Brune     if (dm && dm != ts->dm) {
540d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
550d0b770aSPeter Brune     }
560d0b770aSPeter Brune   }
570d0b770aSPeter Brune   if (Xdot) {
580d0b770aSPeter Brune     if (dm && dm != ts->dm) {
590d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
600d0b770aSPeter Brune     }
610d0b770aSPeter Brune   }
620d0b770aSPeter Brune   PetscFunctionReturn(0);
630d0b770aSPeter Brune }
640d0b770aSPeter Brune 
657445fe48SJed Brown #undef __FUNCT__
667445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
677445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
687445fe48SJed Brown {
697445fe48SJed Brown 
707445fe48SJed Brown   PetscFunctionBegin;
717445fe48SJed Brown   PetscFunctionReturn(0);
727445fe48SJed Brown }
737445fe48SJed Brown 
747445fe48SJed Brown #undef __FUNCT__
757445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
767445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
777445fe48SJed Brown {
787445fe48SJed Brown   TS             ts = (TS)ctx;
797445fe48SJed Brown   PetscErrorCode ierr;
807445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
817445fe48SJed Brown 
827445fe48SJed Brown   PetscFunctionBegin;
837445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
847445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
857445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
867445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
877445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
890d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
900d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
917445fe48SJed Brown   PetscFunctionReturn(0);
927445fe48SJed Brown }
937445fe48SJed Brown 
947445fe48SJed Brown #undef __FUNCT__
95258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
96258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
97258e1594SPeter Brune {
98258e1594SPeter Brune 
99258e1594SPeter Brune   PetscFunctionBegin;
100258e1594SPeter Brune   PetscFunctionReturn(0);
101258e1594SPeter Brune }
102258e1594SPeter Brune 
103258e1594SPeter Brune #undef __FUNCT__
104258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
105258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
106258e1594SPeter Brune {
107258e1594SPeter Brune   TS             ts = (TS)ctx;
108258e1594SPeter Brune   PetscErrorCode ierr;
109258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
110258e1594SPeter Brune 
111258e1594SPeter Brune   PetscFunctionBegin;
112258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
113258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
114258e1594SPeter Brune 
115258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
116258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
117258e1594SPeter Brune 
118258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune 
121258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
122258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
123258e1594SPeter Brune   PetscFunctionReturn(0);
124258e1594SPeter Brune }
125258e1594SPeter Brune 
1263b1890cdSShri Abhyankar #undef __FUNCT__
1273b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1283b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1293b1890cdSShri Abhyankar {
1303b1890cdSShri Abhyankar   PetscErrorCode ierr;
1313b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1323b1890cdSShri Abhyankar 
1333b1890cdSShri Abhyankar   PetscFunctionBegin;
134ce94432eSBarry 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");
1353b1890cdSShri Abhyankar   if (order == th->order) {
1363b1890cdSShri Abhyankar     if (th->endpoint) {
1373b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1383b1890cdSShri Abhyankar     } else {
1393b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1403b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1413b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1423b1890cdSShri Abhyankar     }
1433b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1443b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1453b1890cdSShri Abhyankar   }
1463b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1473b1890cdSShri Abhyankar }
148258e1594SPeter Brune 
149258e1594SPeter Brune #undef __FUNCT__
15024655328SShri #define __FUNCT__ "TSRollBack_Theta"
15124655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15224655328SShri {
15324655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
15424655328SShri   PetscErrorCode ierr;
15524655328SShri 
15624655328SShri   PetscFunctionBegin;
15724655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
15824655328SShri   th->status    = TS_STEP_INCOMPLETE;
15924655328SShri   PetscFunctionReturn(0);
16024655328SShri }
16124655328SShri 
16224655328SShri #undef __FUNCT__
163316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
164193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
165316643e7SJed Brown {
166316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1673b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
1683b1890cdSShri Abhyankar   PetscReal      next_time_step;
1693b1890cdSShri Abhyankar   TSAdapt        adapt;
1704957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
171051f2191SLisandro Dalcin   PetscErrorCode ierr;
172316643e7SJed Brown 
173316643e7SJed Brown   PetscFunctionBegin;
1743b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1753b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
176051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
177b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
178eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
179b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
180b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
181316643e7SJed Brown 
182eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
183eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
184eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
185eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
186eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
187eb284becSJed Brown     }
188ace68cafSJed Brown     if (th->extrapolate) {
189b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
190ace68cafSJed Brown     } else {
1912b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
192ace68cafSJed Brown     }
193eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
194316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
195316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1965ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
197051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
198552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1994957b756SLisandro Dalcin     ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr);
200051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
201051f2191SLisandro Dalcin 
2020298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
203051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2043b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
205552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2063b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2070298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2083b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
209051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
210051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
211051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
212051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
213051f2191SLisandro Dalcin       goto reject_step;
214051f2191SLisandro Dalcin     }
2153b1890cdSShri Abhyankar 
2163b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2172b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
218cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
219316643e7SJed Brown     ts->steps++;
2203b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
221051f2191SLisandro Dalcin     break;
222051f2191SLisandro Dalcin 
223051f2191SLisandro Dalcin reject_step:
224051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
225051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
226051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2273b1890cdSShri Abhyankar     }
228051f2191SLisandro Dalcin     continue;
2293b1890cdSShri Abhyankar   }
230316643e7SJed Brown   PetscFunctionReturn(0);
231316643e7SJed Brown }
232316643e7SJed Brown 
233cd652676SJed Brown #undef __FUNCT__
234cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
235cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
236cd652676SJed Brown {
237cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
2385a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
239cd652676SJed Brown   PetscErrorCode ierr;
240cd652676SJed Brown 
241cd652676SJed Brown   PetscFunctionBegin;
242a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
2435a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
2445a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
245cd652676SJed Brown   PetscFunctionReturn(0);
246cd652676SJed Brown }
247cd652676SJed Brown 
248316643e7SJed Brown /*------------------------------------------------------------*/
249316643e7SJed Brown #undef __FUNCT__
250277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
251277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
252316643e7SJed Brown {
253316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
254316643e7SJed Brown   PetscErrorCode ierr;
255316643e7SJed Brown 
256316643e7SJed Brown   PetscFunctionBegin;
2576bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
2586bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
2593b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
260eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
261277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
262277b19d0SLisandro Dalcin }
263277b19d0SLisandro Dalcin 
264277b19d0SLisandro Dalcin #undef __FUNCT__
265277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
266277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
267277b19d0SLisandro Dalcin {
268277b19d0SLisandro Dalcin   PetscErrorCode ierr;
269277b19d0SLisandro Dalcin 
270277b19d0SLisandro Dalcin   PetscFunctionBegin;
271277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
272277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
273bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
277316643e7SJed Brown   PetscFunctionReturn(0);
278316643e7SJed Brown }
279316643e7SJed Brown 
280316643e7SJed Brown /*
281316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
2822b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
283316643e7SJed Brown */
284316643e7SJed Brown #undef __FUNCT__
2850f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
2860f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
287316643e7SJed Brown {
288316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
289316643e7SJed Brown   PetscErrorCode ierr;
2907445fe48SJed Brown   Vec            X0,Xdot;
2917445fe48SJed Brown   DM             dm,dmsave;
292b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
293316643e7SJed Brown 
294316643e7SJed Brown   PetscFunctionBegin;
2957445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2965a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
2977445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
298b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
2997445fe48SJed Brown 
3007445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
3017445fe48SJed Brown   dmsave = ts->dm;
3027445fe48SJed Brown   ts->dm = dm;
3037445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
3047445fe48SJed Brown   ts->dm = dmsave;
3050d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
306316643e7SJed Brown   PetscFunctionReturn(0);
307316643e7SJed Brown }
308316643e7SJed Brown 
309316643e7SJed Brown #undef __FUNCT__
3100f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
311d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
312316643e7SJed Brown {
313316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
314316643e7SJed Brown   PetscErrorCode ierr;
3157445fe48SJed Brown   Vec            Xdot;
3167445fe48SJed Brown   DM             dm,dmsave;
317b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
318316643e7SJed Brown 
319316643e7SJed Brown   PetscFunctionBegin;
3207445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3217445fe48SJed Brown 
3220f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
3230298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
3247445fe48SJed Brown 
3257445fe48SJed Brown   dmsave = ts->dm;
3267445fe48SJed Brown   ts->dm = dm;
327d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
3287445fe48SJed Brown   ts->dm = dmsave;
3290298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
330316643e7SJed Brown   PetscFunctionReturn(0);
331316643e7SJed Brown }
332316643e7SJed Brown 
333316643e7SJed Brown #undef __FUNCT__
334316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
335316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
336316643e7SJed Brown {
337316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
338316643e7SJed Brown   PetscErrorCode ierr;
3397445fe48SJed Brown   SNES           snes;
340ef749922SLisandro Dalcin   TSAdapt        adapt;
3417445fe48SJed Brown   DM             dm;
342316643e7SJed Brown 
343316643e7SJed Brown   PetscFunctionBegin;
344316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
345316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
3463b1890cdSShri Abhyankar   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
3477445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3487445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
3497445fe48SJed Brown   if (dm) {
3507445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
351258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
3527445fe48SJed Brown   }
3533b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
3543b1890cdSShri Abhyankar   else th->order = 1;
3553b1890cdSShri Abhyankar 
356552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
357ef749922SLisandro Dalcin   if (!th->adapt) {
3583b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
3593b1890cdSShri Abhyankar   }
360316643e7SJed Brown   PetscFunctionReturn(0);
361316643e7SJed Brown }
362316643e7SJed Brown /*------------------------------------------------------------*/
363316643e7SJed Brown 
364316643e7SJed Brown #undef __FUNCT__
365316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
3668c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_Theta(PetscOptions *PetscOptionsObject,TS ts)
367316643e7SJed Brown {
368316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
369316643e7SJed Brown   PetscErrorCode ierr;
370316643e7SJed Brown 
371316643e7SJed Brown   PetscFunctionBegin;
372e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr);
373316643e7SJed Brown   {
3740298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
3750298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
3760298fd71SBarry 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);
3770298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
378d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
379316643e7SJed Brown   }
380316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
381316643e7SJed Brown   PetscFunctionReturn(0);
382316643e7SJed Brown }
383316643e7SJed Brown 
384316643e7SJed Brown #undef __FUNCT__
385316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
386316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
387316643e7SJed Brown {
388316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
389ace3abfcSBarry Smith   PetscBool      iascii;
390316643e7SJed Brown   PetscErrorCode ierr;
391316643e7SJed Brown 
392316643e7SJed Brown   PetscFunctionBegin;
393251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
394316643e7SJed Brown   if (iascii) {
3957c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
396ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
397316643e7SJed Brown   }
398ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
399316643e7SJed Brown   PetscFunctionReturn(0);
400316643e7SJed Brown }
401316643e7SJed Brown 
4020de4c49aSJed Brown #undef __FUNCT__
4030de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
4047087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
4050de4c49aSJed Brown {
4060de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4070de4c49aSJed Brown 
4080de4c49aSJed Brown   PetscFunctionBegin;
4090de4c49aSJed Brown   *theta = th->Theta;
4100de4c49aSJed Brown   PetscFunctionReturn(0);
4110de4c49aSJed Brown }
4120de4c49aSJed Brown 
4130de4c49aSJed Brown #undef __FUNCT__
4140de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
4157087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
4160de4c49aSJed Brown {
4170de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4180de4c49aSJed Brown 
4190de4c49aSJed Brown   PetscFunctionBegin;
4207c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
4210de4c49aSJed Brown   th->Theta = theta;
4220de4c49aSJed Brown   PetscFunctionReturn(0);
4230de4c49aSJed Brown }
424eb284becSJed Brown 
425eb284becSJed Brown #undef __FUNCT__
42678e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
42726f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
42826f2ff8fSLisandro Dalcin {
42926f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
43026f2ff8fSLisandro Dalcin 
43126f2ff8fSLisandro Dalcin   PetscFunctionBegin;
43226f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
43326f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
43426f2ff8fSLisandro Dalcin }
43526f2ff8fSLisandro Dalcin 
43626f2ff8fSLisandro Dalcin #undef __FUNCT__
43726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
438eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
439eb284becSJed Brown {
440eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
441eb284becSJed Brown 
442eb284becSJed Brown   PetscFunctionBegin;
443eb284becSJed Brown   th->endpoint = flg;
444eb284becSJed Brown   PetscFunctionReturn(0);
445eb284becSJed Brown }
4460de4c49aSJed Brown 
447f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
448f9c1d6abSBarry Smith #undef __FUNCT__
449f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
450f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
451f9c1d6abSBarry Smith {
452f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
453f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
4543fd8ae06SJed Brown   const PetscReal one = 1.0;
455f9c1d6abSBarry Smith 
456f9c1d6abSBarry Smith   PetscFunctionBegin;
4573fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
458f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
459f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
460f9c1d6abSBarry Smith   PetscFunctionReturn(0);
461f9c1d6abSBarry Smith }
462f9c1d6abSBarry Smith #endif
463f9c1d6abSBarry Smith 
464*42682096SHong Zhang #undef __FUNCT__
465*42682096SHong Zhang #define __FUNCT__ "TSGetStages_Theta"
466*42682096SHong Zhang static PetscErrorCode  TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y)
467*42682096SHong Zhang {
468*42682096SHong Zhang   TS_Theta     *th = (TS_Theta*)ts->data;
469*42682096SHong Zhang 
470*42682096SHong Zhang   PetscFunctionBegin;
471*42682096SHong Zhang   *ns = (th->endpoint) ? 1 : 0 ;
472*42682096SHong Zhang   if(Y) *Y  = &(th->X);
473*42682096SHong Zhang   PetscFunctionReturn(0);
474*42682096SHong Zhang }
475f9c1d6abSBarry Smith 
476316643e7SJed Brown /* ------------------------------------------------------------ */
477316643e7SJed Brown /*MC
47896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
479316643e7SJed Brown 
480316643e7SJed Brown    Level: beginner
481316643e7SJed Brown 
4824eb428fdSBarry Smith    Options Database:
4830c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
4844eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
4850c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
4864eb428fdSBarry Smith 
487eb284becSJed Brown    Notes:
4880c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
4890c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
4904eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
4914eb428fdSBarry Smith 
4924eb428fdSBarry Smith 
4934eb428fdSBarry Smith 
494eb284becSJed Brown    This method can be applied to DAE.
495eb284becSJed Brown 
496eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
497eb284becSJed Brown 
498eb284becSJed Brown .vb
499eb284becSJed Brown   Theta | Theta
500eb284becSJed Brown   -------------
501eb284becSJed Brown         |  1
502eb284becSJed Brown .ve
503eb284becSJed Brown 
504eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
505eb284becSJed Brown 
506eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
507eb284becSJed Brown 
508eb284becSJed Brown .vb
509eb284becSJed Brown   0 | 0         0
510eb284becSJed Brown   1 | 1-Theta   Theta
511eb284becSJed Brown   -------------------
512eb284becSJed Brown     | 1-Theta   Theta
513eb284becSJed Brown .ve
514eb284becSJed Brown 
515eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
516eb284becSJed Brown 
517eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
518eb284becSJed Brown 
519eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
520eb284becSJed Brown 
5214eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
522eb284becSJed Brown 
523eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
524316643e7SJed Brown 
525316643e7SJed Brown M*/
526316643e7SJed Brown #undef __FUNCT__
527316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
5288cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
529316643e7SJed Brown {
530316643e7SJed Brown   TS_Theta       *th;
531316643e7SJed Brown   PetscErrorCode ierr;
532316643e7SJed Brown 
533316643e7SJed Brown   PetscFunctionBegin;
534277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
535316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
536316643e7SJed Brown   ts->ops->view           = TSView_Theta;
537316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
538316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
539cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
5403b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
54124655328SShri   ts->ops->rollback       = TSRollBack_Theta;
542316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
5430f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
5440f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
545f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
546f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
547f9c1d6abSBarry Smith #endif
548*42682096SHong Zhang   ts->ops->getstages      = TSGetStages_Theta;
549316643e7SJed Brown 
550b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
551316643e7SJed Brown   ts->data = (void*)th;
552316643e7SJed Brown 
5536f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
554316643e7SJed Brown   th->Theta       = 0.5;
5553b1890cdSShri Abhyankar   th->ccfl        = 1.0;
5563b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
557bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
558bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
559bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
560bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
561316643e7SJed Brown   PetscFunctionReturn(0);
562316643e7SJed Brown }
5630de4c49aSJed Brown 
5640de4c49aSJed Brown #undef __FUNCT__
5650de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
5660de4c49aSJed Brown /*@
5670de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
5680de4c49aSJed Brown 
5690de4c49aSJed Brown   Not Collective
5700de4c49aSJed Brown 
5710de4c49aSJed Brown   Input Parameter:
5720de4c49aSJed Brown .  ts - timestepping context
5730de4c49aSJed Brown 
5740de4c49aSJed Brown   Output Parameter:
5750de4c49aSJed Brown .  theta - stage abscissa
5760de4c49aSJed Brown 
5770de4c49aSJed Brown   Note:
5780de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
5790de4c49aSJed Brown 
5800de4c49aSJed Brown   Level: Advanced
5810de4c49aSJed Brown 
5820de4c49aSJed Brown .seealso: TSThetaSetTheta()
5830de4c49aSJed Brown @*/
5847087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
5850de4c49aSJed Brown {
5864ac538c5SBarry Smith   PetscErrorCode ierr;
5870de4c49aSJed Brown 
5880de4c49aSJed Brown   PetscFunctionBegin;
589afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5900de4c49aSJed Brown   PetscValidPointer(theta,2);
5914ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
5920de4c49aSJed Brown   PetscFunctionReturn(0);
5930de4c49aSJed Brown }
5940de4c49aSJed Brown 
5950de4c49aSJed Brown #undef __FUNCT__
5960de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
5970de4c49aSJed Brown /*@
5980de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
5990de4c49aSJed Brown 
6000de4c49aSJed Brown   Not Collective
6010de4c49aSJed Brown 
6020de4c49aSJed Brown   Input Parameter:
6030de4c49aSJed Brown +  ts - timestepping context
6040de4c49aSJed Brown -  theta - stage abscissa
6050de4c49aSJed Brown 
6060de4c49aSJed Brown   Options Database:
6070de4c49aSJed Brown .  -ts_theta_theta <theta>
6080de4c49aSJed Brown 
6090de4c49aSJed Brown   Level: Intermediate
6100de4c49aSJed Brown 
6110de4c49aSJed Brown .seealso: TSThetaGetTheta()
6120de4c49aSJed Brown @*/
6137087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
6140de4c49aSJed Brown {
6154ac538c5SBarry Smith   PetscErrorCode ierr;
6160de4c49aSJed Brown 
6170de4c49aSJed Brown   PetscFunctionBegin;
618afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6194ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
6200de4c49aSJed Brown   PetscFunctionReturn(0);
6210de4c49aSJed Brown }
622f33bbcb6SJed Brown 
623eb284becSJed Brown #undef __FUNCT__
62426f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
62526f2ff8fSLisandro Dalcin /*@
62626f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
62726f2ff8fSLisandro Dalcin 
62826f2ff8fSLisandro Dalcin   Not Collective
62926f2ff8fSLisandro Dalcin 
63026f2ff8fSLisandro Dalcin   Input Parameter:
63126f2ff8fSLisandro Dalcin .  ts - timestepping context
63226f2ff8fSLisandro Dalcin 
63326f2ff8fSLisandro Dalcin   Output Parameter:
63426f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
63526f2ff8fSLisandro Dalcin 
63626f2ff8fSLisandro Dalcin   Level: Advanced
63726f2ff8fSLisandro Dalcin 
63826f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
63926f2ff8fSLisandro Dalcin @*/
64026f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
64126f2ff8fSLisandro Dalcin {
64226f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
64326f2ff8fSLisandro Dalcin 
64426f2ff8fSLisandro Dalcin   PetscFunctionBegin;
64526f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
64626f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
64726f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
64826f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
64926f2ff8fSLisandro Dalcin }
65026f2ff8fSLisandro Dalcin 
65126f2ff8fSLisandro Dalcin #undef __FUNCT__
652eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
653eb284becSJed Brown /*@
654eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
655eb284becSJed Brown 
656eb284becSJed Brown   Not Collective
657eb284becSJed Brown 
658eb284becSJed Brown   Input Parameter:
659eb284becSJed Brown +  ts - timestepping context
660eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
661eb284becSJed Brown 
662eb284becSJed Brown   Options Database:
663eb284becSJed Brown .  -ts_theta_endpoint <flg>
664eb284becSJed Brown 
665eb284becSJed Brown   Level: Intermediate
666eb284becSJed Brown 
667eb284becSJed Brown .seealso: TSTHETA, TSCN
668eb284becSJed Brown @*/
669eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
670eb284becSJed Brown {
671eb284becSJed Brown   PetscErrorCode ierr;
672eb284becSJed Brown 
673eb284becSJed Brown   PetscFunctionBegin;
674eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
675eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
676eb284becSJed Brown   PetscFunctionReturn(0);
677eb284becSJed Brown }
678eb284becSJed Brown 
679f33bbcb6SJed Brown /*
680f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
681f33bbcb6SJed Brown  * The creation functions for these specializations are below.
682f33bbcb6SJed Brown  */
683f33bbcb6SJed Brown 
684f33bbcb6SJed Brown #undef __FUNCT__
685f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
686f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
687f33bbcb6SJed Brown {
688d52bd9f3SBarry Smith   PetscErrorCode ierr;
689d52bd9f3SBarry Smith 
690f33bbcb6SJed Brown   PetscFunctionBegin;
691d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
692f33bbcb6SJed Brown   PetscFunctionReturn(0);
693f33bbcb6SJed Brown }
694f33bbcb6SJed Brown 
695f33bbcb6SJed Brown /*MC
696f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
697f33bbcb6SJed Brown 
698f33bbcb6SJed Brown   Level: beginner
699f33bbcb6SJed Brown 
7004eb428fdSBarry Smith   Notes:
701c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
7024eb428fdSBarry Smith 
7034eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
7044eb428fdSBarry Smith 
705f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
706f33bbcb6SJed Brown 
707f33bbcb6SJed Brown M*/
708f33bbcb6SJed Brown #undef __FUNCT__
709f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
7108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
711f33bbcb6SJed Brown {
712f33bbcb6SJed Brown   PetscErrorCode ierr;
713f33bbcb6SJed Brown 
714f33bbcb6SJed Brown   PetscFunctionBegin;
715f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
716f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
717f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
718f33bbcb6SJed Brown   PetscFunctionReturn(0);
719f33bbcb6SJed Brown }
720f33bbcb6SJed Brown 
721f33bbcb6SJed Brown #undef __FUNCT__
722f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
723f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
724f33bbcb6SJed Brown {
725d52bd9f3SBarry Smith   PetscErrorCode ierr;
726d52bd9f3SBarry Smith 
727f33bbcb6SJed Brown   PetscFunctionBegin;
728d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
729f33bbcb6SJed Brown   PetscFunctionReturn(0);
730f33bbcb6SJed Brown }
731f33bbcb6SJed Brown 
732f33bbcb6SJed Brown /*MC
733f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
734f33bbcb6SJed Brown 
735f33bbcb6SJed Brown   Level: beginner
736f33bbcb6SJed Brown 
737f33bbcb6SJed Brown   Notes:
7387cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
7397cf5af47SJed Brown 
7407cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
741f33bbcb6SJed Brown 
742f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
743f33bbcb6SJed Brown 
744f33bbcb6SJed Brown M*/
745f33bbcb6SJed Brown #undef __FUNCT__
746f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
7478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
748f33bbcb6SJed Brown {
749f33bbcb6SJed Brown   PetscErrorCode ierr;
750f33bbcb6SJed Brown 
751f33bbcb6SJed Brown   PetscFunctionBegin;
752f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
753f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
754eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
755f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
756f33bbcb6SJed Brown   PetscFunctionReturn(0);
757f33bbcb6SJed Brown }
758