xref: /petsc/src/ts/impls/implicit/theta/theta.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1316643e7SJed Brown /*
2316643e7SJed Brown   Code for timestepping with implicit Theta method
3316643e7SJed Brown */
4*b45d2f2cSJed 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) {
277445fe48SJed Brown       ierr = PetscObjectQuery((PetscObject)dm,"TSTheta_X0",(PetscObject*)X0);CHKERRQ(ierr);
287445fe48SJed Brown       if (!*X0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_X0 has not been composed with DM from SNES");
297445fe48SJed Brown     } else *X0 = ts->vec_sol;
307445fe48SJed Brown   }
317445fe48SJed Brown   if (Xdot) {
327445fe48SJed Brown     if (dm && dm != ts->dm) {
337445fe48SJed Brown       ierr = PetscObjectQuery((PetscObject)dm,"TSTheta_Xdot",(PetscObject*)Xdot);CHKERRQ(ierr);
347445fe48SJed Brown       if (!*Xdot) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"TSTheta_Xdot has not been composed with DM from SNES");
357445fe48SJed Brown     } else *Xdot = th->Xdot;
367445fe48SJed Brown   }
377445fe48SJed Brown   PetscFunctionReturn(0);
387445fe48SJed Brown }
397445fe48SJed Brown 
407445fe48SJed Brown #undef __FUNCT__
417445fe48SJed Brown #define __FUNCT__ "DMCoarsenHook_TSTheta"
427445fe48SJed Brown static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx)
437445fe48SJed Brown {
447445fe48SJed Brown   Vec X0,Xdot;
457445fe48SJed Brown   PetscErrorCode ierr;
467445fe48SJed Brown 
477445fe48SJed Brown   PetscFunctionBegin;
487445fe48SJed Brown   ierr = DMCreateGlobalVector(coarse,&X0);CHKERRQ(ierr);
497445fe48SJed Brown   ierr = DMCreateGlobalVector(coarse,&Xdot);CHKERRQ(ierr);
507445fe48SJed Brown   /* Oh noes, this would create a loop because the Vec holds a reference to the DM.
517445fe48SJed Brown      Making a PetscContainer to hold these Vecs would make the following call succeed, but would create a reference loop.
527445fe48SJed Brown      Need to decide on a way to break the reference counting loop.
537445fe48SJed Brown    */
547445fe48SJed Brown   ierr = PetscObjectCompose((PetscObject)coarse,"TSTheta_X0",(PetscObject)X0);CHKERRQ(ierr);
557445fe48SJed Brown   ierr = PetscObjectCompose((PetscObject)coarse,"TSTheta_Xdot",(PetscObject)Xdot);CHKERRQ(ierr);
567445fe48SJed Brown   ierr = VecDestroy(&X0);CHKERRQ(ierr);
577445fe48SJed Brown   ierr = VecDestroy(&Xdot);CHKERRQ(ierr);
587445fe48SJed Brown   PetscFunctionReturn(0);
597445fe48SJed Brown }
607445fe48SJed Brown 
617445fe48SJed Brown #undef __FUNCT__
627445fe48SJed Brown #define __FUNCT__ "DMRestrictHook_TSTheta"
637445fe48SJed Brown static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
647445fe48SJed Brown {
657445fe48SJed Brown   TS ts = (TS)ctx;
667445fe48SJed Brown   PetscErrorCode ierr;
677445fe48SJed Brown   Vec X0,Xdot,X0_c,Xdot_c;
687445fe48SJed Brown 
697445fe48SJed Brown   PetscFunctionBegin;
707445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr);
717445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr);
727445fe48SJed Brown   ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr);
737445fe48SJed Brown   ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr);
747445fe48SJed Brown   ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr);
757445fe48SJed Brown   ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr);
767445fe48SJed Brown   PetscFunctionReturn(0);
777445fe48SJed Brown }
787445fe48SJed Brown 
797445fe48SJed Brown #undef __FUNCT__
80316643e7SJed Brown #define __FUNCT__ "TSStep_Theta"
81193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts)
82316643e7SJed Brown {
83316643e7SJed Brown   TS_Theta            *th = (TS_Theta*)ts->data;
84b70ae86eSJed Brown   PetscInt            its,lits;
85cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
86f1b97656SJed Brown   SNESConvergedReason snesreason;
872b5a38e1SLisandro Dalcin   PetscErrorCode      ierr;
88316643e7SJed Brown 
89316643e7SJed Brown   PetscFunctionBegin;
90cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
91eb284becSJed Brown   th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step;
92316643e7SJed Brown   th->shift = 1./(th->Theta*ts->time_step);
93316643e7SJed Brown 
94eb284becSJed Brown   if (th->endpoint) {           /* This formulation assumes linear time-independent mass matrix */
95eb284becSJed Brown     ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr);
96eb284becSJed Brown     if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);}
97eb284becSJed Brown     ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr);
98eb284becSJed Brown     ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr);
99eb284becSJed Brown   }
100ace68cafSJed Brown   if (th->extrapolate) {
1012b5a38e1SLisandro Dalcin     ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr);
102ace68cafSJed Brown   } else {
1032b5a38e1SLisandro Dalcin     ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
104ace68cafSJed Brown   }
105eb284becSJed Brown   ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr);
106316643e7SJed Brown   ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr);
107316643e7SJed Brown   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
108f1b97656SJed Brown   ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);
109316643e7SJed Brown   ts->nonlinear_its += its; ts->linear_its += lits;
110f1b97656SJed Brown   if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
111f1b97656SJed Brown     ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
112f1b97656SJed 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);
113f1b97656SJed Brown     PetscFunctionReturn(0);
114f1b97656SJed Brown   }
115eb284becSJed Brown   if (th->endpoint) {
116eb284becSJed Brown     ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr);
117eb284becSJed Brown   } else {
1182b5a38e1SLisandro Dalcin     ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr);
1192b5a38e1SLisandro Dalcin     ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr);
120eb284becSJed Brown   }
1212b5a38e1SLisandro Dalcin   ts->ptime += ts->time_step;
122cdbf8f93SLisandro Dalcin   ts->time_step = next_time_step;
123316643e7SJed Brown   ts->steps++;
124316643e7SJed Brown   PetscFunctionReturn(0);
125316643e7SJed Brown }
126316643e7SJed Brown 
127cd652676SJed Brown #undef __FUNCT__
128cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta"
129cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X)
130cd652676SJed Brown {
131cd652676SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
1325a3a76d0SJed Brown   PetscReal      alpha = t - ts->ptime;
133cd652676SJed Brown   PetscErrorCode ierr;
134cd652676SJed Brown 
135cd652676SJed Brown   PetscFunctionBegin;
136a43b19c4SJed Brown   ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr);
1375a3a76d0SJed Brown   if (th->endpoint) alpha *= th->Theta;
1385a3a76d0SJed Brown   ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr);
139cd652676SJed Brown   PetscFunctionReturn(0);
140cd652676SJed Brown }
141cd652676SJed Brown 
142316643e7SJed Brown /*------------------------------------------------------------*/
143316643e7SJed Brown #undef __FUNCT__
144277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta"
145277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts)
146316643e7SJed Brown {
147316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
148316643e7SJed Brown   PetscErrorCode  ierr;
149316643e7SJed Brown 
150316643e7SJed Brown   PetscFunctionBegin;
1516bf464f9SBarry Smith   ierr = VecDestroy(&th->X);CHKERRQ(ierr);
1526bf464f9SBarry Smith   ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr);
153eb284becSJed Brown   ierr = VecDestroy(&th->affine);CHKERRQ(ierr);
154277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
155277b19d0SLisandro Dalcin }
156277b19d0SLisandro Dalcin 
157277b19d0SLisandro Dalcin #undef __FUNCT__
158277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta"
159277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts)
160277b19d0SLisandro Dalcin {
161277b19d0SLisandro Dalcin   PetscErrorCode  ierr;
162277b19d0SLisandro Dalcin 
163277b19d0SLisandro Dalcin   PetscFunctionBegin;
164277b19d0SLisandro Dalcin   ierr = TSReset_Theta(ts);CHKERRQ(ierr);
165277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
166335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
167335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr);
16826f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
169eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr);
170316643e7SJed Brown   PetscFunctionReturn(0);
171316643e7SJed Brown }
172316643e7SJed Brown 
173316643e7SJed Brown /*
174316643e7SJed Brown   This defines the nonlinear equation that is to be solved with SNES
1752b5a38e1SLisandro Dalcin   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
176316643e7SJed Brown */
177316643e7SJed Brown #undef __FUNCT__
1780f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta"
1790f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts)
180316643e7SJed Brown {
181316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
182316643e7SJed Brown   PetscErrorCode ierr;
1837445fe48SJed Brown   Vec            X0,Xdot;
1847445fe48SJed Brown   DM             dm,dmsave;
185316643e7SJed Brown 
186316643e7SJed Brown   PetscFunctionBegin;
1877445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1885a3a76d0SJed Brown   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
1897445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr);
1907445fe48SJed Brown   ierr = VecAXPBYPCZ(Xdot,-th->shift,th->shift,0,X0,x);CHKERRQ(ierr);
1917445fe48SJed Brown 
1927445fe48SJed Brown   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
1937445fe48SJed Brown   dmsave = ts->dm;
1947445fe48SJed Brown   ts->dm = dm;
1957445fe48SJed Brown   ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr);
1967445fe48SJed Brown   ts->dm = dmsave;
197316643e7SJed Brown   PetscFunctionReturn(0);
198316643e7SJed Brown }
199316643e7SJed Brown 
200316643e7SJed Brown #undef __FUNCT__
2010f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta"
2020f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts)
203316643e7SJed Brown {
204316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
205316643e7SJed Brown   PetscErrorCode ierr;
2067445fe48SJed Brown   Vec            Xdot;
2077445fe48SJed Brown   DM             dm,dmsave;
208316643e7SJed Brown 
209316643e7SJed Brown   PetscFunctionBegin;
2107445fe48SJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
2117445fe48SJed Brown 
2120f5c6efeSJed Brown   /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
2137445fe48SJed Brown   ierr = TSThetaGetX0AndXdot(ts,dm,PETSC_NULL,&Xdot);CHKERRQ(ierr);
2147445fe48SJed Brown 
2157445fe48SJed Brown   dmsave = ts->dm;
2167445fe48SJed Brown   ts->dm = dm;
2177445fe48SJed Brown   ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr);
2187445fe48SJed Brown   ts->dm = dmsave;
219316643e7SJed Brown   PetscFunctionReturn(0);
220316643e7SJed Brown }
221316643e7SJed Brown 
222316643e7SJed Brown #undef __FUNCT__
223316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta"
224316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts)
225316643e7SJed Brown {
226316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
227316643e7SJed Brown   PetscErrorCode ierr;
2287445fe48SJed Brown   SNES           snes;
2297445fe48SJed Brown   DM             dm;
230316643e7SJed Brown 
231316643e7SJed Brown   PetscFunctionBegin;
232316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr);
233316643e7SJed Brown   ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr);
2347445fe48SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
2357445fe48SJed Brown   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2367445fe48SJed Brown   if (dm) {
2377445fe48SJed Brown     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr);
2387445fe48SJed Brown   }
239316643e7SJed Brown   PetscFunctionReturn(0);
240316643e7SJed Brown }
241316643e7SJed Brown /*------------------------------------------------------------*/
242316643e7SJed Brown 
243316643e7SJed Brown #undef __FUNCT__
244316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta"
245316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts)
246316643e7SJed Brown {
247316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
248316643e7SJed Brown   PetscErrorCode ierr;
249316643e7SJed Brown 
250316643e7SJed Brown   PetscFunctionBegin;
251d73342a9SJed Brown   ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr);
252316643e7SJed Brown   {
253316643e7SJed Brown     ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr);
254acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr);
255eb284becSJed 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);
256d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
257316643e7SJed Brown   }
258316643e7SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
259316643e7SJed Brown   PetscFunctionReturn(0);
260316643e7SJed Brown }
261316643e7SJed Brown 
262316643e7SJed Brown #undef __FUNCT__
263316643e7SJed Brown #define __FUNCT__ "TSView_Theta"
264316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer)
265316643e7SJed Brown {
266316643e7SJed Brown   TS_Theta       *th = (TS_Theta*)ts->data;
267ace3abfcSBarry Smith   PetscBool       iascii;
268316643e7SJed Brown   PetscErrorCode  ierr;
269316643e7SJed Brown 
270316643e7SJed Brown   PetscFunctionBegin;
2712692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
272316643e7SJed Brown   if (iascii) {
273316643e7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Theta=%G\n",th->Theta);CHKERRQ(ierr);
274ace68cafSJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr);
275316643e7SJed Brown   }
276d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
277316643e7SJed Brown   PetscFunctionReturn(0);
278316643e7SJed Brown }
279316643e7SJed Brown 
2800de4c49aSJed Brown EXTERN_C_BEGIN
2810de4c49aSJed Brown #undef __FUNCT__
2820de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta"
2837087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta_Theta(TS ts,PetscReal *theta)
2840de4c49aSJed Brown {
2850de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2860de4c49aSJed Brown 
2870de4c49aSJed Brown   PetscFunctionBegin;
2880de4c49aSJed Brown   *theta = th->Theta;
2890de4c49aSJed Brown   PetscFunctionReturn(0);
2900de4c49aSJed Brown }
2910de4c49aSJed Brown 
2920de4c49aSJed Brown #undef __FUNCT__
2930de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta"
2947087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta_Theta(TS ts,PetscReal theta)
2950de4c49aSJed Brown {
2960de4c49aSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
2970de4c49aSJed Brown 
2980de4c49aSJed Brown   PetscFunctionBegin;
299e7be1afaSJed Brown   if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta);
3000de4c49aSJed Brown   th->Theta = theta;
3010de4c49aSJed Brown   PetscFunctionReturn(0);
3020de4c49aSJed Brown }
303eb284becSJed Brown 
304eb284becSJed Brown #undef __FUNCT__
30578e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta"
30626f2ff8fSLisandro Dalcin PetscErrorCode  TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint)
30726f2ff8fSLisandro Dalcin {
30826f2ff8fSLisandro Dalcin   TS_Theta *th = (TS_Theta*)ts->data;
30926f2ff8fSLisandro Dalcin 
31026f2ff8fSLisandro Dalcin   PetscFunctionBegin;
31126f2ff8fSLisandro Dalcin   *endpoint = th->endpoint;
31226f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
31326f2ff8fSLisandro Dalcin }
31426f2ff8fSLisandro Dalcin 
31526f2ff8fSLisandro Dalcin #undef __FUNCT__
31626f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta"
317eb284becSJed Brown PetscErrorCode  TSThetaSetEndpoint_Theta(TS ts,PetscBool flg)
318eb284becSJed Brown {
319eb284becSJed Brown   TS_Theta *th = (TS_Theta*)ts->data;
320eb284becSJed Brown 
321eb284becSJed Brown   PetscFunctionBegin;
322eb284becSJed Brown   th->endpoint = flg;
323eb284becSJed Brown   PetscFunctionReturn(0);
324eb284becSJed Brown }
3250de4c49aSJed Brown EXTERN_C_END
3260de4c49aSJed Brown 
327316643e7SJed Brown /* ------------------------------------------------------------ */
328316643e7SJed Brown /*MC
32996f5712cSJed Brown       TSTHETA - DAE solver using the implicit Theta method
330316643e7SJed Brown 
331316643e7SJed Brown    Level: beginner
332316643e7SJed Brown 
333eb284becSJed Brown    Notes:
334eb284becSJed Brown    This method can be applied to DAE.
335eb284becSJed Brown 
336eb284becSJed Brown    This method is cast as a 1-stage implicit Runge-Kutta method.
337eb284becSJed Brown 
338eb284becSJed Brown .vb
339eb284becSJed Brown   Theta | Theta
340eb284becSJed Brown   -------------
341eb284becSJed Brown         |  1
342eb284becSJed Brown .ve
343eb284becSJed Brown 
344eb284becSJed Brown    For the default Theta=0.5, this is also known as the implicit midpoint rule.
345eb284becSJed Brown 
346eb284becSJed Brown    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:
347eb284becSJed Brown 
348eb284becSJed Brown .vb
349eb284becSJed Brown   0 | 0         0
350eb284becSJed Brown   1 | 1-Theta   Theta
351eb284becSJed Brown   -------------------
352eb284becSJed Brown     | 1-Theta   Theta
353eb284becSJed Brown .ve
354eb284becSJed Brown 
355eb284becSJed Brown    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).
356eb284becSJed Brown 
357eb284becSJed Brown    To apply a diagonally implicit RK method to DAE, the stage formula
358eb284becSJed Brown 
359eb284becSJed Brown $  Y_i = X + h sum_j a_ij Y'_j
360eb284becSJed Brown 
361eb284becSJed Brown    is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i)
362eb284becSJed Brown 
363eb284becSJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint()
364316643e7SJed Brown 
365316643e7SJed Brown M*/
366316643e7SJed Brown EXTERN_C_BEGIN
367316643e7SJed Brown #undef __FUNCT__
368316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta"
3697087cfbeSBarry Smith PetscErrorCode  TSCreate_Theta(TS ts)
370316643e7SJed Brown {
371316643e7SJed Brown   TS_Theta       *th;
372316643e7SJed Brown   PetscErrorCode ierr;
373316643e7SJed Brown 
374316643e7SJed Brown   PetscFunctionBegin;
375277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Theta;
376316643e7SJed Brown   ts->ops->destroy        = TSDestroy_Theta;
377316643e7SJed Brown   ts->ops->view           = TSView_Theta;
378316643e7SJed Brown   ts->ops->setup          = TSSetUp_Theta;
379316643e7SJed Brown   ts->ops->step           = TSStep_Theta;
380cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_Theta;
381316643e7SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Theta;
3820f5c6efeSJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
3830f5c6efeSJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
384316643e7SJed Brown 
385316643e7SJed Brown   ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr);
386316643e7SJed Brown   ts->data = (void*)th;
387316643e7SJed Brown 
3886f700aefSJed Brown   th->extrapolate = PETSC_FALSE;
389316643e7SJed Brown   th->Theta       = 0.5;
390316643e7SJed Brown 
3910de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr);
3920de4c49aSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr);
39326f2ff8fSLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr);
394eb284becSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr);
395316643e7SJed Brown   PetscFunctionReturn(0);
396316643e7SJed Brown }
397316643e7SJed Brown EXTERN_C_END
3980de4c49aSJed Brown 
3990de4c49aSJed Brown #undef __FUNCT__
4000de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta"
4010de4c49aSJed Brown /*@
4020de4c49aSJed Brown   TSThetaGetTheta - Get the abscissa of the stage in (0,1].
4030de4c49aSJed Brown 
4040de4c49aSJed Brown   Not Collective
4050de4c49aSJed Brown 
4060de4c49aSJed Brown   Input Parameter:
4070de4c49aSJed Brown .  ts - timestepping context
4080de4c49aSJed Brown 
4090de4c49aSJed Brown   Output Parameter:
4100de4c49aSJed Brown .  theta - stage abscissa
4110de4c49aSJed Brown 
4120de4c49aSJed Brown   Note:
4130de4c49aSJed Brown   Use of this function is normally only required to hack TSTHETA to use a modified integration scheme.
4140de4c49aSJed Brown 
4150de4c49aSJed Brown   Level: Advanced
4160de4c49aSJed Brown 
4170de4c49aSJed Brown .seealso: TSThetaSetTheta()
4180de4c49aSJed Brown @*/
4197087cfbeSBarry Smith PetscErrorCode  TSThetaGetTheta(TS ts,PetscReal *theta)
4200de4c49aSJed Brown {
4214ac538c5SBarry Smith   PetscErrorCode ierr;
4220de4c49aSJed Brown 
4230de4c49aSJed Brown   PetscFunctionBegin;
424afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4250de4c49aSJed Brown   PetscValidPointer(theta,2);
4264ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr);
4270de4c49aSJed Brown   PetscFunctionReturn(0);
4280de4c49aSJed Brown }
4290de4c49aSJed Brown 
4300de4c49aSJed Brown #undef __FUNCT__
4310de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta"
4320de4c49aSJed Brown /*@
4330de4c49aSJed Brown   TSThetaSetTheta - Set the abscissa of the stage in (0,1].
4340de4c49aSJed Brown 
4350de4c49aSJed Brown   Not Collective
4360de4c49aSJed Brown 
4370de4c49aSJed Brown   Input Parameter:
4380de4c49aSJed Brown +  ts - timestepping context
4390de4c49aSJed Brown -  theta - stage abscissa
4400de4c49aSJed Brown 
4410de4c49aSJed Brown   Options Database:
4420de4c49aSJed Brown .  -ts_theta_theta <theta>
4430de4c49aSJed Brown 
4440de4c49aSJed Brown   Level: Intermediate
4450de4c49aSJed Brown 
4460de4c49aSJed Brown .seealso: TSThetaGetTheta()
4470de4c49aSJed Brown @*/
4487087cfbeSBarry Smith PetscErrorCode  TSThetaSetTheta(TS ts,PetscReal theta)
4490de4c49aSJed Brown {
4504ac538c5SBarry Smith   PetscErrorCode ierr;
4510de4c49aSJed Brown 
4520de4c49aSJed Brown   PetscFunctionBegin;
453afb20b64SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4544ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr);
4550de4c49aSJed Brown   PetscFunctionReturn(0);
4560de4c49aSJed Brown }
457f33bbcb6SJed Brown 
458eb284becSJed Brown #undef __FUNCT__
45926f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint"
46026f2ff8fSLisandro Dalcin /*@
46126f2ff8fSLisandro Dalcin   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
46226f2ff8fSLisandro Dalcin 
46326f2ff8fSLisandro Dalcin   Not Collective
46426f2ff8fSLisandro Dalcin 
46526f2ff8fSLisandro Dalcin   Input Parameter:
46626f2ff8fSLisandro Dalcin .  ts - timestepping context
46726f2ff8fSLisandro Dalcin 
46826f2ff8fSLisandro Dalcin   Output Parameter:
46926f2ff8fSLisandro Dalcin .  endpoint - PETSC_TRUE when using the endpoint variant
47026f2ff8fSLisandro Dalcin 
47126f2ff8fSLisandro Dalcin   Level: Advanced
47226f2ff8fSLisandro Dalcin 
47326f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN
47426f2ff8fSLisandro Dalcin @*/
47526f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint)
47626f2ff8fSLisandro Dalcin {
47726f2ff8fSLisandro Dalcin   PetscErrorCode ierr;
47826f2ff8fSLisandro Dalcin 
47926f2ff8fSLisandro Dalcin   PetscFunctionBegin;
48026f2ff8fSLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
48126f2ff8fSLisandro Dalcin   PetscValidPointer(endpoint,2);
48226f2ff8fSLisandro Dalcin   ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr);
48326f2ff8fSLisandro Dalcin   PetscFunctionReturn(0);
48426f2ff8fSLisandro Dalcin }
48526f2ff8fSLisandro Dalcin 
48626f2ff8fSLisandro Dalcin #undef __FUNCT__
487eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint"
488eb284becSJed Brown /*@
489eb284becSJed Brown   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule).
490eb284becSJed Brown 
491eb284becSJed Brown   Not Collective
492eb284becSJed Brown 
493eb284becSJed Brown   Input Parameter:
494eb284becSJed Brown +  ts - timestepping context
495eb284becSJed Brown -  flg - PETSC_TRUE to use the endpoint variant
496eb284becSJed Brown 
497eb284becSJed Brown   Options Database:
498eb284becSJed Brown .  -ts_theta_endpoint <flg>
499eb284becSJed Brown 
500eb284becSJed Brown   Level: Intermediate
501eb284becSJed Brown 
502eb284becSJed Brown .seealso: TSTHETA, TSCN
503eb284becSJed Brown @*/
504eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg)
505eb284becSJed Brown {
506eb284becSJed Brown   PetscErrorCode ierr;
507eb284becSJed Brown 
508eb284becSJed Brown   PetscFunctionBegin;
509eb284becSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
510eb284becSJed Brown   ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
511eb284becSJed Brown   PetscFunctionReturn(0);
512eb284becSJed Brown }
513eb284becSJed Brown 
514f33bbcb6SJed Brown /*
515f33bbcb6SJed Brown  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
516f33bbcb6SJed Brown  * The creation functions for these specializations are below.
517f33bbcb6SJed Brown  */
518f33bbcb6SJed Brown 
519f33bbcb6SJed Brown #undef __FUNCT__
520f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler"
521f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer)
522f33bbcb6SJed Brown {
523d52bd9f3SBarry Smith   PetscErrorCode ierr;
524d52bd9f3SBarry Smith 
525f33bbcb6SJed Brown   PetscFunctionBegin;
526d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
527f33bbcb6SJed Brown   PetscFunctionReturn(0);
528f33bbcb6SJed Brown }
529f33bbcb6SJed Brown 
530f33bbcb6SJed Brown /*MC
531f33bbcb6SJed Brown       TSBEULER - ODE solver using the implicit backward Euler method
532f33bbcb6SJed Brown 
533f33bbcb6SJed Brown   Level: beginner
534f33bbcb6SJed Brown 
535f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA
536f33bbcb6SJed Brown 
537f33bbcb6SJed Brown M*/
538f33bbcb6SJed Brown EXTERN_C_BEGIN
539f33bbcb6SJed Brown #undef __FUNCT__
540f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler"
541f33bbcb6SJed Brown PetscErrorCode  TSCreate_BEuler(TS ts)
542f33bbcb6SJed Brown {
543f33bbcb6SJed Brown   PetscErrorCode ierr;
544f33bbcb6SJed Brown 
545f33bbcb6SJed Brown   PetscFunctionBegin;
546f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
547f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr);
548f33bbcb6SJed Brown   ts->ops->view = TSView_BEuler;
549f33bbcb6SJed Brown   PetscFunctionReturn(0);
550f33bbcb6SJed Brown }
551f33bbcb6SJed Brown EXTERN_C_END
552f33bbcb6SJed Brown 
553f33bbcb6SJed Brown #undef __FUNCT__
554f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN"
555f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer)
556f33bbcb6SJed Brown {
557d52bd9f3SBarry Smith   PetscErrorCode ierr;
558d52bd9f3SBarry Smith 
559f33bbcb6SJed Brown   PetscFunctionBegin;
560d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
561f33bbcb6SJed Brown   PetscFunctionReturn(0);
562f33bbcb6SJed Brown }
563f33bbcb6SJed Brown 
564f33bbcb6SJed Brown /*MC
565f33bbcb6SJed Brown       TSCN - ODE solver using the implicit Crank-Nicolson method.
566f33bbcb6SJed Brown 
567f33bbcb6SJed Brown   Level: beginner
568f33bbcb6SJed Brown 
569f33bbcb6SJed Brown   Notes:
5707cf5af47SJed Brown   TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e.
5717cf5af47SJed Brown 
5727cf5af47SJed Brown $  -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint
573f33bbcb6SJed Brown 
574f33bbcb6SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA
575f33bbcb6SJed Brown 
576f33bbcb6SJed Brown M*/
577f33bbcb6SJed Brown EXTERN_C_BEGIN
578f33bbcb6SJed Brown #undef __FUNCT__
579f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN"
580f33bbcb6SJed Brown PetscErrorCode  TSCreate_CN(TS ts)
581f33bbcb6SJed Brown {
582f33bbcb6SJed Brown   PetscErrorCode ierr;
583f33bbcb6SJed Brown 
584f33bbcb6SJed Brown   PetscFunctionBegin;
585f33bbcb6SJed Brown   ierr = TSCreate_Theta(ts);CHKERRQ(ierr);
586f33bbcb6SJed Brown   ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr);
587eb284becSJed Brown   ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr);
588f33bbcb6SJed Brown   ts->ops->view = TSView_CN;
589f33bbcb6SJed Brown   PetscFunctionReturn(0);
590f33bbcb6SJed Brown }
591f33bbcb6SJed Brown EXTERN_C_END
592