1 /* 2 Code for time stepping with the Runge-Kutta method 3 4 Notes: 5 The general system is written as 6 7 Udot = F(t,U) 8 9 */ 10 11 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 12 #include <petscdm.h> 13 #include <../src/ts/impls/explicit/rk/rk.h> 14 #include <../src/ts/impls/explicit/rk/mrk.h> 15 16 static TSRKType TSRKDefault = TSRK3BS; 17 static PetscBool TSRKRegisterAllCalled; 18 static PetscBool TSRKPackageInitialized; 19 20 static RKTableauLink RKTableauList; 21 22 /*MC 23 TSRK1FE - First order forward Euler scheme. 24 25 This method has one stage. 26 27 Options database: 28 . -ts_rk_type 1fe 29 30 Level: advanced 31 32 .seealso: TSRK, TSRKType, TSRKSetType() 33 M*/ 34 /*MC 35 TSRK2A - Second order RK scheme. 36 37 This method has two stages. 38 39 Options database: 40 . -ts_rk_type 2a 41 42 Level: advanced 43 44 .seealso: TSRK, TSRKType, TSRKSetType() 45 M*/ 46 /*MC 47 TSRK3 - Third order RK scheme. 48 49 This method has three stages. 50 51 Options database: 52 . -ts_rk_type 3 53 54 Level: advanced 55 56 .seealso: TSRK, TSRKType, TSRKSetType() 57 M*/ 58 /*MC 59 TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 60 61 This method has four stages with the First Same As Last (FSAL) property. 62 63 Options database: 64 . -ts_rk_type 3bs 65 66 Level: advanced 67 68 References: https://doi.org/10.1016/0893-9659(89)90079-7 69 70 .seealso: TSRK, TSRKType, TSRKSetType() 71 M*/ 72 /*MC 73 TSRK4 - Fourth order RK scheme. 74 75 This is the classical Runge-Kutta method with four stages. 76 77 Options database: 78 . -ts_rk_type 4 79 80 Level: advanced 81 82 .seealso: TSRK, TSRKType, TSRKSetType() 83 M*/ 84 /*MC 85 TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 86 87 This method has six stages. 88 89 Options database: 90 . -ts_rk_type 5f 91 92 Level: advanced 93 94 .seealso: TSRK, TSRKType, TSRKSetType() 95 M*/ 96 /*MC 97 TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 98 99 This method has seven stages with the First Same As Last (FSAL) property. 100 101 Options database: 102 . -ts_rk_type 5dp 103 104 Level: advanced 105 106 References: https://doi.org/10.1016/0771-050X(80)90013-3 107 108 .seealso: TSRK, TSRKType, TSRKSetType() 109 M*/ 110 /*MC 111 TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 112 113 This method has eight stages with the First Same As Last (FSAL) property. 114 115 Options database: 116 . -ts_rk_type 5bs 117 118 Level: advanced 119 120 References: https://doi.org/10.1016/0898-1221(96)00141-1 121 122 .seealso: TSRK, TSRKType, TSRKSetType() 123 M*/ 124 125 /*@C 126 TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 127 128 Not Collective, but should be called by all processes which will need the schemes to be registered 129 130 Level: advanced 131 132 .keywords: TS, TSRK, register, all 133 134 .seealso: TSRKRegisterDestroy() 135 @*/ 136 PetscErrorCode TSRKRegisterAll(void) 137 { 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 142 TSRKRegisterAllCalled = PETSC_TRUE; 143 144 #define RC PetscRealConstant 145 { 146 const PetscReal 147 A[1][1] = {{0}}, 148 b[1] = {RC(1.0)}; 149 ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 150 } 151 { 152 const PetscReal 153 A[2][2] = {{0,0}, 154 {RC(1.0),0}}, 155 b[2] = {RC(0.5),RC(0.5)}, 156 bembed[2] = {RC(1.0),0}; 157 ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 158 } 159 { 160 const PetscReal 161 A[3][3] = {{0,0,0}, 162 {RC(2.0)/RC(3.0),0,0}, 163 {RC(-1.0)/RC(3.0),RC(1.0),0}}, 164 b[3] = {RC(0.25),RC(0.5),RC(0.25)}; 165 ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 166 } 167 { 168 const PetscReal 169 A[4][4] = {{0,0,0,0}, 170 {RC(1.0)/RC(2.0),0,0,0}, 171 {0,RC(3.0)/RC(4.0),0,0}, 172 {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}}, 173 b[4] = {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}, 174 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)}; 175 ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 176 } 177 { 178 const PetscReal 179 A[4][4] = {{0,0,0,0}, 180 {RC(0.5),0,0,0}, 181 {0,RC(0.5),0,0}, 182 {0,0,RC(1.0),0}}, 183 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)}; 184 ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 185 } 186 { 187 const PetscReal 188 A[6][6] = {{0,0,0,0,0,0}, 189 {RC(0.25),0,0,0,0,0}, 190 {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0}, 191 {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0}, 192 {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0}, 193 {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}}, 194 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)}, 195 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}; 196 ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 197 } 198 { 199 const PetscReal 200 A[7][7] = {{0,0,0,0,0,0,0}, 201 {RC(1.0)/RC(5.0),0,0,0,0,0,0}, 202 {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0}, 203 {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0}, 204 {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}, 205 {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}, 206 {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}}, 207 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}, 208 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)}, 209 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)}, 210 {0,0,0,0,0}, 211 {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)}, 212 {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)}, 213 {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)}, 214 {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)}, 215 {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)}}; 216 217 ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,binterp[0]);CHKERRQ(ierr); 218 } 219 { 220 const PetscReal 221 A[8][8] = {{0,0,0,0,0,0,0,0}, 222 {RC(1.0)/RC(6.0),0,0,0,0,0,0,0}, 223 {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0}, 224 {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0}, 225 {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}, 226 {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}, 227 {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}, 228 {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}}, 229 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}, 230 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)}; 231 ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 232 } 233 #undef RC 234 PetscFunctionReturn(0); 235 } 236 237 /*@C 238 TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 239 240 Not Collective 241 242 Level: advanced 243 244 .keywords: TSRK, register, destroy 245 .seealso: TSRKRegister(), TSRKRegisterAll() 246 @*/ 247 PetscErrorCode TSRKRegisterDestroy(void) 248 { 249 PetscErrorCode ierr; 250 RKTableauLink link; 251 252 PetscFunctionBegin; 253 while ((link = RKTableauList)) { 254 RKTableau t = &link->tab; 255 RKTableauList = link->next; 256 ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 257 ierr = PetscFree (t->bembed); CHKERRQ(ierr); 258 ierr = PetscFree (t->binterp); CHKERRQ(ierr); 259 ierr = PetscFree (t->name); CHKERRQ(ierr); 260 ierr = PetscFree (link); CHKERRQ(ierr); 261 } 262 TSRKRegisterAllCalled = PETSC_FALSE; 263 PetscFunctionReturn(0); 264 } 265 266 /*@C 267 TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 268 from TSInitializePackage(). 269 270 Level: developer 271 272 .keywords: TS, TSRK, initialize, package 273 .seealso: PetscInitialize() 274 @*/ 275 PetscErrorCode TSRKInitializePackage(void) 276 { 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 if (TSRKPackageInitialized) PetscFunctionReturn(0); 281 TSRKPackageInitialized = PETSC_TRUE; 282 ierr = TSRKRegisterAll();CHKERRQ(ierr); 283 ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 /*@C 288 TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 289 called from PetscFinalize(). 290 291 Level: developer 292 293 .keywords: Petsc, destroy, package 294 .seealso: PetscFinalize() 295 @*/ 296 PetscErrorCode TSRKFinalizePackage(void) 297 { 298 PetscErrorCode ierr; 299 300 PetscFunctionBegin; 301 TSRKPackageInitialized = PETSC_FALSE; 302 ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 /*@C 307 TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 308 309 Not Collective, but the same schemes should be registered on all processes on which they will be used 310 311 Input Parameters: 312 + name - identifier for method 313 . order - approximation order of method 314 . s - number of stages, this is the dimension of the matrices below 315 . A - stage coefficients (dimension s*s, row-major) 316 . b - step completion table (dimension s; NULL to use last row of A) 317 . c - abscissa (dimension s; NULL to use row sums of A) 318 . bembed - completion table for embedded method (dimension s; NULL if not available) 319 . p - Order of the interpolation scheme, equal to the number of columns of binterp 320 - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 321 322 Notes: 323 Several RK methods are provided, this function is only needed to create new methods. 324 325 Level: advanced 326 327 .keywords: TS, register 328 329 .seealso: TSRK 330 @*/ 331 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 332 const PetscReal A[],const PetscReal b[],const PetscReal c[], 333 const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 334 { 335 PetscErrorCode ierr; 336 RKTableauLink link; 337 RKTableau t; 338 PetscInt i,j; 339 340 PetscFunctionBegin; 341 PetscValidCharPointer(name,1); 342 PetscValidRealPointer(A,4); 343 if (b) PetscValidRealPointer(b,5); 344 if (c) PetscValidRealPointer(c,6); 345 if (bembed) PetscValidRealPointer(bembed,7); 346 if (binterp || p > 1) PetscValidRealPointer(binterp,9); 347 348 ierr = TSRKInitializePackage();CHKERRQ(ierr); 349 ierr = PetscNew(&link);CHKERRQ(ierr); 350 t = &link->tab; 351 352 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 353 t->order = order; 354 t->s = s; 355 ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 356 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 357 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 358 else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 359 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 360 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 361 t->FSAL = PETSC_TRUE; 362 for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 363 364 if (bembed) { 365 ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 366 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 367 } 368 369 if (!binterp) { p = 1; binterp = t->b; } 370 t->p = p; 371 ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 372 ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 373 374 link->next = RKTableauList; 375 RKTableauList = link; 376 PetscFunctionReturn(0); 377 } 378 379 /* 380 This is for single-step RK method 381 The step completion formula is 382 383 x1 = x0 + h b^T YdotRHS 384 385 This function can be called before or after ts->vec_sol has been updated. 386 Suppose we have a completion formula (b) and an embedded formula (be) of different order. 387 We can write 388 389 x1e = x0 + h be^T YdotRHS 390 = x1 - h b^T YdotRHS + h be^T YdotRHS 391 = x1 + h (be - b)^T YdotRHS 392 393 so we can evaluate the method with different order even after the step has been optimistically completed. 394 */ 395 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 396 { 397 TS_RK *rk = (TS_RK*)ts->data; 398 RKTableau tab = rk->tableau; 399 PetscScalar *w = rk->work; 400 PetscReal h; 401 PetscInt s = tab->s,j; 402 PetscErrorCode ierr; 403 404 PetscFunctionBegin; 405 switch (rk->status) { 406 case TS_STEP_INCOMPLETE: 407 case TS_STEP_PENDING: 408 h = ts->time_step; break; 409 case TS_STEP_COMPLETE: 410 h = ts->ptime - ts->ptime_prev; break; 411 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 412 } 413 if (order == tab->order) { 414 if (rk->status == TS_STEP_INCOMPLETE) { 415 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 416 for (j=0; j<s; j++) w[j] = h*tab->b[j]/rk->dtratio; 417 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 418 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 419 PetscFunctionReturn(0); 420 } else if (order == tab->order-1) { 421 if (!tab->bembed) goto unavailable; 422 if (rk->status == TS_STEP_INCOMPLETE) { /*Complete with the embedded method (be)*/ 423 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 424 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 425 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 426 } else { /*Rollback and re-complete using (be-b) */ 427 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 428 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 429 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 430 } 431 if (done) *done = PETSC_TRUE; 432 PetscFunctionReturn(0); 433 } 434 unavailable: 435 if (done) *done = PETSC_FALSE; 436 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); 437 PetscFunctionReturn(0); 438 } 439 440 static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 441 { 442 TS_RK *rk = (TS_RK*)ts->data; 443 RKTableau tab = rk->tableau; 444 const PetscInt s = tab->s; 445 const PetscReal *b = tab->b,*c = tab->c; 446 Vec *Y = rk->Y; 447 PetscInt i; 448 PetscErrorCode ierr; 449 450 PetscFunctionBegin; 451 /* backup cost integral */ 452 ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 453 for (i=s-1; i>=0; i--) { 454 /* Evolve ts->vec_costintegral to compute integrals */ 455 ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*c[i],Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 456 ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 457 } 458 PetscFunctionReturn(0); 459 } 460 461 static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 462 { 463 TS_RK *rk = (TS_RK*)ts->data; 464 RKTableau tab = rk->tableau; 465 const PetscInt s = tab->s; 466 const PetscReal *b = tab->b,*c = tab->c; 467 Vec *Y = rk->Y; 468 PetscInt i; 469 PetscErrorCode ierr; 470 471 PetscFunctionBegin; 472 for (i=s-1; i>=0; i--) { 473 /* Evolve ts->vec_costintegral to compute integrals */ 474 ierr = TSComputeCostIntegrand(ts,ts->ptime+ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 475 ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 476 } 477 PetscFunctionReturn(0); 478 } 479 480 static PetscErrorCode TSRollBack_RK(TS ts) 481 { 482 TS_RK *rk = (TS_RK*)ts->data; 483 RKTableau tab = rk->tableau; 484 const PetscInt s = tab->s; 485 const PetscReal *b = tab->b; 486 PetscScalar *w = rk->work; 487 Vec *YdotRHS = rk->YdotRHS; 488 PetscInt j; 489 PetscReal h; 490 PetscErrorCode ierr; 491 492 PetscFunctionBegin; 493 switch (rk->status) { 494 case TS_STEP_INCOMPLETE: 495 case TS_STEP_PENDING: 496 h = ts->time_step; break; 497 case TS_STEP_COMPLETE: 498 h = ts->ptime - ts->ptime_prev; break; 499 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 500 } 501 for (j=0; j<s; j++) w[j] = -h*b[j]; 502 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505 506 static PetscErrorCode TSForwardStep_RK(TS ts) 507 { 508 TS_RK *rk = (TS_RK*)ts->data; 509 RKTableau tab = rk->tableau; 510 Mat J,*MatsFwdSensipTemp = rk->MatsFwdSensipTemp; 511 const PetscInt s = tab->s; 512 const PetscReal *A = tab->A,*c = tab->c,*b = tab->b; 513 Vec *Y = rk->Y; 514 PetscInt i,j; 515 PetscReal stage_time,h = ts->time_step; 516 PetscBool zero; 517 PetscErrorCode ierr; 518 519 PetscFunctionBegin; 520 ierr = MatCopy(ts->mat_sensip,rk->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 521 ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 522 523 for (i=0; i<s; i++) { 524 stage_time = ts->ptime + h*c[i]; 525 zero = PETSC_FALSE; 526 if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 527 /* TLM Stage values */ 528 if(!i) { 529 ierr = MatCopy(ts->mat_sensip,rk->MatsFwdStageSensip[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 530 } else if (!zero) { 531 ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 532 for (j=0; j<i; j++) { 533 ierr = MatAXPY(rk->MatsFwdStageSensip[i],h*A[i*s+j],MatsFwdSensipTemp[j],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 534 } 535 ierr = MatAXPY(rk->MatsFwdStageSensip[i],1.,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 536 } else { 537 ierr = MatZeroEntries(rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 538 } 539 540 ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 541 ierr = MatMatMult(J,rk->MatsFwdStageSensip[i],MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatsFwdSensipTemp[i]);CHKERRQ(ierr); 542 if (ts->Jacprhs) { 543 ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 544 if (ts->vecs_sensi2p) { /* TLM used for 2nd-order adjoint */ 545 PetscScalar *xarr; 546 ierr = MatDenseGetColumn(MatsFwdSensipTemp[i],0,&xarr);CHKERRQ(ierr); 547 ierr = VecPlaceArray(rk->VecDeltaFwdSensipCol,xarr);CHKERRQ(ierr); 548 ierr = MatMultAdd(ts->Jacprhs,ts->vec_dir,rk->VecDeltaFwdSensipCol,rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 549 ierr = VecResetArray(rk->VecDeltaFwdSensipCol);CHKERRQ(ierr);CHKERRQ(ierr); 550 ierr = MatDenseRestoreColumn(MatsFwdSensipTemp[i],&xarr);CHKERRQ(ierr); 551 } else { 552 ierr = MatAXPY(MatsFwdSensipTemp[i],1.,ts->Jacprhs,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 553 } 554 } 555 } 556 557 for (i=0; i<s; i++) { 558 ierr = MatAXPY(ts->mat_sensip,h*b[i],rk->MatsFwdSensipTemp[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr); 559 } 560 rk->status = TS_STEP_COMPLETE; 561 PetscFunctionReturn(0); 562 } 563 564 static PetscErrorCode TSForwardGetStages_RK(TS ts,PetscInt *ns,Mat **stagesensip) 565 { 566 TS_RK *rk = (TS_RK*)ts->data; 567 RKTableau tab = rk->tableau; 568 569 PetscFunctionBegin; 570 if (ns) *ns = tab->s; 571 if (stagesensip) *stagesensip = rk->MatsFwdStageSensip; 572 PetscFunctionReturn(0); 573 } 574 575 static PetscErrorCode TSForwardSetUp_RK(TS ts) 576 { 577 TS_RK *rk = (TS_RK*)ts->data; 578 RKTableau tab = rk->tableau; 579 PetscInt i; 580 PetscErrorCode ierr; 581 582 PetscFunctionBegin; 583 /* backup sensitivity results for roll-backs */ 584 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatFwdSensip0);CHKERRQ(ierr); 585 586 ierr = PetscMalloc1(tab->s,&rk->MatsFwdStageSensip);CHKERRQ(ierr); 587 ierr = PetscMalloc1(tab->s,&rk->MatsFwdSensipTemp);CHKERRQ(ierr); 588 for(i=0; i<tab->s; i++) { 589 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 590 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 591 } 592 ierr = VecDuplicate(ts->vec_sol,&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 static PetscErrorCode TSForwardReset_RK(TS ts) 597 { 598 TS_RK *rk = (TS_RK*)ts->data; 599 RKTableau tab = rk->tableau; 600 PetscInt i; 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 ierr = MatDestroy(&rk->MatFwdSensip0);CHKERRQ(ierr); 605 if (rk->MatsFwdStageSensip) { 606 for (i=0; i<tab->s; i++) { 607 ierr = MatDestroy(&rk->MatsFwdStageSensip[i]);CHKERRQ(ierr); 608 } 609 ierr = PetscFree(rk->MatsFwdStageSensip);CHKERRQ(ierr); 610 } 611 if (rk->MatsFwdSensipTemp) { 612 for (i=0; i<tab->s; i++) { 613 ierr = MatDestroy(&rk->MatsFwdSensipTemp[i]);CHKERRQ(ierr); 614 } 615 ierr = PetscFree(rk->MatsFwdSensipTemp);CHKERRQ(ierr); 616 } 617 ierr = VecDestroy(&rk->VecDeltaFwdSensipCol);CHKERRQ(ierr); 618 PetscFunctionReturn(0); 619 } 620 621 static PetscErrorCode TSStep_RK(TS ts) 622 { 623 TS_RK *rk = (TS_RK*)ts->data; 624 RKTableau tab = rk->tableau; 625 const PetscInt s = tab->s; 626 const PetscReal *A = tab->A,*c = tab->c; 627 PetscScalar *w = rk->work; 628 Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 629 PetscBool FSAL = tab->FSAL; 630 TSAdapt adapt; 631 PetscInt i,j; 632 PetscInt rejections = 0; 633 PetscBool stageok,accept = PETSC_TRUE; 634 PetscReal next_time_step = ts->time_step; 635 PetscErrorCode ierr; 636 637 PetscFunctionBegin; 638 if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 639 if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 640 641 rk->status = TS_STEP_INCOMPLETE; 642 while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 643 PetscReal t = ts->ptime; 644 PetscReal h = ts->time_step; 645 for (i=0; i<s; i++) { 646 rk->stage_time = t + h*c[i]; 647 ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 648 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 649 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 650 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 651 ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 652 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 653 ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 654 if (!stageok) goto reject_step; 655 if (FSAL && !i) continue; 656 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 657 } 658 659 rk->status = TS_STEP_INCOMPLETE; 660 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 661 rk->status = TS_STEP_PENDING; 662 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 663 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 664 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 665 ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 666 rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 667 if (!accept) { /* Roll back the current step */ 668 ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 669 ts->time_step = next_time_step; 670 goto reject_step; 671 } 672 673 if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 674 rk->ptime = ts->ptime; 675 rk->time_step = ts->time_step; 676 } 677 678 ts->ptime += ts->time_step; 679 ts->time_step = next_time_step; 680 break; 681 682 reject_step: 683 ts->reject++; accept = PETSC_FALSE; 684 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 685 ts->reason = TS_DIVERGED_STEP_REJECTED; 686 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 687 } 688 } 689 PetscFunctionReturn(0); 690 } 691 692 static PetscErrorCode TSAdjointSetUp_RK(TS ts) 693 { 694 TS_RK *rk = (TS_RK*)ts->data; 695 RKTableau tab = rk->tableau; 696 PetscInt s = tab->s; 697 PetscErrorCode ierr; 698 699 PetscFunctionBegin; 700 if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 701 ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 702 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 703 if(ts->vecs_sensip) { 704 ierr = VecDuplicate(ts->vecs_sensip[0],&rk->VecDeltaMu);CHKERRQ(ierr); 705 } 706 if (ts->vecs_sensi2) { 707 ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 708 ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 709 } 710 if (ts->vecs_sensi2p) { 711 ierr = VecDuplicate(ts->vecs_sensi2p[0],&rk->VecDeltaMu2);CHKERRQ(ierr); 712 } 713 PetscFunctionReturn(0); 714 } 715 716 static PetscErrorCode TSAdjointStep_RK(TS ts) 717 { 718 TS_RK *rk = (TS_RK*)ts->data; 719 RKTableau tab = rk->tableau; 720 Mat J; 721 const PetscInt s = tab->s; 722 const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 723 PetscScalar *w = rk->work,*xarr; 724 Vec *Y = rk->Y,*VecsDeltaLam = rk->VecsDeltaLam,VecDeltaMu = rk->VecDeltaMu,*VecsSensiTemp = rk->VecsSensiTemp; 725 Vec *VecsDeltaLam2 = rk->VecsDeltaLam2,VecDeltaMu2 = rk->VecDeltaMu2,*VecsSensi2Temp = rk->VecsSensi2Temp; 726 PetscInt i,j,nadj; 727 PetscReal t = ts->ptime; 728 PetscReal h = ts->time_step,stage_time; 729 PetscBool zero; 730 PetscErrorCode ierr; 731 732 PetscFunctionBegin; 733 rk->status = TS_STEP_INCOMPLETE; 734 735 for (i=s-1; i>=0; i--) { 736 if (tab->FSAL && i == s-1) { 737 /* VecsDeltaLam[nadj*s+s-1] are initialized with zeros and the values never change.*/ 738 continue; 739 } 740 stage_time = t + h*(1.0-c[i]); 741 zero = PETSC_FALSE; 742 ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 743 ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr); 744 ierr = TSComputeDRDUFunction(ts,stage_time,Y[i],ts->vecs_drdu);CHKERRQ(ierr); 745 if (ts->vecs_sensip) { 746 ierr = TSComputeRHSJacobianP(ts,stage_time,Y[i],ts->Jacprhs);CHKERRQ(ierr); /* get f_p */ 747 ierr = TSComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 748 } 749 750 if (b[i] == 0 && i == s-1) zero = PETSC_TRUE; 751 752 for (j=i+1; j<s; j++) w[j-i-1] = -h*A[j*s+i]; /* coefficients for computing VecsSensiTemp */ 753 754 for (nadj=0; nadj<ts->numcost; nadj++) { 755 /* Stage values of lambda */ 756 if (!zero) { 757 /* b_i*lambda_{n+1} + \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 758 ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); /* VecDeltaLam is an vec array of size s by numcost */ 759 ierr = VecScale(VecsSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 760 ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 761 ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 762 } else { 763 /* \sum_{j=i+1}^s a_{ji}*lambda_{s,j} */ 764 ierr = VecSet(VecsSensiTemp[nadj],0);CHKERRQ(ierr); 765 ierr = VecMAXPY(VecsSensiTemp[nadj],s-i-1,w,&VecsDeltaLam[nadj*s+i+1]);CHKERRQ(ierr); 766 ierr = MatMultTranspose(J,VecsSensiTemp[nadj],VecsDeltaLam[nadj*s+i]);CHKERRQ(ierr); 767 ierr = VecScale(VecsDeltaLam[nadj*s+i],-h);CHKERRQ(ierr); 768 } 769 if (ts->vec_costintegral) { 770 ierr = VecAXPY(VecsDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdu[nadj]);CHKERRQ(ierr); 771 } 772 773 /* Stage values of mu */ 774 if (ts->vecs_sensip) { 775 if (!zero) { 776 ierr = MatMultTranspose(ts->Jacprhs,VecsSensiTemp[nadj],VecDeltaMu);CHKERRQ(ierr); 777 } else { 778 ierr = VecSet(VecDeltaMu,0);CHKERRQ(ierr); 779 } 780 if (ts->vec_costintegral) { 781 ierr = VecAXPY(VecDeltaMu,-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 782 } 783 ierr = VecAXPY(ts->vecs_sensip[nadj],1.,VecDeltaMu);CHKERRQ(ierr); /* update sensip for each stage */ 784 } 785 } 786 787 if (ts->vecs_sensi2 && ts->forward_solve) { /* 2nd-order adjoint, TLM mode has to be turned on */ 788 /* Get w1 at t_{n+1} from TLM matrix */ 789 ierr = MatDenseGetColumn(rk->MatsFwdStageSensip[i],0,&xarr);CHKERRQ(ierr); 790 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 791 /* lambda_s^T F_UU w_1 */ 792 ierr = TSComputeRHSHessianProductFunction1(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_guu);CHKERRQ(ierr); 793 /* lambda_s^T F_UP w_2 */ 794 ierr = TSComputeRHSHessianProductFunction2(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gup);CHKERRQ(ierr); 795 if (ts->vecs_sensi2p) { 796 /* lambda_s^T F_PU w_1 */ 797 ierr = TSComputeRHSHessianProductFunction3(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_sensip_col,ts->vecs_gpu);CHKERRQ(ierr); 798 /* lambda_s^T F_PU w_2 */ 799 ierr = TSComputeRHSHessianProductFunction4(ts,stage_time,Y[i],VecsSensiTemp,ts->vec_dir,ts->vecs_gpp);CHKERRQ(ierr); 800 } 801 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 802 ierr = MatDenseRestoreColumn(rk->MatsFwdStageSensip[i],&xarr);CHKERRQ(ierr); 803 804 for (nadj=0; nadj<ts->numcost; nadj++) { 805 /* Stage values of lambda */ 806 if (!zero) { 807 /* h*J_i^T*(b_i*Lambda_{n+1}+\sum_{j=i+1}^s a_{ji}*Lambda_{s,j} */ 808 ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 809 ierr = VecScale(VecsSensi2Temp[nadj],-h*b[i]);CHKERRQ(ierr); 810 ierr = VecMAXPY(VecsSensi2Temp[nadj],s-i-1,w,&VecsDeltaLam2[nadj*s+i+1]);CHKERRQ(ierr); 811 ierr = MatMultTranspose(J,VecsSensi2Temp[nadj],VecsDeltaLam2[nadj*s+i]);CHKERRQ(ierr); 812 ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],1.,ts->vecs_guu[nadj]);CHKERRQ(ierr); 813 if (ts->vecs_gup) { 814 ierr = VecAXPY(VecsDeltaLam2[nadj*s+i],1.,ts->vecs_gup[nadj]);CHKERRQ(ierr); 815 } 816 } else { 817 ierr = VecSet(VecsDeltaLam2[nadj*s+i],0);CHKERRQ(ierr); 818 } 819 if (ts->vecs_sensi2p) { /* 2nd-order adjoint */ 820 if (!zero) { 821 ierr = MatMultTranspose(ts->Jacprhs,VecsSensi2Temp[nadj],VecDeltaMu2);CHKERRQ(ierr); 822 } else { 823 ierr = VecSet(VecDeltaMu2,0);CHKERRQ(ierr); 824 } 825 if (ts->vecs_gpu) { 826 ierr = VecAXPY(VecDeltaMu2,1,ts->vecs_gpu[nadj]);CHKERRQ(ierr); 827 } 828 if (ts->vecs_gpp) { 829 ierr = VecAXPY(VecDeltaMu2,1,ts->vecs_gpp[nadj]);CHKERRQ(ierr); 830 } 831 ierr = VecAXPY(ts->vecs_sensi2p[nadj],1,VecDeltaMu2);CHKERRQ(ierr); /* update sensi2p for each stage */ 832 } 833 } 834 } 835 } 836 837 for (j=0; j<s; j++) w[j] = 1.0; 838 for (nadj=0; nadj<ts->numcost; nadj++) { /* no need to do this for mu's */ 839 ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecsDeltaLam[nadj*s]);CHKERRQ(ierr); 840 if (ts->vecs_sensi2) { 841 ierr = VecMAXPY(ts->vecs_sensi2[nadj],s,w,&VecsDeltaLam2[nadj*s]);CHKERRQ(ierr); 842 } 843 } 844 rk->status = TS_STEP_COMPLETE; 845 PetscFunctionReturn(0); 846 } 847 848 static PetscErrorCode TSAdjointReset_RK(TS ts) 849 { 850 TS_RK *rk = (TS_RK*)ts->data; 851 RKTableau tab = rk->tableau; 852 PetscErrorCode ierr; 853 854 PetscFunctionBegin; 855 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 856 ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 857 ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 858 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam2);CHKERRQ(ierr); 859 ierr = VecDestroy(&rk->VecDeltaMu2);CHKERRQ(ierr); 860 ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensi2Temp);CHKERRQ(ierr); 861 PetscFunctionReturn(0); 862 } 863 864 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 865 { 866 TS_RK *rk = (TS_RK*)ts->data; 867 PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 868 PetscReal h; 869 PetscReal tt,t; 870 PetscScalar *b; 871 const PetscReal *B = rk->tableau->binterp; 872 PetscErrorCode ierr; 873 874 PetscFunctionBegin; 875 if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 876 877 switch (rk->status) { 878 case TS_STEP_INCOMPLETE: 879 case TS_STEP_PENDING: 880 h = ts->time_step; 881 t = (itime - ts->ptime)/h; 882 break; 883 case TS_STEP_COMPLETE: 884 h = ts->ptime - ts->ptime_prev; 885 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 886 break; 887 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 888 } 889 ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 890 for (i=0; i<s; i++) b[i] = 0; 891 for (j=0,tt=t; j<p; j++,tt*=t) { 892 for (i=0; i<s; i++) { 893 b[i] += h * B[i*p+j] * tt; 894 } 895 } 896 ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 897 ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 898 ierr = PetscFree(b);CHKERRQ(ierr); 899 PetscFunctionReturn(0); 900 } 901 902 /*------------------------------------------------------------*/ 903 904 static PetscErrorCode TSRKTableauReset(TS ts) 905 { 906 TS_RK *rk = (TS_RK*)ts->data; 907 RKTableau tab = rk->tableau; 908 PetscErrorCode ierr; 909 910 PetscFunctionBegin; 911 if (!tab) PetscFunctionReturn(0); 912 ierr = PetscFree(rk->work);CHKERRQ(ierr); 913 ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 914 ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 915 ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecsDeltaLam);CHKERRQ(ierr); 916 ierr = VecDestroyVecs(ts->numcost,&rk->VecsSensiTemp);CHKERRQ(ierr); 917 ierr = VecDestroy(&rk->VecDeltaMu);CHKERRQ(ierr); 918 PetscFunctionReturn(0); 919 } 920 921 static PetscErrorCode TSReset_RK(TS ts) 922 { 923 TS_RK *rk = (TS_RK*)ts->data; 924 PetscErrorCode ierr; 925 926 PetscFunctionBegin; 927 ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 928 ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 929 if (ts->use_splitrhsfunction) { 930 ierr = PetscTryMethod(ts,"TSReset_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 931 } else { 932 ierr = PetscTryMethod(ts,"TSReset_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 933 } 934 PetscFunctionReturn(0); 935 } 936 937 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 938 { 939 PetscFunctionBegin; 940 PetscFunctionReturn(0); 941 } 942 943 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 944 { 945 PetscFunctionBegin; 946 PetscFunctionReturn(0); 947 } 948 949 950 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 951 { 952 PetscFunctionBegin; 953 PetscFunctionReturn(0); 954 } 955 956 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 957 { 958 959 PetscFunctionBegin; 960 PetscFunctionReturn(0); 961 } 962 /* 963 static PetscErrorCode RKSetAdjCoe(RKTableau tab) 964 { 965 PetscReal *A,*b; 966 PetscInt s,i,j; 967 PetscErrorCode ierr; 968 969 PetscFunctionBegin; 970 s = tab->s; 971 ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 972 973 for (i=0; i<s; i++) 974 for (j=0; j<s; j++) { 975 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]; 976 ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 977 } 978 979 for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 980 981 ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 982 ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 983 ierr = PetscFree2(A,b);CHKERRQ(ierr); 984 PetscFunctionReturn(0); 985 } 986 */ 987 988 static PetscErrorCode TSRKTableauSetUp(TS ts) 989 { 990 TS_RK *rk = (TS_RK*)ts->data; 991 RKTableau tab = rk->tableau; 992 PetscErrorCode ierr; 993 994 PetscFunctionBegin; 995 ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 996 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 997 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 998 PetscFunctionReturn(0); 999 } 1000 1001 static PetscErrorCode TSSetUp_RK(TS ts) 1002 { 1003 TS_RK *rk = (TS_RK*)ts->data; 1004 PetscErrorCode ierr; 1005 DM dm; 1006 1007 PetscFunctionBegin; 1008 ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 1009 ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 1010 if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 1011 ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 1012 } 1013 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1014 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1015 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1016 if (ts->use_splitrhsfunction) { 1017 ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateSplit_C",(TS),(ts));CHKERRQ(ierr); 1018 } else { 1019 ierr = PetscTryMethod(ts,"TSSetUp_RK_MultirateNonsplit_C",(TS),(ts));CHKERRQ(ierr); 1020 } 1021 PetscFunctionReturn(0); 1022 } 1023 1024 static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 1025 { 1026 TS_RK *rk = (TS_RK*)ts->data; 1027 PetscErrorCode ierr; 1028 1029 PetscFunctionBegin; 1030 ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 1031 { 1032 RKTableauLink link; 1033 PetscInt count,choice; 1034 PetscBool flg,use_multirate = PETSC_FALSE; 1035 const char **namelist; 1036 1037 for (link=RKTableauList,count=0; link; link=link->next,count++) ; 1038 ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 1039 for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 1040 ierr = PetscOptionsBool("-ts_rk_multirate","Use interpolation-based multirate RK method","TSRKSetMultirate",rk->use_multirate,&use_multirate,&flg);CHKERRQ(ierr); 1041 if (flg) { 1042 ierr = TSRKSetMultirate(ts,use_multirate);CHKERRQ(ierr); 1043 } 1044 ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 1045 if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 1046 ierr = PetscFree(namelist);CHKERRQ(ierr); 1047 } 1048 ierr = PetscOptionsTail();CHKERRQ(ierr); 1049 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Multirate methods options","");CHKERRQ(ierr); 1050 ierr = PetscOptionsInt("-ts_rk_dtratio","time step ratio between slow and fast","",rk->dtratio,&rk->dtratio,NULL);CHKERRQ(ierr); 1051 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1052 PetscFunctionReturn(0); 1053 } 1054 1055 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 1056 { 1057 TS_RK *rk = (TS_RK*)ts->data; 1058 PetscBool iascii; 1059 PetscErrorCode ierr; 1060 1061 PetscFunctionBegin; 1062 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1063 if (iascii) { 1064 RKTableau tab = rk->tableau; 1065 TSRKType rktype; 1066 char buf[512]; 1067 ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 1068 ierr = PetscViewerASCIIPrintf(viewer," RK type %s\n",rktype);CHKERRQ(ierr); 1069 ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 1070 ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 1071 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1072 ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 1073 } 1074 PetscFunctionReturn(0); 1075 } 1076 1077 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 1078 { 1079 PetscErrorCode ierr; 1080 TSAdapt adapt; 1081 1082 PetscFunctionBegin; 1083 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 1084 ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1085 PetscFunctionReturn(0); 1086 } 1087 1088 /*@C 1089 TSRKSetType - Set the type of RK scheme 1090 1091 Logically collective 1092 1093 Input Parameter: 1094 + ts - timestepping context 1095 - rktype - type of RK-scheme 1096 1097 Options Database: 1098 . -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 1099 1100 Level: intermediate 1101 1102 .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 1103 @*/ 1104 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 1105 { 1106 PetscErrorCode ierr; 1107 1108 PetscFunctionBegin; 1109 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1110 PetscValidCharPointer(rktype,2); 1111 ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1112 PetscFunctionReturn(0); 1113 } 1114 1115 /*@C 1116 TSRKGetType - Get the type of RK scheme 1117 1118 Logically collective 1119 1120 Input Parameter: 1121 . ts - timestepping context 1122 1123 Output Parameter: 1124 . rktype - type of RK-scheme 1125 1126 Level: intermediate 1127 1128 .seealso: TSRKGetType() 1129 @*/ 1130 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1131 { 1132 PetscErrorCode ierr; 1133 1134 PetscFunctionBegin; 1135 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1136 ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1137 PetscFunctionReturn(0); 1138 } 1139 1140 static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1141 { 1142 TS_RK *rk = (TS_RK*)ts->data; 1143 1144 PetscFunctionBegin; 1145 *rktype = rk->tableau->name; 1146 PetscFunctionReturn(0); 1147 } 1148 1149 static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1150 { 1151 TS_RK *rk = (TS_RK*)ts->data; 1152 PetscErrorCode ierr; 1153 PetscBool match; 1154 RKTableauLink link; 1155 1156 PetscFunctionBegin; 1157 if (rk->tableau) { 1158 ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1159 if (match) PetscFunctionReturn(0); 1160 } 1161 for (link = RKTableauList; link; link=link->next) { 1162 ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1163 if (match) { 1164 if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1165 rk->tableau = &link->tab; 1166 if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 1167 ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1168 PetscFunctionReturn(0); 1169 } 1170 } 1171 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1172 PetscFunctionReturn(0); 1173 } 1174 1175 static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1176 { 1177 TS_RK *rk = (TS_RK*)ts->data; 1178 1179 PetscFunctionBegin; 1180 if (ns) *ns = rk->tableau->s; 1181 if (Y) *Y = rk->Y; 1182 PetscFunctionReturn(0); 1183 } 1184 1185 static PetscErrorCode TSDestroy_RK(TS ts) 1186 { 1187 PetscErrorCode ierr; 1188 1189 PetscFunctionBegin; 1190 ierr = TSReset_RK(ts);CHKERRQ(ierr); 1191 if (ts->dm) { 1192 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 1193 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 1194 } 1195 ierr = PetscFree(ts->data);CHKERRQ(ierr); 1196 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 1197 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 1198 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",NULL);CHKERRQ(ierr); 1199 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",NULL);CHKERRQ(ierr); 1200 PetscFunctionReturn(0); 1201 } 1202 1203 /*@C 1204 TSRKSetMultirate - Use the interpolation-based multirate RK method 1205 1206 Logically collective 1207 1208 Input Parameter: 1209 + ts - timestepping context 1210 - use_multirate - PETSC_TRUE enables the multirate RK method, sets the basic method to be RK2A and sets the ratio between slow stepsize and fast stepsize to be 2 1211 1212 Options Database: 1213 . -ts_rk_multirate - <true,false> 1214 1215 Notes: 1216 The multirate method requires interpolation. The default interpolation works for 1st- and 2nd- order RK, but not for high-order RKs except TSRK5DP which comes with the interpolation coeffcients (binterp). 1217 1218 Level: intermediate 1219 1220 .seealso: TSRKGetMultirate() 1221 @*/ 1222 PetscErrorCode TSRKSetMultirate(TS ts,PetscBool use_multirate) 1223 { 1224 PetscErrorCode ierr; 1225 1226 PetscFunctionBegin; 1227 ierr = PetscTryMethod(ts,"TSRKSetMultirate_C",(TS,PetscBool),(ts,use_multirate));CHKERRQ(ierr); 1228 PetscFunctionReturn(0); 1229 } 1230 1231 /*@C 1232 TSRKGetMultirate - Gets whether to Use the interpolation-based multirate RK method 1233 1234 Not collective 1235 1236 Input Parameter: 1237 . ts - timestepping context 1238 1239 Output Parameter: 1240 . use_multirate - PETSC_TRUE if the multirate RK method is enabled, PETSC_FALSE otherwise 1241 1242 Level: intermediate 1243 1244 .seealso: TSRKSetMultirate() 1245 @*/ 1246 PetscErrorCode TSRKGetMultirate(TS ts,PetscBool *use_multirate) 1247 { 1248 PetscErrorCode ierr; 1249 1250 PetscFunctionBegin; 1251 ierr = PetscUseMethod(ts,"TSRKGetMultirate_C",(TS,PetscBool*),(ts,use_multirate));CHKERRQ(ierr); 1252 PetscFunctionReturn(0); 1253 } 1254 1255 /*MC 1256 TSRK - ODE and DAE solver using Runge-Kutta schemes 1257 1258 The user should provide the right hand side of the equation 1259 using TSSetRHSFunction(). 1260 1261 Notes: 1262 The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1263 1264 Level: beginner 1265 1266 .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1267 TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister(), TSRKSetMultirate(), TSRKGetMultirate() 1268 1269 M*/ 1270 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1271 { 1272 TS_RK *rk; 1273 PetscErrorCode ierr; 1274 1275 PetscFunctionBegin; 1276 ierr = TSRKInitializePackage();CHKERRQ(ierr); 1277 1278 ts->ops->reset = TSReset_RK; 1279 ts->ops->destroy = TSDestroy_RK; 1280 ts->ops->view = TSView_RK; 1281 ts->ops->load = TSLoad_RK; 1282 ts->ops->setup = TSSetUp_RK; 1283 ts->ops->interpolate = TSInterpolate_RK; 1284 ts->ops->step = TSStep_RK; 1285 ts->ops->evaluatestep = TSEvaluateStep_RK; 1286 ts->ops->rollback = TSRollBack_RK; 1287 ts->ops->setfromoptions = TSSetFromOptions_RK; 1288 ts->ops->getstages = TSGetStages_RK; 1289 1290 ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 1291 ts->ops->adjointsetup = TSAdjointSetUp_RK; 1292 ts->ops->adjointstep = TSAdjointStep_RK; 1293 ts->ops->adjointreset = TSAdjointReset_RK; 1294 1295 ts->ops->forwardintegral = TSForwardCostIntegral_RK; 1296 ts->ops->forwardsetup = TSForwardSetUp_RK; 1297 ts->ops->forwardreset = TSForwardReset_RK; 1298 ts->ops->forwardstep = TSForwardStep_RK; 1299 ts->ops->forwardgetstages= TSForwardGetStages_RK; 1300 1301 ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 1302 ts->data = (void*)rk; 1303 1304 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1305 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1306 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetMultirate_C",TSRKSetMultirate_RK);CHKERRQ(ierr); 1307 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetMultirate_C",TSRKGetMultirate_RK);CHKERRQ(ierr); 1308 1309 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1310 rk->dtratio = 1; 1311 PetscFunctionReturn(0); 1312 } 1313