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_prev = ts->time_step; 54 ts->time_step = next_time_step; 55 ts->steps++; 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "TSInterpolate_Theta" 61 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 62 { 63 TS_Theta *th = (TS_Theta*)ts->data; 64 PetscReal alpha = t - ts->ptime; 65 PetscErrorCode ierr; 66 67 PetscFunctionBegin; 68 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 69 if (th->endpoint) alpha *= th->Theta; 70 ierr = VecWAXPY(X,alpha,th->Xdot,th->X);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 /*------------------------------------------------------------*/ 75 #undef __FUNCT__ 76 #define __FUNCT__ "TSReset_Theta" 77 static PetscErrorCode TSReset_Theta(TS ts) 78 { 79 TS_Theta *th = (TS_Theta*)ts->data; 80 PetscErrorCode ierr; 81 82 PetscFunctionBegin; 83 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 84 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 85 ierr = VecDestroy(&th->affine);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 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","",PETSC_NULL);CHKERRQ(ierr); 99 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_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 TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 214 { 215 TS_Theta *th = (TS_Theta*)ts->data; 216 217 PetscFunctionBegin; 218 th->endpoint = flg; 219 PetscFunctionReturn(0); 220 } 221 EXTERN_C_END 222 223 /* ------------------------------------------------------------ */ 224 /*MC 225 TSTHETA - DAE solver using the implicit Theta method 226 227 Level: beginner 228 229 Notes: 230 This method can be applied to DAE. 231 232 This method is cast as a 1-stage implicit Runge-Kutta method. 233 234 .vb 235 Theta | Theta 236 ------------- 237 | 1 238 .ve 239 240 For the default Theta=0.5, this is also known as the implicit midpoint rule. 241 242 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 243 244 .vb 245 0 | 0 0 246 1 | 1-Theta Theta 247 ------------------- 248 | 1-Theta Theta 249 .ve 250 251 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 252 253 To apply a diagonally implicit RK method to DAE, the stage formula 254 255 $ Y_i = X + h sum_j a_ij Y'_j 256 257 is interpreted as a formula for Y'_i in terms of Y_i and known stuff (Y'_j, j<i) 258 259 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 260 261 M*/ 262 EXTERN_C_BEGIN 263 #undef __FUNCT__ 264 #define __FUNCT__ "TSCreate_Theta" 265 PetscErrorCode TSCreate_Theta(TS ts) 266 { 267 TS_Theta *th; 268 PetscErrorCode ierr; 269 270 PetscFunctionBegin; 271 ts->ops->reset = TSReset_Theta; 272 ts->ops->destroy = TSDestroy_Theta; 273 ts->ops->view = TSView_Theta; 274 ts->ops->setup = TSSetUp_Theta; 275 ts->ops->step = TSStep_Theta; 276 ts->ops->interpolate = TSInterpolate_Theta; 277 ts->ops->setfromoptions = TSSetFromOptions_Theta; 278 ts->ops->snesfunction = SNESTSFormFunction_Theta; 279 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 280 281 ierr = PetscNewLog(ts,TS_Theta,&th);CHKERRQ(ierr); 282 ts->data = (void*)th; 283 284 th->extrapolate = PETSC_FALSE; 285 th->Theta = 0.5; 286 287 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaGetTheta_C","TSThetaGetTheta_Theta",TSThetaGetTheta_Theta);CHKERRQ(ierr); 288 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetTheta_C","TSThetaSetTheta_Theta",TSThetaSetTheta_Theta);CHKERRQ(ierr); 289 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSThetaSetEndpoint_C","TSThetaSetEndpoint_Theta",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 290 PetscFunctionReturn(0); 291 } 292 EXTERN_C_END 293 294 #undef __FUNCT__ 295 #define __FUNCT__ "TSThetaGetTheta" 296 /*@ 297 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 298 299 Not Collective 300 301 Input Parameter: 302 . ts - timestepping context 303 304 Output Parameter: 305 . theta - stage abscissa 306 307 Note: 308 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 309 310 Level: Advanced 311 312 .seealso: TSThetaSetTheta() 313 @*/ 314 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 315 { 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 320 PetscValidPointer(theta,2); 321 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "TSThetaSetTheta" 327 /*@ 328 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 329 330 Not Collective 331 332 Input Parameter: 333 + ts - timestepping context 334 - theta - stage abscissa 335 336 Options Database: 337 . -ts_theta_theta <theta> 338 339 Level: Intermediate 340 341 .seealso: TSThetaGetTheta() 342 @*/ 343 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 344 { 345 PetscErrorCode ierr; 346 347 PetscFunctionBegin; 348 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 349 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 353 #undef __FUNCT__ 354 #define __FUNCT__ "TSThetaSetEndpoint" 355 /*@ 356 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 357 358 Not Collective 359 360 Input Parameter: 361 + ts - timestepping context 362 - flg - PETSC_TRUE to use the endpoint variant 363 364 Options Database: 365 . -ts_theta_endpoint <flg> 366 367 Level: Intermediate 368 369 .seealso: TSTHETA, TSCN 370 @*/ 371 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 372 { 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 377 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 378 PetscFunctionReturn(0); 379 } 380 381 /* 382 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 383 * The creation functions for these specializations are below. 384 */ 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "TSView_BEuler" 388 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 389 { 390 PetscErrorCode ierr; 391 392 PetscFunctionBegin; 393 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 394 PetscFunctionReturn(0); 395 } 396 397 /*MC 398 TSBEULER - ODE solver using the implicit backward Euler method 399 400 Level: beginner 401 402 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 403 404 M*/ 405 EXTERN_C_BEGIN 406 #undef __FUNCT__ 407 #define __FUNCT__ "TSCreate_BEuler" 408 PetscErrorCode TSCreate_BEuler(TS ts) 409 { 410 PetscErrorCode ierr; 411 412 PetscFunctionBegin; 413 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 414 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 415 ts->ops->view = TSView_BEuler; 416 PetscFunctionReturn(0); 417 } 418 EXTERN_C_END 419 420 #undef __FUNCT__ 421 #define __FUNCT__ "TSView_CN" 422 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 423 { 424 PetscErrorCode ierr; 425 426 PetscFunctionBegin; 427 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 431 /*MC 432 TSCN - ODE solver using the implicit Crank-Nicolson method. 433 434 Level: beginner 435 436 Notes: 437 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 438 439 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 440 441 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 442 443 M*/ 444 EXTERN_C_BEGIN 445 #undef __FUNCT__ 446 #define __FUNCT__ "TSCreate_CN" 447 PetscErrorCode TSCreate_CN(TS ts) 448 { 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 453 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 454 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 455 ts->ops->view = TSView_CN; 456 PetscFunctionReturn(0); 457 } 458 EXTERN_C_END 459