1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown */ 4c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 5316643e7SJed Brown 6316643e7SJed Brown typedef struct { 7316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 8eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 9ace3abfcSBarry Smith PetscBool extrapolate; 10eb284becSJed Brown PetscBool endpoint; 11316643e7SJed Brown PetscReal Theta; 12316643e7SJed Brown PetscReal shift; 13316643e7SJed Brown PetscReal stage_time; 14316643e7SJed Brown } TS_Theta; 15316643e7SJed Brown 16316643e7SJed Brown #undef __FUNCT__ 17316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 18193ac0bcSJed Brown static PetscErrorCode TSStep_Theta(TS ts) 19316643e7SJed Brown { 20316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 21b70ae86eSJed Brown PetscInt its,lits; 22cdbf8f93SLisandro Dalcin PetscReal next_time_step; 232b5a38e1SLisandro Dalcin PetscErrorCode ierr; 24316643e7SJed Brown 25316643e7SJed Brown PetscFunctionBegin; 26cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 27eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 28316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 29316643e7SJed Brown 30eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 31eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 32eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 33eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 34eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 35eb284becSJed Brown } 36ace68cafSJed Brown if (th->extrapolate) { 372b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 38ace68cafSJed Brown } else { 392b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 40ace68cafSJed Brown } 41eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 42316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 43316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 44316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 452b5a38e1SLisandro Dalcin 46eb284becSJed Brown if (th->endpoint) { 47eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 48eb284becSJed Brown } else { 492b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 502b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 51eb284becSJed Brown } 522b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 53cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 54316643e7SJed Brown ts->steps++; 55316643e7SJed Brown PetscFunctionReturn(0); 56316643e7SJed Brown } 57316643e7SJed Brown 58cd652676SJed Brown #undef __FUNCT__ 59cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 60cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 61cd652676SJed Brown { 62cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 635a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 64cd652676SJed Brown PetscErrorCode ierr; 65cd652676SJed Brown 66cd652676SJed Brown PetscFunctionBegin; 67a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 685a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 695a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 70cd652676SJed Brown PetscFunctionReturn(0); 71cd652676SJed Brown } 72cd652676SJed Brown 73316643e7SJed Brown /*------------------------------------------------------------*/ 74316643e7SJed Brown #undef __FUNCT__ 75277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 76277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 77316643e7SJed Brown { 78316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 79316643e7SJed Brown PetscErrorCode ierr; 80316643e7SJed Brown 81316643e7SJed Brown PetscFunctionBegin; 826bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 836bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 84eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 85277b19d0SLisandro Dalcin PetscFunctionReturn(0); 86277b19d0SLisandro Dalcin } 87277b19d0SLisandro Dalcin 88277b19d0SLisandro Dalcin #undef __FUNCT__ 89277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 90277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 91277b19d0SLisandro Dalcin { 92277b19d0SLisandro Dalcin PetscErrorCode ierr; 93277b19d0SLisandro Dalcin 94277b19d0SLisandro Dalcin PetscFunctionBegin; 95277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 96277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 97335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 98335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 99*26f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 100eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 101316643e7SJed Brown PetscFunctionReturn(0); 102316643e7SJed Brown } 103316643e7SJed Brown 104316643e7SJed Brown /* 105316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1062b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 107316643e7SJed Brown */ 108316643e7SJed Brown #undef __FUNCT__ 1090f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1100f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 111316643e7SJed Brown { 112316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 113316643e7SJed Brown PetscErrorCode ierr; 114316643e7SJed Brown 115316643e7SJed Brown PetscFunctionBegin; 1165a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 1172b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 118214bc6a2SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 119316643e7SJed Brown PetscFunctionReturn(0); 120316643e7SJed Brown } 121316643e7SJed Brown 122316643e7SJed Brown #undef __FUNCT__ 1230f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1240f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 125316643e7SJed Brown { 126316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 127316643e7SJed Brown PetscErrorCode ierr; 128316643e7SJed Brown 129316643e7SJed Brown PetscFunctionBegin; 1300f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 131214bc6a2SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 132316643e7SJed Brown PetscFunctionReturn(0); 133316643e7SJed Brown } 134316643e7SJed Brown 135316643e7SJed Brown 136316643e7SJed Brown #undef __FUNCT__ 137316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 138316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 139316643e7SJed Brown { 140316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 141316643e7SJed Brown PetscErrorCode ierr; 142316643e7SJed Brown 143316643e7SJed Brown PetscFunctionBegin; 144316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 145316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 146316643e7SJed Brown PetscFunctionReturn(0); 147316643e7SJed Brown } 148316643e7SJed Brown /*------------------------------------------------------------*/ 149316643e7SJed Brown 150316643e7SJed Brown #undef __FUNCT__ 151316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 152316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 153316643e7SJed Brown { 154316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 155316643e7SJed Brown PetscErrorCode ierr; 156316643e7SJed Brown 157316643e7SJed Brown PetscFunctionBegin; 158d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 159316643e7SJed Brown { 160316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 161acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 162eb284becSJed 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); 163d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 164316643e7SJed Brown } 165316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 166316643e7SJed Brown PetscFunctionReturn(0); 167316643e7SJed Brown } 168316643e7SJed Brown 169316643e7SJed Brown #undef __FUNCT__ 170316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 171316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 172316643e7SJed Brown { 173316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 174ace3abfcSBarry Smith PetscBool iascii; 175316643e7SJed Brown PetscErrorCode ierr; 176316643e7SJed Brown 177316643e7SJed Brown PetscFunctionBegin; 1782692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 179316643e7SJed Brown if (iascii) { 180316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 181ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 182316643e7SJed Brown } 183d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 184316643e7SJed Brown PetscFunctionReturn(0); 185316643e7SJed Brown } 186316643e7SJed Brown 1870de4c49aSJed Brown EXTERN_C_BEGIN 1880de4c49aSJed Brown #undef __FUNCT__ 1890de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 1907087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1910de4c49aSJed Brown { 1920de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1930de4c49aSJed Brown 1940de4c49aSJed Brown PetscFunctionBegin; 1950de4c49aSJed Brown *theta = th->Theta; 1960de4c49aSJed Brown PetscFunctionReturn(0); 1970de4c49aSJed Brown } 1980de4c49aSJed Brown 1990de4c49aSJed Brown #undef __FUNCT__ 2000de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 2017087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 2020de4c49aSJed Brown { 2030de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2040de4c49aSJed Brown 2050de4c49aSJed Brown PetscFunctionBegin; 206e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 2070de4c49aSJed Brown th->Theta = theta; 2080de4c49aSJed Brown PetscFunctionReturn(0); 2090de4c49aSJed Brown } 210eb284becSJed Brown 211eb284becSJed Brown #undef __FUNCT__ 212eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint_Theta" 213*26f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 214*26f2ff8fSLisandro Dalcin { 215*26f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 216*26f2ff8fSLisandro Dalcin 217*26f2ff8fSLisandro Dalcin PetscFunctionBegin; 218*26f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 219*26f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 220*26f2ff8fSLisandro Dalcin } 221*26f2ff8fSLisandro Dalcin 222*26f2ff8fSLisandro Dalcin #undef __FUNCT__ 223*26f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 224eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 225eb284becSJed Brown { 226eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 227eb284becSJed Brown 228eb284becSJed Brown PetscFunctionBegin; 229eb284becSJed Brown th->endpoint = flg; 230eb284becSJed Brown PetscFunctionReturn(0); 231eb284becSJed Brown } 2320de4c49aSJed Brown EXTERN_C_END 2330de4c49aSJed Brown 234316643e7SJed Brown /* ------------------------------------------------------------ */ 235316643e7SJed Brown /*MC 23696f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 237316643e7SJed Brown 238316643e7SJed Brown Level: beginner 239316643e7SJed Brown 240eb284becSJed Brown Notes: 241eb284becSJed Brown This method can be applied to DAE. 242eb284becSJed Brown 243eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 244eb284becSJed Brown 245eb284becSJed Brown .vb 246eb284becSJed Brown Theta | Theta 247eb284becSJed Brown ------------- 248eb284becSJed Brown | 1 249eb284becSJed Brown .ve 250eb284becSJed Brown 251eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 252eb284becSJed Brown 253eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 254eb284becSJed Brown 255eb284becSJed Brown .vb 256eb284becSJed Brown 0 | 0 0 257eb284becSJed Brown 1 | 1-Theta Theta 258eb284becSJed Brown ------------------- 259eb284becSJed Brown | 1-Theta Theta 260eb284becSJed Brown .ve 261eb284becSJed Brown 262eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 263eb284becSJed Brown 264eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 265eb284becSJed Brown 266eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 267eb284becSJed Brown 268eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 269eb284becSJed Brown 270eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 271316643e7SJed Brown 272316643e7SJed Brown M*/ 273316643e7SJed Brown EXTERN_C_BEGIN 274316643e7SJed Brown #undef __FUNCT__ 275316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2767087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 277316643e7SJed Brown { 278316643e7SJed Brown TS_Theta *th; 279316643e7SJed Brown PetscErrorCode ierr; 280316643e7SJed Brown 281316643e7SJed Brown PetscFunctionBegin; 282277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 283316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 284316643e7SJed Brown ts->ops->view = TSView_Theta; 285316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 286316643e7SJed Brown ts->ops->step = TSStep_Theta; 287cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 288316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2890f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2900f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 291316643e7SJed Brown 292316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 293316643e7SJed Brown ts->data = (void*)th; 294316643e7SJed Brown 2956f700aefSJed Brown th->extrapolate = PETSC_FALSE; 296316643e7SJed Brown th->Theta = 0.5; 297316643e7SJed Brown 2980de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 2990de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 300*26f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 301eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 302316643e7SJed Brown PetscFunctionReturn(0); 303316643e7SJed Brown } 304316643e7SJed Brown EXTERN_C_END 3050de4c49aSJed Brown 3060de4c49aSJed Brown #undef __FUNCT__ 3070de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 3080de4c49aSJed Brown /*@ 3090de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 3100de4c49aSJed Brown 3110de4c49aSJed Brown Not Collective 3120de4c49aSJed Brown 3130de4c49aSJed Brown Input Parameter: 3140de4c49aSJed Brown . ts - timestepping context 3150de4c49aSJed Brown 3160de4c49aSJed Brown Output Parameter: 3170de4c49aSJed Brown . theta - stage abscissa 3180de4c49aSJed Brown 3190de4c49aSJed Brown Note: 3200de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 3210de4c49aSJed Brown 3220de4c49aSJed Brown Level: Advanced 3230de4c49aSJed Brown 3240de4c49aSJed Brown .seealso: TSThetaSetTheta() 3250de4c49aSJed Brown @*/ 3267087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 3270de4c49aSJed Brown { 3284ac538c5SBarry Smith PetscErrorCode ierr; 3290de4c49aSJed Brown 3300de4c49aSJed Brown PetscFunctionBegin; 331afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3320de4c49aSJed Brown PetscValidPointer(theta,2); 3334ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 3340de4c49aSJed Brown PetscFunctionReturn(0); 3350de4c49aSJed Brown } 3360de4c49aSJed Brown 3370de4c49aSJed Brown #undef __FUNCT__ 3380de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 3390de4c49aSJed Brown /*@ 3400de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 3410de4c49aSJed Brown 3420de4c49aSJed Brown Not Collective 3430de4c49aSJed Brown 3440de4c49aSJed Brown Input Parameter: 3450de4c49aSJed Brown + ts - timestepping context 3460de4c49aSJed Brown - theta - stage abscissa 3470de4c49aSJed Brown 3480de4c49aSJed Brown Options Database: 3490de4c49aSJed Brown . -ts_theta_theta <theta> 3500de4c49aSJed Brown 3510de4c49aSJed Brown Level: Intermediate 3520de4c49aSJed Brown 3530de4c49aSJed Brown .seealso: TSThetaGetTheta() 3540de4c49aSJed Brown @*/ 3557087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 3560de4c49aSJed Brown { 3574ac538c5SBarry Smith PetscErrorCode ierr; 3580de4c49aSJed Brown 3590de4c49aSJed Brown PetscFunctionBegin; 360afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3614ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 3620de4c49aSJed Brown PetscFunctionReturn(0); 3630de4c49aSJed Brown } 364f33bbcb6SJed Brown 365eb284becSJed Brown #undef __FUNCT__ 366*26f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 367*26f2ff8fSLisandro Dalcin /*@ 368*26f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 369*26f2ff8fSLisandro Dalcin 370*26f2ff8fSLisandro Dalcin Not Collective 371*26f2ff8fSLisandro Dalcin 372*26f2ff8fSLisandro Dalcin Input Parameter: 373*26f2ff8fSLisandro Dalcin . ts - timestepping context 374*26f2ff8fSLisandro Dalcin 375*26f2ff8fSLisandro Dalcin Output Parameter: 376*26f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 377*26f2ff8fSLisandro Dalcin 378*26f2ff8fSLisandro Dalcin Level: Advanced 379*26f2ff8fSLisandro Dalcin 380*26f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 381*26f2ff8fSLisandro Dalcin @*/ 382*26f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 383*26f2ff8fSLisandro Dalcin { 384*26f2ff8fSLisandro Dalcin PetscErrorCode ierr; 385*26f2ff8fSLisandro Dalcin 386*26f2ff8fSLisandro Dalcin PetscFunctionBegin; 387*26f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 388*26f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 389*26f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 390*26f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 391*26f2ff8fSLisandro Dalcin } 392*26f2ff8fSLisandro Dalcin 393*26f2ff8fSLisandro Dalcin #undef __FUNCT__ 394eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 395eb284becSJed Brown /*@ 396eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 397eb284becSJed Brown 398eb284becSJed Brown Not Collective 399eb284becSJed Brown 400eb284becSJed Brown Input Parameter: 401eb284becSJed Brown + ts - timestepping context 402eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 403eb284becSJed Brown 404eb284becSJed Brown Options Database: 405eb284becSJed Brown . -ts_theta_endpoint <flg> 406eb284becSJed Brown 407eb284becSJed Brown Level: Intermediate 408eb284becSJed Brown 409eb284becSJed Brown .seealso: TSTHETA, TSCN 410eb284becSJed Brown @*/ 411eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 412eb284becSJed Brown { 413eb284becSJed Brown PetscErrorCode ierr; 414eb284becSJed Brown 415eb284becSJed Brown PetscFunctionBegin; 416eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 417eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 418eb284becSJed Brown PetscFunctionReturn(0); 419eb284becSJed Brown } 420eb284becSJed Brown 421f33bbcb6SJed Brown /* 422f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 423f33bbcb6SJed Brown * The creation functions for these specializations are below. 424f33bbcb6SJed Brown */ 425f33bbcb6SJed Brown 426f33bbcb6SJed Brown #undef __FUNCT__ 427f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 428f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 429f33bbcb6SJed Brown { 430d52bd9f3SBarry Smith PetscErrorCode ierr; 431d52bd9f3SBarry Smith 432f33bbcb6SJed Brown PetscFunctionBegin; 433d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 434f33bbcb6SJed Brown PetscFunctionReturn(0); 435f33bbcb6SJed Brown } 436f33bbcb6SJed Brown 437f33bbcb6SJed Brown /*MC 438f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 439f33bbcb6SJed Brown 440f33bbcb6SJed Brown Level: beginner 441f33bbcb6SJed Brown 442f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 443f33bbcb6SJed Brown 444f33bbcb6SJed Brown M*/ 445f33bbcb6SJed Brown EXTERN_C_BEGIN 446f33bbcb6SJed Brown #undef __FUNCT__ 447f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 448f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 449f33bbcb6SJed Brown { 450f33bbcb6SJed Brown PetscErrorCode ierr; 451f33bbcb6SJed Brown 452f33bbcb6SJed Brown PetscFunctionBegin; 453f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 454f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 455f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 456f33bbcb6SJed Brown PetscFunctionReturn(0); 457f33bbcb6SJed Brown } 458f33bbcb6SJed Brown EXTERN_C_END 459f33bbcb6SJed Brown 460f33bbcb6SJed Brown #undef __FUNCT__ 461f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 462f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 463f33bbcb6SJed Brown { 464d52bd9f3SBarry Smith PetscErrorCode ierr; 465d52bd9f3SBarry Smith 466f33bbcb6SJed Brown PetscFunctionBegin; 467d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 468f33bbcb6SJed Brown PetscFunctionReturn(0); 469f33bbcb6SJed Brown } 470f33bbcb6SJed Brown 471f33bbcb6SJed Brown /*MC 472f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 473f33bbcb6SJed Brown 474f33bbcb6SJed Brown Level: beginner 475f33bbcb6SJed Brown 476f33bbcb6SJed Brown Notes: 4777cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 4787cf5af47SJed Brown 4797cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 480f33bbcb6SJed Brown 481f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 482f33bbcb6SJed Brown 483f33bbcb6SJed Brown M*/ 484f33bbcb6SJed Brown EXTERN_C_BEGIN 485f33bbcb6SJed Brown #undef __FUNCT__ 486f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 487f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 488f33bbcb6SJed Brown { 489f33bbcb6SJed Brown PetscErrorCode ierr; 490f33bbcb6SJed Brown 491f33bbcb6SJed Brown PetscFunctionBegin; 492f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 493f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 494eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 495f33bbcb6SJed Brown ts->ops->view = TSView_CN; 496f33bbcb6SJed Brown PetscFunctionReturn(0); 497f33bbcb6SJed Brown } 498f33bbcb6SJed Brown EXTERN_C_END 499