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 */ 8*eb284becSJed Brown Vec affine; /* Affine vector needed for residual at beginning of step */ 9ace3abfcSBarry Smith PetscBool extrapolate; 10*eb284becSJed 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; 222b5a38e1SLisandro Dalcin PetscErrorCode ierr; 23316643e7SJed Brown 24316643e7SJed Brown PetscFunctionBegin; 25193ac0bcSJed Brown ts->time_step = ts->next_time_step; 26*eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 27316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 28316643e7SJed Brown 29*eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 30*eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 31*eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 32*eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 33*eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 34*eb284becSJed Brown } 35ace68cafSJed Brown if (th->extrapolate) { 362b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 37ace68cafSJed Brown } else { 382b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 39ace68cafSJed Brown } 40*eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 41316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 42316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 43316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 442b5a38e1SLisandro Dalcin 45*eb284becSJed Brown if (th->endpoint) { 46*eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 47*eb284becSJed Brown } else { 482b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 492b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 50*eb284becSJed Brown } 512b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 52193ac0bcSJed Brown ts->next_time_step = ts->time_step; 53316643e7SJed Brown ts->steps++; 54316643e7SJed Brown PetscFunctionReturn(0); 55316643e7SJed Brown } 56316643e7SJed Brown 57cd652676SJed Brown #undef __FUNCT__ 58cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 59cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 60cd652676SJed Brown { 61cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 62cd652676SJed Brown PetscErrorCode ierr; 63cd652676SJed Brown 64cd652676SJed Brown PetscFunctionBegin; 65cd652676SJed Brown ierr = VecWAXPY(X,t-ts->ptime,ts->vec_sol,th->Xdot);CHKERRQ(ierr); 66cd652676SJed Brown PetscFunctionReturn(0); 67cd652676SJed Brown } 68cd652676SJed Brown 69316643e7SJed Brown /*------------------------------------------------------------*/ 70316643e7SJed Brown #undef __FUNCT__ 71277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 72277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 73316643e7SJed Brown { 74316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 75316643e7SJed Brown PetscErrorCode ierr; 76316643e7SJed Brown 77316643e7SJed Brown PetscFunctionBegin; 786bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 796bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 80*eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 81277b19d0SLisandro Dalcin PetscFunctionReturn(0); 82277b19d0SLisandro Dalcin } 83277b19d0SLisandro Dalcin 84277b19d0SLisandro Dalcin #undef __FUNCT__ 85277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 86277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 87277b19d0SLisandro Dalcin { 88277b19d0SLisandro Dalcin PetscErrorCode ierr; 89277b19d0SLisandro Dalcin 90277b19d0SLisandro Dalcin PetscFunctionBegin; 91277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 92277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 93335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 94335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 95*eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 96316643e7SJed Brown PetscFunctionReturn(0); 97316643e7SJed Brown } 98316643e7SJed Brown 99316643e7SJed Brown /* 100316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1012b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 102316643e7SJed Brown */ 103316643e7SJed Brown #undef __FUNCT__ 1040f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1050f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 106316643e7SJed Brown { 107316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 108316643e7SJed Brown PetscErrorCode ierr; 109316643e7SJed Brown 110316643e7SJed Brown PetscFunctionBegin; 1112b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 112214bc6a2SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 113316643e7SJed Brown PetscFunctionReturn(0); 114316643e7SJed Brown } 115316643e7SJed Brown 116316643e7SJed Brown #undef __FUNCT__ 1170f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1180f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 119316643e7SJed Brown { 120316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 121316643e7SJed Brown PetscErrorCode ierr; 122316643e7SJed Brown 123316643e7SJed Brown PetscFunctionBegin; 1240f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 125214bc6a2SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 126316643e7SJed Brown PetscFunctionReturn(0); 127316643e7SJed Brown } 128316643e7SJed Brown 129316643e7SJed Brown 130316643e7SJed Brown #undef __FUNCT__ 131316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 132316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 133316643e7SJed Brown { 134316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 135316643e7SJed Brown PetscErrorCode ierr; 136316643e7SJed Brown 137316643e7SJed Brown PetscFunctionBegin; 138316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 139316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 140316643e7SJed Brown PetscFunctionReturn(0); 141316643e7SJed Brown } 142316643e7SJed Brown /*------------------------------------------------------------*/ 143316643e7SJed Brown 144316643e7SJed Brown #undef __FUNCT__ 145316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 146316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 147316643e7SJed Brown { 148316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 149316643e7SJed Brown PetscErrorCode ierr; 150316643e7SJed Brown 151316643e7SJed Brown PetscFunctionBegin; 152d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 153316643e7SJed Brown { 154316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 155acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 156*eb284becSJed 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); 157316643e7SJed Brown } 158316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 159316643e7SJed Brown PetscFunctionReturn(0); 160316643e7SJed Brown } 161316643e7SJed Brown 162316643e7SJed Brown #undef __FUNCT__ 163316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 164316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 165316643e7SJed Brown { 166316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 167ace3abfcSBarry Smith PetscBool iascii; 168316643e7SJed Brown PetscErrorCode ierr; 169316643e7SJed Brown 170316643e7SJed Brown PetscFunctionBegin; 1712692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 172316643e7SJed Brown if (iascii) { 173316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 174ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 175316643e7SJed Brown } 176316643e7SJed Brown PetscFunctionReturn(0); 177316643e7SJed Brown } 178316643e7SJed Brown 1790de4c49aSJed Brown EXTERN_C_BEGIN 1800de4c49aSJed Brown #undef __FUNCT__ 1810de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 1827087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1830de4c49aSJed Brown { 1840de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1850de4c49aSJed Brown 1860de4c49aSJed Brown PetscFunctionBegin; 1870de4c49aSJed Brown *theta = th->Theta; 1880de4c49aSJed Brown PetscFunctionReturn(0); 1890de4c49aSJed Brown } 1900de4c49aSJed Brown 1910de4c49aSJed Brown #undef __FUNCT__ 1920de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 1937087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 1940de4c49aSJed Brown { 1950de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1960de4c49aSJed Brown 1970de4c49aSJed Brown PetscFunctionBegin; 198e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 1990de4c49aSJed Brown th->Theta = theta; 2000de4c49aSJed Brown PetscFunctionReturn(0); 2010de4c49aSJed Brown } 202*eb284becSJed Brown 203*eb284becSJed Brown #undef __FUNCT__ 204*eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint_Theta" 205*eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 206*eb284becSJed Brown { 207*eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 208*eb284becSJed Brown 209*eb284becSJed Brown PetscFunctionBegin; 210*eb284becSJed Brown th->endpoint = flg; 211*eb284becSJed Brown PetscFunctionReturn(0); 212*eb284becSJed Brown } 2130de4c49aSJed Brown EXTERN_C_END 2140de4c49aSJed Brown 215316643e7SJed Brown /* ------------------------------------------------------------ */ 216316643e7SJed Brown /*MC 21796f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 218316643e7SJed Brown 219316643e7SJed Brown Level: beginner 220316643e7SJed Brown 221*eb284becSJed Brown Notes: 222*eb284becSJed Brown This method can be applied to DAE. 223*eb284becSJed Brown 224*eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 225*eb284becSJed Brown 226*eb284becSJed Brown .vb 227*eb284becSJed Brown Theta | Theta 228*eb284becSJed Brown ------------- 229*eb284becSJed Brown | 1 230*eb284becSJed Brown .ve 231*eb284becSJed Brown 232*eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 233*eb284becSJed Brown 234*eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 235*eb284becSJed Brown 236*eb284becSJed Brown .vb 237*eb284becSJed Brown 0 | 0 0 238*eb284becSJed Brown 1 | 1-Theta Theta 239*eb284becSJed Brown ------------------- 240*eb284becSJed Brown | 1-Theta Theta 241*eb284becSJed Brown .ve 242*eb284becSJed Brown 243*eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 244*eb284becSJed Brown 245*eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 246*eb284becSJed Brown 247*eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 248*eb284becSJed Brown 249*eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 250*eb284becSJed Brown 251*eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 252316643e7SJed Brown 253316643e7SJed Brown M*/ 254316643e7SJed Brown EXTERN_C_BEGIN 255316643e7SJed Brown #undef __FUNCT__ 256316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2577087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 258316643e7SJed Brown { 259316643e7SJed Brown TS_Theta *th; 260316643e7SJed Brown PetscErrorCode ierr; 261316643e7SJed Brown 262316643e7SJed Brown PetscFunctionBegin; 263277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 264316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 265316643e7SJed Brown ts->ops->view = TSView_Theta; 266316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 267316643e7SJed Brown ts->ops->step = TSStep_Theta; 268cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 269316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2700f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2710f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 272316643e7SJed Brown 273316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 274316643e7SJed Brown ts->data = (void*)th; 275316643e7SJed Brown 2766f700aefSJed Brown th->extrapolate = PETSC_FALSE; 277316643e7SJed Brown th->Theta = 0.5; 278316643e7SJed Brown 2790de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 2800de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 281*eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 282316643e7SJed Brown PetscFunctionReturn(0); 283316643e7SJed Brown } 284316643e7SJed Brown EXTERN_C_END 2850de4c49aSJed Brown 2860de4c49aSJed Brown #undef __FUNCT__ 2870de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 2880de4c49aSJed Brown /*@ 2890de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 2900de4c49aSJed Brown 2910de4c49aSJed Brown Not Collective 2920de4c49aSJed Brown 2930de4c49aSJed Brown Input Parameter: 2940de4c49aSJed Brown . ts - timestepping context 2950de4c49aSJed Brown 2960de4c49aSJed Brown Output Parameter: 2970de4c49aSJed Brown . theta - stage abscissa 2980de4c49aSJed Brown 2990de4c49aSJed Brown Note: 3000de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 3010de4c49aSJed Brown 3020de4c49aSJed Brown Level: Advanced 3030de4c49aSJed Brown 3040de4c49aSJed Brown .seealso: TSThetaSetTheta() 3050de4c49aSJed Brown @*/ 3067087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 3070de4c49aSJed Brown { 3084ac538c5SBarry Smith PetscErrorCode ierr; 3090de4c49aSJed Brown 3100de4c49aSJed Brown PetscFunctionBegin; 311afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3120de4c49aSJed Brown PetscValidPointer(theta,2); 3134ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 3140de4c49aSJed Brown PetscFunctionReturn(0); 3150de4c49aSJed Brown } 3160de4c49aSJed Brown 3170de4c49aSJed Brown #undef __FUNCT__ 3180de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 3190de4c49aSJed Brown /*@ 3200de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 3210de4c49aSJed Brown 3220de4c49aSJed Brown Not Collective 3230de4c49aSJed Brown 3240de4c49aSJed Brown Input Parameter: 3250de4c49aSJed Brown + ts - timestepping context 3260de4c49aSJed Brown - theta - stage abscissa 3270de4c49aSJed Brown 3280de4c49aSJed Brown Options Database: 3290de4c49aSJed Brown . -ts_theta_theta <theta> 3300de4c49aSJed Brown 3310de4c49aSJed Brown Level: Intermediate 3320de4c49aSJed Brown 3330de4c49aSJed Brown .seealso: TSThetaGetTheta() 3340de4c49aSJed Brown @*/ 3357087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 3360de4c49aSJed Brown { 3374ac538c5SBarry Smith PetscErrorCode ierr; 3380de4c49aSJed Brown 3390de4c49aSJed Brown PetscFunctionBegin; 340afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3414ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 3420de4c49aSJed Brown PetscFunctionReturn(0); 3430de4c49aSJed Brown } 344f33bbcb6SJed Brown 345*eb284becSJed Brown #undef __FUNCT__ 346*eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 347*eb284becSJed Brown /*@ 348*eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 349*eb284becSJed Brown 350*eb284becSJed Brown Not Collective 351*eb284becSJed Brown 352*eb284becSJed Brown Input Parameter: 353*eb284becSJed Brown + ts - timestepping context 354*eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 355*eb284becSJed Brown 356*eb284becSJed Brown Options Database: 357*eb284becSJed Brown . -ts_theta_endpoint <flg> 358*eb284becSJed Brown 359*eb284becSJed Brown Level: Intermediate 360*eb284becSJed Brown 361*eb284becSJed Brown .seealso: TSTHETA, TSCN 362*eb284becSJed Brown @*/ 363*eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 364*eb284becSJed Brown { 365*eb284becSJed Brown PetscErrorCode ierr; 366*eb284becSJed Brown 367*eb284becSJed Brown PetscFunctionBegin; 368*eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 369*eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 370*eb284becSJed Brown PetscFunctionReturn(0); 371*eb284becSJed Brown } 372*eb284becSJed Brown 373f33bbcb6SJed Brown /* 374f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 375f33bbcb6SJed Brown * The creation functions for these specializations are below. 376f33bbcb6SJed Brown */ 377f33bbcb6SJed Brown 378f33bbcb6SJed Brown #undef __FUNCT__ 379f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 380f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 381f33bbcb6SJed Brown { 382f33bbcb6SJed Brown PetscFunctionBegin; 383f33bbcb6SJed Brown PetscFunctionReturn(0); 384f33bbcb6SJed Brown } 385f33bbcb6SJed Brown 386f33bbcb6SJed Brown /*MC 387f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 388f33bbcb6SJed Brown 389f33bbcb6SJed Brown Level: beginner 390f33bbcb6SJed Brown 391f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 392f33bbcb6SJed Brown 393f33bbcb6SJed Brown M*/ 394f33bbcb6SJed Brown EXTERN_C_BEGIN 395f33bbcb6SJed Brown #undef __FUNCT__ 396f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 397f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 398f33bbcb6SJed Brown { 399f33bbcb6SJed Brown PetscErrorCode ierr; 400f33bbcb6SJed Brown 401f33bbcb6SJed Brown PetscFunctionBegin; 402f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 403f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 404f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 405f33bbcb6SJed Brown PetscFunctionReturn(0); 406f33bbcb6SJed Brown } 407f33bbcb6SJed Brown EXTERN_C_END 408f33bbcb6SJed Brown 409f33bbcb6SJed Brown #undef __FUNCT__ 410f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 411f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 412f33bbcb6SJed Brown { 413f33bbcb6SJed Brown PetscFunctionBegin; 414f33bbcb6SJed Brown PetscFunctionReturn(0); 415f33bbcb6SJed Brown } 416f33bbcb6SJed Brown 417f33bbcb6SJed Brown /*MC 418f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 419f33bbcb6SJed Brown 420f33bbcb6SJed Brown Level: beginner 421f33bbcb6SJed Brown 422f33bbcb6SJed Brown Notes: 423f33bbcb6SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 424f33bbcb6SJed Brown 425f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 426f33bbcb6SJed Brown 427f33bbcb6SJed Brown M*/ 428f33bbcb6SJed Brown EXTERN_C_BEGIN 429f33bbcb6SJed Brown #undef __FUNCT__ 430f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 431f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 432f33bbcb6SJed Brown { 433f33bbcb6SJed Brown PetscErrorCode ierr; 434f33bbcb6SJed Brown 435f33bbcb6SJed Brown PetscFunctionBegin; 436f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 437f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 438*eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 439f33bbcb6SJed Brown ts->ops->view = TSView_CN; 440f33bbcb6SJed Brown PetscFunctionReturn(0); 441f33bbcb6SJed Brown } 442f33bbcb6SJed Brown EXTERN_C_END 443