1 /* 2 Code for timestepping with implicit Theta method 3 */ 4 #include <private/tsimpl.h> /*I "petscts.h" I*/ 5 6 typedef struct { 7 Vec X,Xdot; /* Storage for one stage */ 8 Vec affine; /* Affine vector needed for residual at beginning of step */ 9 PetscBool extrapolate; 10 PetscBool endpoint; 11 PetscReal Theta; 12 PetscReal shift; 13 PetscReal stage_time; 14 } TS_Theta; 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "TSStep_Theta" 18 static PetscErrorCode TSStep_Theta(TS ts) 19 { 20 TS_Theta *th = (TS_Theta*)ts->data; 21 PetscInt its,lits; 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 ts->time_step = ts->next_time_step; 26 th->stage_time = ts->ptime + (th->endpoint ? 1. : th->Theta)*ts->time_step; 27 th->shift = 1./(th->Theta*ts->time_step); 28 29 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 30 ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 31 if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 32 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 33 ierr = VecScale(th->affine,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 34 } 35 if (th->extrapolate) { 36 ierr = VecWAXPY(th->X,1./th->shift,th->Xdot,ts->vec_sol);CHKERRQ(ierr); 37 } else { 38 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 39 } 40 ierr = SNESSolve(ts->snes,th->affine,th->X);CHKERRQ(ierr); 41 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 42 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 43 ts->nonlinear_its += its; ts->linear_its += lits; 44 45 if (th->endpoint) { 46 ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 47 } else { 48 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,th->X);CHKERRQ(ierr); 49 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 50 } 51 ts->ptime += ts->time_step; 52 ts->next_time_step = ts->time_step; 53 ts->steps++; 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "TSInterpolate_Theta" 59 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 60 { 61 TS_Theta *th = (TS_Theta*)ts->data; 62 PetscReal alpha = t - ts->ptime; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 67 if (th->endpoint) alpha *= th->Theta; 68 ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 /*------------------------------------------------------------*/ 73 #undef __FUNCT__ 74 #define __FUNCT__ "TSReset_Theta" 75 static PetscErrorCode TSReset_Theta(TS ts) 76 { 77 TS_Theta *th = (TS_Theta*)ts->data; 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 82 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 83 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "TSDestroy_Theta" 89 static PetscErrorCode TSDestroy_Theta(TS ts) 90 { 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 95 ierr = PetscFree(ts->data);CHKERRQ(ierr); 96 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 97 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 98 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","",PETSC_NULL);CHKERRQ(ierr); 99 PetscFunctionReturn(0); 100 } 101 102 /* 103 This defines the nonlinear equation that is to be solved with SNES 104 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 105 */ 106 #undef __FUNCT__ 107 #define __FUNCT__ "SNESTSFormFunction_Theta" 108 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 109 { 110 TS_Theta *th = (TS_Theta*)ts->data; 111 PetscErrorCode ierr; 112 113 PetscFunctionBegin; 114 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 115 ierr = VecAXPBYPCZ(th->Xdot,-th->shift,th->shift,0,ts->vec_sol,x);CHKERRQ(ierr); 116 ierr = TSComputeIFunction(ts,th->stage_time,x,th->Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 #undef __FUNCT__ 121 #define __FUNCT__ "SNESTSFormJacobian_Theta" 122 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *str,TS ts) 123 { 124 TS_Theta *th = (TS_Theta*)ts->data; 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 /* th->Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 129 ierr = TSComputeIJacobian(ts,th->stage_time,x,th->Xdot,th->shift,A,B,str,PETSC_FALSE);CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "TSSetUp_Theta" 136 static PetscErrorCode TSSetUp_Theta(TS ts) 137 { 138 TS_Theta *th = (TS_Theta*)ts->data; 139 PetscErrorCode ierr; 140 141 PetscFunctionBegin; 142 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 143 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 144 PetscFunctionReturn(0); 145 } 146 /*------------------------------------------------------------*/ 147 148 #undef __FUNCT__ 149 #define __FUNCT__ "TSSetFromOptions_Theta" 150 static PetscErrorCode TSSetFromOptions_Theta(TS ts) 151 { 152 TS_Theta *th = (TS_Theta*)ts->data; 153 PetscErrorCode ierr; 154 155 PetscFunctionBegin; 156 ierr = PetscOptionsHead("Theta ODE solver options");CHKERRQ(ierr); 157 { 158 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,PETSC_NULL);CHKERRQ(ierr); 159 ierr = PetscOptionsBool("-ts_theta_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,PETSC_NULL);CHKERRQ(ierr); 160 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 } 162 ierr = PetscOptionsTail();CHKERRQ(ierr); 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "TSView_Theta" 168 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 169 { 170 TS_Theta *th = (TS_Theta*)ts->data; 171 PetscBool iascii; 172 PetscErrorCode ierr; 173 174 PetscFunctionBegin; 175 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 176 if (iascii) { 177 ierr = PetscViewerASCIIPrintf(viewer," Theta=%G\n",th->Theta);CHKERRQ(ierr); 178 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate?"yes":"no");CHKERRQ(ierr); 179 } 180 PetscFunctionReturn(0); 181 } 182 183 EXTERN_C_BEGIN 184 #undef __FUNCT__ 185 #define __FUNCT__ "TSThetaGetTheta_Theta" 186 PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 187 { 188 TS_Theta *th = (TS_Theta*)ts->data; 189 190 PetscFunctionBegin; 191 *theta = th->Theta; 192 PetscFunctionReturn(0); 193 } 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "TSThetaSetTheta_Theta" 197 PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 198 { 199 TS_Theta *th = (TS_Theta*)ts->data; 200 201 PetscFunctionBegin; 202 if (theta <= 0 || 1 < theta) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Theta %G not in range (0,1]",theta); 203 th->Theta = theta; 204 PetscFunctionReturn(0); 205 } 206 207 #undef __FUNCT__ 208 #define __FUNCT__ "TSThetaSetEndpoint_Theta" 209 PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 210 { 211 TS_Theta *th = (TS_Theta*)ts->data; 212 213 PetscFunctionBegin; 214 th->endpoint = flg; 215 PetscFunctionReturn(0); 216 } 217 EXTERN_C_END 218 219 /* ------------------------------------------------------------ */ 220 /*MC 221 TSTHETA - DAE solver using the implicit Theta method 222 223 Level: beginner 224 225 Notes: 226 This method can be applied to DAE. 227 228 This method is cast as a 1-stage implicit Runge-Kutta method. 229 230 .vb 231 Theta | Theta 232 ------------- 233 | 1 234 .ve 235 236 For the default Theta=0.5, this is also known as the implicit midpoint rule. 237 238 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 239 240 .vb 241 0 | 0 0 242 1 | 1-Theta Theta 243 ------------------- 244 | 1-Theta Theta 245 .ve 246 247 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 248 249 To apply a diagonally implicit RK method to DAE, the stage formula 250 251 $ Y_i = X + h sum_j a_ij Y'_j 252 253 is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 254 255 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 256 257 M*/ 258 EXTERN_C_BEGIN 259 #undef __FUNCT__ 260 #define __FUNCT__ "TSCreate_Theta" 261 PetscErrorCode TSCreate_Theta(TS ts) 262 { 263 TS_Theta *th; 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 ts->ops->reset = TSReset_Theta; 268 ts->ops->destroy = TSDestroy_Theta; 269 ts->ops->view = TSView_Theta; 270 ts->ops->setup = TSSetUp_Theta; 271 ts->ops->step = TSStep_Theta; 272 ts->ops->interpolate = TSInterpolate_Theta; 273 ts->ops->setfromoptions = TSSetFromOptions_Theta; 274 ts->ops->snesfunction = SNESTSFormFunction_Theta; 275 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 276 277 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 278 ts->data = (void*)th; 279 280 th->extrapolate = PETSC_FALSE; 281 th->Theta = 0.5; 282 283 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 284 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 285 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 EXTERN_C_END 289 290 #undef __FUNCT__ 291 #define __FUNCT__ "TSThetaGetTheta" 292 /*@ 293 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 294 295 Not Collective 296 297 Input Parameter: 298 . ts - timestepping context 299 300 Output Parameter: 301 . theta - stage abscissa 302 303 Note: 304 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 305 306 Level: Advanced 307 308 .seealso: TSThetaSetTheta() 309 @*/ 310 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 311 { 312 PetscErrorCode ierr; 313 314 PetscFunctionBegin; 315 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 316 PetscValidPointer(theta,2); 317 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 318 PetscFunctionReturn(0); 319 } 320 321 #undef __FUNCT__ 322 #define __FUNCT__ "TSThetaSetTheta" 323 /*@ 324 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 325 326 Not Collective 327 328 Input Parameter: 329 + ts - timestepping context 330 - theta - stage abscissa 331 332 Options Database: 333 . -ts_theta_theta <theta> 334 335 Level: Intermediate 336 337 .seealso: TSThetaGetTheta() 338 @*/ 339 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 340 { 341 PetscErrorCode ierr; 342 343 PetscFunctionBegin; 344 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 345 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "TSThetaSetEndpoint" 351 /*@ 352 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 353 354 Not Collective 355 356 Input Parameter: 357 + ts - timestepping context 358 - flg - PETSC_TRUE to use the endpoint variant 359 360 Options Database: 361 . -ts_theta_endpoint <flg> 362 363 Level: Intermediate 364 365 .seealso: TSTHETA, TSCN 366 @*/ 367 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 368 { 369 PetscErrorCode ierr; 370 371 PetscFunctionBegin; 372 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 373 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 374 PetscFunctionReturn(0); 375 } 376 377 /* 378 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 379 * The creation functions for these specializations are below. 380 */ 381 382 #undef __FUNCT__ 383 #define __FUNCT__ "TSView_BEuler" 384 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 385 { 386 PetscFunctionBegin; 387 PetscFunctionReturn(0); 388 } 389 390 /*MC 391 TSBEULER - ODE solver using the implicit backward Euler method 392 393 Level: beginner 394 395 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 396 397 M*/ 398 EXTERN_C_BEGIN 399 #undef __FUNCT__ 400 #define __FUNCT__ "TSCreate_BEuler" 401 PetscErrorCode TSCreate_BEuler(TS ts) 402 { 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 407 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 408 ts->ops->view = TSView_BEuler; 409 PetscFunctionReturn(0); 410 } 411 EXTERN_C_END 412 413 #undef __FUNCT__ 414 #define __FUNCT__ "TSView_CN" 415 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 416 { 417 PetscFunctionBegin; 418 PetscFunctionReturn(0); 419 } 420 421 /*MC 422 TSCN - ODE solver using the implicit Crank-Nicolson method. 423 424 Level: beginner 425 426 Notes: 427 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 428 429 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 430 431 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 432 433 M*/ 434 EXTERN_C_BEGIN 435 #undef __FUNCT__ 436 #define __FUNCT__ "TSCreate_CN" 437 PetscErrorCode TSCreate_CN(TS ts) 438 { 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 443 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 444 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 445 ts->ops->view = TSView_CN; 446 PetscFunctionReturn(0); 447 } 448 EXTERN_C_END 449