xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision 4eb428fdc6f36238aef3665003db33fd4cbae8e2)
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;
93b70ae86eSJed Brown   PetscInt            its,lits;
94cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
95f1b97656SJed Brown   SNESConvergedReason snesreason;
962b5a38e1SLisandro Dalcin   PetscErrorCode      ierr;
97316643e7SJed Brown 
98316643e7SJed Brown   PetscFunctionBegin;
99cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
100eb284becSJed Brown   th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
101316643e7SJed Brown   th->shift = 1./(th->Theta*ts->time_step);
102316643e7SJed Brown 
103eb284becSJed Brown   if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
104eb284becSJed Brown     ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
105eb284becSJed Brown     if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
106eb284becSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
107eb284becSJed Brown     ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
108eb284becSJed Brown   }
109ace68cafSJed Brown   if (th->extrapolate) {
1102b5a38e1SLisandro Dalcin     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
111ace68cafSJed Brown   } else {
1122b5a38e1SLisandro Dalcin     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
113ace68cafSJed Brown   }
114eb284becSJed Brown   ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
115316643e7SJed Brown   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
116316643e7SJed Brown   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
117f1b97656SJed Brown   ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
1185ef26d82SJed Brown   ts->snes_its += its; ts->ksp_its += lits;
119f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
120f1b97656SJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
121f1b97656SJed 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);
122f1b97656SJed Brown     PetscFunctionReturn(0);
123f1b97656SJed Brown   }
124eb284becSJed Brown   if (th->endpoint) {
125eb284becSJed Brown     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
126eb284becSJed Brown   } else {
1272b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
1282b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
129eb284becSJed Brown   }
1302b5a38e1SLisandro Dalcin   ts->ptime += ts->time_step;
131cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
132316643e7SJed Brown   ts->steps++;
133316643e7SJed Brown   PetscFunctionReturn(0);
134316643e7SJed Brown }
135316643e7SJed Brown 
136cd652676SJed Brown #undef __FUNCT__
137cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
138cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
139cd652676SJed Brown {
140cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1415a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
142cd652676SJed Brown   PetscErrorCode ierr;
143cd652676SJed Brown 
144cd652676SJed Brown   PetscFunctionBegin;
145a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
1465a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
1475a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
148cd652676SJed Brown   PetscFunctionReturn(0);
149cd652676SJed Brown }
150cd652676SJed Brown 
151316643e7SJed Brown /*------------------------------------------------------------*/
152316643e7SJed Brown #undef __FUNCT__
153277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
154277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
155316643e7SJed Brown {
156316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
157316643e7SJed Brown   PetscErrorCode  ierr;
158316643e7SJed Brown 
159316643e7SJed Brown   PetscFunctionBegin;
1606bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
1616bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
162eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
163277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
164277b19d0SLisandro Dalcin }
165277b19d0SLisandro Dalcin 
166277b19d0SLisandro Dalcin #undef __FUNCT__
167277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
168277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
169277b19d0SLisandro Dalcin {
170277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
171277b19d0SLisandro Dalcin 
172277b19d0SLisandro Dalcin   PetscFunctionBegin;
173277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
174277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
175335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
176335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
17726f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
178eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
179316643e7SJed Brown   PetscFunctionReturn(0);
180316643e7SJed Brown }
181316643e7SJed Brown 
182316643e7SJed Brown /*
183316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
1842b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
185316643e7SJed Brown */
186316643e7SJed Brown #undef __FUNCT__
1870f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
1880f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
189316643e7SJed Brown {
190316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
191316643e7SJed Brown   PetscErrorCode ierr;
1927445fe48SJed Brown   Vec            X0,Xdot;
1937445fe48SJed Brown   DM             dm,dmsave;
194316643e7SJed Brown 
195316643e7SJed Brown   PetscFunctionBegin;
1967445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1975a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
1987445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
1997445fe48SJed Brown   ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr);
2007445fe48SJed Brown 
2017445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
2027445fe48SJed Brown   dmsave = ts->dm;
2037445fe48SJed Brown   ts->dm = dm;
2047445fe48SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
2057445fe48SJed Brown   ts->dm = dmsave;
2060d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
207316643e7SJed Brown   PetscFunctionReturn(0);
208316643e7SJed Brown }
209316643e7SJed Brown 
210316643e7SJed Brown #undef __FUNCT__
2110f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
2120f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
213316643e7SJed Brown {
214316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
215316643e7SJed Brown   PetscErrorCode ierr;
2167445fe48SJed Brown   Vec            Xdot;
2177445fe48SJed Brown   DM             dm,dmsave;
218316643e7SJed Brown 
219316643e7SJed Brown   PetscFunctionBegin;
2207445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2217445fe48SJed Brown 
2220f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
2237445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
2247445fe48SJed Brown 
2257445fe48SJed Brown   dmsave = ts->dm;
2267445fe48SJed Brown   ts->dm = dm;
2277445fe48SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
2287445fe48SJed Brown   ts->dm = dmsave;
2290d0b770aSPeter Brune   ierr = TSThetaRestoreX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
230316643e7SJed Brown   PetscFunctionReturn(0);
231316643e7SJed Brown }
232316643e7SJed Brown 
233316643e7SJed Brown #undef __FUNCT__
234316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
235316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
236316643e7SJed Brown {
237316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
238316643e7SJed Brown   PetscErrorCode ierr;
2397445fe48SJed Brown   SNES           snes;
2407445fe48SJed Brown   DM             dm;
241316643e7SJed Brown 
242316643e7SJed Brown   PetscFunctionBegin;
243316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
244316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
2457445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2467445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2477445fe48SJed Brown   if (dm) {
2487445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
2497445fe48SJed Brown   }
250316643e7SJed Brown   PetscFunctionReturn(0);
251316643e7SJed Brown }
252316643e7SJed Brown /*------------------------------------------------------------*/
253316643e7SJed Brown 
254316643e7SJed Brown #undef __FUNCT__
255316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
256316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
257316643e7SJed Brown {
258316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
259316643e7SJed Brown   PetscErrorCode ierr;
260316643e7SJed Brown 
261316643e7SJed Brown   PetscFunctionBegin;
262d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
263316643e7SJed Brown   {
264316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
265acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
266eb284becSJed 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);
267d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
268316643e7SJed Brown   }
269316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
270316643e7SJed Brown   PetscFunctionReturn(0);
271316643e7SJed Brown }
272316643e7SJed Brown 
273316643e7SJed Brown #undef __FUNCT__
274316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
275316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
276316643e7SJed Brown {
277316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
278ace3abfcSBarry Smith   PetscBool       iascii;
279316643e7SJed Brown   PetscErrorCode  ierr;
280316643e7SJed Brown 
281316643e7SJed Brown   PetscFunctionBegin;
282251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
283316643e7SJed Brown   if (iascii) {
284316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
285ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
286316643e7SJed Brown   }
287d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
288316643e7SJed Brown   PetscFunctionReturn(0);
289316643e7SJed Brown }
290316643e7SJed Brown 
2910de4c49aSJed Brown EXTERN_C_BEGIN
2920de4c49aSJed Brown #undef __FUNCT__
2930de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
2947087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
2950de4c49aSJed Brown {
2960de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2970de4c49aSJed Brown 
2980de4c49aSJed Brown   PetscFunctionBegin;
2990de4c49aSJed Brown   *theta = th->Theta;
3000de4c49aSJed Brown   PetscFunctionReturn(0);
3010de4c49aSJed Brown }
3020de4c49aSJed Brown 
3030de4c49aSJed Brown #undef __FUNCT__
3040de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
3057087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
3060de4c49aSJed Brown {
3070de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
3080de4c49aSJed Brown 
3090de4c49aSJed Brown   PetscFunctionBegin;
310e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
3110de4c49aSJed Brown   th->Theta = theta;
3120de4c49aSJed Brown   PetscFunctionReturn(0);
3130de4c49aSJed Brown }
314eb284becSJed Brown 
315eb284becSJed Brown #undef __FUNCT__
31678e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
31726f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
31826f2ff8fSLisandro Dalcin {
31926f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
32026f2ff8fSLisandro Dalcin 
32126f2ff8fSLisandro Dalcin   PetscFunctionBegin;
32226f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
32326f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
32426f2ff8fSLisandro Dalcin }
32526f2ff8fSLisandro Dalcin 
32626f2ff8fSLisandro Dalcin #undef __FUNCT__
32726f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
328eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
329eb284becSJed Brown {
330eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
331eb284becSJed Brown 
332eb284becSJed Brown   PetscFunctionBegin;
333eb284becSJed Brown   th->endpoint = flg;
334eb284becSJed Brown   PetscFunctionReturn(0);
335eb284becSJed Brown }
3360de4c49aSJed Brown EXTERN_C_END
3370de4c49aSJed Brown 
338316643e7SJed Brown /* ------------------------------------------------------------ */
339316643e7SJed Brown /*MC
34096f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
341316643e7SJed Brown 
342316643e7SJed Brown    Level: beginner
343316643e7SJed Brown 
344*4eb428fdSBarry Smith    Options Database:
345*4eb428fdSBarry Smith       -ts_theta_theta <Theta> - Location of stage (0<Theta<=1);  Theta = 1.0 (
346*4eb428fdSBarry Smith       -ts_theta_extrapolate <flg> Extrapolate stage solution from previous solution (sometimes unstable)
347*4eb428fdSBarry Smith       -ts_theta_endpoint <flag> - Use the endpoint instead of midpoint form of the Theta method
348*4eb428fdSBarry Smith 
349eb284becSJed Brown    Notes:
350*4eb428fdSBarry Smith $  -ts_type theta -ts_theta 1.0 corresponds to backward Euler (TSBEULER)
351*4eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
352*4eb428fdSBarry Smith 
353*4eb428fdSBarry Smith 
354*4eb428fdSBarry Smith 
355eb284becSJed Brown    This method can be applied to DAE.
356eb284becSJed Brown 
357eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
358eb284becSJed Brown 
359eb284becSJed Brown .vb
360eb284becSJed Brown   Theta | Theta
361eb284becSJed Brown   -------------
362eb284becSJed Brown         |  1
363eb284becSJed Brown .ve
364eb284becSJed Brown 
365eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
366eb284becSJed Brown 
367eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
368eb284becSJed Brown 
369eb284becSJed Brown .vb
370eb284becSJed Brown   0 | 0         0
371eb284becSJed Brown   1 | 1-Theta   Theta
372eb284becSJed Brown   -------------------
373eb284becSJed Brown     | 1-Theta   Theta
374eb284becSJed Brown .ve
375eb284becSJed Brown 
376eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
377eb284becSJed Brown 
378eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
379eb284becSJed Brown 
380eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
381eb284becSJed Brown 
382*4eb428fdSBarry Smith    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)
383eb284becSJed Brown 
384eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
385316643e7SJed Brown 
386316643e7SJed Brown M*/
387316643e7SJed Brown EXTERN_C_BEGIN
388316643e7SJed Brown #undef __FUNCT__
389316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
3907087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
391316643e7SJed Brown {
392316643e7SJed Brown   TS_Theta       *th;
393316643e7SJed Brown   PetscErrorCode ierr;
394316643e7SJed Brown 
395316643e7SJed Brown   PetscFunctionBegin;
396277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
397316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
398316643e7SJed Brown   ts->ops->view           = TSView_Theta;
399316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
400316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
401cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
402316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
4030f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
4040f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
405316643e7SJed Brown 
406316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
407316643e7SJed Brown   ts->data = (void*)th;
408316643e7SJed Brown 
4096f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
410316643e7SJed Brown   th->Theta       = 0.5;
411316643e7SJed Brown 
4120de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
4130de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
41426f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
415eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
416316643e7SJed Brown   PetscFunctionReturn(0);
417316643e7SJed Brown }
418316643e7SJed Brown EXTERN_C_END
4190de4c49aSJed Brown 
4200de4c49aSJed Brown #undef __FUNCT__
4210de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
4220de4c49aSJed Brown /*@
4230de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
4240de4c49aSJed Brown 
4250de4c49aSJed Brown   Not Collective
4260de4c49aSJed Brown 
4270de4c49aSJed Brown   Input Parameter:
4280de4c49aSJed Brown .  ts - timestepping context
4290de4c49aSJed Brown 
4300de4c49aSJed Brown   Output Parameter:
4310de4c49aSJed Brown .  theta - stage abscissa
4320de4c49aSJed Brown 
4330de4c49aSJed Brown   Note:
4340de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
4350de4c49aSJed Brown 
4360de4c49aSJed Brown   Level: Advanced
4370de4c49aSJed Brown 
4380de4c49aSJed Brown .seealso: TSThetaSetTheta()
4390de4c49aSJed Brown @*/
4407087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
4410de4c49aSJed Brown {
4424ac538c5SBarry Smith   PetscErrorCode ierr;
4430de4c49aSJed Brown 
4440de4c49aSJed Brown   PetscFunctionBegin;
445afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4460de4c49aSJed Brown   PetscValidPointer(theta,2);
4474ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
4480de4c49aSJed Brown   PetscFunctionReturn(0);
4490de4c49aSJed Brown }
4500de4c49aSJed Brown 
4510de4c49aSJed Brown #undef __FUNCT__
4520de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
4530de4c49aSJed Brown /*@
4540de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
4550de4c49aSJed Brown 
4560de4c49aSJed Brown   Not Collective
4570de4c49aSJed Brown 
4580de4c49aSJed Brown   Input Parameter:
4590de4c49aSJed Brown +  ts - timestepping context
4600de4c49aSJed Brown -  theta - stage abscissa
4610de4c49aSJed Brown 
4620de4c49aSJed Brown   Options Database:
4630de4c49aSJed Brown .  -ts_theta_theta <theta>
4640de4c49aSJed Brown 
4650de4c49aSJed Brown   Level: Intermediate
4660de4c49aSJed Brown 
4670de4c49aSJed Brown .seealso: TSThetaGetTheta()
4680de4c49aSJed Brown @*/
4697087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
4700de4c49aSJed Brown {
4714ac538c5SBarry Smith   PetscErrorCode ierr;
4720de4c49aSJed Brown 
4730de4c49aSJed Brown   PetscFunctionBegin;
474afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4754ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
4760de4c49aSJed Brown   PetscFunctionReturn(0);
4770de4c49aSJed Brown }
478f33bbcb6SJed Brown 
479eb284becSJed Brown #undef __FUNCT__
48026f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
48126f2ff8fSLisandro Dalcin /*@
48226f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
48326f2ff8fSLisandro Dalcin 
48426f2ff8fSLisandro Dalcin   Not Collective
48526f2ff8fSLisandro Dalcin 
48626f2ff8fSLisandro Dalcin   Input Parameter:
48726f2ff8fSLisandro Dalcin .  ts - timestepping context
48826f2ff8fSLisandro Dalcin 
48926f2ff8fSLisandro Dalcin   Output Parameter:
49026f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
49126f2ff8fSLisandro Dalcin 
49226f2ff8fSLisandro Dalcin   Level: Advanced
49326f2ff8fSLisandro Dalcin 
49426f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
49526f2ff8fSLisandro Dalcin @*/
49626f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
49726f2ff8fSLisandro Dalcin {
49826f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
49926f2ff8fSLisandro Dalcin 
50026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
50126f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
50226f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
50326f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
50426f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
50526f2ff8fSLisandro Dalcin }
50626f2ff8fSLisandro Dalcin 
50726f2ff8fSLisandro Dalcin #undef __FUNCT__
508eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
509eb284becSJed Brown /*@
510eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
511eb284becSJed Brown 
512eb284becSJed Brown   Not Collective
513eb284becSJed Brown 
514eb284becSJed Brown   Input Parameter:
515eb284becSJed Brown +  ts - timestepping context
516eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
517eb284becSJed Brown 
518eb284becSJed Brown   Options Database:
519eb284becSJed Brown .  -ts_theta_endpoint <flg>
520eb284becSJed Brown 
521eb284becSJed Brown   Level: Intermediate
522eb284becSJed Brown 
523eb284becSJed Brown .seealso: TSTHETA, TSCN
524eb284becSJed Brown @*/
525eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
526eb284becSJed Brown {
527eb284becSJed Brown   PetscErrorCode ierr;
528eb284becSJed Brown 
529eb284becSJed Brown   PetscFunctionBegin;
530eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
531eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
532eb284becSJed Brown   PetscFunctionReturn(0);
533eb284becSJed Brown }
534eb284becSJed Brown 
535f33bbcb6SJed Brown /*
536f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
537f33bbcb6SJed Brown  * The creation functions for these specializations are below.
538f33bbcb6SJed Brown  */
539f33bbcb6SJed Brown 
540f33bbcb6SJed Brown #undef __FUNCT__
541f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
542f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
543f33bbcb6SJed Brown {
544d52bd9f3SBarry Smith   PetscErrorCode ierr;
545d52bd9f3SBarry Smith 
546f33bbcb6SJed Brown   PetscFunctionBegin;
547d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
548f33bbcb6SJed Brown   PetscFunctionReturn(0);
549f33bbcb6SJed Brown }
550f33bbcb6SJed Brown 
551f33bbcb6SJed Brown /*MC
552f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
553f33bbcb6SJed Brown 
554f33bbcb6SJed Brown   Level: beginner
555f33bbcb6SJed Brown 
556*4eb428fdSBarry Smith   Notes:
557*4eb428fdSBarry Smith   TSCN is equivalent to TSTHETA with Theta=1.0
558*4eb428fdSBarry Smith 
559*4eb428fdSBarry Smith $  -ts_type theta -ts_theta_theta 1.
560*4eb428fdSBarry Smith 
561f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
562f33bbcb6SJed Brown 
563f33bbcb6SJed Brown M*/
564f33bbcb6SJed Brown EXTERN_C_BEGIN
565f33bbcb6SJed Brown #undef __FUNCT__
566f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
567f33bbcb6SJed Brown PetscErrorCode  TSCreate_BEuler(TS ts)
568f33bbcb6SJed Brown {
569f33bbcb6SJed Brown   PetscErrorCode ierr;
570f33bbcb6SJed Brown 
571f33bbcb6SJed Brown   PetscFunctionBegin;
572f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
573f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
574f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
575f33bbcb6SJed Brown   PetscFunctionReturn(0);
576f33bbcb6SJed Brown }
577f33bbcb6SJed Brown EXTERN_C_END
578f33bbcb6SJed Brown 
579f33bbcb6SJed Brown #undef __FUNCT__
580f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
581f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
582f33bbcb6SJed Brown {
583d52bd9f3SBarry Smith   PetscErrorCode ierr;
584d52bd9f3SBarry Smith 
585f33bbcb6SJed Brown   PetscFunctionBegin;
586d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
587f33bbcb6SJed Brown   PetscFunctionReturn(0);
588f33bbcb6SJed Brown }
589f33bbcb6SJed Brown 
590f33bbcb6SJed Brown /*MC
591f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
592f33bbcb6SJed Brown 
593f33bbcb6SJed Brown   Level: beginner
594f33bbcb6SJed Brown 
595f33bbcb6SJed Brown   Notes:
5967cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
5977cf5af47SJed Brown 
5987cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
599f33bbcb6SJed Brown 
600f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
601f33bbcb6SJed Brown 
602f33bbcb6SJed Brown M*/
603f33bbcb6SJed Brown EXTERN_C_BEGIN
604f33bbcb6SJed Brown #undef __FUNCT__
605f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
606f33bbcb6SJed Brown PetscErrorCode  TSCreate_CN(TS ts)
607f33bbcb6SJed Brown {
608f33bbcb6SJed Brown   PetscErrorCode ierr;
609f33bbcb6SJed Brown 
610f33bbcb6SJed Brown   PetscFunctionBegin;
611f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
612f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
613eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
614f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
615f33bbcb6SJed Brown   PetscFunctionReturn(0);
616f33bbcb6SJed Brown }
617f33bbcb6SJed Brown EXTERN_C_END
618