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; 23*f1b97656SJed Brown SNESConvergedReason snesreason; 242b5a38e1SLisandro Dalcin PetscErrorCode ierr; 25316643e7SJed Brown 26316643e7SJed Brown PetscFunctionBegin; 27cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 28eb284becSJed Brown th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 29316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 30316643e7SJed Brown 31eb284becSJed Brown if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 32eb284becSJed Brown ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 33eb284becSJed Brown if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 34eb284becSJed Brown ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 35eb284becSJed Brown ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 36eb284becSJed Brown } 37ace68cafSJed Brown if (th->extrapolate) { 382b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 39ace68cafSJed Brown } else { 402b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 41ace68cafSJed Brown } 42eb284becSJed Brown ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 43316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 44316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 45*f1b97656SJed Brown ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 46316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 47*f1b97656SJed Brown if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 48*f1b97656SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 49*f1b97656SJed 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); 50*f1b97656SJed Brown PetscFunctionReturn(0); 51*f1b97656SJed Brown } 52eb284becSJed Brown if (th->endpoint) { 53eb284becSJed Brown ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 54eb284becSJed Brown } else { 552b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 562b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 57eb284becSJed Brown } 582b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 59cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 60316643e7SJed Brown ts->steps++; 61316643e7SJed Brown PetscFunctionReturn(0); 62316643e7SJed Brown } 63316643e7SJed Brown 64cd652676SJed Brown #undef __FUNCT__ 65cd652676SJed Brown #define __FUNCT__ "TSInterpolate_Theta" 66cd652676SJed Brown static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 67cd652676SJed Brown { 68cd652676SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 695a3a76d0SJed Brown PetscReal alpha = t - ts->ptime; 70cd652676SJed Brown PetscErrorCode ierr; 71cd652676SJed Brown 72cd652676SJed Brown PetscFunctionBegin; 73a43b19c4SJed Brown ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 745a3a76d0SJed Brown if (th->endpoint) alpha *= th->Theta; 755a3a76d0SJed Brown ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 76cd652676SJed Brown PetscFunctionReturn(0); 77cd652676SJed Brown } 78cd652676SJed Brown 79316643e7SJed Brown /*------------------------------------------------------------*/ 80316643e7SJed Brown #undef __FUNCT__ 81277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 82277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 83316643e7SJed Brown { 84316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 85316643e7SJed Brown PetscErrorCode ierr; 86316643e7SJed Brown 87316643e7SJed Brown PetscFunctionBegin; 886bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 896bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 90eb284becSJed Brown ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 91277b19d0SLisandro Dalcin PetscFunctionReturn(0); 92277b19d0SLisandro Dalcin } 93277b19d0SLisandro Dalcin 94277b19d0SLisandro Dalcin #undef __FUNCT__ 95277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_Theta" 96277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_Theta(TS ts) 97277b19d0SLisandro Dalcin { 98277b19d0SLisandro Dalcin PetscErrorCode ierr; 99277b19d0SLisandro Dalcin 100277b19d0SLisandro Dalcin PetscFunctionBegin; 101277b19d0SLisandro Dalcin ierr = TSReset_Theta(ts);CHKERRQ(ierr); 102277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 103335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 104335f802eSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 10526f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 106eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 107316643e7SJed Brown PetscFunctionReturn(0); 108316643e7SJed Brown } 109316643e7SJed Brown 110316643e7SJed Brown /* 111316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1122b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 113316643e7SJed Brown */ 114316643e7SJed Brown #undef __FUNCT__ 1150f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1160f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 117316643e7SJed Brown { 118316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 119316643e7SJed Brown PetscErrorCode ierr; 120316643e7SJed Brown 121316643e7SJed Brown PetscFunctionBegin; 1225a3a76d0SJed Brown /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 1232b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 124214bc6a2SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 125316643e7SJed Brown PetscFunctionReturn(0); 126316643e7SJed Brown } 127316643e7SJed Brown 128316643e7SJed Brown #undef __FUNCT__ 1290f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1300f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 131316643e7SJed Brown { 132316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 133316643e7SJed Brown PetscErrorCode ierr; 134316643e7SJed Brown 135316643e7SJed Brown PetscFunctionBegin; 1360f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 137214bc6a2SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 138316643e7SJed Brown PetscFunctionReturn(0); 139316643e7SJed Brown } 140316643e7SJed Brown 141316643e7SJed Brown 142316643e7SJed Brown #undef __FUNCT__ 143316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 144316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 145316643e7SJed Brown { 146316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 147316643e7SJed Brown PetscErrorCode ierr; 148316643e7SJed Brown 149316643e7SJed Brown PetscFunctionBegin; 150316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 151316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 152316643e7SJed Brown PetscFunctionReturn(0); 153316643e7SJed Brown } 154316643e7SJed Brown /*------------------------------------------------------------*/ 155316643e7SJed Brown 156316643e7SJed Brown #undef __FUNCT__ 157316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 158316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 159316643e7SJed Brown { 160316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 161316643e7SJed Brown PetscErrorCode ierr; 162316643e7SJed Brown 163316643e7SJed Brown PetscFunctionBegin; 164d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 165316643e7SJed Brown { 166316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 167acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 168eb284becSJed 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); 169d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 170316643e7SJed Brown } 171316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 172316643e7SJed Brown PetscFunctionReturn(0); 173316643e7SJed Brown } 174316643e7SJed Brown 175316643e7SJed Brown #undef __FUNCT__ 176316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 177316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 178316643e7SJed Brown { 179316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 180ace3abfcSBarry Smith PetscBool iascii; 181316643e7SJed Brown PetscErrorCode ierr; 182316643e7SJed Brown 183316643e7SJed Brown PetscFunctionBegin; 1842692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 185316643e7SJed Brown if (iascii) { 186316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 187ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 188316643e7SJed Brown } 189d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 190316643e7SJed Brown PetscFunctionReturn(0); 191316643e7SJed Brown } 192316643e7SJed Brown 1930de4c49aSJed Brown EXTERN_C_BEGIN 1940de4c49aSJed Brown #undef __FUNCT__ 1950de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 1967087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 1970de4c49aSJed Brown { 1980de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 1990de4c49aSJed Brown 2000de4c49aSJed Brown PetscFunctionBegin; 2010de4c49aSJed Brown *theta = th->Theta; 2020de4c49aSJed Brown PetscFunctionReturn(0); 2030de4c49aSJed Brown } 2040de4c49aSJed Brown 2050de4c49aSJed Brown #undef __FUNCT__ 2060de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 2077087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 2080de4c49aSJed Brown { 2090de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2100de4c49aSJed Brown 2110de4c49aSJed Brown PetscFunctionBegin; 212e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 2130de4c49aSJed Brown th->Theta = theta; 2140de4c49aSJed Brown PetscFunctionReturn(0); 2150de4c49aSJed Brown } 216eb284becSJed Brown 217eb284becSJed Brown #undef __FUNCT__ 21878e224dfSJed Brown #define __FUNCT__ "TSThetaGetEndpoint_Theta" 21926f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 22026f2ff8fSLisandro Dalcin { 22126f2ff8fSLisandro Dalcin TS_Theta *th = (TS_Theta*)ts->data; 22226f2ff8fSLisandro Dalcin 22326f2ff8fSLisandro Dalcin PetscFunctionBegin; 22426f2ff8fSLisandro Dalcin *endpoint = th->endpoint; 22526f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 22626f2ff8fSLisandro Dalcin } 22726f2ff8fSLisandro Dalcin 22826f2ff8fSLisandro Dalcin #undef __FUNCT__ 22926f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaSetEndpoint_Theta" 230eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 231eb284becSJed Brown { 232eb284becSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 233eb284becSJed Brown 234eb284becSJed Brown PetscFunctionBegin; 235eb284becSJed Brown th->endpoint = flg; 236eb284becSJed Brown PetscFunctionReturn(0); 237eb284becSJed Brown } 2380de4c49aSJed Brown EXTERN_C_END 2390de4c49aSJed Brown 240316643e7SJed Brown /* ------------------------------------------------------------ */ 241316643e7SJed Brown /*MC 24296f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 243316643e7SJed Brown 244316643e7SJed Brown Level: beginner 245316643e7SJed Brown 246eb284becSJed Brown Notes: 247eb284becSJed Brown This method can be applied to DAE. 248eb284becSJed Brown 249eb284becSJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 250eb284becSJed Brown 251eb284becSJed Brown .vb 252eb284becSJed Brown Theta | Theta 253eb284becSJed Brown ------------- 254eb284becSJed Brown | 1 255eb284becSJed Brown .ve 256eb284becSJed Brown 257eb284becSJed Brown For the default Theta=0.5, this is also known as the implicit midpoint rule. 258eb284becSJed Brown 259eb284becSJed Brown When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 260eb284becSJed Brown 261eb284becSJed Brown .vb 262eb284becSJed Brown 0 | 0 0 263eb284becSJed Brown 1 | 1-Theta Theta 264eb284becSJed Brown ------------------- 265eb284becSJed Brown | 1-Theta Theta 266eb284becSJed Brown .ve 267eb284becSJed Brown 268eb284becSJed Brown For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 269eb284becSJed Brown 270eb284becSJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 271eb284becSJed Brown 272eb284becSJed Brown $ Y_i = X + h sum_j a_ij Y'_j 273eb284becSJed Brown 274eb284becSJed Brown is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 275eb284becSJed Brown 276eb284becSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 277316643e7SJed Brown 278316643e7SJed Brown M*/ 279316643e7SJed Brown EXTERN_C_BEGIN 280316643e7SJed Brown #undef __FUNCT__ 281316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2827087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 283316643e7SJed Brown { 284316643e7SJed Brown TS_Theta *th; 285316643e7SJed Brown PetscErrorCode ierr; 286316643e7SJed Brown 287316643e7SJed Brown PetscFunctionBegin; 288277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 289316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 290316643e7SJed Brown ts->ops->view = TSView_Theta; 291316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 292316643e7SJed Brown ts->ops->step = TSStep_Theta; 293cd652676SJed Brown ts->ops->interpolate = TSInterpolate_Theta; 294316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2950f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2960f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 297316643e7SJed Brown 298316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 299316643e7SJed Brown ts->data = (void*)th; 300316643e7SJed Brown 3016f700aefSJed Brown th->extrapolate = PETSC_FALSE; 302316643e7SJed Brown th->Theta = 0.5; 303316643e7SJed Brown 3040de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 3050de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 30626f2ff8fSLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetEndpoint_C","TSThetaGetEndpoint_Theta",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 307eb284becSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 308316643e7SJed Brown PetscFunctionReturn(0); 309316643e7SJed Brown } 310316643e7SJed Brown EXTERN_C_END 3110de4c49aSJed Brown 3120de4c49aSJed Brown #undef __FUNCT__ 3130de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 3140de4c49aSJed Brown /*@ 3150de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 3160de4c49aSJed Brown 3170de4c49aSJed Brown Not Collective 3180de4c49aSJed Brown 3190de4c49aSJed Brown Input Parameter: 3200de4c49aSJed Brown . ts - timestepping context 3210de4c49aSJed Brown 3220de4c49aSJed Brown Output Parameter: 3230de4c49aSJed Brown . theta - stage abscissa 3240de4c49aSJed Brown 3250de4c49aSJed Brown Note: 3260de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 3270de4c49aSJed Brown 3280de4c49aSJed Brown Level: Advanced 3290de4c49aSJed Brown 3300de4c49aSJed Brown .seealso: TSThetaSetTheta() 3310de4c49aSJed Brown @*/ 3327087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 3330de4c49aSJed Brown { 3344ac538c5SBarry Smith PetscErrorCode ierr; 3350de4c49aSJed Brown 3360de4c49aSJed Brown PetscFunctionBegin; 337afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3380de4c49aSJed Brown PetscValidPointer(theta,2); 3394ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 3400de4c49aSJed Brown PetscFunctionReturn(0); 3410de4c49aSJed Brown } 3420de4c49aSJed Brown 3430de4c49aSJed Brown #undef __FUNCT__ 3440de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 3450de4c49aSJed Brown /*@ 3460de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 3470de4c49aSJed Brown 3480de4c49aSJed Brown Not Collective 3490de4c49aSJed Brown 3500de4c49aSJed Brown Input Parameter: 3510de4c49aSJed Brown + ts - timestepping context 3520de4c49aSJed Brown - theta - stage abscissa 3530de4c49aSJed Brown 3540de4c49aSJed Brown Options Database: 3550de4c49aSJed Brown . -ts_theta_theta <theta> 3560de4c49aSJed Brown 3570de4c49aSJed Brown Level: Intermediate 3580de4c49aSJed Brown 3590de4c49aSJed Brown .seealso: TSThetaGetTheta() 3600de4c49aSJed Brown @*/ 3617087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 3620de4c49aSJed Brown { 3634ac538c5SBarry Smith PetscErrorCode ierr; 3640de4c49aSJed Brown 3650de4c49aSJed Brown PetscFunctionBegin; 366afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3674ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 3680de4c49aSJed Brown PetscFunctionReturn(0); 3690de4c49aSJed Brown } 370f33bbcb6SJed Brown 371eb284becSJed Brown #undef __FUNCT__ 37226f2ff8fSLisandro Dalcin #define __FUNCT__ "TSThetaGetEndpoint" 37326f2ff8fSLisandro Dalcin /*@ 37426f2ff8fSLisandro Dalcin TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 37526f2ff8fSLisandro Dalcin 37626f2ff8fSLisandro Dalcin Not Collective 37726f2ff8fSLisandro Dalcin 37826f2ff8fSLisandro Dalcin Input Parameter: 37926f2ff8fSLisandro Dalcin . ts - timestepping context 38026f2ff8fSLisandro Dalcin 38126f2ff8fSLisandro Dalcin Output Parameter: 38226f2ff8fSLisandro Dalcin . endpoint - PETSC_TRUE when using the endpoint variant 38326f2ff8fSLisandro Dalcin 38426f2ff8fSLisandro Dalcin Level: Advanced 38526f2ff8fSLisandro Dalcin 38626f2ff8fSLisandro Dalcin .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 38726f2ff8fSLisandro Dalcin @*/ 38826f2ff8fSLisandro Dalcin PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 38926f2ff8fSLisandro Dalcin { 39026f2ff8fSLisandro Dalcin PetscErrorCode ierr; 39126f2ff8fSLisandro Dalcin 39226f2ff8fSLisandro Dalcin PetscFunctionBegin; 39326f2ff8fSLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 39426f2ff8fSLisandro Dalcin PetscValidPointer(endpoint,2); 39526f2ff8fSLisandro Dalcin ierr = PetscTryMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 39626f2ff8fSLisandro Dalcin PetscFunctionReturn(0); 39726f2ff8fSLisandro Dalcin } 39826f2ff8fSLisandro Dalcin 39926f2ff8fSLisandro Dalcin #undef __FUNCT__ 400eb284becSJed Brown #define __FUNCT__ "TSThetaSetEndpoint" 401eb284becSJed Brown /*@ 402eb284becSJed Brown TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 403eb284becSJed Brown 404eb284becSJed Brown Not Collective 405eb284becSJed Brown 406eb284becSJed Brown Input Parameter: 407eb284becSJed Brown + ts - timestepping context 408eb284becSJed Brown - flg - PETSC_TRUE to use the endpoint variant 409eb284becSJed Brown 410eb284becSJed Brown Options Database: 411eb284becSJed Brown . -ts_theta_endpoint <flg> 412eb284becSJed Brown 413eb284becSJed Brown Level: Intermediate 414eb284becSJed Brown 415eb284becSJed Brown .seealso: TSTHETA, TSCN 416eb284becSJed Brown @*/ 417eb284becSJed Brown PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 418eb284becSJed Brown { 419eb284becSJed Brown PetscErrorCode ierr; 420eb284becSJed Brown 421eb284becSJed Brown PetscFunctionBegin; 422eb284becSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 423eb284becSJed Brown ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 424eb284becSJed Brown PetscFunctionReturn(0); 425eb284becSJed Brown } 426eb284becSJed Brown 427f33bbcb6SJed Brown /* 428f33bbcb6SJed Brown * TSBEULER and TSCN are straightforward specializations of TSTHETA. 429f33bbcb6SJed Brown * The creation functions for these specializations are below. 430f33bbcb6SJed Brown */ 431f33bbcb6SJed Brown 432f33bbcb6SJed Brown #undef __FUNCT__ 433f33bbcb6SJed Brown #define __FUNCT__ "TSView_BEuler" 434f33bbcb6SJed Brown static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 435f33bbcb6SJed Brown { 436d52bd9f3SBarry Smith PetscErrorCode ierr; 437d52bd9f3SBarry Smith 438f33bbcb6SJed Brown PetscFunctionBegin; 439d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 440f33bbcb6SJed Brown PetscFunctionReturn(0); 441f33bbcb6SJed Brown } 442f33bbcb6SJed Brown 443f33bbcb6SJed Brown /*MC 444f33bbcb6SJed Brown TSBEULER - ODE solver using the implicit backward Euler method 445f33bbcb6SJed Brown 446f33bbcb6SJed Brown Level: beginner 447f33bbcb6SJed Brown 448f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 449f33bbcb6SJed Brown 450f33bbcb6SJed Brown M*/ 451f33bbcb6SJed Brown EXTERN_C_BEGIN 452f33bbcb6SJed Brown #undef __FUNCT__ 453f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_BEuler" 454f33bbcb6SJed Brown PetscErrorCode TSCreate_BEuler(TS ts) 455f33bbcb6SJed Brown { 456f33bbcb6SJed Brown PetscErrorCode ierr; 457f33bbcb6SJed Brown 458f33bbcb6SJed Brown PetscFunctionBegin; 459f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 460f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 461f33bbcb6SJed Brown ts->ops->view = TSView_BEuler; 462f33bbcb6SJed Brown PetscFunctionReturn(0); 463f33bbcb6SJed Brown } 464f33bbcb6SJed Brown EXTERN_C_END 465f33bbcb6SJed Brown 466f33bbcb6SJed Brown #undef __FUNCT__ 467f33bbcb6SJed Brown #define __FUNCT__ "TSView_CN" 468f33bbcb6SJed Brown static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 469f33bbcb6SJed Brown { 470d52bd9f3SBarry Smith PetscErrorCode ierr; 471d52bd9f3SBarry Smith 472f33bbcb6SJed Brown PetscFunctionBegin; 473d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 474f33bbcb6SJed Brown PetscFunctionReturn(0); 475f33bbcb6SJed Brown } 476f33bbcb6SJed Brown 477f33bbcb6SJed Brown /*MC 478f33bbcb6SJed Brown TSCN - ODE solver using the implicit Crank-Nicolson method. 479f33bbcb6SJed Brown 480f33bbcb6SJed Brown Level: beginner 481f33bbcb6SJed Brown 482f33bbcb6SJed Brown Notes: 4837cf5af47SJed Brown TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 4847cf5af47SJed Brown 4857cf5af47SJed Brown $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 486f33bbcb6SJed Brown 487f33bbcb6SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 488f33bbcb6SJed Brown 489f33bbcb6SJed Brown M*/ 490f33bbcb6SJed Brown EXTERN_C_BEGIN 491f33bbcb6SJed Brown #undef __FUNCT__ 492f33bbcb6SJed Brown #define __FUNCT__ "TSCreate_CN" 493f33bbcb6SJed Brown PetscErrorCode TSCreate_CN(TS ts) 494f33bbcb6SJed Brown { 495f33bbcb6SJed Brown PetscErrorCode ierr; 496f33bbcb6SJed Brown 497f33bbcb6SJed Brown PetscFunctionBegin; 498f33bbcb6SJed Brown ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 499f33bbcb6SJed Brown ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 500eb284becSJed Brown ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 501f33bbcb6SJed Brown ts->ops->view = TSView_CN; 502f33bbcb6SJed Brown PetscFunctionReturn(0); 503f33bbcb6SJed Brown } 504f33bbcb6SJed Brown EXTERN_C_END 505