1316643e7SJed Brown /* 2316643e7SJed Brown Code for timestepping with implicit Theta method 3316643e7SJed Brown 4316643e7SJed Brown Notes: 5316643e7SJed Brown This method can be applied to DAE. 6316643e7SJed Brown 7316643e7SJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 8316643e7SJed Brown 9316643e7SJed Brown Theta | Theta 10316643e7SJed Brown ------------- 11316643e7SJed Brown | 1 12316643e7SJed Brown 13316643e7SJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 14316643e7SJed Brown 15316643e7SJed Brown X_i = x + h sum_j a_ij X'_j 16316643e7SJed Brown 17316643e7SJed Brown is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i) 18316643e7SJed Brown */ 19c6db04a5SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 20316643e7SJed Brown 21316643e7SJed Brown typedef struct { 22316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 23ace3abfcSBarry Smith PetscBool extrapolate; 24316643e7SJed Brown PetscReal Theta; 25316643e7SJed Brown PetscReal shift; 26316643e7SJed Brown PetscReal stage_time; 27316643e7SJed Brown } TS_Theta; 28316643e7SJed Brown 29316643e7SJed Brown #undef __FUNCT__ 30316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 31316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime) 32316643e7SJed Brown { 33316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 34186e87acSLisandro Dalcin PetscInt i,its,lits; 352b5a38e1SLisandro Dalcin PetscErrorCode ierr; 36316643e7SJed Brown 37316643e7SJed Brown PetscFunctionBegin; 38316643e7SJed Brown *steps = -ts->steps; 392b5a38e1SLisandro Dalcin *ptime = ts->ptime; 402b5a38e1SLisandro Dalcin 412b5a38e1SLisandro Dalcin ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 42316643e7SJed Brown 43186e87acSLisandro Dalcin for (i=0; i<ts->max_steps; i++) { 44316643e7SJed Brown if (ts->ptime + ts->time_step > ts->max_time) break; 453f2090d5SJed Brown ierr = TSPreStep(ts);CHKERRQ(ierr); 462b5a38e1SLisandro Dalcin 47316643e7SJed Brown th->stage_time = ts->ptime + th->Theta*ts->time_step; 48316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 49316643e7SJed Brown 50ace68cafSJed Brown if (th->extrapolate) { 512b5a38e1SLisandro Dalcin ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 52ace68cafSJed Brown } else { 532b5a38e1SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 54ace68cafSJed Brown } 55316643e7SJed Brown ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr); 56316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 57316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 58316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 592b5a38e1SLisandro Dalcin 602b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 612b5a38e1SLisandro Dalcin ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 622b5a38e1SLisandro Dalcin ts->ptime += ts->time_step; 63316643e7SJed Brown ts->steps++; 642b5a38e1SLisandro Dalcin 653f2090d5SJed Brown ierr = TSPostStep(ts);CHKERRQ(ierr); 662b5a38e1SLisandro Dalcin ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 67316643e7SJed Brown } 68316643e7SJed Brown 69316643e7SJed Brown *steps += ts->steps; 70316643e7SJed Brown *ptime = ts->ptime; 71316643e7SJed Brown PetscFunctionReturn(0); 72316643e7SJed Brown } 73316643e7SJed Brown 74316643e7SJed Brown /*------------------------------------------------------------*/ 75316643e7SJed Brown #undef __FUNCT__ 76277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Theta" 77277b19d0SLisandro Dalcin static PetscErrorCode TSReset_Theta(TS ts) 78316643e7SJed Brown { 79316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 80316643e7SJed Brown PetscErrorCode ierr; 81316643e7SJed Brown 82316643e7SJed Brown PetscFunctionBegin; 83*6bf464f9SBarry Smith ierr = VecDestroy(&th->X);CHKERRQ(ierr); 84*6bf464f9SBarry Smith ierr = VecDestroy(&th->Xdot);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); 97316643e7SJed Brown PetscFunctionReturn(0); 98316643e7SJed Brown } 99316643e7SJed Brown 100316643e7SJed Brown /* 101316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 1022b5a38e1SLisandro Dalcin G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 103316643e7SJed Brown */ 104316643e7SJed Brown #undef __FUNCT__ 1050f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormFunction_Theta" 1060f5c6efeSJed Brown static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 107316643e7SJed Brown { 108316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 109316643e7SJed Brown PetscErrorCode ierr; 110316643e7SJed Brown 111316643e7SJed Brown PetscFunctionBegin; 1122b5a38e1SLisandro Dalcin ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 113316643e7SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr); 114316643e7SJed Brown PetscFunctionReturn(0); 115316643e7SJed Brown } 116316643e7SJed Brown 117316643e7SJed Brown #undef __FUNCT__ 1180f5c6efeSJed Brown #define __FUNCT__ "SNESTSFormJacobian_Theta" 1190f5c6efeSJed Brown static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 120316643e7SJed Brown { 121316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 122316643e7SJed Brown PetscErrorCode ierr; 123316643e7SJed Brown 124316643e7SJed Brown PetscFunctionBegin; 1250f5c6efeSJed Brown /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 126316643e7SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr); 127316643e7SJed Brown PetscFunctionReturn(0); 128316643e7SJed Brown } 129316643e7SJed Brown 130316643e7SJed Brown 131316643e7SJed Brown #undef __FUNCT__ 132316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 133316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 134316643e7SJed Brown { 135316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 136316643e7SJed Brown PetscErrorCode ierr; 137a6e0575dSJed Brown Vec res; 138316643e7SJed Brown 139316643e7SJed Brown PetscFunctionBegin; 14017186662SBarry Smith if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 141316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 142316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 143a6e0575dSJed Brown ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr); 144a6e0575dSJed Brown ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 145*6bf464f9SBarry Smith ierr = VecDestroy(&res);CHKERRQ(ierr); 1467c0b301bSJed Brown /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 1477c0b301bSJed Brown replace A and we don't want to mess with that. With -snes_mf, A and B will be replaced as well as the function and 1487c0b301bSJed Brown context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 1497c0b301bSJed Brown the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 1507c0b301bSJed Brown snes->jacP should be the TS. */ 1517c0b301bSJed Brown { 1527c0b301bSJed Brown Mat A,B; 1537c0b301bSJed Brown PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 1547c0b301bSJed Brown void *ctx; 1557c0b301bSJed Brown ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 1560f5c6efeSJed Brown ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 1577c0b301bSJed Brown } 158316643e7SJed Brown PetscFunctionReturn(0); 159316643e7SJed Brown } 160316643e7SJed Brown /*------------------------------------------------------------*/ 161316643e7SJed Brown 162316643e7SJed Brown #undef __FUNCT__ 163316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 164316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 165316643e7SJed Brown { 166316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 167316643e7SJed Brown PetscErrorCode ierr; 168316643e7SJed Brown 169316643e7SJed Brown PetscFunctionBegin; 170d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 171316643e7SJed Brown { 172316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 173acfcf0e5SJed Brown ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 174316643e7SJed Brown } 175316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 176316643e7SJed Brown PetscFunctionReturn(0); 177316643e7SJed Brown } 178316643e7SJed Brown 179316643e7SJed Brown #undef __FUNCT__ 180316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 181316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 182316643e7SJed Brown { 183316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 184ace3abfcSBarry Smith PetscBool iascii; 185316643e7SJed Brown PetscErrorCode ierr; 186316643e7SJed Brown 187316643e7SJed Brown PetscFunctionBegin; 1882692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 189316643e7SJed Brown if (iascii) { 190316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 191ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 192316643e7SJed Brown } else { 193e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name); 194316643e7SJed Brown } 195316643e7SJed Brown PetscFunctionReturn(0); 196316643e7SJed Brown } 197316643e7SJed Brown 1980de4c49aSJed Brown EXTERN_C_BEGIN 1990de4c49aSJed Brown #undef __FUNCT__ 2000de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta_Theta" 2017087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 2020de4c49aSJed Brown { 2030de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2040de4c49aSJed Brown 2050de4c49aSJed Brown PetscFunctionBegin; 2060de4c49aSJed Brown *theta = th->Theta; 2070de4c49aSJed Brown PetscFunctionReturn(0); 2080de4c49aSJed Brown } 2090de4c49aSJed Brown 2100de4c49aSJed Brown #undef __FUNCT__ 2110de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta_Theta" 2127087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 2130de4c49aSJed Brown { 2140de4c49aSJed Brown TS_Theta *th = (TS_Theta*)ts->data; 2150de4c49aSJed Brown 2160de4c49aSJed Brown PetscFunctionBegin; 217e7be1afaSJed Brown if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 2180de4c49aSJed Brown th->Theta = theta; 2190de4c49aSJed Brown PetscFunctionReturn(0); 2200de4c49aSJed Brown } 2210de4c49aSJed Brown EXTERN_C_END 2220de4c49aSJed Brown 223316643e7SJed Brown /* ------------------------------------------------------------ */ 224316643e7SJed Brown /*MC 22596f5712cSJed Brown TSTHETA - DAE solver using the implicit Theta method 226316643e7SJed Brown 227316643e7SJed Brown Level: beginner 228316643e7SJed Brown 229316643e7SJed Brown .seealso: TSCreate(), TS, TSSetType() 230316643e7SJed Brown 231316643e7SJed Brown M*/ 232316643e7SJed Brown EXTERN_C_BEGIN 233316643e7SJed Brown #undef __FUNCT__ 234316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 2357087cfbeSBarry Smith PetscErrorCode TSCreate_Theta(TS ts) 236316643e7SJed Brown { 237316643e7SJed Brown TS_Theta *th; 238316643e7SJed Brown PetscErrorCode ierr; 239316643e7SJed Brown 240316643e7SJed Brown PetscFunctionBegin; 241277b19d0SLisandro Dalcin ts->ops->reset = TSReset_Theta; 242316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 243316643e7SJed Brown ts->ops->view = TSView_Theta; 244316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 245316643e7SJed Brown ts->ops->step = TSStep_Theta; 246316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 2470f5c6efeSJed Brown ts->ops->snesfunction = SNESTSFormFunction_Theta; 2480f5c6efeSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 249316643e7SJed Brown 2502b5a38e1SLisandro Dalcin ts->problem_type = TS_NONLINEAR; 251316643e7SJed Brown ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 252316643e7SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 253316643e7SJed Brown 254316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 255316643e7SJed Brown ts->data = (void*)th; 256316643e7SJed Brown 2576f700aefSJed Brown th->extrapolate = PETSC_FALSE; 258316643e7SJed Brown th->Theta = 0.5; 259316643e7SJed Brown 2600de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 2610de4c49aSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 262316643e7SJed Brown PetscFunctionReturn(0); 263316643e7SJed Brown } 264316643e7SJed Brown EXTERN_C_END 2650de4c49aSJed Brown 2660de4c49aSJed Brown #undef __FUNCT__ 2670de4c49aSJed Brown #define __FUNCT__ "TSThetaGetTheta" 2680de4c49aSJed Brown /*@ 2690de4c49aSJed Brown TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 2700de4c49aSJed Brown 2710de4c49aSJed Brown Not Collective 2720de4c49aSJed Brown 2730de4c49aSJed Brown Input Parameter: 2740de4c49aSJed Brown . ts - timestepping context 2750de4c49aSJed Brown 2760de4c49aSJed Brown Output Parameter: 2770de4c49aSJed Brown . theta - stage abscissa 2780de4c49aSJed Brown 2790de4c49aSJed Brown Note: 2800de4c49aSJed Brown Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 2810de4c49aSJed Brown 2820de4c49aSJed Brown Level: Advanced 2830de4c49aSJed Brown 2840de4c49aSJed Brown .seealso: TSThetaSetTheta() 2850de4c49aSJed Brown @*/ 2867087cfbeSBarry Smith PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 2870de4c49aSJed Brown { 2884ac538c5SBarry Smith PetscErrorCode ierr; 2890de4c49aSJed Brown 2900de4c49aSJed Brown PetscFunctionBegin; 291afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2920de4c49aSJed Brown PetscValidPointer(theta,2); 2934ac538c5SBarry Smith ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 2940de4c49aSJed Brown PetscFunctionReturn(0); 2950de4c49aSJed Brown } 2960de4c49aSJed Brown 2970de4c49aSJed Brown #undef __FUNCT__ 2980de4c49aSJed Brown #define __FUNCT__ "TSThetaSetTheta" 2990de4c49aSJed Brown /*@ 3000de4c49aSJed Brown TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 3010de4c49aSJed Brown 3020de4c49aSJed Brown Not Collective 3030de4c49aSJed Brown 3040de4c49aSJed Brown Input Parameter: 3050de4c49aSJed Brown + ts - timestepping context 3060de4c49aSJed Brown - theta - stage abscissa 3070de4c49aSJed Brown 3080de4c49aSJed Brown Options Database: 3090de4c49aSJed Brown . -ts_theta_theta <theta> 3100de4c49aSJed Brown 3110de4c49aSJed Brown Level: Intermediate 3120de4c49aSJed Brown 3130de4c49aSJed Brown .seealso: TSThetaGetTheta() 3140de4c49aSJed Brown @*/ 3157087cfbeSBarry Smith PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 3160de4c49aSJed Brown { 3174ac538c5SBarry Smith PetscErrorCode ierr; 3180de4c49aSJed Brown 3190de4c49aSJed Brown PetscFunctionBegin; 320afb20b64SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3214ac538c5SBarry Smith ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 3220de4c49aSJed Brown PetscFunctionReturn(0); 3230de4c49aSJed Brown } 324