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