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; 222b5a38e1SLisandro Dalcin PetscErrorCode ierr; 23316643e7SJed Brown 24316643e7SJed Brown PetscFunctionBegin; 25193ac0bcSJed Brown ts->time_step = ts->next_time_step; 26eb284becSJed 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 29eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 30eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 31eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 32eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 33eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 34eb284becSJed 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 } 40eb284becSJed 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 45eb284becSJed Brown if (th->endpoint) { 46eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 47eb284becSJed 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); 50eb284becSJed 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; 625a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 63cd652676SJed Brown PetscErrorCode ierr; 64cd652676SJed Brown 65cd652676SJed Brown PetscFunctionBegin; 66a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 675a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 685a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 69cd652676SJed Brown PetscFunctionReturn(0); 70cd652676SJed Brown } 71cd652676SJed Brown 72316643e7SJed Brown /*------------------------------------------------------------*/ 73316643e7SJed Brown #undef __FUNCT__ 74277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 75277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 76316643e7SJed Brown { 77316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 78316643e7SJed Brown PetscErrorCode ierr; 79316643e7SJed Brown 80316643e7SJed Brown PetscFunctionBegin; 816bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 826bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 83eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 84277b19d0SLisandro Dalcin PetscFunctionReturn(0); 85277b19d0SLisandro Dalcin } 86277b19d0SLisandro Dalcin 87277b19d0SLisandro Dalcin #undef __FUNCT__ 88277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 89277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 90277b19d0SLisandro Dalcin { 91277b19d0SLisandro Dalcin PetscErrorCode ierr; 92277b19d0SLisandro Dalcin 93277b19d0SLisandro Dalcin PetscFunctionBegin; 94277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 95277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 96335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 97335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 98eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 99316643e7SJed Brown PetscFunctionReturn(0); 100316643e7SJed Brown } 101316643e7SJed Brown 102316643e7SJed Brown /* 103316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1042b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 105316643e7SJed Brown */ 106316643e7SJed Brown #undef __FUNCT__ 1070f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1080f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 109316643e7SJed Brown { 110316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 111316643e7SJed Brown PetscErrorCode ierr; 112316643e7SJed Brown 113316643e7SJed Brown PetscFunctionBegin; 1145a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 1152b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 116214bc6a2SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 117316643e7SJed Brown PetscFunctionReturn(0); 118316643e7SJed Brown } 119316643e7SJed Brown 120316643e7SJed Brown #undef __FUNCT__ 1210f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1220f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 123316643e7SJed Brown { 124316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 125316643e7SJed Brown PetscErrorCode ierr; 126316643e7SJed Brown 127316643e7SJed Brown PetscFunctionBegin; 1280f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 129214bc6a2SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 130316643e7SJed Brown PetscFunctionReturn(0); 131316643e7SJed Brown } 132316643e7SJed Brown 133316643e7SJed Brown 134316643e7SJed Brown #undef __FUNCT__ 135316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 136316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 137316643e7SJed Brown { 138316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 139316643e7SJed Brown PetscErrorCode ierr; 140316643e7SJed Brown 141316643e7SJed Brown PetscFunctionBegin; 142316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 143316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 144316643e7SJed Brown PetscFunctionReturn(0); 145316643e7SJed Brown } 146316643e7SJed Brown /*------------------------------------------------------------*/ 147316643e7SJed Brown 148316643e7SJed Brown #undef __FUNCT__ 149316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 150316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 151316643e7SJed Brown { 152316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 153316643e7SJed Brown PetscErrorCode ierr; 154316643e7SJed Brown 155316643e7SJed Brown PetscFunctionBegin; 156d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 157316643e7SJed Brown { 158316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 159acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 160eb284becSJed 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); 161*d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 162316643e7SJed Brown } 163316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 164316643e7SJed Brown PetscFunctionReturn(0); 165316643e7SJed Brown } 166316643e7SJed Brown 167316643e7SJed Brown #undef __FUNCT__ 168316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 169316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 170316643e7SJed Brown { 171316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 172ace3abfcSBarry Smith PetscBool iascii; 173316643e7SJed Brown PetscErrorCode ierr; 174316643e7SJed Brown 175316643e7SJed Brown PetscFunctionBegin; 1762692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 177316643e7SJed Brown if (iascii) { 178316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 179ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 180316643e7SJed Brown } 181*d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 182316643e7SJed Brown PetscFunctionReturn(0); 183316643e7SJed Brown } 184316643e7SJed Brown 1850de4c49aSJed Brown EXTERN_C_BEGIN 1860de4c49aSJed Brown #undef __FUNCT__ 1870de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 1887087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1890de4c49aSJed Brown { 1900de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1910de4c49aSJed Brown 1920de4c49aSJed Brown PetscFunctionBegin; 1930de4c49aSJed Brown *theta = th->Theta; 1940de4c49aSJed Brown PetscFunctionReturn(0); 1950de4c49aSJed Brown } 1960de4c49aSJed Brown 1970de4c49aSJed Brown #undef __FUNCT__ 1980de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 1997087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 2000de4c49aSJed Brown { 2010de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2020de4c49aSJed Brown 2030de4c49aSJed Brown PetscFunctionBegin; 204e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 2050de4c49aSJed Brown th->Theta = theta; 2060de4c49aSJed Brown PetscFunctionReturn(0); 2070de4c49aSJed Brown } 208eb284becSJed Brown 209eb284becSJed Brown #undef __FUNCT__ 210eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint_Theta" 211eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 212eb284becSJed Brown { 213eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 214eb284becSJed Brown 215eb284becSJed Brown PetscFunctionBegin; 216eb284becSJed Brown th->endpoint = flg; 217eb284becSJed Brown PetscFunctionReturn(0); 218eb284becSJed Brown } 2190de4c49aSJed Brown EXTERN_C_END 2200de4c49aSJed Brown 221316643e7SJed Brown /* ------------------------------------------------------------ */ 222316643e7SJed Brown /*MC 22396f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 224316643e7SJed Brown 225316643e7SJed Brown Level: beginner 226316643e7SJed Brown 227eb284becSJed Brown Notes: 228eb284becSJed Brown This method can be applied to DAE. 229eb284becSJed Brown 230eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 231eb284becSJed Brown 232eb284becSJed Brown .vb 233eb284becSJed Brown Theta | Theta 234eb284becSJed Brown ------------- 235eb284becSJed Brown | 1 236eb284becSJed Brown .ve 237eb284becSJed Brown 238eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 239eb284becSJed Brown 240eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 241eb284becSJed Brown 242eb284becSJed Brown .vb 243eb284becSJed Brown 0 | 0 0 244eb284becSJed Brown 1 | 1-Theta Theta 245eb284becSJed Brown ------------------- 246eb284becSJed Brown | 1-Theta Theta 247eb284becSJed Brown .ve 248eb284becSJed Brown 249eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 250eb284becSJed Brown 251eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 252eb284becSJed Brown 253eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 254eb284becSJed Brown 255eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 256eb284becSJed Brown 257eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 258316643e7SJed Brown 259316643e7SJed Brown M*/ 260316643e7SJed Brown EXTERN_C_BEGIN 261316643e7SJed Brown #undef __FUNCT__ 262316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2637087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 264316643e7SJed Brown { 265316643e7SJed Brown TS_Theta *th; 266316643e7SJed Brown PetscErrorCode ierr; 267316643e7SJed Brown 268316643e7SJed Brown PetscFunctionBegin; 269277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 270316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 271316643e7SJed Brown ts->ops->view = TSView_Theta; 272316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 273316643e7SJed Brown ts->ops->step = TSStep_Theta; 274cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 275316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2760f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2770f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 278316643e7SJed Brown 279316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 280316643e7SJed Brown ts->data = (void*)th; 281316643e7SJed Brown 2826f700aefSJed Brown th->extrapolate = PETSC_FALSE; 283316643e7SJed Brown th->Theta = 0.5; 284316643e7SJed Brown 2850de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 2860de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 287eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 288316643e7SJed Brown PetscFunctionReturn(0); 289316643e7SJed Brown } 290316643e7SJed Brown EXTERN_C_END 2910de4c49aSJed Brown 2920de4c49aSJed Brown #undef __FUNCT__ 2930de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 2940de4c49aSJed Brown /*@ 2950de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 2960de4c49aSJed Brown 2970de4c49aSJed Brown Not Collective 2980de4c49aSJed Brown 2990de4c49aSJed Brown Input Parameter: 3000de4c49aSJed Brown . ts - timestepping context 3010de4c49aSJed Brown 3020de4c49aSJed Brown Output Parameter: 3030de4c49aSJed Brown . theta - stage abscissa 3040de4c49aSJed Brown 3050de4c49aSJed Brown Note: 3060de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 3070de4c49aSJed Brown 3080de4c49aSJed Brown Level: Advanced 3090de4c49aSJed Brown 3100de4c49aSJed Brown .seealso: TSThetaSetTheta() 3110de4c49aSJed Brown @*/ 3127087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 3130de4c49aSJed Brown { 3144ac538c5SBarry Smith PetscErrorCode ierr; 3150de4c49aSJed Brown 3160de4c49aSJed Brown PetscFunctionBegin; 317afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3180de4c49aSJed Brown PetscValidPointer(theta,2); 3194ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 3200de4c49aSJed Brown PetscFunctionReturn(0); 3210de4c49aSJed Brown } 3220de4c49aSJed Brown 3230de4c49aSJed Brown #undef __FUNCT__ 3240de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 3250de4c49aSJed Brown /*@ 3260de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 3270de4c49aSJed Brown 3280de4c49aSJed Brown Not Collective 3290de4c49aSJed Brown 3300de4c49aSJed Brown Input Parameter: 3310de4c49aSJed Brown + ts - timestepping context 3320de4c49aSJed Brown - theta - stage abscissa 3330de4c49aSJed Brown 3340de4c49aSJed Brown Options Database: 3350de4c49aSJed Brown . -ts_theta_theta <theta> 3360de4c49aSJed Brown 3370de4c49aSJed Brown Level: Intermediate 3380de4c49aSJed Brown 3390de4c49aSJed Brown .seealso: TSThetaGetTheta() 3400de4c49aSJed Brown @*/ 3417087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 3420de4c49aSJed Brown { 3434ac538c5SBarry Smith PetscErrorCode ierr; 3440de4c49aSJed Brown 3450de4c49aSJed Brown PetscFunctionBegin; 346afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3474ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 3480de4c49aSJed Brown PetscFunctionReturn(0); 3490de4c49aSJed Brown } 350f33bbcb6SJed Brown 351eb284becSJed Brown #undef __FUNCT__ 352eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 353eb284becSJed Brown /*@ 354eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 355eb284becSJed Brown 356eb284becSJed Brown Not Collective 357eb284becSJed Brown 358eb284becSJed Brown Input Parameter: 359eb284becSJed Brown + ts - timestepping context 360eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 361eb284becSJed Brown 362eb284becSJed Brown Options Database: 363eb284becSJed Brown . -ts_theta_endpoint <flg> 364eb284becSJed Brown 365eb284becSJed Brown Level: Intermediate 366eb284becSJed Brown 367eb284becSJed Brown .seealso: TSTHETA, TSCN 368eb284becSJed Brown @*/ 369eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 370eb284becSJed Brown { 371eb284becSJed Brown PetscErrorCode ierr; 372eb284becSJed Brown 373eb284becSJed Brown PetscFunctionBegin; 374eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 375eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 376eb284becSJed Brown PetscFunctionReturn(0); 377eb284becSJed Brown } 378eb284becSJed Brown 379f33bbcb6SJed Brown /* 380f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 381f33bbcb6SJed Brown * The creation functions for these specializations are below. 382f33bbcb6SJed Brown */ 383f33bbcb6SJed Brown 384f33bbcb6SJed Brown #undef __FUNCT__ 385f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 386f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 387f33bbcb6SJed Brown { 388*d52bd9f3SBarry Smith PetscErrorCode ierr; 389*d52bd9f3SBarry Smith 390f33bbcb6SJed Brown PetscFunctionBegin; 391*d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 392f33bbcb6SJed Brown PetscFunctionReturn(0); 393f33bbcb6SJed Brown } 394f33bbcb6SJed Brown 395f33bbcb6SJed Brown /*MC 396f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 397f33bbcb6SJed Brown 398f33bbcb6SJed Brown Level: beginner 399f33bbcb6SJed Brown 400f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 401f33bbcb6SJed Brown 402f33bbcb6SJed Brown M*/ 403f33bbcb6SJed Brown EXTERN_C_BEGIN 404f33bbcb6SJed Brown #undef __FUNCT__ 405f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 406f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 407f33bbcb6SJed Brown { 408f33bbcb6SJed Brown PetscErrorCode ierr; 409f33bbcb6SJed Brown 410f33bbcb6SJed Brown PetscFunctionBegin; 411f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 412f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 413f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 414f33bbcb6SJed Brown PetscFunctionReturn(0); 415f33bbcb6SJed Brown } 416f33bbcb6SJed Brown EXTERN_C_END 417f33bbcb6SJed Brown 418f33bbcb6SJed Brown #undef __FUNCT__ 419f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 420f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 421f33bbcb6SJed Brown { 422*d52bd9f3SBarry Smith PetscErrorCode ierr; 423*d52bd9f3SBarry Smith 424f33bbcb6SJed Brown PetscFunctionBegin; 425*d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 426f33bbcb6SJed Brown PetscFunctionReturn(0); 427f33bbcb6SJed Brown } 428f33bbcb6SJed Brown 429f33bbcb6SJed Brown /*MC 430f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 431f33bbcb6SJed Brown 432f33bbcb6SJed Brown Level: beginner 433f33bbcb6SJed Brown 434f33bbcb6SJed Brown Notes: 4357cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 4367cf5af47SJed Brown 4377cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 438f33bbcb6SJed Brown 439f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 440f33bbcb6SJed Brown 441f33bbcb6SJed Brown M*/ 442f33bbcb6SJed Brown EXTERN_C_BEGIN 443f33bbcb6SJed Brown #undef __FUNCT__ 444f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 445f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 446f33bbcb6SJed Brown { 447f33bbcb6SJed Brown PetscErrorCode ierr; 448f33bbcb6SJed Brown 449f33bbcb6SJed Brown PetscFunctionBegin; 450f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 451f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 452eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 453f33bbcb6SJed Brown ts->ops->view = TSView_CN; 454f33bbcb6SJed Brown PetscFunctionReturn(0); 455f33bbcb6SJed Brown } 456f33bbcb6SJed Brown EXTERN_C_END 457