1 /* 2 Code for time stepping with the Runge-Kutta method 3 4 Notes: 5 The general system is written as 6 7 Udot = F(t,U) 8 9 */ 10 11 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12 #include <petscdm.h> 13 #include <../src/ts/impls/explicit/rk/rk.h> 14 #include <../src/ts/impls/explicit/rk/mrk.h> 15 16 static TSRKType TSRKDefault = TSRK3BS; 17 static PetscBool TSRKRegisterAllCalled; 18 static PetscBool TSRKPackageInitialized; 19 20 /*MC 21 TSRK1FE - First order forward Euler scheme. 22 23 This method has one stage. 24 25 Options database: 26 . -ts_rk_type 1fe 27 28 Level: advanced 29 30 .seealso: TSRK, TSRKType, TSRKSetType() 31 M*/ 32 /*MC 33 TSRK2A - Second order RK scheme. 34 35 This method has two stages. 36 37 Options database: 38 . -ts_rk_type 2a 39 40 Level: advanced 41 42 .seealso: TSRK, TSRKType, TSRKSetType() 43 M*/ 44 /*MC 45 TSRK3 - Third order RK scheme. 46 47 This method has three stages. 48 49 Options database: 50 . -ts_rk_type 3 51 52 Level: advanced 53 54 .seealso: TSRK, TSRKType, TSRKSetType() 55 M*/ 56 /*MC 57 TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 58 59 This method has four stages with the First Same As Last (FSAL) property. 60 61 Options database: 62 . -ts_rk_type 3bs 63 64 Level: advanced 65 66 References: https://doi.org/10.1016/0893-9659(89)90079-7 67 68 .seealso: TSRK, TSRKType, TSRKSetType() 69 M*/ 70 /*MC 71 TSRK4 - Fourth order RK scheme. 72 73 This is the classical Runge-Kutta method with four stages. 74 75 Options database: 76 . -ts_rk_type 4 77 78 Level: advanced 79 80 .seealso: TSRK, TSRKType, TSRKSetType() 81 M*/ 82 /*MC 83 TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 84 85 This method has six stages. 86 87 Options database: 88 . -ts_rk_type 5f 89 90 Level: advanced 91 92 .seealso: TSRK, TSRKType, TSRKSetType() 93 M*/ 94 /*MC 95 TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 96 97 This method has seven stages with the First Same As Last (FSAL) property. 98 99 Options database: 100 . -ts_rk_type 5dp 101 102 Level: advanced 103 104 References: https://doi.org/10.1016/0771-050X(80)90013-3 105 106 .seealso: TSRK, TSRKType, TSRKSetType() 107 M*/ 108 /*MC 109 TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 110 111 This method has eight stages with the First Same As Last (FSAL) property. 112 113 Options database: 114 . -ts_rk_type 5bs 115 116 Level: advanced 117 118 References: https://doi.org/10.1016/0898-1221(96)00141-1 119 120 .seealso: TSRK, TSRKType, TSRKSetType() 121 M*/ 122 123 /*@C 124 TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 125 126 Not Collective, but should be called by all processes which will need the schemes to be registered 127 128 Level: advanced 129 130 .keywords: TS, TSRK, register, all 131 132 .seealso: TSRKRegisterDestroy() 133 @*/ 134 PetscErrorCode TSRKRegisterAll(void) 135 { 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 140 TSRKRegisterAllCalled = PETSC_TRUE; 141 142 #define RC PetscRealConstant 143 { 144 const PetscReal 145 A[1][1] = {{0}}, 146 b[1] = {RC(1.0)}; 147 ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 148 } 149 { 150 const PetscReal 151 A[2][2] = {{0,0}, 152 {RC(1.0),0}}, 153 b[2] = {RC(0.5),RC(0.5)}, 154 bembed[2] = {RC(1.0),0}; 155 ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 156 } 157 { 158 const PetscReal 159 A[3][3] = {{0,0,0}, 160 {RC(2.0)/RC(3.0),0,0}, 161 {RC(-1.0)/RC(3.0),RC(1.0),0}}, 162 b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 163 ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 164 } 165 { 166 const PetscReal 167 A[4][4] = {{0,0,0,0}, 168 {RC(1.0)/RC(2.0),0,0,0}, 169 {0,RC(3.0)/RC(4.0),0,0}, 170 {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 171 b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 172 bembed[4] = {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)}; 173 ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 174 } 175 { 176 const PetscReal 177 A[4][4] = {{0,0,0,0}, 178 {RC(0.5),0,0,0}, 179 {0,RC(0.5),0,0}, 180 {0,0,RC(1.0),0}}, 181 b[4] = {RC(1.0)/RC(6.0),RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),RC(1.0)/RC(6.0)}; 182 ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 183 } 184 { 185 const PetscReal 186 A[6][6] = {{0,0,0,0,0,0}, 187 {RC(0.25),0,0,0,0,0}, 188 {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 189 {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 190 {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 191 {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}}, 192 b[6] = {RC(16.0)/RC(135.0),0,RC(6656.0)/RC(12825.0),RC(28561.0)/RC(56430.0),RC(-9.0)/RC(50.0),RC(2.0)/RC(55.0)}, 193 bembed[6] = {RC(25.0)/RC(216.0),0,RC(1408.0)/RC(2565.0),RC(2197.0)/RC(4104.0),RC(-1.0)/RC(5.0),0}; 194 ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 195 } 196 { 197 const PetscReal 198 A[7][7] = {{0,0,0,0,0,0,0}, 199 {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 200 {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 201 {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 202 {RC(19372.0)/RC(6561.0),RC(-25360.0)/RC(2187.0),RC(64448.0)/RC(6561.0),RC(-212.0)/RC(729.0),0,0,0}, 203 {RC(9017.0)/RC(3168.0),RC(-355.0)/RC(33.0),RC(46732.0)/RC(5247.0),RC(49.0)/RC(176.0),RC(-5103.0)/RC(18656.0),0,0}, 204 {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}}, 205 b[7] = {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}, 206 bembed[7] = {RC(5179.0)/RC(57600.0),0,RC(7571.0)/RC(16695.0),RC(393.0)/RC(640.0),RC(-92097.0)/RC(339200.0),RC(187.0)/RC(2100.0),RC(1.0)/RC(40.0)}, 207 binterp[7][5] = {{RC(1.0),RC(-4034104133.0)/RC(1410260304.0),RC(105330401.0)/RC(33982176.0),RC(-13107642775.0)/RC(11282082432.0),RC(6542295.0)/RC(470086768.0)}, 208 {0,0,0,0,0}, 209 {0,RC(132343189600.0)/RC(32700410799.0),RC(-833316000.0)/RC(131326951.0),RC(91412856700.0)/RC(32700410799.0),RC(-523383600.0)/RC(10900136933.0)}, 210 {0,RC(-115792950.0)/RC(29380423.0),RC(185270875.0)/RC(16991088.0),RC(-12653452475.0)/RC(1880347072.0),RC(98134425.0)/RC(235043384.0)}, 211 {0,RC(70805911779.0)/RC(24914598704.0),RC(-4531260609.0)/RC(600351776.0),RC(988140236175.0)/RC(199316789632.0),RC(-14307999165.0)/RC(24914598704.0)}, 212 {0,RC(-331320693.0)/RC(205662961.0),RC(31361737.0)/RC(7433601.0),RC(-2426908385.0)/RC(822651844.0),RC(97305120.0)/RC(205662961.0)}, 213 {0,RC(44764047.0)/RC(29380423.0),RC(-1532549.0)/RC(353981.0),RC(90730570.0)/RC(29380423.0),RC(-8293050.0)/RC(29380423.0)}}; 214 215 ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 216 } 217 { 218 const PetscReal 219 A[8][8] = {{0,0,0,0,0,0,0,0}, 220 {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 221 {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 222 {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 223 {RC(68.0)/RC(297.0),RC(-4.0)/RC(11.0),RC(42.0)/RC(143.0),RC(1960.0)/RC(3861.0),0,0,0,0}, 224 {RC(597.0)/RC(22528.0),RC(81.0)/RC(352.0),RC(63099.0)/RC(585728.0),RC(58653.0)/RC(366080.0),RC(4617.0)/RC(20480.0),0,0,0}, 225 {RC(174197.0)/RC(959244.0),RC(-30942.0)/RC(79937.0),RC(8152137.0)/RC(19744439.0),RC(666106.0)/RC(1039181.0),RC(-29421.0)/RC(29068.0),RC(482048.0)/RC(414219.0),0,0}, 226 {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}}, 227 b[8] = {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}, 228 bembed[8] = {RC(2479.0)/RC(34992.0),0,RC(123.0)/RC(416.0),RC(612941.0)/RC(3411720.0),RC(43.0)/RC(1440.0),RC(2272.0)/RC(6561.0),RC(79937.0)/RC(1113912.0),RC(3293.0)/RC(556956.0)}; 229 ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 230 } 231 #undef RC 232 PetscFunctionReturn(0); 233 } 234 235 /*@C 236 TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 237 238 Not Collective 239 240 Level: advanced 241 242 .keywords: TSRK, register, destroy 243 .seealso: TSRKRegister(), TSRKRegisterAll() 244 @*/ 245 PetscErrorCode TSRKRegisterDestroy(void) 246 { 247 PetscErrorCode ierr; 248 RKTableauLink link; 249 250 PetscFunctionBegin; 251 while ((link = RKTableauList)) { 252 RKTableau t = &link->tab; 253 RKTableauList = link->next; 254 ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 255 ierr = PetscFree (t->bembed); CHKERRQ(ierr); 256 ierr = PetscFree (t->binterp); CHKERRQ(ierr); 257 ierr = PetscFree (t->name); CHKERRQ(ierr); 258 ierr = PetscFree (link); CHKERRQ(ierr); 259 } 260 TSRKRegisterAllCalled = PETSC_FALSE; 261 PetscFunctionReturn(0); 262 } 263 264 /*@C 265 TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 266 from TSInitializePackage(). 267 268 Level: developer 269 270 .keywords: TS, TSRK, initialize, package 271 .seealso: PetscInitialize() 272 @*/ 273 PetscErrorCode TSRKInitializePackage(void) 274 { 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 if (TSRKPackageInitialized) PetscFunctionReturn(0); 279 TSRKPackageInitialized = PETSC_TRUE; 280 ierr = TSRKRegisterAll();CHKERRQ(ierr); 281 ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 /*@C 286 TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 287 called from PetscFinalize(). 288 289 Level: developer 290 291 .keywords: Petsc, destroy, package 292 .seealso: PetscFinalize() 293 @*/ 294 PetscErrorCode TSRKFinalizePackage(void) 295 { 296 PetscErrorCode ierr; 297 298 PetscFunctionBegin; 299 TSRKPackageInitialized = PETSC_FALSE; 300 ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 301 PetscFunctionReturn(0); 302 } 303 304 /*@C 305 TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 306 307 Not Collective, but the same schemes should be registered on all processes on which they will be used 308 309 Input Parameters: 310 + name - identifier for method 311 . order - approximation order of method 312 . s - number of stages, this is the dimension of the matrices below 313 . A - stage coefficients (dimension s*s, row-major) 314 . b - step completion table (dimension s; NULL to use last row of A) 315 . c - abscissa (dimension s; NULL to use row sums of A) 316 . bembed - completion table for embedded method (dimension s; NULL if not available) 317 . p - Order of the interpolation scheme, equal to the number of columns of binterp 318 - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 319 320 Notes: 321 Several RK methods are provided, this function is only needed to create new methods. 322 323 Level: advanced 324 325 .keywords: TS, register 326 327 .seealso: TSRK 328 @*/ 329 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 330 const PetscReal A[],const PetscReal b[],const PetscReal c[], 331 const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 332 { 333 PetscErrorCode ierr; 334 RKTableauLink link; 335 RKTableau t; 336 PetscInt i,j; 337 338 PetscFunctionBegin; 339 PetscValidCharPointer(name,1); 340 PetscValidRealPointer(A,4); 341 if (b) PetscValidRealPointer(b,5); 342 if (c) PetscValidRealPointer(c,6); 343 if (bembed) PetscValidRealPointer(bembed,7); 344 if (binterp || p > 1) PetscValidRealPointer(binterp,9); 345 346 ierr = TSRKInitializePackage();CHKERRQ(ierr); 347 ierr = PetscNew(&link);CHKERRQ(ierr); 348 t = &link->tab; 349 350 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 351 t->order = order; 352 t->s = s; 353 ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 354 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 355 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 356 else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 357 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 358 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 359 t->FSAL = PETSC_TRUE; 360 for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 361 362 if (bembed) { 363 ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 364 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 365 } 366 367 if (!binterp) { p = 1; binterp = t->b; } 368 t->p = p; 369 ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 370 ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 371 372 link->next = RKTableauList; 373 RKTableauList = link; 374 PetscFunctionReturn(0); 375 } 376 377 /* 378 This is for single-step RK method 379 The step completion formula is 380 381 x1 = x0 + h b^T YdotRHS 382 383 This function can be called before or after ts->vec_sol has been updated. 384 Suppose we have a completion formula (b) and an embedded formula (be) of different order. 385 We can write 386 387 x1e = x0 + h be^T YdotRHS 388 = x1 - h b^T YdotRHS + h be^T YdotRHS 389 = x1 + h (be - b)^T YdotRHS 390 391 so we can evaluate the method with different order even after the step has been optimistically completed. 392 */ 393 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 394 { 395 TS_RK *rk = (TS_RK*)ts->data; 396 RKTableau tab = rk->tableau; 397 PetscScalar *w = rk->work; 398 PetscReal h; 399 PetscInt s = tab->s,j; 400 PetscErrorCode ierr; 401 402 PetscFunctionBegin; 403 switch (rk->status) { 404 case TS_STEP_INCOMPLETE: 405 case TS_STEP_PENDING: 406 h = ts->time_step; break; 407 case TS_STEP_COMPLETE: 408 h = ts->ptime - ts->ptime_prev; break; 409 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 410 } 411 if (order == tab->order) { 412 if (rk->status == TS_STEP_INCOMPLETE) { 413 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 414 for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 415 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 416 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 417 PetscFunctionReturn(0); 418 } else if (order == tab->order-1) { 419 if (!tab->bembed) goto unavailable; 420 if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 421 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 422 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 423 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 424 } else { /*Rollback and re-complete using (be-b) */ 425 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 426 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 427 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 428 } 429 if (done) *done = PETSC_TRUE; 430 PetscFunctionReturn(0); 431 } 432 unavailable: 433 if (done) *done = PETSC_FALSE; 434 else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 435 PetscFunctionReturn(0); 436 } 437 438 static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 439 { 440 TS_RK *rk = (TS_RK*)ts->data; 441 RKTableau tab = rk->tableau; 442 const PetscInt s = tab->s; 443 const PetscReal *b = tab->b,*c = tab->c; 444 Vec *Y = rk->Y; 445 PetscInt i; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 /* backup cost integral */ 450 ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 451 for (i=s-1; i>=0; i--) { 452 /* Evolve ts->vec_costintegral to compute integrals */ 453 ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 454 ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 455 } 456 PetscFunctionReturn(0); 457 } 458 459 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 460 { 461 TS_RK *rk = (TS_RK*)ts->data; 462 RKTableau tab = rk->tableau; 463 const PetscInt s = tab->s; 464 const PetscReal *b = tab->b,*c = tab->c; 465 Vec *Y = rk->Y; 466 PetscInt i; 467 PetscErrorCode ierr; 468 469 PetscFunctionBegin; 470 for (i=s-1; i>=0; i--) { 471 /* Evolve ts->vec_costintegral to compute integrals */ 472 ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 473 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 474 } 475 PetscFunctionReturn(0); 476 } 477 478 static PetscErrorCode TSRollBack_RK(TS ts) 479 { 480 TS_RK *rk = (TS_RK*)ts->data; 481 RKTableau tab = rk->tableau; 482 const PetscInt s = tab->s; 483 const PetscReal *b = tab->b; 484 PetscScalar *w = rk->work; 485 Vec *YdotRHS = rk->YdotRHS; 486 PetscInt j; 487 PetscReal h; 488 PetscErrorCode ierr; 489 490 PetscFunctionBegin; 491 switch (rk->status) { 492 case TS_STEP_INCOMPLETE: 493 case TS_STEP_PENDING: 494 h = ts->time_step; break; 495 case TS_STEP_COMPLETE: 496 h = ts->ptime - ts->ptime_prev; break; 497 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 498 } 499 for (j=0; j<s; j++) w[j] = -h*b[j]; 500 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 501 PetscFunctionReturn(0); 502 } 503 504 static PetscErrorCode TSStep_RK(TS ts) 505 { 506 TS_RK *rk = (TS_RK*)ts->data; 507 RKTableau tab = rk->tableau; 508 const PetscInt s = tab->s; 509 const PetscReal *A = tab->A,*c = tab->c; 510 PetscScalar *w = rk->work; 511 Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 512 PetscBool FSAL = tab->FSAL; 513 TSAdapt adapt; 514 PetscInt i,j; 515 PetscInt rejections = 0; 516 PetscBool stageok,accept = PETSC_TRUE; 517 PetscReal next_time_step = ts->time_step; 518 PetscErrorCode ierr; 519 520 PetscFunctionBegin; 521 if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 522 if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 523 524 rk->status = TS_STEP_INCOMPLETE; 525 while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 526 PetscReal t = ts->ptime; 527 PetscReal h = ts->time_step; 528 for (i=0; i<s; i++) { 529 rk->stage_time = t + h*c[i]; 530 ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 531 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 532 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 533 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 534 ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 535 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 536 ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 537 if (!stageok) goto reject_step; 538 if (FSAL && !i) continue; 539 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 540 } 541 542 rk->status = TS_STEP_INCOMPLETE; 543 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 544 rk->status = TS_STEP_PENDING; 545 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 546 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 547 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 548 ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 549 rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 550 if (!accept) { /* Roll back the current step */ 551 ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 552 ts->time_step = next_time_step; 553 goto reject_step; 554 } 555 556 if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 557 rk->ptime = ts->ptime; 558 rk->time_step = ts->time_step; 559 } 560 561 ts->ptime += ts->time_step; 562 ts->time_step = next_time_step; 563 break; 564 565 reject_step: 566 ts->reject++; accept = PETSC_FALSE; 567 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 568 ts->reason = TS_DIVERGED_STEP_REJECTED; 569 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 570 } 571 } 572 PetscFunctionReturn(0); 573 } 574 575 static PetscErrorCode TSAdjointSetUp_RK(TS ts) 576 { 577 TS_RK *rk = (TS_RK*)ts->data; 578 RKTableau tab = rk->tableau; 579 PetscInt s = tab->s; 580 PetscErrorCode ierr; 581 582 PetscFunctionBegin; 583 if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 584 ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 585 if(ts->vecs_sensip) { 586 ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 587 } 588 PetscFunctionReturn(0); 589 } 590 591 static PetscErrorCode TSAdjointStep_RK(TS ts) 592 { 593 TS_RK *rk = (TS_RK*)ts->data; 594 RKTableau tab = rk->tableau; 595 const PetscInt s = tab->s; 596 const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 597 PetscScalar *w = rk->work; 598 Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu; 599 PetscInt i,j,nadj; 600 PetscReal t = ts->ptime; 601 PetscErrorCode ierr; 602 PetscReal h = ts->time_step; 603 604 PetscFunctionBegin; 605 rk->status = TS_STEP_INCOMPLETE; 606 for (i=s-1; i>=0; i--) { 607 Mat J; 608 PetscReal stage_time = t + h*(1.0-c[i]); 609 PetscBool zero = PETSC_FALSE; 610 611 ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 612 ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 613 if (ts->vec_costintegral) { 614 ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 615 } 616 /* Stage values of mu */ 617 if (ts->vecs_sensip) { 618 ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 619 if (ts->vec_costintegral) { 620 ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 621 } 622 } 623 624 if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 625 for (nadj=0; nadj<ts->numcost; nadj++) { 626 DM dm; 627 Vec VecSensiTemp; 628 629 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 630 ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 631 /* Stage values of lambda */ 632 if (!zero) { 633 ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr); 634 ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr); 635 for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; 636 ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 637 ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 638 } else { 639 ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr); 640 } 641 if (ts->vec_costintegral) { 642 ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 643 } 644 645 /* Stage values of mu */ 646 if (ts->vecs_sensip) { 647 if (!zero) { 648 ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 649 } else { 650 ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr); 651 } 652 if (ts->vec_costintegral) { 653 ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 654 } 655 } 656 ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 657 } 658 } 659 660 for (j=0; j<s; j++) w[j] = 1.0; 661 for (nadj=0; nadj<ts->numcost; nadj++) { 662 ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 663 if (ts->vecs_sensip) { 664 ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 665 } 666 } 667 rk->status = TS_STEP_COMPLETE; 668 PetscFunctionReturn(0); 669 } 670 671 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 672 { 673 TS_RK *rk = (TS_RK*)ts->data; 674 PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 675 PetscReal h; 676 PetscReal tt,t; 677 PetscScalar *b; 678 const PetscReal *B = rk->tableau->binterp; 679 PetscErrorCode ierr; 680 681 PetscFunctionBegin; 682 if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 683 684 switch (rk->status) { 685 case TS_STEP_INCOMPLETE: 686 case TS_STEP_PENDING: 687 h = ts->time_step; 688 t = (itime - ts->ptime)/h; 689 break; 690 case TS_STEP_COMPLETE: 691 h = ts->ptime - ts->ptime_prev; 692 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 693 break; 694 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 695 } 696 ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 697 for (i=0; i<s; i++) b[i] = 0; 698 for (j=0,tt=t; j<p; j++,tt*=t) { 699 for (i=0; i<s; i++) { 700 b[i] += h * B[i*p+j] * tt; 701 } 702 } 703 ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 704 ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 705 ierr = PetscFree(b);CHKERRQ(ierr); 706 PetscFunctionReturn(0); 707 } 708 709 /*------------------------------------------------------------*/ 710 711 static PetscErrorCode TSRKTableauReset(TS ts) 712 { 713 TS_RK *rk = (TS_RK*)ts->data; 714 RKTableau tab = rk->tableau; 715 PetscErrorCode ierr; 716 717 PetscFunctionBegin; 718 if (!tab) PetscFunctionReturn(0); 719 ierr = PetscFree(rk->work);CHKERRQ(ierr); 720 ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 721 ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 722 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 723 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 724 PetscFunctionReturn(0); 725 } 726 727 static PetscErrorCode TSReset_RK(TS ts) 728 { 729 TS_RK *rk = (TS_RK*)ts->data; 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 734 ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 735 if (ts->use_splitrhsfunction) { 736 ierr = PetscTryMethod(ts,"TSReset_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 737 } else { 738 ierr = PetscTryMethod(ts,"TSReset_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 739 } 740 PetscFunctionReturn(0); 741 } 742 743 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 744 { 745 PetscFunctionBegin; 746 PetscFunctionReturn(0); 747 } 748 749 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 750 { 751 PetscFunctionBegin; 752 PetscFunctionReturn(0); 753 } 754 755 756 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 757 { 758 PetscFunctionBegin; 759 PetscFunctionReturn(0); 760 } 761 762 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 763 { 764 765 PetscFunctionBegin; 766 PetscFunctionReturn(0); 767 } 768 /* 769 static PetscErrorCode RKSetAdjCoe(RKTableau tab) 770 { 771 PetscReal *A,*b; 772 PetscInt s,i,j; 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 s = tab->s; 777 ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 778 779 for (i=0; i<s; i++) 780 for (j=0; j<s; j++) { 781 A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i]; 782 ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 783 } 784 785 for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 786 787 ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 788 ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 789 ierr = PetscFree2(A,b);CHKERRQ(ierr); 790 PetscFunctionReturn(0); 791 } 792 */ 793 794 static PetscErrorCode TSRKTableauSetUp(TS ts) 795 { 796 TS_RK *rk = (TS_RK*)ts->data; 797 RKTableau tab = rk->tableau; 798 PetscErrorCode ierr; 799 800 PetscFunctionBegin; 801 ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 802 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 803 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 804 PetscFunctionReturn(0); 805 } 806 807 static PetscErrorCode TSSetUp_RK(TS ts) 808 { 809 TS_RK *rk = (TS_RK*)ts->data; 810 PetscErrorCode ierr; 811 DM dm; 812 813 PetscFunctionBegin; 814 ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 815 ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 816 if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 817 ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 818 } 819 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 820 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 821 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 822 if (ts->use_splitrhsfunction) { 823 ierr = PetscTryMethod(ts,"TSSetUp_MRKSPLIT_C",(TS),(ts));CHKERRQ(ierr); 824 } else { 825 ierr = PetscTryMethod(ts,"TSSetUp_MRKNONSPLIT_C",(TS),(ts));CHKERRQ(ierr); 826 } 827 PetscFunctionReturn(0); 828 } 829 830 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 831 { 832 TS_RK *rk = (TS_RK*)ts->data; 833 PetscErrorCode ierr; 834 835 PetscFunctionBegin; 836 ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 837 { 838 RKTableauLink link; 839 PetscInt count,choice; 840 PetscBool flg,use_multirate = PETSC_FALSE; 841 const char **namelist; 842 843 for (link=RKTableauList,count=0; link; link=link->next,count++) ; 844 ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 845 for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 846 ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 847 if (flg) { 848 ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 849 } 850 ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 851 if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 852 ierr = PetscFree(namelist);CHKERRQ(ierr); 853 } 854 ierr = PetscOptionsTail();CHKERRQ(ierr); 855 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 856 ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 857 ierr = PetscOptionsEnd();CHKERRQ(ierr); 858 PetscFunctionReturn(0); 859 } 860 861 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 862 { 863 TS_RK *rk = (TS_RK*)ts->data; 864 PetscBool iascii; 865 PetscErrorCode ierr; 866 867 PetscFunctionBegin; 868 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 869 if (iascii) { 870 RKTableau tab = rk->tableau; 871 TSRKType rktype; 872 char buf[512]; 873 ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 874 ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 875 ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 876 ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 877 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 878 ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 879 } 880 PetscFunctionReturn(0); 881 } 882 883 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 884 { 885 PetscErrorCode ierr; 886 TSAdapt adapt; 887 888 PetscFunctionBegin; 889 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 890 ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 891 PetscFunctionReturn(0); 892 } 893 894 /*@C 895 TSRKSetType - Set the type of RK scheme 896 897 Logically collective 898 899 Input Parameter: 900 + ts - timestepping context 901 - rktype - type of RK-scheme 902 903 Options Database: 904 . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 905 906 Level: intermediate 907 908 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 909 @*/ 910 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 911 { 912 PetscErrorCode ierr; 913 914 PetscFunctionBegin; 915 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 916 PetscValidCharPointer(rktype,2); 917 ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 918 PetscFunctionReturn(0); 919 } 920 921 /*@C 922 TSRKGetType - Get the type of RK scheme 923 924 Logically collective 925 926 Input Parameter: 927 . ts - timestepping context 928 929 Output Parameter: 930 . rktype - type of RK-scheme 931 932 Level: intermediate 933 934 .seealso: TSRKGetType() 935 @*/ 936 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 937 { 938 PetscErrorCode ierr; 939 940 PetscFunctionBegin; 941 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 942 ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 947 { 948 TS_RK *rk = (TS_RK*)ts->data; 949 950 PetscFunctionBegin; 951 *rktype = rk->tableau->name; 952 PetscFunctionReturn(0); 953 } 954 955 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 956 { 957 TS_RK *rk = (TS_RK*)ts->data; 958 PetscErrorCode ierr; 959 PetscBool match; 960 RKTableauLink link; 961 962 PetscFunctionBegin; 963 if (rk->tableau) { 964 ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 965 if (match) PetscFunctionReturn(0); 966 } 967 for (link = RKTableauList; link; link=link->next) { 968 ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 969 if (match) { 970 if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 971 rk->tableau = &link->tab; 972 if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 973 ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 974 PetscFunctionReturn(0); 975 } 976 } 977 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 978 PetscFunctionReturn(0); 979 } 980 981 static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 982 { 983 TS_RK *rk = (TS_RK*)ts->data; 984 985 PetscFunctionBegin; 986 if (ns) *ns = rk->tableau->s; 987 if (Y) *Y = rk->Y; 988 PetscFunctionReturn(0); 989 } 990 991 static PetscErrorCode TSDestroy_RK(TS ts) 992 { 993 PetscErrorCode ierr; 994 995 PetscFunctionBegin; 996 ierr = TSReset_RK(ts);CHKERRQ(ierr); 997 if (ts->dm) { 998 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 999 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1000 } 1001 ierr = PetscFree(ts->data);CHKERRQ(ierr); 1002 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1003 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 1004 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 1005 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1006 PetscFunctionReturn(0); 1007 } 1008 1009 /*@C 1010 TSRKSetMultirate - Use the interpolation-based multirate RK method 1011 1012 Logically collective 1013 1014 Input Parameter: 1015 + ts - timestepping context 1016 - use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2 1017 1018 Options Database: 1019 . -ts_rk_multirate - <true,false> 1020 1021 Notes: 1022 The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp). 1023 1024 Level: intermediate 1025 1026 .seealso: TSRKGetMultirate() 1027 @*/ 1028 PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate) 1029 { 1030 PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate)); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /*@C 1035 TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 1036 1037 Not collective 1038 1039 Input Parameter: 1040 . ts - timestepping context 1041 1042 Output Parameter: 1043 . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 1044 1045 Level: intermediate 1046 1047 .seealso: TSRKSetMultirate() 1048 @*/ 1049 PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate) 1050 { 1051 PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate)); 1052 PetscFunctionReturn(0); 1053 } 1054 1055 /*MC 1056 TSRK - ODE and DAE solver using Runge-Kutta schemes 1057 1058 The user should provide the right hand side of the equation 1059 using TSSetRHSFunction(). 1060 1061 Notes: 1062 The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1063 1064 Level: beginner 1065 1066 .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1067 TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1068 1069 M*/ 1070 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1071 { 1072 TS_RK *rk; 1073 PetscErrorCode ierr; 1074 1075 PetscFunctionBegin; 1076 ierr = TSRKInitializePackage();CHKERRQ(ierr); 1077 1078 ts->ops->reset = TSReset_RK; 1079 ts->ops->destroy = TSDestroy_RK; 1080 ts->ops->view = TSView_RK; 1081 ts->ops->load = TSLoad_RK; 1082 ts->ops->setup = TSSetUp_RK; 1083 ts->ops->adjointsetup = TSAdjointSetUp_RK; 1084 ts->ops->interpolate = TSInterpolate_RK; 1085 ts->ops->step = TSStep_RK; 1086 ts->ops->evaluatestep = TSEvaluateStep_RK; 1087 ts->ops->rollback = TSRollBack_RK; 1088 ts->ops->setfromoptions = TSSetFromOptions_RK; 1089 ts->ops->getstages = TSGetStages_RK; 1090 ts->ops->adjointstep = TSAdjointStep_RK; 1091 1092 ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 1093 ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1094 1095 ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 1096 ts->data = (void*)rk; 1097 1098 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1099 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1100 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr); 1101 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr); 1102 1103 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1104 rk->dtratio = 1; 1105 PetscFunctionReturn(0); 1106 } 1107