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