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