xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 4e17fbf5bafa47185df4021e6c4f5fb738d34ecb)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
57445fe48SJed Brown #include <petscsnesfas.h>
6316643e7SJed Brown 
7316643e7SJed Brown typedef struct {
8316643e7SJed Brown   Vec       X,Xdot;                   /* Storage for one stage */
9eb284becSJed Brown   Vec       affine;                   /* Affine vector needed for residual at beginning of step */
10ace3abfcSBarry Smith   PetscBool extrapolate;
11eb284becSJed Brown   PetscBool endpoint;
12316643e7SJed Brown   PetscReal Theta;
13316643e7SJed Brown   PetscReal shift;
14316643e7SJed Brown   PetscReal stage_time;
15316643e7SJed Brown } TS_Theta;
16316643e7SJed Brown 
17316643e7SJed Brown #undef __FUNCT__
187445fe48SJed Brown #define __FUNCT__ "TSThetaGetX0AndXdot"
197445fe48SJed Brown static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
207445fe48SJed Brown {
217445fe48SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
227445fe48SJed Brown   PetscErrorCode ierr;
237445fe48SJed Brown 
247445fe48SJed Brown   PetscFunctionBegin;
257445fe48SJed Brown   if (X0) {
267445fe48SJed Brown     if (dm && dm != ts->dm) {
270d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
287445fe48SJed Brown     } else *X0 = ts->vec_sol;
297445fe48SJed Brown   }
307445fe48SJed Brown   if (Xdot) {
317445fe48SJed Brown     if (dm && dm != ts->dm) {
320d0b770aSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
337445fe48SJed Brown     } else *Xdot = th->Xdot;
347445fe48SJed Brown   }
357445fe48SJed Brown   PetscFunctionReturn(0);
367445fe48SJed Brown }
377445fe48SJed Brown 
380d0b770aSPeter Brune 
390d0b770aSPeter Brune #undef __FUNCT__
400d0b770aSPeter Brune #define __FUNCT__ "TSThetaRestoreX0AndXdot"
410d0b770aSPeter Brune static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot)
420d0b770aSPeter Brune {
430d0b770aSPeter Brune   PetscErrorCode ierr;
440d0b770aSPeter Brune 
450d0b770aSPeter Brune   PetscFunctionBegin;
460d0b770aSPeter Brune   if (X0) {
470d0b770aSPeter Brune     if (dm && dm != ts->dm) {
480d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr);
490d0b770aSPeter Brune     }
500d0b770aSPeter Brune   }
510d0b770aSPeter Brune   if (Xdot) {
520d0b770aSPeter Brune     if (dm && dm != ts->dm) {
530d0b770aSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr);
540d0b770aSPeter Brune     }
550d0b770aSPeter Brune   }
560d0b770aSPeter Brune   PetscFunctionReturn(0);
570d0b770aSPeter Brune }
580d0b770aSPeter Brune 
597445fe48SJed Brown #undef __FUNCT__
607445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
617445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
627445fe48SJed Brown {
637445fe48SJed Brown 
647445fe48SJed Brown   PetscFunctionBegin;
657445fe48SJed Brown   PetscFunctionReturn(0);
667445fe48SJed Brown }
677445fe48SJed Brown 
687445fe48SJed Brown #undef __FUNCT__
697445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
707445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
717445fe48SJed Brown {
727445fe48SJed Brown   TS ts = (TS)ctx;
737445fe48SJed Brown   PetscErrorCode ierr;
747445fe48SJed Brown   Vec X0,Xdot,X0_c,Xdot_c;
757445fe48SJed Brown 
767445fe48SJed Brown   PetscFunctionBegin;
777445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
787445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
797445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
807445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
817445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
827445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
830d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
840d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
857445fe48SJed Brown   PetscFunctionReturn(0);
867445fe48SJed Brown }
877445fe48SJed Brown 
887445fe48SJed Brown #undef __FUNCT__
89316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
90193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
91316643e7SJed Brown {
92316643e7SJed Brown   TS_Theta            *th = (TS_Theta*)ts->data;
93*4e17fbf5SBarry Smith   PetscInt            its,lits,reject;
94cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
95f1b97656SJed Brown   SNESConvergedReason snesreason;
962b5a38e1SLisandro Dalcin   PetscErrorCode      ierr;
97316643e7SJed Brown 
98316643e7SJed Brown   PetscFunctionBegin;
99*4e17fbf5SBarry Smith   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
100cdbf8f93SLisandro Dalcin     next_time_step = ts->time_step;
101eb284becSJed Brown     th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
102316643e7SJed Brown     th->shift = 1./(th->Theta*ts->time_step);
103b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
104b8123daeSJed Brown     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
105316643e7SJed Brown 
106eb284becSJed Brown     if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
107eb284becSJed Brown       ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
108eb284becSJed Brown       if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
109eb284becSJed Brown       ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
110eb284becSJed Brown       ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
111eb284becSJed Brown     }
112ace68cafSJed Brown     if (th->extrapolate) {
1132b5a38e1SLisandro Dalcin       ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
114ace68cafSJed Brown     } else {
1152b5a38e1SLisandro Dalcin       ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
116ace68cafSJed Brown     }
117eb284becSJed Brown     ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
118316643e7SJed Brown     ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
119316643e7SJed Brown     ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
120f1b97656SJed Brown     ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
1215ef26d82SJed Brown     ts->snes_its += its; ts->ksp_its += lits;
122*4e17fbf5SBarry Smith     if (snesreason > 0) break;
123*4e17fbf5SBarry Smith     ierr = PetscInfo3(ts,"Step=%D, Cutting time-step from %g to %g\n",ts->steps,(double)ts->time_step,(double).5*ts->time_step);CHKERRQ(ierr);
124*4e17fbf5SBarry Smith     ts->time_step = .5*ts->time_step;
125*4e17fbf5SBarry Smith   }
126f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
127f1b97656SJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
128f1b97656SJed Brown     ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
129f1b97656SJed Brown     PetscFunctionReturn(0);
130f1b97656SJed Brown   }
131eb284becSJed Brown   if (th->endpoint) {
132eb284becSJed Brown     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
133eb284becSJed Brown   } else {
1342b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
1352b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
136eb284becSJed Brown   }
1372b5a38e1SLisandro Dalcin   ts->ptime += ts->time_step;
138cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
139316643e7SJed Brown   ts->steps++;
140316643e7SJed Brown   PetscFunctionReturn(0);
141316643e7SJed Brown }
142316643e7SJed Brown 
143cd652676SJed Brown #undef __FUNCT__
144cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
145cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
146cd652676SJed Brown {
147cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1485a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
149cd652676SJed Brown   PetscErrorCode ierr;
150cd652676SJed Brown 
151cd652676SJed Brown   PetscFunctionBegin;
152a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
1535a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
1545a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
155cd652676SJed Brown   PetscFunctionReturn(0);
156cd652676SJed Brown }
157cd652676SJed Brown 
158316643e7SJed Brown /*------------------------------------------------------------*/
159316643e7SJed Brown #undef __FUNCT__
160277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
161277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
162316643e7SJed Brown {
163316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
164316643e7SJed Brown   PetscErrorCode  ierr;
165316643e7SJed Brown 
166316643e7SJed Brown   PetscFunctionBegin;
1676bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
1686bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
169eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
170277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
171277b19d0SLisandro Dalcin }
172277b19d0SLisandro Dalcin 
173277b19d0SLisandro Dalcin #undef __FUNCT__
174277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
175277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
176277b19d0SLisandro Dalcin {
177277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
178277b19d0SLisandro Dalcin 
179277b19d0SLisandro Dalcin   PetscFunctionBegin;
180277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
181277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
182335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
183335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
18426f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
185eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
186316643e7SJed Brown   PetscFunctionReturn(0);
187316643e7SJed Brown }
188316643e7SJed Brown 
189316643e7SJed Brown /*
190316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
1912b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
192316643e7SJed Brown */
193316643e7SJed Brown #undef __FUNCT__
1940f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
1950f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
196316643e7SJed Brown {
197316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
198316643e7SJed Brown   PetscErrorCode ierr;
1997445fe48SJed Brown   Vec            X0,Xdot;
2007445fe48SJed Brown   DM             dm,dmsave;
201316643e7SJed Brown 
202316643e7SJed Brown   PetscFunctionBegin;
2037445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2045a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
2057445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
2067445fe48SJed Brown   ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr);
2077445fe48SJed Brown 
2087445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
2097445fe48SJed Brown   dmsave = ts->dm;
2107445fe48SJed Brown   ts->dm = dm;
2117445fe48SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
2127445fe48SJed Brown   ts->dm = dmsave;
2130d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
214316643e7SJed Brown   PetscFunctionReturn(0);
215316643e7SJed Brown }
216316643e7SJed Brown 
217316643e7SJed Brown #undef __FUNCT__
2180f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
2190f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
220316643e7SJed Brown {
221316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
222316643e7SJed Brown   PetscErrorCode ierr;
2237445fe48SJed Brown   Vec            Xdot;
2247445fe48SJed Brown   DM             dm,dmsave;
225316643e7SJed Brown 
226316643e7SJed Brown   PetscFunctionBegin;
2277445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2287445fe48SJed Brown 
2290f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
2307445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
2317445fe48SJed Brown 
2327445fe48SJed Brown   dmsave = ts->dm;
2337445fe48SJed Brown   ts->dm = dm;
2347445fe48SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
2357445fe48SJed Brown   ts->dm = dmsave;
2360d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
237316643e7SJed Brown   PetscFunctionReturn(0);
238316643e7SJed Brown }
239316643e7SJed Brown 
240316643e7SJed Brown #undef __FUNCT__
241316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
242316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
243316643e7SJed Brown {
244316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
245316643e7SJed Brown   PetscErrorCode ierr;
2467445fe48SJed Brown   SNES           snes;
2477445fe48SJed Brown   DM             dm;
248316643e7SJed Brown 
249316643e7SJed Brown   PetscFunctionBegin;
250316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
251316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
2527445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2537445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2547445fe48SJed Brown   if (dm) {
2557445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
2567445fe48SJed Brown   }
257316643e7SJed Brown   PetscFunctionReturn(0);
258316643e7SJed Brown }
259316643e7SJed Brown /*------------------------------------------------------------*/
260316643e7SJed Brown 
261316643e7SJed Brown #undef __FUNCT__
262316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
263316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
264316643e7SJed Brown {
265316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
266316643e7SJed Brown   PetscErrorCode ierr;
267316643e7SJed Brown 
268316643e7SJed Brown   PetscFunctionBegin;
269d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
270316643e7SJed Brown   {
271316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
272acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
273eb284becSJed Brown     ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,PETSC_NULL);CHKERRQ(ierr);
274d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
275316643e7SJed Brown   }
276316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
277316643e7SJed Brown   PetscFunctionReturn(0);
278316643e7SJed Brown }
279316643e7SJed Brown 
280316643e7SJed Brown #undef __FUNCT__
281316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
282316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
283316643e7SJed Brown {
284316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
285ace3abfcSBarry Smith   PetscBool       iascii;
286316643e7SJed Brown   PetscErrorCode  ierr;
287316643e7SJed Brown 
288316643e7SJed Brown   PetscFunctionBegin;
289251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
290316643e7SJed Brown   if (iascii) {
291316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
292ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
293316643e7SJed Brown   }
294d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
295316643e7SJed Brown   PetscFunctionReturn(0);
296316643e7SJed Brown }
297316643e7SJed Brown 
2980de4c49aSJed Brown EXTERN_C_BEGIN
2990de4c49aSJed Brown #undef __FUNCT__
3000de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
3017087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
3020de4c49aSJed Brown {
3030de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
3040de4c49aSJed Brown 
3050de4c49aSJed Brown   PetscFunctionBegin;
3060de4c49aSJed Brown   *theta = th->Theta;
3070de4c49aSJed Brown   PetscFunctionReturn(0);
3080de4c49aSJed Brown }
3090de4c49aSJed Brown 
3100de4c49aSJed Brown #undef __FUNCT__
3110de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
3127087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
3130de4c49aSJed Brown {
3140de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
3150de4c49aSJed Brown 
3160de4c49aSJed Brown   PetscFunctionBegin;
317e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
3180de4c49aSJed Brown   th->Theta = theta;
3190de4c49aSJed Brown   PetscFunctionReturn(0);
3200de4c49aSJed Brown }
321eb284becSJed Brown 
322eb284becSJed Brown #undef __FUNCT__
32378e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
32426f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
32526f2ff8fSLisandro Dalcin {
32626f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
32726f2ff8fSLisandro Dalcin 
32826f2ff8fSLisandro Dalcin   PetscFunctionBegin;
32926f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
33026f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
33126f2ff8fSLisandro Dalcin }
33226f2ff8fSLisandro Dalcin 
33326f2ff8fSLisandro Dalcin #undef __FUNCT__
33426f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
335eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
336eb284becSJed Brown {
337eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
338eb284becSJed Brown 
339eb284becSJed Brown   PetscFunctionBegin;
340eb284becSJed Brown   th->endpoint = flg;
341eb284becSJed Brown   PetscFunctionReturn(0);
342eb284becSJed Brown }
3430de4c49aSJed Brown EXTERN_C_END
3440de4c49aSJed Brown 
345316643e7SJed Brown /* ------------------------------------------------------------ */
346316643e7SJed Brown /*MC
34796f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
348316643e7SJed Brown 
349316643e7SJed Brown    Level: beginner
350316643e7SJed Brown 
3514eb428fdSBarry Smith    Options Database:
3524eb428fdSBarry Smith       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1);  Theta = 1.0 (
3534eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
3544eb428fdSBarry Smith       -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method
3554eb428fdSBarry Smith 
356eb284becSJed Brown    Notes:
3574eb428fdSBarry Smith $  -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER)
3584eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
3594eb428fdSBarry Smith 
3604eb428fdSBarry Smith 
3614eb428fdSBarry Smith 
362eb284becSJed Brown    This method can be applied to DAE.
363eb284becSJed Brown 
364eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
365eb284becSJed Brown 
366eb284becSJed Brown .vb
367eb284becSJed Brown   Theta | Theta
368eb284becSJed Brown   -------------
369eb284becSJed Brown         |  1
370eb284becSJed Brown .ve
371eb284becSJed Brown 
372eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
373eb284becSJed Brown 
374eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
375eb284becSJed Brown 
376eb284becSJed Brown .vb
377eb284becSJed Brown   0 | 0         0
378eb284becSJed Brown   1 | 1-Theta   Theta
379eb284becSJed Brown   -------------------
380eb284becSJed Brown     | 1-Theta   Theta
381eb284becSJed Brown .ve
382eb284becSJed Brown 
383eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
384eb284becSJed Brown 
385eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
386eb284becSJed Brown 
387eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
388eb284becSJed Brown 
3894eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
390eb284becSJed Brown 
391eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
392316643e7SJed Brown 
393316643e7SJed Brown M*/
394316643e7SJed Brown EXTERN_C_BEGIN
395316643e7SJed Brown #undef __FUNCT__
396316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
3977087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
398316643e7SJed Brown {
399316643e7SJed Brown   TS_Theta       *th;
400316643e7SJed Brown   PetscErrorCode ierr;
401316643e7SJed Brown 
402316643e7SJed Brown   PetscFunctionBegin;
403277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
404316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
405316643e7SJed Brown   ts->ops->view           = TSView_Theta;
406316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
407316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
408cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
409316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
4100f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
4110f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
412316643e7SJed Brown 
413316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
414316643e7SJed Brown   ts->data = (void*)th;
415316643e7SJed Brown 
4166f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
417316643e7SJed Brown   th->Theta       = 0.5;
418316643e7SJed Brown 
4190de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
4200de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
42126f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
422eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
423316643e7SJed Brown   PetscFunctionReturn(0);
424316643e7SJed Brown }
425316643e7SJed Brown EXTERN_C_END
4260de4c49aSJed Brown 
4270de4c49aSJed Brown #undef __FUNCT__
4280de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
4290de4c49aSJed Brown /*@
4300de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
4310de4c49aSJed Brown 
4320de4c49aSJed Brown   Not Collective
4330de4c49aSJed Brown 
4340de4c49aSJed Brown   Input Parameter:
4350de4c49aSJed Brown .  ts - timestepping context
4360de4c49aSJed Brown 
4370de4c49aSJed Brown   Output Parameter:
4380de4c49aSJed Brown .  theta - stage abscissa
4390de4c49aSJed Brown 
4400de4c49aSJed Brown   Note:
4410de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
4420de4c49aSJed Brown 
4430de4c49aSJed Brown   Level: Advanced
4440de4c49aSJed Brown 
4450de4c49aSJed Brown .seealso: TSThetaSetTheta()
4460de4c49aSJed Brown @*/
4477087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
4480de4c49aSJed Brown {
4494ac538c5SBarry Smith   PetscErrorCode ierr;
4500de4c49aSJed Brown 
4510de4c49aSJed Brown   PetscFunctionBegin;
452afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4530de4c49aSJed Brown   PetscValidPointer(theta,2);
4544ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
4550de4c49aSJed Brown   PetscFunctionReturn(0);
4560de4c49aSJed Brown }
4570de4c49aSJed Brown 
4580de4c49aSJed Brown #undef __FUNCT__
4590de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
4600de4c49aSJed Brown /*@
4610de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
4620de4c49aSJed Brown 
4630de4c49aSJed Brown   Not Collective
4640de4c49aSJed Brown 
4650de4c49aSJed Brown   Input Parameter:
4660de4c49aSJed Brown +  ts - timestepping context
4670de4c49aSJed Brown -  theta - stage abscissa
4680de4c49aSJed Brown 
4690de4c49aSJed Brown   Options Database:
4700de4c49aSJed Brown .  -ts_theta_theta <theta>
4710de4c49aSJed Brown 
4720de4c49aSJed Brown   Level: Intermediate
4730de4c49aSJed Brown 
4740de4c49aSJed Brown .seealso: TSThetaGetTheta()
4750de4c49aSJed Brown @*/
4767087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
4770de4c49aSJed Brown {
4784ac538c5SBarry Smith   PetscErrorCode ierr;
4790de4c49aSJed Brown 
4800de4c49aSJed Brown   PetscFunctionBegin;
481afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4824ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
4830de4c49aSJed Brown   PetscFunctionReturn(0);
4840de4c49aSJed Brown }
485f33bbcb6SJed Brown 
486eb284becSJed Brown #undef __FUNCT__
48726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
48826f2ff8fSLisandro Dalcin /*@
48926f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
49026f2ff8fSLisandro Dalcin 
49126f2ff8fSLisandro Dalcin   Not Collective
49226f2ff8fSLisandro Dalcin 
49326f2ff8fSLisandro Dalcin   Input Parameter:
49426f2ff8fSLisandro Dalcin .  ts - timestepping context
49526f2ff8fSLisandro Dalcin 
49626f2ff8fSLisandro Dalcin   Output Parameter:
49726f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
49826f2ff8fSLisandro Dalcin 
49926f2ff8fSLisandro Dalcin   Level: Advanced
50026f2ff8fSLisandro Dalcin 
50126f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
50226f2ff8fSLisandro Dalcin @*/
50326f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
50426f2ff8fSLisandro Dalcin {
50526f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
50626f2ff8fSLisandro Dalcin 
50726f2ff8fSLisandro Dalcin   PetscFunctionBegin;
50826f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
50926f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
51026f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
51126f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
51226f2ff8fSLisandro Dalcin }
51326f2ff8fSLisandro Dalcin 
51426f2ff8fSLisandro Dalcin #undef __FUNCT__
515eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
516eb284becSJed Brown /*@
517eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
518eb284becSJed Brown 
519eb284becSJed Brown   Not Collective
520eb284becSJed Brown 
521eb284becSJed Brown   Input Parameter:
522eb284becSJed Brown +  ts - timestepping context
523eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
524eb284becSJed Brown 
525eb284becSJed Brown   Options Database:
526eb284becSJed Brown .  -ts_theta_endpoint <flg>
527eb284becSJed Brown 
528eb284becSJed Brown   Level: Intermediate
529eb284becSJed Brown 
530eb284becSJed Brown .seealso: TSTHETA, TSCN
531eb284becSJed Brown @*/
532eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
533eb284becSJed Brown {
534eb284becSJed Brown   PetscErrorCode ierr;
535eb284becSJed Brown 
536eb284becSJed Brown   PetscFunctionBegin;
537eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
538eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
539eb284becSJed Brown   PetscFunctionReturn(0);
540eb284becSJed Brown }
541eb284becSJed Brown 
542f33bbcb6SJed Brown /*
543f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
544f33bbcb6SJed Brown  * The creation functions for these specializations are below.
545f33bbcb6SJed Brown  */
546f33bbcb6SJed Brown 
547f33bbcb6SJed Brown #undef __FUNCT__
548f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
549f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
550f33bbcb6SJed Brown {
551d52bd9f3SBarry Smith   PetscErrorCode ierr;
552d52bd9f3SBarry Smith 
553f33bbcb6SJed Brown   PetscFunctionBegin;
554d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
555f33bbcb6SJed Brown   PetscFunctionReturn(0);
556f33bbcb6SJed Brown }
557f33bbcb6SJed Brown 
558f33bbcb6SJed Brown /*MC
559f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
560f33bbcb6SJed Brown 
561f33bbcb6SJed Brown   Level: beginner
562f33bbcb6SJed Brown 
5634eb428fdSBarry Smith   Notes:
5644eb428fdSBarry Smith   TSCN is equivalent to TSTHETA with Theta=1.0
5654eb428fdSBarry Smith 
5664eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
5674eb428fdSBarry Smith 
568f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
569f33bbcb6SJed Brown 
570f33bbcb6SJed Brown M*/
571f33bbcb6SJed Brown EXTERN_C_BEGIN
572f33bbcb6SJed Brown #undef __FUNCT__
573f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
574f33bbcb6SJed Brown PetscErrorCode  TSCreate_BEuler(TS ts)
575f33bbcb6SJed Brown {
576f33bbcb6SJed Brown   PetscErrorCode ierr;
577f33bbcb6SJed Brown 
578f33bbcb6SJed Brown   PetscFunctionBegin;
579f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
580f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
581f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
582f33bbcb6SJed Brown   PetscFunctionReturn(0);
583f33bbcb6SJed Brown }
584f33bbcb6SJed Brown EXTERN_C_END
585f33bbcb6SJed Brown 
586f33bbcb6SJed Brown #undef __FUNCT__
587f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
588f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
589f33bbcb6SJed Brown {
590d52bd9f3SBarry Smith   PetscErrorCode ierr;
591d52bd9f3SBarry Smith 
592f33bbcb6SJed Brown   PetscFunctionBegin;
593d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
594f33bbcb6SJed Brown   PetscFunctionReturn(0);
595f33bbcb6SJed Brown }
596f33bbcb6SJed Brown 
597f33bbcb6SJed Brown /*MC
598f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
599f33bbcb6SJed Brown 
600f33bbcb6SJed Brown   Level: beginner
601f33bbcb6SJed Brown 
602f33bbcb6SJed Brown   Notes:
6037cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
6047cf5af47SJed Brown 
6057cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
606f33bbcb6SJed Brown 
607f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
608f33bbcb6SJed Brown 
609f33bbcb6SJed Brown M*/
610f33bbcb6SJed Brown EXTERN_C_BEGIN
611f33bbcb6SJed Brown #undef __FUNCT__
612f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
613f33bbcb6SJed Brown PetscErrorCode  TSCreate_CN(TS ts)
614f33bbcb6SJed Brown {
615f33bbcb6SJed Brown   PetscErrorCode ierr;
616f33bbcb6SJed Brown 
617f33bbcb6SJed Brown   PetscFunctionBegin;
618f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
619f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
620eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
621f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
622f33bbcb6SJed Brown   PetscFunctionReturn(0);
623f33bbcb6SJed Brown }
624f33bbcb6SJed Brown EXTERN_C_END
625