xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision ac75fa188fe85ab501cb360879ce65a617d54b59)
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;
170f1b97656SJed Brown   SNESConvergedReason snesreason;
1712b5a38e1SLisandro Dalcin   PetscErrorCode      ierr;
1723b1890cdSShri Abhyankar   TSAdapt             adapt;
1733b1890cdSShri Abhyankar   PetscBool           accept;
174316643e7SJed Brown 
175316643e7SJed Brown   PetscFunctionBegin;
17624655328SShri   next_time_step = ts->time_step;
1773b1890cdSShri Abhyankar   th->status = TS_STEP_INCOMPLETE;
1783b1890cdSShri Abhyankar   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
1793b1890cdSShri Abhyankar   for (reject=0; reject<ts->max_reject && !ts->reason && th->status != TS_STEP_COMPLETE; reject++,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);
199f1b97656SJed Brown     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
2009be3e283SDebojyoti Ghosh     ierr = TSPostStage(ts,th->stage_time,0,&(th->X));CHKERRQ(ierr);
2015ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
202552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2033b1890cdSShri Abhyankar     ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
2043b1890cdSShri Abhyankar     if (!accept) continue;
2050298fd71SBarry Smith     ierr = TSEvaluateStep(ts,th->order,ts->vec_sol,NULL);CHKERRQ(ierr);
2063b1890cdSShri Abhyankar     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
207552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
2083b1890cdSShri Abhyankar     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
2090298fd71SBarry Smith     ierr = TSAdaptCandidateAdd(adapt,NULL,th->order,1,th->ccfl,1.0,PETSC_TRUE);CHKERRQ(ierr);
2103b1890cdSShri Abhyankar     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
2113b1890cdSShri Abhyankar 
2123b1890cdSShri Abhyankar     if (accept) {
2133b1890cdSShri Abhyankar       /* ignore next_scheme for now */
2142b5a38e1SLisandro Dalcin       ts->ptime    += ts->time_step;
215cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
216316643e7SJed Brown       ts->steps++;
2173b1890cdSShri Abhyankar       th->status = TS_STEP_COMPLETE;
2183b1890cdSShri Abhyankar     } else {                    /* Roll back the current step */
21924655328SShri       ts->ptime += next_time_step; /* This will be undone in rollback */
220ec5563edSShri       th->status = TS_STEP_INCOMPLETE;
22124655328SShri       ierr = TSRollBack(ts);CHKERRQ(ierr);
2223b1890cdSShri Abhyankar     }
2233b1890cdSShri Abhyankar   }
224316643e7SJed Brown   PetscFunctionReturn(0);
225316643e7SJed Brown }
226316643e7SJed Brown 
227cd652676SJed Brown #undef __FUNCT__
228cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
229cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
230cd652676SJed Brown {
231cd652676SJed Brown   TS_Theta       *th   = (TS_Theta*)ts->data;
2325a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
233cd652676SJed Brown   PetscErrorCode ierr;
234cd652676SJed Brown 
235cd652676SJed Brown   PetscFunctionBegin;
236a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
2375a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
2385a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
239cd652676SJed Brown   PetscFunctionReturn(0);
240cd652676SJed Brown }
241cd652676SJed Brown 
242316643e7SJed Brown /*------------------------------------------------------------*/
243316643e7SJed Brown #undef __FUNCT__
244277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
245277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
246316643e7SJed Brown {
247316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
248316643e7SJed Brown   PetscErrorCode ierr;
249316643e7SJed Brown 
250316643e7SJed Brown   PetscFunctionBegin;
2516bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
2526bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
2533b1890cdSShri Abhyankar   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
254eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
255277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
256277b19d0SLisandro Dalcin }
257277b19d0SLisandro Dalcin 
258277b19d0SLisandro Dalcin #undef __FUNCT__
259277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
260277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
261277b19d0SLisandro Dalcin {
262277b19d0SLisandro Dalcin   PetscErrorCode ierr;
263277b19d0SLisandro Dalcin 
264277b19d0SLisandro Dalcin   PetscFunctionBegin;
265277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
266277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
267bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr);
268bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr);
269bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr);
270bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr);
271316643e7SJed Brown   PetscFunctionReturn(0);
272316643e7SJed Brown }
273316643e7SJed Brown 
274316643e7SJed Brown /*
275316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
2762b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
277316643e7SJed Brown */
278316643e7SJed Brown #undef __FUNCT__
2790f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
2800f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
281316643e7SJed Brown {
282316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
283316643e7SJed Brown   PetscErrorCode ierr;
2847445fe48SJed Brown   Vec            X0,Xdot;
2857445fe48SJed Brown   DM             dm,dmsave;
286b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
287316643e7SJed Brown 
288316643e7SJed Brown   PetscFunctionBegin;
2897445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2905a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
2917445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
292b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr);
2937445fe48SJed Brown 
2947445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
2957445fe48SJed Brown   dmsave = ts->dm;
2967445fe48SJed Brown   ts->dm = dm;
2977445fe48SJed Brown   ierr   = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
2987445fe48SJed Brown   ts->dm = dmsave;
2990d0b770aSPeter Brune   ierr   = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
300316643e7SJed Brown   PetscFunctionReturn(0);
301316643e7SJed Brown }
302316643e7SJed Brown 
303316643e7SJed Brown #undef __FUNCT__
3040f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
305d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts)
306316643e7SJed Brown {
307316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
308316643e7SJed Brown   PetscErrorCode ierr;
3097445fe48SJed Brown   Vec            Xdot;
3107445fe48SJed Brown   DM             dm,dmsave;
311b296d7d5SJed Brown   PetscReal      shift = 1./(th->Theta*ts->time_step);
312316643e7SJed Brown 
313316643e7SJed Brown   PetscFunctionBegin;
3147445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
3157445fe48SJed Brown 
3160f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
3170298fd71SBarry Smith   ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
3187445fe48SJed Brown 
3197445fe48SJed Brown   dmsave = ts->dm;
3207445fe48SJed Brown   ts->dm = dm;
321d1e9a80fSBarry Smith   ierr   = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr);
3227445fe48SJed Brown   ts->dm = dmsave;
3230298fd71SBarry Smith   ierr   = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr);
324316643e7SJed Brown   PetscFunctionReturn(0);
325316643e7SJed Brown }
326316643e7SJed Brown 
327316643e7SJed Brown #undef __FUNCT__
328316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
329316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
330316643e7SJed Brown {
331316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
332316643e7SJed Brown   PetscErrorCode ierr;
3337445fe48SJed Brown   SNES           snes;
334ef749922SLisandro Dalcin   TSAdapt        adapt;
3357445fe48SJed Brown   DM             dm;
336316643e7SJed Brown 
337316643e7SJed Brown   PetscFunctionBegin;
338316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
339316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
3403b1890cdSShri Abhyankar   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
3417445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
3427445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
3437445fe48SJed Brown   if (dm) {
3447445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
345258e1594SPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr);
3467445fe48SJed Brown   }
3473b1890cdSShri Abhyankar   if (th->Theta == 0.5 && th->endpoint) th->order = 2;
3483b1890cdSShri Abhyankar   else th->order = 1;
3493b1890cdSShri Abhyankar 
350552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
351ef749922SLisandro Dalcin   if (!th->adapt) {
3523b1890cdSShri Abhyankar     ierr = TSAdaptSetType(adapt,TSADAPTNONE);CHKERRQ(ierr);
3533b1890cdSShri Abhyankar   }
354316643e7SJed Brown   PetscFunctionReturn(0);
355316643e7SJed Brown }
356316643e7SJed Brown /*------------------------------------------------------------*/
357316643e7SJed Brown 
358316643e7SJed Brown #undef __FUNCT__
359316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
360316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
361316643e7SJed Brown {
362316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
363316643e7SJed Brown   PetscErrorCode ierr;
364316643e7SJed Brown 
365316643e7SJed Brown   PetscFunctionBegin;
366d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
367316643e7SJed Brown   {
3680298fd71SBarry Smith     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr);
3690298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr);
3700298fd71SBarry 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);
3710298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_theta_adapt","Use time-step adaptivity with the Theta method","",th->adapt,&th->adapt,NULL);CHKERRQ(ierr);
372d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
373316643e7SJed Brown   }
374316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
375316643e7SJed Brown   PetscFunctionReturn(0);
376316643e7SJed Brown }
377316643e7SJed Brown 
378316643e7SJed Brown #undef __FUNCT__
379316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
380316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
381316643e7SJed Brown {
382316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
383ace3abfcSBarry Smith   PetscBool      iascii;
384316643e7SJed Brown   PetscErrorCode ierr;
385316643e7SJed Brown 
386316643e7SJed Brown   PetscFunctionBegin;
387251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
388316643e7SJed Brown   if (iascii) {
3897c8652ddSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%g\n",(double)th->Theta);CHKERRQ(ierr);
390ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr);
391316643e7SJed Brown   }
392*ac75fa18SLisandro Dalcin   if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
393316643e7SJed Brown   PetscFunctionReturn(0);
394316643e7SJed Brown }
395316643e7SJed Brown 
3960de4c49aSJed Brown #undef __FUNCT__
3970de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
3987087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
3990de4c49aSJed Brown {
4000de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4010de4c49aSJed Brown 
4020de4c49aSJed Brown   PetscFunctionBegin;
4030de4c49aSJed Brown   *theta = th->Theta;
4040de4c49aSJed Brown   PetscFunctionReturn(0);
4050de4c49aSJed Brown }
4060de4c49aSJed Brown 
4070de4c49aSJed Brown #undef __FUNCT__
4080de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
4097087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
4100de4c49aSJed Brown {
4110de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
4120de4c49aSJed Brown 
4130de4c49aSJed Brown   PetscFunctionBegin;
4147c8652ddSBarry Smith   if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta);
4150de4c49aSJed Brown   th->Theta = theta;
4160de4c49aSJed Brown   PetscFunctionReturn(0);
4170de4c49aSJed Brown }
418eb284becSJed Brown 
419eb284becSJed Brown #undef __FUNCT__
42078e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
42126f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
42226f2ff8fSLisandro Dalcin {
42326f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
42426f2ff8fSLisandro Dalcin 
42526f2ff8fSLisandro Dalcin   PetscFunctionBegin;
42626f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
42726f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
42826f2ff8fSLisandro Dalcin }
42926f2ff8fSLisandro Dalcin 
43026f2ff8fSLisandro Dalcin #undef __FUNCT__
43126f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
432eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
433eb284becSJed Brown {
434eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
435eb284becSJed Brown 
436eb284becSJed Brown   PetscFunctionBegin;
437eb284becSJed Brown   th->endpoint = flg;
438eb284becSJed Brown   PetscFunctionReturn(0);
439eb284becSJed Brown }
4400de4c49aSJed Brown 
441f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
442f9c1d6abSBarry Smith #undef __FUNCT__
443f9c1d6abSBarry Smith #define __FUNCT__ "TSComputeLinearStability_Theta"
444f9c1d6abSBarry Smith static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi)
445f9c1d6abSBarry Smith {
446f9c1d6abSBarry Smith   PetscComplex z   = xr + xi*PETSC_i,f;
447f9c1d6abSBarry Smith   TS_Theta     *th = (TS_Theta*)ts->data;
4483fd8ae06SJed Brown   const PetscReal one = 1.0;
449f9c1d6abSBarry Smith 
450f9c1d6abSBarry Smith   PetscFunctionBegin;
4513fd8ae06SJed Brown   f   = (one + (one - th->Theta)*z)/(one - th->Theta*z);
452f9c1d6abSBarry Smith   *yr = PetscRealPartComplex(f);
453f9c1d6abSBarry Smith   *yi = PetscImaginaryPartComplex(f);
454f9c1d6abSBarry Smith   PetscFunctionReturn(0);
455f9c1d6abSBarry Smith }
456f9c1d6abSBarry Smith #endif
457f9c1d6abSBarry Smith 
458f9c1d6abSBarry Smith 
459316643e7SJed Brown /* ------------------------------------------------------------ */
460316643e7SJed Brown /*MC
46196f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
462316643e7SJed Brown 
463316643e7SJed Brown    Level: beginner
464316643e7SJed Brown 
4654eb428fdSBarry Smith    Options Database:
4660c3ba866SJed Brown       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1)
4674eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
4680c3ba866SJed Brown       -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
4694eb428fdSBarry Smith 
470eb284becSJed Brown    Notes:
4710c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
4720c3ba866SJed Brown $  -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
4734eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
4744eb428fdSBarry Smith 
4754eb428fdSBarry Smith 
4764eb428fdSBarry Smith 
477eb284becSJed Brown    This method can be applied to DAE.
478eb284becSJed Brown 
479eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
480eb284becSJed Brown 
481eb284becSJed Brown .vb
482eb284becSJed Brown   Theta | Theta
483eb284becSJed Brown   -------------
484eb284becSJed Brown         |  1
485eb284becSJed Brown .ve
486eb284becSJed Brown 
487eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
488eb284becSJed Brown 
489eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
490eb284becSJed Brown 
491eb284becSJed Brown .vb
492eb284becSJed Brown   0 | 0         0
493eb284becSJed Brown   1 | 1-Theta   Theta
494eb284becSJed Brown   -------------------
495eb284becSJed Brown     | 1-Theta   Theta
496eb284becSJed Brown .ve
497eb284becSJed Brown 
498eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
499eb284becSJed Brown 
500eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
501eb284becSJed Brown 
502eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
503eb284becSJed Brown 
5044eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
505eb284becSJed Brown 
506eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
507316643e7SJed Brown 
508316643e7SJed Brown M*/
509316643e7SJed Brown #undef __FUNCT__
510316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
5118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
512316643e7SJed Brown {
513316643e7SJed Brown   TS_Theta       *th;
514316643e7SJed Brown   PetscErrorCode ierr;
515316643e7SJed Brown 
516316643e7SJed Brown   PetscFunctionBegin;
517277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
518316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
519316643e7SJed Brown   ts->ops->view           = TSView_Theta;
520316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
521316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
522cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
5233b1890cdSShri Abhyankar   ts->ops->evaluatestep   = TSEvaluateStep_Theta;
52424655328SShri   ts->ops->rollback       = TSRollBack_Theta;
525316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
5260f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
5270f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
528f9c1d6abSBarry Smith #if defined(PETSC_HAVE_COMPLEX)
529f9c1d6abSBarry Smith   ts->ops->linearstability = TSComputeLinearStability_Theta;
530f9c1d6abSBarry Smith #endif
531316643e7SJed Brown 
532b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
533316643e7SJed Brown   ts->data = (void*)th;
534316643e7SJed Brown 
5356f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
536316643e7SJed Brown   th->Theta       = 0.5;
5373b1890cdSShri Abhyankar   th->ccfl        = 1.0;
5383b1890cdSShri Abhyankar   th->adapt       = PETSC_FALSE;
539bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr);
540bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr);
541bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
542bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
543316643e7SJed Brown   PetscFunctionReturn(0);
544316643e7SJed Brown }
5450de4c49aSJed Brown 
5460de4c49aSJed Brown #undef __FUNCT__
5470de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
5480de4c49aSJed Brown /*@
5490de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
5500de4c49aSJed Brown 
5510de4c49aSJed Brown   Not Collective
5520de4c49aSJed Brown 
5530de4c49aSJed Brown   Input Parameter:
5540de4c49aSJed Brown .  ts - timestepping context
5550de4c49aSJed Brown 
5560de4c49aSJed Brown   Output Parameter:
5570de4c49aSJed Brown .  theta - stage abscissa
5580de4c49aSJed Brown 
5590de4c49aSJed Brown   Note:
5600de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
5610de4c49aSJed Brown 
5620de4c49aSJed Brown   Level: Advanced
5630de4c49aSJed Brown 
5640de4c49aSJed Brown .seealso: TSThetaSetTheta()
5650de4c49aSJed Brown @*/
5667087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
5670de4c49aSJed Brown {
5684ac538c5SBarry Smith   PetscErrorCode ierr;
5690de4c49aSJed Brown 
5700de4c49aSJed Brown   PetscFunctionBegin;
571afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5720de4c49aSJed Brown   PetscValidPointer(theta,2);
5734ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
5740de4c49aSJed Brown   PetscFunctionReturn(0);
5750de4c49aSJed Brown }
5760de4c49aSJed Brown 
5770de4c49aSJed Brown #undef __FUNCT__
5780de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
5790de4c49aSJed Brown /*@
5800de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
5810de4c49aSJed Brown 
5820de4c49aSJed Brown   Not Collective
5830de4c49aSJed Brown 
5840de4c49aSJed Brown   Input Parameter:
5850de4c49aSJed Brown +  ts - timestepping context
5860de4c49aSJed Brown -  theta - stage abscissa
5870de4c49aSJed Brown 
5880de4c49aSJed Brown   Options Database:
5890de4c49aSJed Brown .  -ts_theta_theta <theta>
5900de4c49aSJed Brown 
5910de4c49aSJed Brown   Level: Intermediate
5920de4c49aSJed Brown 
5930de4c49aSJed Brown .seealso: TSThetaGetTheta()
5940de4c49aSJed Brown @*/
5957087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
5960de4c49aSJed Brown {
5974ac538c5SBarry Smith   PetscErrorCode ierr;
5980de4c49aSJed Brown 
5990de4c49aSJed Brown   PetscFunctionBegin;
600afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
6014ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
6020de4c49aSJed Brown   PetscFunctionReturn(0);
6030de4c49aSJed Brown }
604f33bbcb6SJed Brown 
605eb284becSJed Brown #undef __FUNCT__
60626f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
60726f2ff8fSLisandro Dalcin /*@
60826f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
60926f2ff8fSLisandro Dalcin 
61026f2ff8fSLisandro Dalcin   Not Collective
61126f2ff8fSLisandro Dalcin 
61226f2ff8fSLisandro Dalcin   Input Parameter:
61326f2ff8fSLisandro Dalcin .  ts - timestepping context
61426f2ff8fSLisandro Dalcin 
61526f2ff8fSLisandro Dalcin   Output Parameter:
61626f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
61726f2ff8fSLisandro Dalcin 
61826f2ff8fSLisandro Dalcin   Level: Advanced
61926f2ff8fSLisandro Dalcin 
62026f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
62126f2ff8fSLisandro Dalcin @*/
62226f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
62326f2ff8fSLisandro Dalcin {
62426f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
62526f2ff8fSLisandro Dalcin 
62626f2ff8fSLisandro Dalcin   PetscFunctionBegin;
62726f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
62826f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
62926f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
63026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
63126f2ff8fSLisandro Dalcin }
63226f2ff8fSLisandro Dalcin 
63326f2ff8fSLisandro Dalcin #undef __FUNCT__
634eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
635eb284becSJed Brown /*@
636eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
637eb284becSJed Brown 
638eb284becSJed Brown   Not Collective
639eb284becSJed Brown 
640eb284becSJed Brown   Input Parameter:
641eb284becSJed Brown +  ts - timestepping context
642eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
643eb284becSJed Brown 
644eb284becSJed Brown   Options Database:
645eb284becSJed Brown .  -ts_theta_endpoint <flg>
646eb284becSJed Brown 
647eb284becSJed Brown   Level: Intermediate
648eb284becSJed Brown 
649eb284becSJed Brown .seealso: TSTHETA, TSCN
650eb284becSJed Brown @*/
651eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
652eb284becSJed Brown {
653eb284becSJed Brown   PetscErrorCode ierr;
654eb284becSJed Brown 
655eb284becSJed Brown   PetscFunctionBegin;
656eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
657eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
658eb284becSJed Brown   PetscFunctionReturn(0);
659eb284becSJed Brown }
660eb284becSJed Brown 
661f33bbcb6SJed Brown /*
662f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
663f33bbcb6SJed Brown  * The creation functions for these specializations are below.
664f33bbcb6SJed Brown  */
665f33bbcb6SJed Brown 
666f33bbcb6SJed Brown #undef __FUNCT__
667f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
668f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
669f33bbcb6SJed Brown {
670d52bd9f3SBarry Smith   PetscErrorCode ierr;
671d52bd9f3SBarry Smith 
672f33bbcb6SJed Brown   PetscFunctionBegin;
673d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
674f33bbcb6SJed Brown   PetscFunctionReturn(0);
675f33bbcb6SJed Brown }
676f33bbcb6SJed Brown 
677f33bbcb6SJed Brown /*MC
678f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
679f33bbcb6SJed Brown 
680f33bbcb6SJed Brown   Level: beginner
681f33bbcb6SJed Brown 
6824eb428fdSBarry Smith   Notes:
683c7eb6c99SShri Abhyankar   TSBEULER is equivalent to TSTHETA with Theta=1.0
6844eb428fdSBarry Smith 
6854eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
6864eb428fdSBarry Smith 
687f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
688f33bbcb6SJed Brown 
689f33bbcb6SJed Brown M*/
690f33bbcb6SJed Brown #undef __FUNCT__
691f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
6928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
693f33bbcb6SJed Brown {
694f33bbcb6SJed Brown   PetscErrorCode ierr;
695f33bbcb6SJed Brown 
696f33bbcb6SJed Brown   PetscFunctionBegin;
697f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
698f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
699f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
700f33bbcb6SJed Brown   PetscFunctionReturn(0);
701f33bbcb6SJed Brown }
702f33bbcb6SJed Brown 
703f33bbcb6SJed Brown #undef __FUNCT__
704f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
705f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
706f33bbcb6SJed Brown {
707d52bd9f3SBarry Smith   PetscErrorCode ierr;
708d52bd9f3SBarry Smith 
709f33bbcb6SJed Brown   PetscFunctionBegin;
710d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
711f33bbcb6SJed Brown   PetscFunctionReturn(0);
712f33bbcb6SJed Brown }
713f33bbcb6SJed Brown 
714f33bbcb6SJed Brown /*MC
715f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
716f33bbcb6SJed Brown 
717f33bbcb6SJed Brown   Level: beginner
718f33bbcb6SJed Brown 
719f33bbcb6SJed Brown   Notes:
7207cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
7217cf5af47SJed Brown 
7227cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
723f33bbcb6SJed Brown 
724f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
725f33bbcb6SJed Brown 
726f33bbcb6SJed Brown M*/
727f33bbcb6SJed Brown #undef __FUNCT__
728f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
7298cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
730f33bbcb6SJed Brown {
731f33bbcb6SJed Brown   PetscErrorCode ierr;
732f33bbcb6SJed Brown 
733f33bbcb6SJed Brown   PetscFunctionBegin;
734f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
735f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
736eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
737f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
738f33bbcb6SJed Brown   PetscFunctionReturn(0);
739f33bbcb6SJed Brown }
740