1 2 /* 3 Code for timestepping with implicit Theta method 4 5 Notes: 6 This method can be applied to DAE. 7 8 This method is cast as a 1-stage implicit Runge-Kutta method. 9 10 Theta | Theta 11 ------------- 12 | 1 13 14 To apply a diagonally implicit RK method to DAE, the stage formula 15 16 X_i = x + h sum_j a_ij X'_j 17 18 is interpreted as a formula for X'_i in terms of X_i and known stuff (X'_j, j<i) 19 */ 20 #include <private/tsimpl.h> /*I "petscts.h" I*/ 21 22 typedef struct { 23 Vec X,Xdot; /* Storage for one stage */ 24 PetscBool extrapolate; 25 PetscReal Theta; 26 PetscReal shift; 27 PetscReal stage_time; 28 } TS_Theta; 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "TSStep_Theta" 32 static PetscErrorCode TSStep_Theta(TS ts,PetscInt *steps,PetscReal *ptime) 33 { 34 TS_Theta *th = (TS_Theta*)ts->data; 35 PetscInt i,max_steps = ts->max_steps,its,lits; 36 PetscErrorCode ierr; 37 38 PetscFunctionBegin; 39 *steps = -ts->steps; 40 *ptime = ts->ptime; 41 42 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 43 44 for (i=0; i<max_steps; i++) { 45 if (ts->ptime + ts->time_step > ts->max_time) break; 46 ierr = TSPreStep(ts);CHKERRQ(ierr); 47 48 th->stage_time = ts->ptime + th->Theta*ts->time_step; 49 th->shift = 1./(th->Theta*ts->time_step); 50 51 if (th->extrapolate) { 52 ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 53 } else { 54 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 55 } 56 ierr = SNESSolve(ts->snes,PETSC_NULL,th->X);CHKERRQ(ierr); 57 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 58 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 59 ts->nonlinear_its += its; ts->linear_its += lits; 60 61 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 62 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 63 ts->ptime += ts->time_step; 64 ts->steps++; 65 66 ierr = TSPostStep(ts);CHKERRQ(ierr); 67 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 68 } 69 70 *steps += ts->steps; 71 *ptime = ts->ptime; 72 PetscFunctionReturn(0); 73 } 74 75 /*------------------------------------------------------------*/ 76 #undef __FUNCT__ 77 #define __FUNCT__ "TSReset_Theta" 78 static PetscErrorCode TSReset_Theta(TS ts) 79 { 80 TS_Theta *th = (TS_Theta*)ts->data; 81 PetscErrorCode ierr; 82 83 PetscFunctionBegin; 84 if (th->X) {ierr = VecDestroy(th->X);CHKERRQ(ierr);} 85 if (th->Xdot) {ierr = VecDestroy(th->Xdot);CHKERRQ(ierr);} 86 PetscFunctionReturn(0); 87 } 88 89 #undef __FUNCT__ 90 #define __FUNCT__ "TSDestroy_Theta" 91 static PetscErrorCode TSDestroy_Theta(TS ts) 92 { 93 PetscErrorCode ierr; 94 95 PetscFunctionBegin; 96 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 97 ierr = PetscFree(ts->data);CHKERRQ(ierr); 98 PetscFunctionReturn(0); 99 } 100 101 /* 102 This defines the nonlinear equation that is to be solved with SNES 103 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 104 */ 105 #undef __FUNCT__ 106 #define __FUNCT__ "SNESTSFormFunction_Theta" 107 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 108 { 109 TS_Theta *th = (TS_Theta*)ts->data; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 114 ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "SNESTSFormJacobian_Theta" 120 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 121 { 122 TS_Theta *th = (TS_Theta*)ts->data; 123 PetscErrorCode ierr; 124 125 PetscFunctionBegin; 126 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 127 ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 132 #undef __FUNCT__ 133 #define __FUNCT__ "TSSetUp_Theta" 134 static PetscErrorCode TSSetUp_Theta(TS ts) 135 { 136 TS_Theta *th = (TS_Theta*)ts->data; 137 PetscErrorCode ierr; 138 Vec res; 139 140 PetscFunctionBegin; 141 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 142 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 143 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 144 ierr = VecDuplicate(ts->vec_sol,&res);CHKERRQ(ierr); 145 ierr = SNESSetFunction(ts->snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 146 ierr = VecDestroy(res);CHKERRQ(ierr); 147 /* This is nasty. SNESSetFromOptions() is usually called in TSSetFromOptions(). With -snes_mf_operator, it will 148 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 149 context. Note that SNESSetFunction() normally has not been called before SNESSetFromOptions(), so when -snes_mf sets 150 the Jacobian user context to snes->funP, it will actually be NULL. This is not a problem because both snes->funP and 151 snes->jacP should be the TS. */ 152 { 153 Mat A,B; 154 PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 155 void *ctx; 156 ierr = SNESGetJacobian(ts->snes,&A,&B,&func,&ctx);CHKERRQ(ierr); 157 ierr = SNESSetJacobian(ts->snes,A?A:ts->A,B?B:ts->B,func?func:&SNESTSFormJacobian,ctx?ctx:ts);CHKERRQ(ierr); 158 } 159 PetscFunctionReturn(0); 160 } 161 /*------------------------------------------------------------*/ 162 163 #undef __FUNCT__ 164 #define __FUNCT__ "TSSetFromOptions_Theta" 165 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 166 { 167 TS_Theta *th = (TS_Theta*)ts->data; 168 PetscErrorCode ierr; 169 170 PetscFunctionBegin; 171 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 172 { 173 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 174 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 175 } 176 ierr = PetscOptionsTail();CHKERRQ(ierr); 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "TSView_Theta" 182 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 183 { 184 TS_Theta *th = (TS_Theta*)ts->data; 185 PetscBool iascii; 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 190 if (iascii) { 191 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 192 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 193 } else { 194 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_Theta",((PetscObject)viewer)->type_name); 195 } 196 PetscFunctionReturn(0); 197 } 198 199 EXTERN_C_BEGIN 200 #undef __FUNCT__ 201 #define __FUNCT__ "TSThetaGetTheta_Theta" 202 PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 203 { 204 TS_Theta *th = (TS_Theta*)ts->data; 205 206 PetscFunctionBegin; 207 *theta = th->Theta; 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNCT__ 212 #define __FUNCT__ "TSThetaSetTheta_Theta" 213 PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 214 { 215 TS_Theta *th = (TS_Theta*)ts->data; 216 217 PetscFunctionBegin; 218 if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 219 th->Theta = theta; 220 PetscFunctionReturn(0); 221 } 222 EXTERN_C_END 223 224 /* ------------------------------------------------------------ */ 225 /*MC 226 TSTHETA - DAE solver using the implicit Theta method 227 228 Level: beginner 229 230 .seealso: TSCreate(), TS, TSSetType() 231 232 M*/ 233 EXTERN_C_BEGIN 234 #undef __FUNCT__ 235 #define __FUNCT__ "TSCreate_Theta" 236 PetscErrorCode TSCreate_Theta(TS ts) 237 { 238 TS_Theta *th; 239 PetscErrorCode ierr; 240 241 PetscFunctionBegin; 242 ts->ops->reset = TSReset_Theta; 243 ts->ops->destroy = TSDestroy_Theta; 244 ts->ops->view = TSView_Theta; 245 ts->ops->setup = TSSetUp_Theta; 246 ts->ops->step = TSStep_Theta; 247 ts->ops->setfromoptions = TSSetFromOptions_Theta; 248 ts->ops->snesfunction = SNESTSFormFunction_Theta; 249 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 250 251 ts->problem_type = TS_NONLINEAR; 252 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 253 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 254 255 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 256 ts->data = (void*)th; 257 258 th->extrapolate = PETSC_FALSE; 259 th->Theta = 0.5; 260 261 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 262 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 EXTERN_C_END 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "TSThetaGetTheta" 269 /*@ 270 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 271 272 Not Collective 273 274 Input Parameter: 275 . ts - timestepping context 276 277 Output Parameter: 278 . theta - stage abscissa 279 280 Note: 281 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 282 283 Level: Advanced 284 285 .seealso: TSThetaSetTheta() 286 @*/ 287 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 288 { 289 PetscErrorCode ierr; 290 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 293 PetscValidPointer(theta,2); 294 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 #undef __FUNCT__ 299 #define __FUNCT__ "TSThetaSetTheta" 300 /*@ 301 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 302 303 Not Collective 304 305 Input Parameter: 306 + ts - timestepping context 307 - theta - stage abscissa 308 309 Options Database: 310 . -ts_theta_theta <theta> 311 312 Level: Intermediate 313 314 .seealso: TSThetaGetTheta() 315 @*/ 316 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 317 { 318 PetscErrorCode ierr; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 322 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 323 PetscFunctionReturn(0); 324 } 325