1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown 4316643e7SJed Brown Notes: 5316643e7SJed Brown This method can be applied to DAE. 6316643e7SJed Brown 7316643e7SJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 8316643e7SJed Brown 9316643e7SJed Brown Theta | Theta 10316643e7SJed Brown ------------- 11316643e7SJed Brown | 1 12316643e7SJed Brown 13316643e7SJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 14316643e7SJed Brown 15316643e7SJed Brown X_i = x + h sum_j a_ij X'_j 16316643e7SJed Brown 17316643e7SJed Brown is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i) 18316643e7SJed Brown */ 19c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 20316643e7SJed Brown 21316643e7SJed Brown typedef struct { 22316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 23ace3abfcSBarry Smith PetscBool extrapolate; 24316643e7SJed Brown PetscReal Theta; 25316643e7SJed Brown PetscReal shift; 26316643e7SJed Brown PetscReal stage_time; 27316643e7SJed Brown } TS_Theta; 28316643e7SJed Brown 29316643e7SJed Brown #undef __FUNCT__ 30316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 31193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 32316643e7SJed Brown { 33316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 34b70ae86eSJed Brown PetscInt its,lits; 352b5a38e1SLisandro Dalcin PetscErrorCode ierr; 36316643e7SJed Brown 37316643e7SJed Brown PetscFunctionBegin; 38193ac0bcSJed Brown ts->time_step = ts->next_time_step; 39316643e7SJed Brown th->stage_time = ts->ptime + th->Theta*ts->time_step; 40316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 41316643e7SJed Brown 42ace68cafSJed Brown if (th->extrapolate) { 432b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 44ace68cafSJed Brown } else { 452b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 46ace68cafSJed Brown } 47316643e7SJed Brown ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr); 48316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 49316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 50316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 512b5a38e1SLisandro Dalcin 522b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 532b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 542b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 55193ac0bcSJed Brown ts->next_time_step = ts->time_step; 56316643e7SJed Brown ts->steps++; 57316643e7SJed Brown PetscFunctionReturn(0); 58316643e7SJed Brown } 59316643e7SJed Brown 60316643e7SJed Brown /*------------------------------------------------------------*/ 61316643e7SJed Brown #undef __FUNCT__ 62277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 63277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 64316643e7SJed Brown { 65316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 66316643e7SJed Brown PetscErrorCode ierr; 67316643e7SJed Brown 68316643e7SJed Brown PetscFunctionBegin; 696bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 706bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 71277b19d0SLisandro Dalcin PetscFunctionReturn(0); 72277b19d0SLisandro Dalcin } 73277b19d0SLisandro Dalcin 74277b19d0SLisandro Dalcin #undef __FUNCT__ 75277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 76277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 77277b19d0SLisandro Dalcin { 78277b19d0SLisandro Dalcin PetscErrorCode ierr; 79277b19d0SLisandro Dalcin 80277b19d0SLisandro Dalcin PetscFunctionBegin; 81277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 82277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 83*335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 84*335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 85316643e7SJed Brown PetscFunctionReturn(0); 86316643e7SJed Brown } 87316643e7SJed Brown 88316643e7SJed Brown /* 89316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 902b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 91316643e7SJed Brown */ 92316643e7SJed Brown #undef __FUNCT__ 930f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 940f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 95316643e7SJed Brown { 96316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 97316643e7SJed Brown PetscErrorCode ierr; 98316643e7SJed Brown 99316643e7SJed Brown PetscFunctionBegin; 1002b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 101214bc6a2SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 102316643e7SJed Brown PetscFunctionReturn(0); 103316643e7SJed Brown } 104316643e7SJed Brown 105316643e7SJed Brown #undef __FUNCT__ 1060f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1070f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 108316643e7SJed Brown { 109316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 110316643e7SJed Brown PetscErrorCode ierr; 111316643e7SJed Brown 112316643e7SJed Brown PetscFunctionBegin; 1130f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 114214bc6a2SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 115316643e7SJed Brown PetscFunctionReturn(0); 116316643e7SJed Brown } 117316643e7SJed Brown 118316643e7SJed Brown 119316643e7SJed Brown #undef __FUNCT__ 120316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 121316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 122316643e7SJed Brown { 123316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 124316643e7SJed Brown PetscErrorCode ierr; 125316643e7SJed Brown 126316643e7SJed Brown PetscFunctionBegin; 127316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 128316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 129316643e7SJed Brown PetscFunctionReturn(0); 130316643e7SJed Brown } 131316643e7SJed Brown /*------------------------------------------------------------*/ 132316643e7SJed Brown 133316643e7SJed Brown #undef __FUNCT__ 134316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 135316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 136316643e7SJed Brown { 137316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 138316643e7SJed Brown PetscErrorCode ierr; 139316643e7SJed Brown 140316643e7SJed Brown PetscFunctionBegin; 141d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 142316643e7SJed Brown { 143316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 144acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 145316643e7SJed Brown } 146316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 147316643e7SJed Brown PetscFunctionReturn(0); 148316643e7SJed Brown } 149316643e7SJed Brown 150316643e7SJed Brown #undef __FUNCT__ 151316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 152316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 153316643e7SJed Brown { 154316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 155ace3abfcSBarry Smith PetscBool iascii; 156316643e7SJed Brown PetscErrorCode ierr; 157316643e7SJed Brown 158316643e7SJed Brown PetscFunctionBegin; 1592692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 160316643e7SJed Brown if (iascii) { 161316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 162ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 163316643e7SJed Brown } else { 164e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name); 165316643e7SJed Brown } 166316643e7SJed Brown PetscFunctionReturn(0); 167316643e7SJed Brown } 168316643e7SJed Brown 1690de4c49aSJed Brown EXTERN_C_BEGIN 1700de4c49aSJed Brown #undef __FUNCT__ 1710de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 1727087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1730de4c49aSJed Brown { 1740de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1750de4c49aSJed Brown 1760de4c49aSJed Brown PetscFunctionBegin; 1770de4c49aSJed Brown *theta = th->Theta; 1780de4c49aSJed Brown PetscFunctionReturn(0); 1790de4c49aSJed Brown } 1800de4c49aSJed Brown 1810de4c49aSJed Brown #undef __FUNCT__ 1820de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 1837087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 1840de4c49aSJed Brown { 1850de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1860de4c49aSJed Brown 1870de4c49aSJed Brown PetscFunctionBegin; 188e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 1890de4c49aSJed Brown th->Theta = theta; 1900de4c49aSJed Brown PetscFunctionReturn(0); 1910de4c49aSJed Brown } 1920de4c49aSJed Brown EXTERN_C_END 1930de4c49aSJed Brown 194316643e7SJed Brown /* ------------------------------------------------------------ */ 195316643e7SJed Brown /*MC 19696f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 197316643e7SJed Brown 198316643e7SJed Brown Level: beginner 199316643e7SJed Brown 200316643e7SJed Brown .seealso: TSCreate(), TS, TSSetType() 201316643e7SJed Brown 202316643e7SJed Brown M*/ 203316643e7SJed Brown EXTERN_C_BEGIN 204316643e7SJed Brown #undef __FUNCT__ 205316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2067087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 207316643e7SJed Brown { 208316643e7SJed Brown TS_Theta *th; 209316643e7SJed Brown PetscErrorCode ierr; 210316643e7SJed Brown 211316643e7SJed Brown PetscFunctionBegin; 212277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 213316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 214316643e7SJed Brown ts->ops->view = TSView_Theta; 215316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 216316643e7SJed Brown ts->ops->step = TSStep_Theta; 217316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2180f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2190f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 220316643e7SJed Brown 221316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 222316643e7SJed Brown ts->data = (void*)th; 223316643e7SJed Brown 2246f700aefSJed Brown th->extrapolate = PETSC_FALSE; 225316643e7SJed Brown th->Theta = 0.5; 226316643e7SJed Brown 2270de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 2280de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 229316643e7SJed Brown PetscFunctionReturn(0); 230316643e7SJed Brown } 231316643e7SJed Brown EXTERN_C_END 2320de4c49aSJed Brown 2330de4c49aSJed Brown #undef __FUNCT__ 2340de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 2350de4c49aSJed Brown /*@ 2360de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 2370de4c49aSJed Brown 2380de4c49aSJed Brown Not Collective 2390de4c49aSJed Brown 2400de4c49aSJed Brown Input Parameter: 2410de4c49aSJed Brown . ts - timestepping context 2420de4c49aSJed Brown 2430de4c49aSJed Brown Output Parameter: 2440de4c49aSJed Brown . theta - stage abscissa 2450de4c49aSJed Brown 2460de4c49aSJed Brown Note: 2470de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 2480de4c49aSJed Brown 2490de4c49aSJed Brown Level: Advanced 2500de4c49aSJed Brown 2510de4c49aSJed Brown .seealso: TSThetaSetTheta() 2520de4c49aSJed Brown @*/ 2537087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 2540de4c49aSJed Brown { 2554ac538c5SBarry Smith PetscErrorCode ierr; 2560de4c49aSJed Brown 2570de4c49aSJed Brown PetscFunctionBegin; 258afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2590de4c49aSJed Brown PetscValidPointer(theta,2); 2604ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 2610de4c49aSJed Brown PetscFunctionReturn(0); 2620de4c49aSJed Brown } 2630de4c49aSJed Brown 2640de4c49aSJed Brown #undef __FUNCT__ 2650de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 2660de4c49aSJed Brown /*@ 2670de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 2680de4c49aSJed Brown 2690de4c49aSJed Brown Not Collective 2700de4c49aSJed Brown 2710de4c49aSJed Brown Input Parameter: 2720de4c49aSJed Brown + ts - timestepping context 2730de4c49aSJed Brown - theta - stage abscissa 2740de4c49aSJed Brown 2750de4c49aSJed Brown Options Database: 2760de4c49aSJed Brown . -ts_theta_theta <theta> 2770de4c49aSJed Brown 2780de4c49aSJed Brown Level: Intermediate 2790de4c49aSJed Brown 2800de4c49aSJed Brown .seealso: TSThetaGetTheta() 2810de4c49aSJed Brown @*/ 2827087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 2830de4c49aSJed Brown { 2844ac538c5SBarry Smith PetscErrorCode ierr; 2850de4c49aSJed Brown 2860de4c49aSJed Brown PetscFunctionBegin; 287afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2884ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 2890de4c49aSJed Brown PetscFunctionReturn(0); 2900de4c49aSJed Brown } 291f33bbcb6SJed Brown 292f33bbcb6SJed Brown /* 293f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 294f33bbcb6SJed Brown * The creation functions for these specializations are below. 295f33bbcb6SJed Brown */ 296f33bbcb6SJed Brown 297f33bbcb6SJed Brown #undef __FUNCT__ 298f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 299f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 300f33bbcb6SJed Brown { 301f33bbcb6SJed Brown PetscFunctionBegin; 302f33bbcb6SJed Brown PetscFunctionReturn(0); 303f33bbcb6SJed Brown } 304f33bbcb6SJed Brown 305f33bbcb6SJed Brown /*MC 306f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 307f33bbcb6SJed Brown 308f33bbcb6SJed Brown Level: beginner 309f33bbcb6SJed Brown 310f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 311f33bbcb6SJed Brown 312f33bbcb6SJed Brown M*/ 313f33bbcb6SJed Brown EXTERN_C_BEGIN 314f33bbcb6SJed Brown #undef __FUNCT__ 315f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 316f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 317f33bbcb6SJed Brown { 318f33bbcb6SJed Brown PetscErrorCode ierr; 319f33bbcb6SJed Brown 320f33bbcb6SJed Brown PetscFunctionBegin; 321f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 322f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 323f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 324f33bbcb6SJed Brown PetscFunctionReturn(0); 325f33bbcb6SJed Brown } 326f33bbcb6SJed Brown EXTERN_C_END 327f33bbcb6SJed Brown 328f33bbcb6SJed Brown #undef __FUNCT__ 329f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 330f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 331f33bbcb6SJed Brown { 332f33bbcb6SJed Brown PetscFunctionBegin; 333f33bbcb6SJed Brown PetscFunctionReturn(0); 334f33bbcb6SJed Brown } 335f33bbcb6SJed Brown 336f33bbcb6SJed Brown /*MC 337f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 338f33bbcb6SJed Brown 339f33bbcb6SJed Brown Level: beginner 340f33bbcb6SJed Brown 341f33bbcb6SJed Brown Notes: 342f33bbcb6SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 343f33bbcb6SJed Brown 344f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 345f33bbcb6SJed Brown 346f33bbcb6SJed Brown M*/ 347f33bbcb6SJed Brown EXTERN_C_BEGIN 348f33bbcb6SJed Brown #undef __FUNCT__ 349f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 350f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 351f33bbcb6SJed Brown { 352f33bbcb6SJed Brown PetscErrorCode ierr; 353f33bbcb6SJed Brown 354f33bbcb6SJed Brown PetscFunctionBegin; 355f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 356f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 357f33bbcb6SJed Brown ts->ops->view = TSView_CN; 358f33bbcb6SJed Brown PetscFunctionReturn(0); 359f33bbcb6SJed Brown } 360f33bbcb6SJed Brown EXTERN_C_END 361