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 = TSRKInitializePackage();CHKERRQ(ierr); 370 ierr = PetscNew(&link);CHKERRQ(ierr); 371 t = &link->tab; 372 373 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 374 t->order = order; 375 t->s = s; 376 ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 377 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 378 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 379 else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 380 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 381 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 382 t->FSAL = PETSC_TRUE; 383 for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 384 385 if (bembed) { 386 ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 387 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 388 } 389 390 if (!binterp) { p = 1; binterp = t->b; } 391 t->p = p; 392 ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 393 ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 394 395 link->next = RKTableauList; 396 RKTableauList = link; 397 PetscFunctionReturn(0); 398 } 399 400 /* 401 The step completion formula is 402 403 x1 = x0 + h b^T YdotRHS 404 405 This function can be called before or after ts->vec_sol has been updated. 406 Suppose we have a completion formula (b) and an embedded formula (be) of different order. 407 We can write 408 409 x1e = x0 + h be^T YdotRHS 410 = x1 - h b^T YdotRHS + h be^T YdotRHS 411 = x1 + h (be - b)^T YdotRHS 412 413 so we can evaluate the method with different order even after the step has been optimistically completed. 414 */ 415 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 416 { 417 TS_RK *rk = (TS_RK*)ts->data; 418 RKTableau tab = rk->tableau; 419 PetscScalar *w = rk->work; 420 PetscReal h; 421 PetscInt s = tab->s,j; 422 PetscErrorCode ierr; 423 424 PetscFunctionBegin; 425 switch (rk->status) { 426 case TS_STEP_INCOMPLETE: 427 case TS_STEP_PENDING: 428 h = ts->time_step; break; 429 case TS_STEP_COMPLETE: 430 h = ts->ptime - ts->ptime_prev; break; 431 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 432 } 433 if (order == tab->order) { 434 if (rk->status == TS_STEP_INCOMPLETE) { 435 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 436 for (j=0; j<s; j++) w[j] = h*tab->b[j]; 437 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 438 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 439 PetscFunctionReturn(0); 440 } else if (order == tab->order-1) { 441 if (!tab->bembed) goto unavailable; 442 if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 443 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 444 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 445 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 446 } else { /* Rollback and re-complete using (be-b) */ 447 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 448 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 449 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 450 } 451 if (done) *done = PETSC_TRUE; 452 PetscFunctionReturn(0); 453 } 454 unavailable: 455 if (done) *done = PETSC_FALSE; 456 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); 457 PetscFunctionReturn(0); 458 } 459 460 static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 461 { 462 TS_RK *rk = (TS_RK*)ts->data; 463 RKTableau tab = rk->tableau; 464 const PetscInt s = tab->s; 465 const PetscReal *b = tab->b,*c = tab->c; 466 Vec *Y = rk->Y; 467 PetscInt i; 468 PetscErrorCode ierr; 469 470 PetscFunctionBegin; 471 /* backup cost integral */ 472 ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 473 for (i=s-1; i>=0; i--) { 474 /* Evolve ts->vec_costintegral to compute integrals */ 475 ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 476 ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 477 } 478 PetscFunctionReturn(0); 479 } 480 481 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 482 { 483 TS_RK *rk = (TS_RK*)ts->data; 484 RKTableau tab = rk->tableau; 485 const PetscInt s = tab->s; 486 const PetscReal *b = tab->b,*c = tab->c; 487 Vec *Y = rk->Y; 488 PetscInt i; 489 PetscErrorCode ierr; 490 491 PetscFunctionBegin; 492 for (i=s-1; i>=0; i--) { 493 /* Evolve ts->vec_costintegral to compute integrals */ 494 ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 495 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 496 } 497 PetscFunctionReturn(0); 498 } 499 500 static PetscErrorCode TSRollBack_RK(TS ts) 501 { 502 TS_RK *rk = (TS_RK*)ts->data; 503 RKTableau tab = rk->tableau; 504 const PetscInt s = tab->s; 505 const PetscReal *b = tab->b; 506 PetscScalar *w = rk->work; 507 Vec *YdotRHS = rk->YdotRHS; 508 PetscInt j; 509 PetscReal h; 510 PetscErrorCode ierr; 511 512 PetscFunctionBegin; 513 switch (rk->status) { 514 case TS_STEP_INCOMPLETE: 515 case TS_STEP_PENDING: 516 h = ts->time_step; break; 517 case TS_STEP_COMPLETE: 518 h = ts->ptime - ts->ptime_prev; break; 519 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 520 } 521 for (j=0; j<s; j++) w[j] = -h*b[j]; 522 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 523 PetscFunctionReturn(0); 524 } 525 526 static PetscErrorCode TSStep_RK(TS ts) 527 { 528 TS_RK *rk = (TS_RK*)ts->data; 529 RKTableau tab = rk->tableau; 530 const PetscInt s = tab->s; 531 const PetscReal *A = tab->A,*c = tab->c; 532 PetscScalar *w = rk->work; 533 Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 534 PetscBool FSAL = tab->FSAL; 535 TSAdapt adapt; 536 PetscInt i,j; 537 PetscInt rejections = 0; 538 PetscBool stageok,accept = PETSC_TRUE; 539 PetscReal next_time_step = ts->time_step; 540 PetscErrorCode ierr; 541 542 PetscFunctionBegin; 543 if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 544 if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 545 546 rk->status = TS_STEP_INCOMPLETE; 547 while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 548 PetscReal t = ts->ptime; 549 PetscReal h = ts->time_step; 550 for (i=0; i<s; i++) { 551 rk->stage_time = t + h*c[i]; 552 ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 553 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 554 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 555 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 556 ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 557 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 558 ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 559 if (!stageok) goto reject_step; 560 if (FSAL && !i) continue; 561 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 562 } 563 564 rk->status = TS_STEP_INCOMPLETE; 565 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 566 rk->status = TS_STEP_PENDING; 567 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 568 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 569 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 570 ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 571 rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 572 if (!accept) { /* Roll back the current step */ 573 ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 574 ts->time_step = next_time_step; 575 goto reject_step; 576 } 577 578 if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 579 rk->ptime = ts->ptime; 580 rk->time_step = ts->time_step; 581 } 582 583 ts->ptime += ts->time_step; 584 ts->time_step = next_time_step; 585 break; 586 587 reject_step: 588 ts->reject++; accept = PETSC_FALSE; 589 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 590 ts->reason = TS_DIVERGED_STEP_REJECTED; 591 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 592 } 593 } 594 PetscFunctionReturn(0); 595 } 596 597 static PetscErrorCode TSAdjointSetUp_RK(TS ts) 598 { 599 TS_RK *rk = (TS_RK*)ts->data; 600 RKTableau tab = rk->tableau; 601 PetscInt s = tab->s; 602 PetscErrorCode ierr; 603 604 PetscFunctionBegin; 605 if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 606 ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 607 if(ts->vecs_sensip) { 608 ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 609 } 610 PetscFunctionReturn(0); 611 } 612 613 static PetscErrorCode TSAdjointStep_RK(TS ts) 614 { 615 TS_RK *rk = (TS_RK*)ts->data; 616 RKTableau tab = rk->tableau; 617 const PetscInt s = tab->s; 618 const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 619 PetscScalar *w = rk->work; 620 Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu; 621 PetscInt i,j,nadj; 622 PetscReal t = ts->ptime; 623 PetscErrorCode ierr; 624 PetscReal h = ts->time_step; 625 626 PetscFunctionBegin; 627 rk->status = TS_STEP_INCOMPLETE; 628 for (i=s-1; i>=0; i--) { 629 Mat J; 630 PetscReal stage_time = t + h*(1.0-c[i]); 631 PetscBool zero = PETSC_FALSE; 632 633 ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 634 ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 635 if (ts->vec_costintegral) { 636 ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 637 } 638 /* Stage values of mu */ 639 if (ts->vecs_sensip) { 640 ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 641 if (ts->vec_costintegral) { 642 ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 643 } 644 } 645 646 if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 647 for (nadj=0; nadj<ts->numcost; nadj++) { 648 DM dm; 649 Vec VecSensiTemp; 650 651 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 652 ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 653 /* Stage values of lambda */ 654 if (!zero) { 655 ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr); 656 ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr); 657 for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; 658 ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 659 ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 660 } else { 661 ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr); 662 } 663 if (ts->vec_costintegral) { 664 ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 665 } 666 667 /* Stage values of mu */ 668 if (ts->vecs_sensip) { 669 if (!zero) { 670 ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 671 } else { 672 ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr); 673 } 674 if (ts->vec_costintegral) { 675 ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 676 } 677 } 678 ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr); 679 } 680 } 681 682 for (j=0; j<s; j++) w[j] = 1.0; 683 for (nadj=0; nadj<ts->numcost; nadj++) { 684 ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 685 if (ts->vecs_sensip) { 686 ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 687 } 688 } 689 rk->status = TS_STEP_COMPLETE; 690 PetscFunctionReturn(0); 691 } 692 693 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 694 { 695 TS_RK *rk = (TS_RK*)ts->data; 696 PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 697 PetscReal h; 698 PetscReal tt,t; 699 PetscScalar *b; 700 const PetscReal *B = rk->tableau->binterp; 701 PetscErrorCode ierr; 702 703 PetscFunctionBegin; 704 if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 705 706 switch (rk->status) { 707 case TS_STEP_INCOMPLETE: 708 case TS_STEP_PENDING: 709 h = ts->time_step; 710 t = (itime - ts->ptime)/h; 711 break; 712 case TS_STEP_COMPLETE: 713 h = ts->ptime - ts->ptime_prev; 714 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 715 break; 716 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 717 } 718 ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 719 for (i=0; i<s; i++) b[i] = 0; 720 for (j=0,tt=t; j<p; j++,tt*=t) { 721 for (i=0; i<s; i++) { 722 b[i] += h * B[i*p+j] * tt; 723 } 724 } 725 726 ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 727 ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 728 729 ierr = PetscFree(b);CHKERRQ(ierr); 730 PetscFunctionReturn(0); 731 } 732 733 /*------------------------------------------------------------*/ 734 735 static PetscErrorCode TSRKTableauReset(TS ts) 736 { 737 TS_RK *rk = (TS_RK*)ts->data; 738 RKTableau tab = rk->tableau; 739 PetscErrorCode ierr; 740 741 PetscFunctionBegin; 742 if (!tab) PetscFunctionReturn(0); 743 ierr = PetscFree(rk->work);CHKERRQ(ierr); 744 ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 745 ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 746 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 747 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 748 PetscFunctionReturn(0); 749 } 750 751 static PetscErrorCode TSReset_RK(TS ts) 752 { 753 TS_RK *rk = (TS_RK*)ts->data; 754 PetscErrorCode ierr; 755 756 PetscFunctionBegin; 757 ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 758 ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 759 PetscFunctionReturn(0); 760 } 761 762 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 763 { 764 PetscFunctionBegin; 765 PetscFunctionReturn(0); 766 } 767 768 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 769 { 770 PetscFunctionBegin; 771 PetscFunctionReturn(0); 772 } 773 774 775 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 776 { 777 PetscFunctionBegin; 778 PetscFunctionReturn(0); 779 } 780 781 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 782 { 783 784 PetscFunctionBegin; 785 PetscFunctionReturn(0); 786 } 787 /* 788 static PetscErrorCode RKSetAdjCoe(RKTableau tab) 789 { 790 PetscReal *A,*b; 791 PetscInt s,i,j; 792 PetscErrorCode ierr; 793 794 PetscFunctionBegin; 795 s = tab->s; 796 ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 797 798 for (i=0; i<s; i++) 799 for (j=0; j<s; j++) { 800 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]; 801 ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 802 } 803 804 for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 805 806 ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 807 ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 808 ierr = PetscFree2(A,b);CHKERRQ(ierr); 809 PetscFunctionReturn(0); 810 } 811 */ 812 813 static PetscErrorCode TSRKTableauSetUp(TS ts) 814 { 815 TS_RK *rk = (TS_RK*)ts->data; 816 RKTableau tab = rk->tableau; 817 PetscErrorCode ierr; 818 819 PetscFunctionBegin; 820 ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 821 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 822 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 823 PetscFunctionReturn(0); 824 } 825 826 827 static PetscErrorCode TSSetUp_RK(TS ts) 828 { 829 TS_RK *rk = (TS_RK*)ts->data; 830 PetscErrorCode ierr; 831 DM dm; 832 833 PetscFunctionBegin; 834 ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 835 ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 836 if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 837 ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 838 } 839 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 840 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 841 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845 846 /*------------------------------------------------------------*/ 847 848 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 849 { 850 TS_RK *rk = (TS_RK*)ts->data; 851 PetscErrorCode ierr; 852 853 PetscFunctionBegin; 854 ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 855 { 856 RKTableauLink link; 857 PetscInt count,choice; 858 PetscBool flg; 859 const char **namelist; 860 for (link=RKTableauList,count=0; link; link=link->next,count++) ; 861 ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 862 for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 863 ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 864 if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 865 ierr = PetscFree(namelist);CHKERRQ(ierr); 866 } 867 ierr = PetscOptionsTail();CHKERRQ(ierr); 868 PetscFunctionReturn(0); 869 } 870 871 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 872 { 873 TS_RK *rk = (TS_RK*)ts->data; 874 PetscBool iascii; 875 PetscErrorCode ierr; 876 877 PetscFunctionBegin; 878 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 879 if (iascii) { 880 RKTableau tab = rk->tableau; 881 TSRKType rktype; 882 char buf[512]; 883 ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 884 ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 885 ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 886 ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 887 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 888 ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 889 } 890 PetscFunctionReturn(0); 891 } 892 893 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 894 { 895 PetscErrorCode ierr; 896 TSAdapt adapt; 897 898 PetscFunctionBegin; 899 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 900 ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 901 PetscFunctionReturn(0); 902 } 903 904 /*@C 905 TSRKSetType - Set the type of RK scheme 906 907 Logically collective 908 909 Input Parameter: 910 + ts - timestepping context 911 - rktype - type of RK-scheme 912 913 Options Database: 914 . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 915 916 Level: intermediate 917 918 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 919 @*/ 920 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 921 { 922 PetscErrorCode ierr; 923 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 926 PetscValidCharPointer(rktype,2); 927 ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 928 PetscFunctionReturn(0); 929 } 930 931 /*@C 932 TSRKGetType - Get the type of RK scheme 933 934 Logically collective 935 936 Input Parameter: 937 . ts - timestepping context 938 939 Output Parameter: 940 . rktype - type of RK-scheme 941 942 Level: intermediate 943 944 .seealso: TSRKGetType() 945 @*/ 946 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 947 { 948 PetscErrorCode ierr; 949 950 PetscFunctionBegin; 951 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 952 ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 953 PetscFunctionReturn(0); 954 } 955 956 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 957 { 958 TS_RK *rk = (TS_RK*)ts->data; 959 960 PetscFunctionBegin; 961 *rktype = rk->tableau->name; 962 PetscFunctionReturn(0); 963 } 964 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 965 { 966 TS_RK *rk = (TS_RK*)ts->data; 967 PetscErrorCode ierr; 968 PetscBool match; 969 RKTableauLink link; 970 971 PetscFunctionBegin; 972 if (rk->tableau) { 973 ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 974 if (match) PetscFunctionReturn(0); 975 } 976 for (link = RKTableauList; link; link=link->next) { 977 ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 978 if (match) { 979 if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 980 rk->tableau = &link->tab; 981 if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 982 ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 983 PetscFunctionReturn(0); 984 } 985 } 986 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 987 PetscFunctionReturn(0); 988 } 989 990 static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 991 { 992 TS_RK *rk = (TS_RK*)ts->data; 993 994 PetscFunctionBegin; 995 if (ns) *ns = rk->tableau->s; 996 if (Y) *Y = rk->Y; 997 PetscFunctionReturn(0); 998 } 999 1000 static PetscErrorCode TSDestroy_RK(TS ts) 1001 { 1002 PetscErrorCode ierr; 1003 1004 PetscFunctionBegin; 1005 ierr = TSReset_RK(ts);CHKERRQ(ierr); 1006 if (ts->dm) { 1007 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1008 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1009 } 1010 ierr = PetscFree(ts->data);CHKERRQ(ierr); 1011 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1012 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 1013 PetscFunctionReturn(0); 1014 } 1015 1016 /* ------------------------------------------------------------ */ 1017 /*MC 1018 TSRK - ODE and DAE solver using Runge-Kutta schemes 1019 1020 The user should provide the right hand side of the equation 1021 using TSSetRHSFunction(). 1022 1023 Notes: 1024 The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1025 1026 Level: beginner 1027 1028 .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1029 TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1030 1031 M*/ 1032 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1033 { 1034 TS_RK *rk; 1035 PetscErrorCode ierr; 1036 1037 PetscFunctionBegin; 1038 ierr = TSRKInitializePackage();CHKERRQ(ierr); 1039 1040 ts->ops->reset = TSReset_RK; 1041 ts->ops->destroy = TSDestroy_RK; 1042 ts->ops->view = TSView_RK; 1043 ts->ops->load = TSLoad_RK; 1044 ts->ops->setup = TSSetUp_RK; 1045 ts->ops->adjointsetup = TSAdjointSetUp_RK; 1046 ts->ops->step = TSStep_RK; 1047 ts->ops->interpolate = TSInterpolate_RK; 1048 ts->ops->evaluatestep = TSEvaluateStep_RK; 1049 ts->ops->rollback = TSRollBack_RK; 1050 ts->ops->setfromoptions = TSSetFromOptions_RK; 1051 ts->ops->getstages = TSGetStages_RK; 1052 ts->ops->adjointstep = TSAdjointStep_RK; 1053 1054 ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 1055 ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1056 1057 ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 1058 ts->data = (void*)rk; 1059 1060 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1061 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1062 1063 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1064 PetscFunctionReturn(0); 1065 } 1066