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