xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 051f21914e098ffbbbb8cc4b9383afc3070b8f37)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4f9c1d6abSBarry Smith #define PETSC_DESIRE_COMPLEX
5b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
6a055c8caSBarry Smith #include <petscsnes.h>
71e25c274SJed Brown #include <petscdm.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 */
13ace3abfcSBarry Smith   PetscBool    extrapolate;
14eb284becSJed Brown   PetscBool    endpoint;
15316643e7SJed Brown   PetscReal    Theta;
16316643e7SJed Brown   PetscReal    stage_time;
173b1890cdSShri Abhyankar   TSStepStatus status;
183b1890cdSShri Abhyankar   char         *name;
193b1890cdSShri Abhyankar   PetscInt     order;
203b1890cdSShri Abhyankar   PetscReal    ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
213b1890cdSShri Abhyankar   PetscBool    adapt;  /* use time-step adaptivity ? */
22316643e7SJed Brown } TS_Theta;
23316643e7SJed Brown 
24316643e7SJed Brown #undef __FUNCT__
257445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
267445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
277445fe48SJed Brown {
287445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
297445fe48SJed Brown   PetscErrorCode ierr;
307445fe48SJed Brown 
317445fe48SJed Brown   PetscFunctionBegin;
327445fe48SJed Brown   if (X0) {
337445fe48SJed Brown     if (dm && dm != ts->dm) {
340d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
357445fe48SJed Brown     } else *X0 = ts->vec_sol;
367445fe48SJed Brown   }
377445fe48SJed Brown   if (Xdot) {
387445fe48SJed Brown     if (dm && dm != ts->dm) {
390d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
407445fe48SJed Brown     } else *Xdot = th->Xdot;
417445fe48SJed Brown   }
427445fe48SJed Brown   PetscFunctionReturn(0);
437445fe48SJed Brown }
447445fe48SJed Brown 
450d0b770aSPeter Brune 
460d0b770aSPeter Brune #undef __FUNCT__
470d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
480d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
490d0b770aSPeter Brune {
500d0b770aSPeter Brune   PetscErrorCode ierr;
510d0b770aSPeter Brune 
520d0b770aSPeter Brune   PetscFunctionBegin;
530d0b770aSPeter Brune   if (X0) {
540d0b770aSPeter Brune     if (dm && dm != ts->dm) {
550d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
560d0b770aSPeter Brune     }
570d0b770aSPeter Brune   }
580d0b770aSPeter Brune   if (Xdot) {
590d0b770aSPeter Brune     if (dm && dm != ts->dm) {
600d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
610d0b770aSPeter Brune     }
620d0b770aSPeter Brune   }
630d0b770aSPeter Brune   PetscFunctionReturn(0);
640d0b770aSPeter Brune }
650d0b770aSPeter Brune 
667445fe48SJed Brown #undef __FUNCT__
677445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
687445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
697445fe48SJed Brown {
707445fe48SJed Brown 
717445fe48SJed Brown   PetscFunctionBegin;
727445fe48SJed Brown   PetscFunctionReturn(0);
737445fe48SJed Brown }
747445fe48SJed Brown 
757445fe48SJed Brown #undef __FUNCT__
767445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
777445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
787445fe48SJed Brown {
797445fe48SJed Brown   TS             ts = (TS)ctx;
807445fe48SJed Brown   PetscErrorCode ierr;
817445fe48SJed Brown   Vec            X0,Xdot,X0_c,Xdot_c;
827445fe48SJed Brown 
837445fe48SJed Brown   PetscFunctionBegin;
847445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
857445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
867445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
877445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
887445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
897445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
900d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
910d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
927445fe48SJed Brown   PetscFunctionReturn(0);
937445fe48SJed Brown }
947445fe48SJed Brown 
957445fe48SJed Brown #undef __FUNCT__
96258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_TSTheta"
97258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx)
98258e1594SPeter Brune {
99258e1594SPeter Brune 
100258e1594SPeter Brune   PetscFunctionBegin;
101258e1594SPeter Brune   PetscFunctionReturn(0);
102258e1594SPeter Brune }
103258e1594SPeter Brune 
104258e1594SPeter Brune #undef __FUNCT__
105258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSTheta"
106258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
107258e1594SPeter Brune {
108258e1594SPeter Brune   TS             ts = (TS)ctx;
109258e1594SPeter Brune   PetscErrorCode ierr;
110258e1594SPeter Brune   Vec            X0,Xdot,X0_sub,Xdot_sub;
111258e1594SPeter Brune 
112258e1594SPeter Brune   PetscFunctionBegin;
113258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
114258e1594SPeter Brune   ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
115258e1594SPeter Brune 
116258e1594SPeter Brune   ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
117258e1594SPeter Brune   ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118258e1594SPeter Brune 
119258e1594SPeter Brune   ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120258e1594SPeter Brune   ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121258e1594SPeter Brune 
122258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
123258e1594SPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr);
124258e1594SPeter Brune   PetscFunctionReturn(0);
125258e1594SPeter Brune }
126258e1594SPeter Brune 
1273b1890cdSShri Abhyankar #undef __FUNCT__
1283b1890cdSShri Abhyankar #define __FUNCT__ "TSEvaluateStep_Theta"
1293b1890cdSShri Abhyankar static PetscErrorCode TSEvaluateStep_Theta(TS ts,PetscInt order,Vec U,PetscBool *done)
1303b1890cdSShri Abhyankar {
1313b1890cdSShri Abhyankar   PetscErrorCode ierr;
1323b1890cdSShri Abhyankar   TS_Theta       *th = (TS_Theta*)ts->data;
1333b1890cdSShri Abhyankar 
1343b1890cdSShri Abhyankar   PetscFunctionBegin;
135ce94432eSBarry 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");
1363b1890cdSShri Abhyankar   if (order == th->order) {
1373b1890cdSShri Abhyankar     if (th->endpoint) {
1383b1890cdSShri Abhyankar       ierr = VecCopy(th->X,U);CHKERRQ(ierr);
1393b1890cdSShri Abhyankar     } else {
1403b1890cdSShri Abhyankar       PetscReal shift = 1./(th->Theta*ts->time_step);
1413b1890cdSShri Abhyankar       ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,U,th->X);CHKERRQ(ierr);
1423b1890cdSShri Abhyankar       ierr = VecAXPY(U,ts->time_step,th->Xdot);CHKERRQ(ierr);
1433b1890cdSShri Abhyankar     }
1443b1890cdSShri Abhyankar   } else if (order == th->order-1 && order) {
1453b1890cdSShri Abhyankar     ierr = VecWAXPY(U,ts->time_step,th->Xdot,th->X0);CHKERRQ(ierr);
1463b1890cdSShri Abhyankar   }
1473b1890cdSShri Abhyankar   PetscFunctionReturn(0);
1483b1890cdSShri Abhyankar }
149258e1594SPeter Brune 
150258e1594SPeter Brune #undef __FUNCT__
15124655328SShri #define __FUNCT__ "TSRollBack_Theta"
15224655328SShri static PetscErrorCode TSRollBack_Theta(TS ts)
15324655328SShri {
15424655328SShri   TS_Theta       *th = (TS_Theta*)ts->data;
15524655328SShri   PetscErrorCode ierr;
15624655328SShri 
15724655328SShri   PetscFunctionBegin;
15824655328SShri   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
15924655328SShri   th->status    = TS_STEP_INCOMPLETE;
16024655328SShri   PetscFunctionReturn(0);
16124655328SShri }
16224655328SShri 
16324655328SShri #undef __FUNCT__
164316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
165193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
166316643e7SJed Brown {
167316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1683b1890cdSShri Abhyankar   PetscInt       its,lits,reject,next_scheme;
1693b1890cdSShri Abhyankar   PetscReal      next_time_step;
1703b1890cdSShri Abhyankar   TSAdapt        adapt;
1714957b756SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
172*051f2191SLisandro Dalcin   PetscErrorCode ierr;
173316643e7SJed Brown 
174316643e7SJed Brown   PetscFunctionBegin;
1753b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1763b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
177*051f2191SLisandro Dalcin   for (reject=0; !ts->reason && th->status != TS_STEP_COMPLETE; ts->reject++) {
178b296d7d5SJed Brown     PetscReal shift = 1./(th->Theta*ts->time_step);
179eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
180b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
181b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
182316643e7SJed Brown 
183eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
184eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
185eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
186eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
187eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
188eb284becSJed Brown     }
189ace68cafSJed Brown     if (th->extrapolate) {
190b296d7d5SJed Brown       ierr = VecWAXPY(th->X,1./shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
191ace68cafSJed Brown     } else {
1922b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
193ace68cafSJed Brown     }
194eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
195316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
196316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
1975ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
198*051f2191SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
199552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2004957b756SLisandro Dalcin     ierr = TSAdaptCheckStage(adapt,ts,&stageok);CHKERRQ(ierr);
201*051f2191SLisandro Dalcin     if (!stageok) {accept = PETSC_FALSE; goto reject_step;}
202*051f2191SLisandro Dalcin 
2030298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
204*051f2191SLisandro Dalcin     th->status = TS_STEP_PENDING;
2053b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
206552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2073b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2080298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2093b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
210*051f2191SLisandro Dalcin     if (!accept) {           /* Roll back the current step */
211*051f2191SLisandro Dalcin       ts->ptime += next_time_step; /* This will be undone in rollback */
212*051f2191SLisandro Dalcin       th->status = TS_STEP_INCOMPLETE;
213*051f2191SLisandro Dalcin       ierr = TSRollBack(ts);CHKERRQ(ierr);
214*051f2191SLisandro Dalcin       goto reject_step;
215*051f2191SLisandro Dalcin     }
2163b1890cdSShri Abhyankar 
2173b1890cdSShri Abhyankar     /* ignore next_scheme for now */
2182b5a38e1SLisandro Dalcin     ts->ptime    += ts->time_step;
219cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
220316643e7SJed Brown     ts->steps++;
2213b1890cdSShri Abhyankar     th->status = TS_STEP_COMPLETE;
222*051f2191SLisandro Dalcin     break;
223*051f2191SLisandro Dalcin 
224*051f2191SLisandro Dalcin reject_step:
225*051f2191SLisandro Dalcin     if (!ts->reason && ++reject > ts->max_reject && ts->max_reject >= 0) {
226*051f2191SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
227*051f2191SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,reject);CHKERRQ(ierr);
2283b1890cdSShri Abhyankar     }
229*051f2191SLisandro Dalcin     continue;
2303b1890cdSShri Abhyankar   }
231316643e7SJed Brown   PetscFunctionReturn(0);
232316643e7SJed Brown }
233316643e7SJed Brown 
234cd652676SJed Brown #undef __FUNCT__
235cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
236cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
237cd652676SJed Brown {
238cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
2395a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
240cd652676SJed Brown   PetscErrorCode ierr;
241cd652676SJed Brown 
242cd652676SJed Brown   PetscFunctionBegin;
243a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
2445a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
2455a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
246cd652676SJed Brown   PetscFunctionReturn(0);
247cd652676SJed Brown }
248cd652676SJed Brown 
249316643e7SJed Brown /*------------------------------------------------------------*/
250316643e7SJed Brown #undef __FUNCT__
251277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
252277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
253316643e7SJed Brown {
254316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
255316643e7SJed Brown   PetscErrorCode ierr;
256316643e7SJed Brown 
257316643e7SJed Brown   PetscFunctionBegin;
2586bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
2596bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
2603b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
261eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
262277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
263277b19d0SLisandro Dalcin }
264277b19d0SLisandro Dalcin 
265277b19d0SLisandro Dalcin #undef __FUNCT__
266277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
267277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
268277b19d0SLisandro Dalcin {
269277b19d0SLisandro Dalcin   PetscErrorCode ierr;
270277b19d0SLisandro Dalcin 
271277b19d0SLisandro Dalcin   PetscFunctionBegin;
272277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
273277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
274bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
275bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
276bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
277bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
278316643e7SJed Brown   PetscFunctionReturn(0);
279316643e7SJed Brown }
280316643e7SJed Brown 
281316643e7SJed Brown /*
282316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
2832b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
284316643e7SJed Brown */
285316643e7SJed Brown #undef __FUNCT__
2860f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
2870f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
288316643e7SJed Brown {
289316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
290316643e7SJed Brown   PetscErrorCode ierr;
2917445fe48SJed Brown   Vec            X0,Xdot;
2927445fe48SJed Brown   DM             dm,dmsave;
293b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
294316643e7SJed Brown 
295316643e7SJed Brown   PetscFunctionBegin;
2967445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2975a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
2987445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
299b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
3007445fe48SJed Brown 
3017445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
3027445fe48SJed Brown   dmsave = ts->dm;
3037445fe48SJed Brown   ts->dm = dm;
3047445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
3057445fe48SJed Brown   ts->dm = dmsave;
3060d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
307316643e7SJed Brown   PetscFunctionReturn(0);
308316643e7SJed Brown }
309316643e7SJed Brown 
310316643e7SJed Brown #undef __FUNCT__
3110f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
312d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
313316643e7SJed Brown {
314316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
315316643e7SJed Brown   PetscErrorCode ierr;
3167445fe48SJed Brown   Vec            Xdot;
3177445fe48SJed Brown   DM             dm,dmsave;
318b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
319316643e7SJed Brown 
320316643e7SJed Brown   PetscFunctionBegin;
3217445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3227445fe48SJed Brown 
3230f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
3240298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
3257445fe48SJed Brown 
3267445fe48SJed Brown   dmsave = ts->dm;
3277445fe48SJed Brown   ts->dm = dm;
328d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
3297445fe48SJed Brown   ts->dm = dmsave;
3300298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
331316643e7SJed Brown   PetscFunctionReturn(0);
332316643e7SJed Brown }
333316643e7SJed Brown 
334316643e7SJed Brown #undef __FUNCT__
335316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
336316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
337316643e7SJed Brown {
338316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
339316643e7SJed Brown   PetscErrorCode ierr;
3407445fe48SJed Brown   SNES           snes;
341ef749922SLisandro Dalcin   TSAdapt        adapt;
3427445fe48SJed Brown   DM             dm;
343316643e7SJed Brown 
344316643e7SJed Brown   PetscFunctionBegin;
345316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
346316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
3473b1890cdSShri Abhyankar   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
3487445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3497445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
3507445fe48SJed Brown   if (dm) {
3517445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
352258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
3537445fe48SJed Brown   }
3543b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
3553b1890cdSShri Abhyankar   else th->order = 1;
3563b1890cdSShri Abhyankar 
357552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
358ef749922SLisandro Dalcin   if (!th->adapt) {
3593b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
3603b1890cdSShri Abhyankar   }
361316643e7SJed Brown   PetscFunctionReturn(0);
362316643e7SJed Brown }
363316643e7SJed Brown /*------------------------------------------------------------*/
364316643e7SJed Brown 
365316643e7SJed Brown #undef __FUNCT__
366316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
367316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
368316643e7SJed Brown {
369316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
370316643e7SJed Brown   PetscErrorCode ierr;
371316643e7SJed Brown 
372316643e7SJed Brown   PetscFunctionBegin;
373d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
374316643e7SJed Brown   {
3750298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
3760298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
3770298fd71SBarry 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);
3780298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
379d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
380316643e7SJed Brown   }
381316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
382316643e7SJed Brown   PetscFunctionReturn(0);
383316643e7SJed Brown }
384316643e7SJed Brown 
385316643e7SJed Brown #undef __FUNCT__
386316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
387316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
388316643e7SJed Brown {
389316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
390ace3abfcSBarry Smith   PetscBool      iascii;
391316643e7SJed Brown   PetscErrorCode ierr;
392316643e7SJed Brown 
393316643e7SJed Brown   PetscFunctionBegin;
394251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
395316643e7SJed Brown   if (iascii) {
3967c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
397ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
398316643e7SJed Brown   }
399ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
400316643e7SJed Brown   PetscFunctionReturn(0);
401316643e7SJed Brown }
402316643e7SJed Brown 
4030de4c49aSJed Brown #undef __FUNCT__
4040de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
4057087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
4060de4c49aSJed Brown {
4070de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4080de4c49aSJed Brown 
4090de4c49aSJed Brown   PetscFunctionBegin;
4100de4c49aSJed Brown   *theta = th->Theta;
4110de4c49aSJed Brown   PetscFunctionReturn(0);
4120de4c49aSJed Brown }
4130de4c49aSJed Brown 
4140de4c49aSJed Brown #undef __FUNCT__
4150de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
4167087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
4170de4c49aSJed Brown {
4180de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4190de4c49aSJed Brown 
4200de4c49aSJed Brown   PetscFunctionBegin;
4217c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
4220de4c49aSJed Brown   th->Theta = theta;
4230de4c49aSJed Brown   PetscFunctionReturn(0);
4240de4c49aSJed Brown }
425eb284becSJed Brown 
426eb284becSJed Brown #undef __FUNCT__
42778e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
42826f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
42926f2ff8fSLisandro Dalcin {
43026f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
43126f2ff8fSLisandro Dalcin 
43226f2ff8fSLisandro Dalcin   PetscFunctionBegin;
43326f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
43426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
43526f2ff8fSLisandro Dalcin }
43626f2ff8fSLisandro Dalcin 
43726f2ff8fSLisandro Dalcin #undef __FUNCT__
43826f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
439eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
440eb284becSJed Brown {
441eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
442eb284becSJed Brown 
443eb284becSJed Brown   PetscFunctionBegin;
444eb284becSJed Brown   th->endpoint = flg;
445eb284becSJed Brown   PetscFunctionReturn(0);
446eb284becSJed Brown }
4470de4c49aSJed Brown 
448f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
449f9c1d6abSBarry Smith #undef __FUNCT__
450f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
451f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
452f9c1d6abSBarry Smith {
453f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
454f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
4553fd8ae06SJed Brown   const PetscReal one = 1.0;
456f9c1d6abSBarry Smith 
457f9c1d6abSBarry Smith   PetscFunctionBegin;
4583fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
459f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
460f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
461f9c1d6abSBarry Smith   PetscFunctionReturn(0);
462f9c1d6abSBarry Smith }
463f9c1d6abSBarry Smith #endif
464f9c1d6abSBarry Smith 
465f9c1d6abSBarry Smith 
466316643e7SJed Brown /* ------------------------------------------------------------ */
467316643e7SJed Brown /*MC
46896f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
469316643e7SJed Brown 
470316643e7SJed Brown    Level: beginner
471316643e7SJed Brown 
4724eb428fdSBarry Smith    Options Database:
4730c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
4744eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
4750c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
4764eb428fdSBarry Smith 
477eb284becSJed Brown    Notes:
4780c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
4790c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
4804eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
4814eb428fdSBarry Smith 
4824eb428fdSBarry Smith 
4834eb428fdSBarry Smith 
484eb284becSJed Brown    This method can be applied to DAE.
485eb284becSJed Brown 
486eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
487eb284becSJed Brown 
488eb284becSJed Brown .vb
489eb284becSJed Brown   Theta | Theta
490eb284becSJed Brown   -------------
491eb284becSJed Brown         |  1
492eb284becSJed Brown .ve
493eb284becSJed Brown 
494eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
495eb284becSJed Brown 
496eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
497eb284becSJed Brown 
498eb284becSJed Brown .vb
499eb284becSJed Brown   0 | 0         0
500eb284becSJed Brown   1 | 1-Theta   Theta
501eb284becSJed Brown   -------------------
502eb284becSJed Brown     | 1-Theta   Theta
503eb284becSJed Brown .ve
504eb284becSJed Brown 
505eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
506eb284becSJed Brown 
507eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
508eb284becSJed Brown 
509eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
510eb284becSJed Brown 
5114eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
512eb284becSJed Brown 
513eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
514316643e7SJed Brown 
515316643e7SJed Brown M*/
516316643e7SJed Brown #undef __FUNCT__
517316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
5188cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
519316643e7SJed Brown {
520316643e7SJed Brown   TS_Theta       *th;
521316643e7SJed Brown   PetscErrorCode ierr;
522316643e7SJed Brown 
523316643e7SJed Brown   PetscFunctionBegin;
524277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
525316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
526316643e7SJed Brown   ts->ops->view           = TSView_Theta;
527316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
528316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
529cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
5303b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
53124655328SShri   ts->ops->rollback       = TSRollBack_Theta;
532316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
5330f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
5340f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
535f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
536f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
537f9c1d6abSBarry Smith #endif
538316643e7SJed Brown 
539b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
540316643e7SJed Brown   ts->data = (void*)th;
541316643e7SJed Brown 
5426f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
543316643e7SJed Brown   th->Theta       = 0.5;
5443b1890cdSShri Abhyankar   th->ccfl        = 1.0;
5453b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
546bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
547bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
548bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
549bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
550316643e7SJed Brown   PetscFunctionReturn(0);
551316643e7SJed Brown }
5520de4c49aSJed Brown 
5530de4c49aSJed Brown #undef __FUNCT__
5540de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
5550de4c49aSJed Brown /*@
5560de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
5570de4c49aSJed Brown 
5580de4c49aSJed Brown   Not Collective
5590de4c49aSJed Brown 
5600de4c49aSJed Brown   Input Parameter:
5610de4c49aSJed Brown .  ts - timestepping context
5620de4c49aSJed Brown 
5630de4c49aSJed Brown   Output Parameter:
5640de4c49aSJed Brown .  theta - stage abscissa
5650de4c49aSJed Brown 
5660de4c49aSJed Brown   Note:
5670de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
5680de4c49aSJed Brown 
5690de4c49aSJed Brown   Level: Advanced
5700de4c49aSJed Brown 
5710de4c49aSJed Brown .seealso: TSThetaSetTheta()
5720de4c49aSJed Brown @*/
5737087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
5740de4c49aSJed Brown {
5754ac538c5SBarry Smith   PetscErrorCode ierr;
5760de4c49aSJed Brown 
5770de4c49aSJed Brown   PetscFunctionBegin;
578afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5790de4c49aSJed Brown   PetscValidPointer(theta,2);
5804ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
5810de4c49aSJed Brown   PetscFunctionReturn(0);
5820de4c49aSJed Brown }
5830de4c49aSJed Brown 
5840de4c49aSJed Brown #undef __FUNCT__
5850de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
5860de4c49aSJed Brown /*@
5870de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
5880de4c49aSJed Brown 
5890de4c49aSJed Brown   Not Collective
5900de4c49aSJed Brown 
5910de4c49aSJed Brown   Input Parameter:
5920de4c49aSJed Brown +  ts - timestepping context
5930de4c49aSJed Brown -  theta - stage abscissa
5940de4c49aSJed Brown 
5950de4c49aSJed Brown   Options Database:
5960de4c49aSJed Brown .  -ts_theta_theta <theta>
5970de4c49aSJed Brown 
5980de4c49aSJed Brown   Level: Intermediate
5990de4c49aSJed Brown 
6000de4c49aSJed Brown .seealso: TSThetaGetTheta()
6010de4c49aSJed Brown @*/
6027087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
6030de4c49aSJed Brown {
6044ac538c5SBarry Smith   PetscErrorCode ierr;
6050de4c49aSJed Brown 
6060de4c49aSJed Brown   PetscFunctionBegin;
607afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6084ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
6090de4c49aSJed Brown   PetscFunctionReturn(0);
6100de4c49aSJed Brown }
611f33bbcb6SJed Brown 
612eb284becSJed Brown #undef __FUNCT__
61326f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
61426f2ff8fSLisandro Dalcin /*@
61526f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
61626f2ff8fSLisandro Dalcin 
61726f2ff8fSLisandro Dalcin   Not Collective
61826f2ff8fSLisandro Dalcin 
61926f2ff8fSLisandro Dalcin   Input Parameter:
62026f2ff8fSLisandro Dalcin .  ts - timestepping context
62126f2ff8fSLisandro Dalcin 
62226f2ff8fSLisandro Dalcin   Output Parameter:
62326f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
62426f2ff8fSLisandro Dalcin 
62526f2ff8fSLisandro Dalcin   Level: Advanced
62626f2ff8fSLisandro Dalcin 
62726f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
62826f2ff8fSLisandro Dalcin @*/
62926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
63026f2ff8fSLisandro Dalcin {
63126f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
63226f2ff8fSLisandro Dalcin 
63326f2ff8fSLisandro Dalcin   PetscFunctionBegin;
63426f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
63526f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
63626f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
63726f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
63826f2ff8fSLisandro Dalcin }
63926f2ff8fSLisandro Dalcin 
64026f2ff8fSLisandro Dalcin #undef __FUNCT__
641eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
642eb284becSJed Brown /*@
643eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
644eb284becSJed Brown 
645eb284becSJed Brown   Not Collective
646eb284becSJed Brown 
647eb284becSJed Brown   Input Parameter:
648eb284becSJed Brown +  ts - timestepping context
649eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
650eb284becSJed Brown 
651eb284becSJed Brown   Options Database:
652eb284becSJed Brown .  -ts_theta_endpoint <flg>
653eb284becSJed Brown 
654eb284becSJed Brown   Level: Intermediate
655eb284becSJed Brown 
656eb284becSJed Brown .seealso: TSTHETA, TSCN
657eb284becSJed Brown @*/
658eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
659eb284becSJed Brown {
660eb284becSJed Brown   PetscErrorCode ierr;
661eb284becSJed Brown 
662eb284becSJed Brown   PetscFunctionBegin;
663eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
664eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
665eb284becSJed Brown   PetscFunctionReturn(0);
666eb284becSJed Brown }
667eb284becSJed Brown 
668f33bbcb6SJed Brown /*
669f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
670f33bbcb6SJed Brown  * The creation functions for these specializations are below.
671f33bbcb6SJed Brown  */
672f33bbcb6SJed Brown 
673f33bbcb6SJed Brown #undef __FUNCT__
674f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
675f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
676f33bbcb6SJed Brown {
677d52bd9f3SBarry Smith   PetscErrorCode ierr;
678d52bd9f3SBarry Smith 
679f33bbcb6SJed Brown   PetscFunctionBegin;
680d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
681f33bbcb6SJed Brown   PetscFunctionReturn(0);
682f33bbcb6SJed Brown }
683f33bbcb6SJed Brown 
684f33bbcb6SJed Brown /*MC
685f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
686f33bbcb6SJed Brown 
687f33bbcb6SJed Brown   Level: beginner
688f33bbcb6SJed Brown 
6894eb428fdSBarry Smith   Notes:
690c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
6914eb428fdSBarry Smith 
6924eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
6934eb428fdSBarry Smith 
694f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
695f33bbcb6SJed Brown 
696f33bbcb6SJed Brown M*/
697f33bbcb6SJed Brown #undef __FUNCT__
698f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
6998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
700f33bbcb6SJed Brown {
701f33bbcb6SJed Brown   PetscErrorCode ierr;
702f33bbcb6SJed Brown 
703f33bbcb6SJed Brown   PetscFunctionBegin;
704f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
705f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
706f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
707f33bbcb6SJed Brown   PetscFunctionReturn(0);
708f33bbcb6SJed Brown }
709f33bbcb6SJed Brown 
710f33bbcb6SJed Brown #undef __FUNCT__
711f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
712f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
713f33bbcb6SJed Brown {
714d52bd9f3SBarry Smith   PetscErrorCode ierr;
715d52bd9f3SBarry Smith 
716f33bbcb6SJed Brown   PetscFunctionBegin;
717d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
718f33bbcb6SJed Brown   PetscFunctionReturn(0);
719f33bbcb6SJed Brown }
720f33bbcb6SJed Brown 
721f33bbcb6SJed Brown /*MC
722f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
723f33bbcb6SJed Brown 
724f33bbcb6SJed Brown   Level: beginner
725f33bbcb6SJed Brown 
726f33bbcb6SJed Brown   Notes:
7277cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
7287cf5af47SJed Brown 
7297cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
730f33bbcb6SJed Brown 
731f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
732f33bbcb6SJed Brown 
733f33bbcb6SJed Brown M*/
734f33bbcb6SJed Brown #undef __FUNCT__
735f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
7368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
737f33bbcb6SJed Brown {
738f33bbcb6SJed Brown   PetscErrorCode ierr;
739f33bbcb6SJed Brown 
740f33bbcb6SJed Brown   PetscFunctionBegin;
741f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
742f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
743eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
744f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
745f33bbcb6SJed Brown   PetscFunctionReturn(0);
746f33bbcb6SJed Brown }
747