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