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