1316643e7SJed Brown #define PETSCTS_DLL 2316643e7SJed Brown 3316643e7SJed Brown /* 4316643e7SJed Brown Code for timestepping with implicit Theta method 5316643e7SJed Brown 6316643e7SJed Brown Notes: 7316643e7SJed Brown This method can be applied to DAE. 8316643e7SJed Brown 9316643e7SJed Brown This method is cast as a 1-stage implicit Runge-Kutta method. 10316643e7SJed Brown 11316643e7SJed Brown Theta | Theta 12316643e7SJed Brown ------------- 13316643e7SJed Brown | 1 14316643e7SJed Brown 15316643e7SJed Brown To apply a diagonally implicit RK method to DAE, the stage formula 16316643e7SJed Brown 17316643e7SJed Brown X_i = x + h sum_j a_ij X'_j 18316643e7SJed Brown 19316643e7SJed Brown is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i) 20316643e7SJed Brown */ 21316643e7SJed Brown #include "private/tsimpl.h" /*I "petscts.h" I*/ 22316643e7SJed Brown 23316643e7SJed Brown typedef struct { 24316643e7SJed Brown Vec Xold; 25316643e7SJed Brown Vec X,Xdot; /* Storage for one stage */ 26316643e7SJed Brown Vec res; /* DAE residuals */ 27ace68cafSJed Brown PetscTruth extrapolate; 28316643e7SJed Brown PetscReal Theta; 29316643e7SJed Brown PetscReal shift; 30316643e7SJed Brown PetscReal stage_time; 31316643e7SJed Brown } TS_Theta; 32316643e7SJed Brown 33316643e7SJed Brown #undef __FUNCT__ 34316643e7SJed Brown #define __FUNCT__ "TSStep_Theta" 35316643e7SJed Brown static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime) 36316643e7SJed Brown { 37316643e7SJed Brown Vec sol = ts->vec_sol; 38316643e7SJed Brown PetscErrorCode ierr; 39316643e7SJed Brown PetscInt i,max_steps = ts->max_steps,its,lits; 40316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 41316643e7SJed Brown 42316643e7SJed Brown PetscFunctionBegin; 43316643e7SJed Brown *steps = -ts->steps; 44316643e7SJed Brown ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 45316643e7SJed Brown 46316643e7SJed Brown for (i=0; i<max_steps; i++) { 47316643e7SJed Brown if (ts->ptime + ts->time_step > ts->max_time) break; 48316643e7SJed Brown th->stage_time = ts->ptime + th->Theta*ts->time_step; 49316643e7SJed Brown th->shift = 1./(th->Theta*ts->time_step); 50316643e7SJed Brown ts->ptime += ts->time_step; 51316643e7SJed Brown 52316643e7SJed Brown ierr = VecCopy(sol,th->Xold);CHKERRQ(ierr); /* Used within function evalutaion */ 53ace68cafSJed Brown if (th->extrapolate) { 54316643e7SJed Brown ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,sol);CHKERRQ(ierr); 55ace68cafSJed Brown } else { 56ace68cafSJed Brown ierr = VecCopy(sol,th->X);CHKERRQ(ierr); 57ace68cafSJed Brown } 58316643e7SJed Brown ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr); 59316643e7SJed Brown ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 60316643e7SJed Brown ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 61316643e7SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 62316643e7SJed Brown ierr = VecAXPY(sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 63316643e7SJed Brown ts->steps++; 64316643e7SJed Brown ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 65316643e7SJed Brown } 66316643e7SJed Brown 67316643e7SJed Brown *steps += ts->steps; 68316643e7SJed Brown *ptime = ts->ptime; 69316643e7SJed Brown PetscFunctionReturn(0); 70316643e7SJed Brown } 71316643e7SJed Brown 72316643e7SJed Brown /*------------------------------------------------------------*/ 73316643e7SJed Brown #undef __FUNCT__ 74316643e7SJed Brown #define __FUNCT__ "TSDestroy_Theta" 75316643e7SJed Brown static PetscErrorCode TSDestroy_Theta(TS ts) 76316643e7SJed Brown { 77316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 78316643e7SJed Brown PetscErrorCode ierr; 79316643e7SJed Brown 80316643e7SJed Brown PetscFunctionBegin; 81316643e7SJed Brown ierr = VecDestroy(th->Xold);CHKERRQ(ierr); 82316643e7SJed Brown ierr = VecDestroy(th->X);CHKERRQ(ierr); 83316643e7SJed Brown ierr = VecDestroy(th->Xdot);CHKERRQ(ierr); 84316643e7SJed Brown ierr = VecDestroy(th->res);CHKERRQ(ierr); 85316643e7SJed Brown ierr = PetscFree(th);CHKERRQ(ierr); 86316643e7SJed Brown PetscFunctionReturn(0); 87316643e7SJed Brown } 88316643e7SJed Brown 89316643e7SJed Brown /* 90316643e7SJed Brown This defines the nonlinear equation that is to be solved with SNES 91316643e7SJed Brown G(U) = F[t0+T*dt, U, (U-U0)*shift] = 0 92316643e7SJed Brown */ 93316643e7SJed Brown #undef __FUNCT__ 94316643e7SJed Brown #define __FUNCT__ "TSThetaFunction" 95316643e7SJed Brown static PetscErrorCode TSThetaFunction(SNES snes,Vec x,Vec y,void *ctx) 96316643e7SJed Brown { 97316643e7SJed Brown TS ts = (TS)ctx; 98316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 99316643e7SJed Brown PetscErrorCode ierr; 100316643e7SJed Brown 101316643e7SJed Brown PetscFunctionBegin; 102316643e7SJed Brown ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,th->Xold,x);CHKERRQ(ierr); 103316643e7SJed Brown ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr); 104316643e7SJed Brown PetscFunctionReturn(0); 105316643e7SJed Brown } 106316643e7SJed Brown 107316643e7SJed Brown #undef __FUNCT__ 108316643e7SJed Brown #define __FUNCT__ "TSThetaJacobian" 109316643e7SJed Brown static PetscErrorCode TSThetaJacobian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,void *ctx) 110316643e7SJed Brown { 111316643e7SJed Brown TS ts = (TS)ctx; 112316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 113316643e7SJed Brown PetscErrorCode ierr; 114316643e7SJed Brown 115316643e7SJed Brown PetscFunctionBegin; 116316643e7SJed Brown /* th->Xdot will have already been computed in TSThetaFunction */ 117316643e7SJed Brown ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr); 118316643e7SJed Brown PetscFunctionReturn(0); 119316643e7SJed Brown } 120316643e7SJed Brown 121316643e7SJed Brown 122316643e7SJed Brown #undef __FUNCT__ 123316643e7SJed Brown #define __FUNCT__ "TSSetUp_Theta" 124316643e7SJed Brown static PetscErrorCode TSSetUp_Theta(TS ts) 125316643e7SJed Brown { 126316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 127316643e7SJed Brown PetscErrorCode ierr; 128316643e7SJed Brown 129316643e7SJed Brown PetscFunctionBegin; 130316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xold);CHKERRQ(ierr); 131316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 132316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 133316643e7SJed Brown ierr = VecDuplicate(ts->vec_sol,&th->res);CHKERRQ(ierr); 134316643e7SJed Brown ierr = SNESSetFunction(ts->snes,th->res,TSThetaFunction,ts);CHKERRQ(ierr); 1357c0b301bSJed Brown /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 1367c0b301bSJed 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 1377c0b301bSJed Brown context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 1387c0b301bSJed Brown the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 1397c0b301bSJed Brown snes->jacP should be the TS. */ 1407c0b301bSJed Brown { 1417c0b301bSJed Brown Mat A,B; 1427c0b301bSJed Brown PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 1437c0b301bSJed Brown void *ctx; 1447c0b301bSJed Brown ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 1457c0b301bSJed Brown ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&TSThetaJacobian,ctx?ctx:ts);CHKERRQ(ierr); 1467c0b301bSJed Brown } 147316643e7SJed Brown PetscFunctionReturn(0); 148316643e7SJed Brown } 149316643e7SJed Brown /*------------------------------------------------------------*/ 150316643e7SJed Brown 151316643e7SJed Brown #undef __FUNCT__ 152316643e7SJed Brown #define __FUNCT__ "TSSetFromOptions_Theta" 153316643e7SJed Brown static PetscErrorCode TSSetFromOptions_Theta(TS ts) 154316643e7SJed Brown { 155316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 156316643e7SJed Brown PetscErrorCode ierr; 157316643e7SJed Brown 158316643e7SJed Brown PetscFunctionBegin; 159*d73342a9SJed Brown ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 160316643e7SJed Brown { 161316643e7SJed Brown ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 162ace68cafSJed Brown ierr = PetscOptionsTruth("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 163316643e7SJed Brown } 164316643e7SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 165316643e7SJed Brown PetscFunctionReturn(0); 166316643e7SJed Brown } 167316643e7SJed Brown 168316643e7SJed Brown #undef __FUNCT__ 169316643e7SJed Brown #define __FUNCT__ "TSView_Theta" 170316643e7SJed Brown static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 171316643e7SJed Brown { 172316643e7SJed Brown TS_Theta *th = (TS_Theta*)ts->data; 173316643e7SJed Brown PetscTruth iascii; 174316643e7SJed Brown PetscErrorCode ierr; 175316643e7SJed Brown 176316643e7SJed Brown PetscFunctionBegin; 177316643e7SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 178316643e7SJed Brown if (iascii) { 179316643e7SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 180ace68cafSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 181316643e7SJed Brown } else { 182316643e7SJed Brown SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name); 183316643e7SJed Brown } 184316643e7SJed Brown PetscFunctionReturn(0); 185316643e7SJed Brown } 186316643e7SJed Brown 187316643e7SJed Brown /* ------------------------------------------------------------ */ 188316643e7SJed Brown /*MC 1897c0b301bSJed Brown TS_THETA - DAE solver using the implicit Theta method 190316643e7SJed Brown 191316643e7SJed Brown Level: beginner 192316643e7SJed Brown 193316643e7SJed Brown .seealso: TSCreate(), TS, TSSetType() 194316643e7SJed Brown 195316643e7SJed Brown M*/ 196316643e7SJed Brown EXTERN_C_BEGIN 197316643e7SJed Brown #undef __FUNCT__ 198316643e7SJed Brown #define __FUNCT__ "TSCreate_Theta" 199316643e7SJed Brown PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Theta(TS ts) 200316643e7SJed Brown { 201316643e7SJed Brown TS_Theta *th; 202316643e7SJed Brown PetscErrorCode ierr; 203316643e7SJed Brown 204316643e7SJed Brown PetscFunctionBegin; 205316643e7SJed Brown ts->ops->destroy = TSDestroy_Theta; 206316643e7SJed Brown ts->ops->view = TSView_Theta; 207316643e7SJed Brown ts->ops->setup = TSSetUp_Theta; 208316643e7SJed Brown ts->ops->step = TSStep_Theta; 209316643e7SJed Brown ts->ops->setfromoptions = TSSetFromOptions_Theta; 210316643e7SJed Brown 211316643e7SJed Brown ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 212316643e7SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 213316643e7SJed Brown 214316643e7SJed Brown ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 215316643e7SJed Brown ts->data = (void*)th; 216316643e7SJed Brown 217ace68cafSJed Brown th->extrapolate = PETSC_TRUE; 218316643e7SJed Brown th->Theta = 0.5; 219316643e7SJed Brown 220316643e7SJed Brown PetscFunctionReturn(0); 221316643e7SJed Brown } 222316643e7SJed Brown EXTERN_C_END 223