1 /* 2 Code for timestepping with Runge-Kutta method 3 4 Notes: 5 The general system is written as 6 7 F(t,U,Udot) = G(t,U) 8 9 where F represents the stiff part of the physics and G represents the non-stiff part. 10 11 */ 12 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 13 #include <petscdm.h> 14 15 static TSRKType TSRKDefault = TSRK3; 16 static PetscBool TSRKRegisterAllCalled; 17 static PetscBool TSRKPackageInitialized; 18 static PetscInt explicit_stage_time_id; 19 20 typedef struct _RKTableau *RKTableau; 21 struct _RKTableau { 22 char *name; 23 PetscInt order; /* Classical approximation order of the method */ 24 PetscInt s; /* Number of stages */ 25 PetscInt pinterp; /* Interpolation order */ 26 PetscReal *A,*b,*c; /* Tableau */ 27 PetscReal *bembed; /* Embedded formula of order one less (order-1) */ 28 PetscReal *binterp; /* Dense output formula */ 29 PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 30 }; 31 typedef struct _RKTableauLink *RKTableauLink; 32 struct _RKTableauLink { 33 struct _RKTableau tab; 34 RKTableauLink next; 35 }; 36 static RKTableauLink RKTableauList; 37 38 typedef struct { 39 RKTableau tableau; 40 Vec *Y; /* States computed during the step */ 41 Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 42 PetscScalar *work; /* Scalar work */ 43 PetscReal stage_time; 44 TSStepStatus status; 45 } TS_RK; 46 47 /*MC 48 TSRK1 - First order forward Euler scheme. 49 50 This method has one stage. 51 52 Level: advanced 53 54 .seealso: TSRK 55 M*/ 56 /*MC 57 TSRK2 - Second order RK scheme. 58 59 This method has two stages. 60 61 Level: advanced 62 63 .seealso: TSRK 64 M*/ 65 /*MC 66 TSRK3 - Third order RK scheme. 67 68 This method has three stages. 69 70 Level: advanced 71 72 .seealso: TSRK 73 M*/ 74 /*MC 75 TSRK4 - Fourth order RK scheme. 76 77 This method has four stages. 78 79 Level: advanced 80 81 .seealso: TSRK 82 M*/ 83 /*MC 84 TSRK5F - Fifth order Fehlberg RK scheme with 4th order embedded method. 85 86 This method has six stages. 87 88 Level: advanced 89 90 .seealso: TSRK 91 M*/ 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "TSRKRegisterAll" 95 /*@C 96 TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 97 98 Not Collective, but should be called by all processes which will need the schemes to be registered 99 100 Level: advanced 101 102 .keywords: TS, TSRK, register, all 103 104 .seealso: TSRKRegisterDestroy() 105 @*/ 106 PetscErrorCode TSRKRegisterAll(void) 107 { 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 112 TSRKRegisterAllCalled = PETSC_TRUE; 113 114 { 115 const PetscReal 116 A[1][1] = {{0.0}}, 117 b[1] = {1.0}; 118 ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr); 119 } 120 { 121 const PetscReal 122 A[2][2] = {{0.0,0.0}, 123 {1.0,0.0}}, 124 b[2] = {0.5,0.5}, 125 bembed[2] = {1.0,0}; 126 ierr = TSRKRegister(TSRK2,2,2,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr); 127 } 128 { 129 const PetscReal 130 A[3][3] = {{0,0,0}, 131 {2.0/3.0,0,0}, 132 {-1.0/3.0,1.0,0}}, 133 b[3] = {0.25,0.5,0.25}; 134 ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,2,b);CHKERRQ(ierr); 135 } 136 { 137 const PetscReal 138 A[4][4] = {{0,0,0,0}, 139 {0.5,0,0,0}, 140 {0,0.5,0,0}, 141 {0,0,1.0,0}}, 142 b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 143 ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,2,b);CHKERRQ(ierr); 144 } 145 { 146 const PetscReal 147 A[6][6] = {{0,0,0,0,0,0}, 148 {0.25,0,0,0,0,0}, 149 {3.0/32.0,9.0/32.0,0,0,0,0}, 150 {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 151 {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 152 {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 153 b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 154 bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 155 ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr); 156 } 157 PetscFunctionReturn(0); 158 } 159 160 #undef __FUNCT__ 161 #define __FUNCT__ "TSRKRegisterDestroy" 162 /*@C 163 TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 164 165 Not Collective 166 167 Level: advanced 168 169 .keywords: TSRK, register, destroy 170 .seealso: TSRKRegister(), TSRKRegisterAll() 171 @*/ 172 PetscErrorCode TSRKRegisterDestroy(void) 173 { 174 PetscErrorCode ierr; 175 RKTableauLink link; 176 177 PetscFunctionBegin; 178 while ((link = RKTableauList)) { 179 RKTableau t = &link->tab; 180 RKTableauList = link->next; 181 ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 182 ierr = PetscFree (t->bembed); CHKERRQ(ierr); 183 ierr = PetscFree (t->binterp); CHKERRQ(ierr); 184 ierr = PetscFree (t->name); CHKERRQ(ierr); 185 ierr = PetscFree (link); CHKERRQ(ierr); 186 } 187 TSRKRegisterAllCalled = PETSC_FALSE; 188 PetscFunctionReturn(0); 189 } 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "TSRKInitializePackage" 193 /*@C 194 TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 195 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 196 when using static libraries. 197 198 Level: developer 199 200 .keywords: TS, TSRK, initialize, package 201 .seealso: PetscInitialize() 202 @*/ 203 PetscErrorCode TSRKInitializePackage(void) 204 { 205 PetscErrorCode ierr; 206 207 PetscFunctionBegin; 208 if (TSRKPackageInitialized) PetscFunctionReturn(0); 209 TSRKPackageInitialized = PETSC_TRUE; 210 ierr = TSRKRegisterAll();CHKERRQ(ierr); 211 ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr); 212 ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 213 PetscFunctionReturn(0); 214 } 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "TSRKFinalizePackage" 218 /*@C 219 TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 220 called from PetscFinalize(). 221 222 Level: developer 223 224 .keywords: Petsc, destroy, package 225 .seealso: PetscFinalize() 226 @*/ 227 PetscErrorCode TSRKFinalizePackage(void) 228 { 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 TSRKPackageInitialized = PETSC_FALSE; 233 ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 234 PetscFunctionReturn(0); 235 } 236 237 #undef __FUNCT__ 238 #define __FUNCT__ "TSRKRegister" 239 /*@C 240 TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 241 242 Not Collective, but the same schemes should be registered on all processes on which they will be used 243 244 Input Parameters: 245 + name - identifier for method 246 . order - approximation order of method 247 . s - number of stages, this is the dimension of the matrices below 248 . A - stage coefficients (dimension s*s, row-major) 249 . b - step completion table (dimension s; NULL to use last row of A) 250 . c - abscissa (dimension s; NULL to use row sums of A) 251 . bembed - completion table for embedded method (dimension s; NULL if not available) 252 . pinterp - Order of the interpolation scheme, equal to the number of columns of binterp 253 - binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt) 254 255 Notes: 256 Several RK methods are provided, this function is only needed to create new methods. 257 258 Level: advanced 259 260 .keywords: TS, register 261 262 .seealso: TSRK 263 @*/ 264 PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 265 const PetscReal A[],const PetscReal b[],const PetscReal c[], 266 const PetscReal bembed[], 267 PetscInt pinterp,const PetscReal binterp[]) 268 { 269 PetscErrorCode ierr; 270 RKTableauLink link; 271 RKTableau t; 272 PetscInt i,j; 273 274 PetscFunctionBegin; 275 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 276 ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 277 t = &link->tab; 278 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 279 t->order = order; 280 t->s = s; 281 ierr = PetscMalloc3(s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 282 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 283 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 284 else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 285 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 286 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 287 288 if (bembed) { 289 ierr = PetscMalloc(s*sizeof(PetscReal),&t->bembed);CHKERRQ(ierr); 290 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 291 } 292 293 t->pinterp = pinterp; 294 ierr = PetscMalloc(s*pinterp*sizeof(PetscReal),&t->binterp);CHKERRQ(ierr); 295 ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 296 link->next = RKTableauList; 297 RKTableauList = link; 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "TSEvaluateStep_RK" 303 /* 304 The step completion formula is 305 306 x1 = x0 + h b^T YdotRHS 307 308 This function can be called before or after ts->vec_sol has been updated. 309 Suppose we have a completion formula (b) and an embedded formula (be) of different order. 310 We can write 311 312 x1e = x0 + h be^T YdotRHS 313 = x1 - h b^T YdotRHS + h be^T YdotRHS 314 = x1 + h (be - b)^T YdotRHS 315 316 so we can evaluate the method with different order even after the step has been optimistically completed. 317 */ 318 static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 319 { 320 TS_RK *rk = (TS_RK*)ts->data; 321 RKTableau tab = rk->tableau; 322 PetscScalar *w = rk->work; 323 PetscReal h; 324 PetscInt s = tab->s,j; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 switch (rk->status) { 329 case TS_STEP_INCOMPLETE: 330 case TS_STEP_PENDING: 331 h = ts->time_step; break; 332 case TS_STEP_COMPLETE: 333 h = ts->time_step_prev; break; 334 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 335 } 336 if (order == tab->order) { 337 if (rk->status == TS_STEP_INCOMPLETE) { 338 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 339 for (j=0; j<s; j++) w[j] = h*tab->b[j]; 340 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 341 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 342 if (done) *done = PETSC_TRUE; 343 PetscFunctionReturn(0); 344 } else if (order == tab->order-1) { 345 if (!tab->bembed) goto unavailable; 346 if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 347 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 348 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 349 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 350 } else { /* Rollback and re-complete using (be-b) */ 351 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 352 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 353 ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 354 } 355 if (done) *done = PETSC_TRUE; 356 PetscFunctionReturn(0); 357 } 358 unavailable: 359 if (done) *done = PETSC_FALSE; 360 else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 361 PetscFunctionReturn(0); 362 } 363 364 #undef __FUNCT__ 365 #define __FUNCT__ "TSStep_RK" 366 static PetscErrorCode TSStep_RK(TS ts) 367 { 368 TS_RK *rk = (TS_RK*)ts->data; 369 RKTableau tab = rk->tableau; 370 const PetscInt s = tab->s; 371 const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 372 PetscScalar *w = rk->work; 373 Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 374 TSAdapt adapt; 375 PetscInt i,j,reject,next_scheme; 376 PetscReal next_time_step; 377 PetscReal t; 378 PetscBool accept; 379 PetscErrorCode ierr; 380 381 PetscFunctionBegin; 382 383 next_time_step = ts->time_step; 384 t = ts->ptime; 385 accept = PETSC_TRUE; 386 rk->status = TS_STEP_INCOMPLETE; 387 388 389 for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 390 PetscReal h = ts->time_step; 391 ierr = TSPreStep(ts);CHKERRQ(ierr); 392 for (i=0; i<s; i++) { 393 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 394 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 395 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 396 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 397 ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 398 if (!accept) goto reject_step; 399 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 400 } 401 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 402 rk->status = TS_STEP_PENDING; 403 404 /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 405 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 406 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 407 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 408 ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 409 if (accept) { 410 /* ignore next_scheme for now */ 411 ts->ptime += ts->time_step; 412 ts->time_step = next_time_step; 413 ts->steps++; 414 rk->status = TS_STEP_COMPLETE; 415 ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 416 417 break; 418 } else { /* Roll back the current step */ 419 for (j=0; j<s; j++) w[j] = -h*b[j]; 420 ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr); 421 ts->time_step = next_time_step; 422 rk->status = TS_STEP_INCOMPLETE; 423 } 424 reject_step: continue; 425 } 426 if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "TSInterpolate_RK" 432 static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 433 { 434 TS_RK *rk = (TS_RK*)ts->data; 435 PetscInt s = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j; 436 PetscReal h; 437 PetscReal tt,t; 438 PetscScalar *b; 439 const PetscReal *B = rk->tableau->binterp; 440 PetscErrorCode ierr; 441 442 PetscFunctionBegin; 443 if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 444 switch (rk->status) { 445 case TS_STEP_INCOMPLETE: 446 case TS_STEP_PENDING: 447 h = ts->time_step; 448 t = (itime - ts->ptime)/h; 449 break; 450 case TS_STEP_COMPLETE: 451 h = ts->time_step_prev; 452 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 453 break; 454 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 455 } 456 ierr = PetscMalloc(s*sizeof(PetscScalar),&b);CHKERRQ(ierr); 457 for (i=0; i<s; i++) b[i] = 0; 458 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 459 for (i=0; i<s; i++) { 460 b[i] += h * B[i*pinterp+j] * tt; 461 } 462 } 463 ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 464 ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 465 ierr = PetscFree(b);CHKERRQ(ierr); 466 PetscFunctionReturn(0); 467 } 468 469 /*------------------------------------------------------------*/ 470 #undef __FUNCT__ 471 #define __FUNCT__ "TSReset_RK" 472 static PetscErrorCode TSReset_RK(TS ts) 473 { 474 TS_RK *rk = (TS_RK*)ts->data; 475 PetscInt s; 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 if (!rk->tableau) PetscFunctionReturn(0); 480 s = rk->tableau->s; 481 ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr); 482 ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr); 483 ierr = PetscFree(rk->work);CHKERRQ(ierr); 484 PetscFunctionReturn(0); 485 } 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "TSDestroy_RK" 489 static PetscErrorCode TSDestroy_RK(TS ts) 490 { 491 PetscErrorCode ierr; 492 493 PetscFunctionBegin; 494 ierr = TSReset_RK(ts);CHKERRQ(ierr); 495 ierr = PetscFree(ts->data);CHKERRQ(ierr); 496 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 497 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 498 PetscFunctionReturn(0); 499 } 500 501 502 #undef __FUNCT__ 503 #define __FUNCT__ "DMCoarsenHook_TSRK" 504 static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 505 { 506 PetscFunctionBegin; 507 PetscFunctionReturn(0); 508 } 509 510 #undef __FUNCT__ 511 #define __FUNCT__ "DMRestrictHook_TSRK" 512 static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 513 { 514 PetscFunctionBegin; 515 PetscFunctionReturn(0); 516 } 517 518 519 #undef __FUNCT__ 520 #define __FUNCT__ "DMSubDomainHook_TSRK" 521 static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 522 { 523 PetscFunctionBegin; 524 PetscFunctionReturn(0); 525 } 526 527 #undef __FUNCT__ 528 #define __FUNCT__ "DMSubDomainRestrictHook_TSRK" 529 static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 530 { 531 532 PetscFunctionBegin; 533 PetscFunctionReturn(0); 534 } 535 536 #undef __FUNCT__ 537 #define __FUNCT__ "TSSetUp_RK" 538 static PetscErrorCode TSSetUp_RK(TS ts) 539 { 540 TS_RK *rk = (TS_RK*)ts->data; 541 RKTableau tab; 542 PetscInt s; 543 PetscErrorCode ierr; 544 DM dm; 545 546 PetscFunctionBegin; 547 if (!rk->tableau) { 548 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 549 } 550 tab = rk->tableau; 551 s = tab->s; 552 ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr); 553 ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr); 554 ierr = PetscMalloc(s*sizeof(rk->work[0]),&rk->work);CHKERRQ(ierr); 555 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 556 if (dm) { 557 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 558 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 559 } 560 PetscFunctionReturn(0); 561 } 562 /*------------------------------------------------------------*/ 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "TSSetFromOptions_RK" 566 static PetscErrorCode TSSetFromOptions_RK(TS ts) 567 { 568 PetscErrorCode ierr; 569 char rktype[256]; 570 571 PetscFunctionBegin; 572 ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr); 573 { 574 RKTableauLink link; 575 PetscInt count,choice; 576 PetscBool flg; 577 const char **namelist; 578 ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr); 579 for (link=RKTableauList,count=0; link; link=link->next,count++) ; 580 ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 581 for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 582 ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr); 583 ierr = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr); 584 ierr = PetscFree(namelist);CHKERRQ(ierr); 585 } 586 ierr = PetscOptionsTail();CHKERRQ(ierr); 587 PetscFunctionReturn(0); 588 } 589 590 #undef __FUNCT__ 591 #define __FUNCT__ "PetscFormatRealArray" 592 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 593 { 594 PetscErrorCode ierr; 595 PetscInt i; 596 size_t left,count; 597 char *p; 598 599 PetscFunctionBegin; 600 for (i=0,p=buf,left=len; i<n; i++) { 601 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 602 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 603 left -= count; 604 p += count; 605 *p++ = ' '; 606 } 607 p[i ? 0 : -1] = 0; 608 PetscFunctionReturn(0); 609 } 610 611 #undef __FUNCT__ 612 #define __FUNCT__ "TSView_RK" 613 static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 614 { 615 TS_RK *rk = (TS_RK*)ts->data; 616 RKTableau tab = rk->tableau; 617 PetscBool iascii; 618 PetscErrorCode ierr; 619 TSAdapt adapt; 620 621 PetscFunctionBegin; 622 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 623 if (iascii) { 624 TSRKType rktype; 625 char buf[512]; 626 ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 627 ierr = PetscViewerASCIIPrintf(viewer," RK %s\n",rktype);CHKERRQ(ierr); 628 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 629 ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 630 } 631 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 632 ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr); 633 PetscFunctionReturn(0); 634 } 635 636 #undef __FUNCT__ 637 #define __FUNCT__ "TSLoad_RK" 638 static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 639 { 640 PetscErrorCode ierr; 641 TSAdapt tsadapt; 642 643 PetscFunctionBegin; 644 ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr); 645 ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr); 646 PetscFunctionReturn(0); 647 } 648 649 #undef __FUNCT__ 650 #define __FUNCT__ "TSRKSetType" 651 /*@C 652 TSRKSetType - Set the type of RK scheme 653 654 Logically collective 655 656 Input Parameter: 657 + ts - timestepping context 658 - rktype - type of RK-scheme 659 660 Level: intermediate 661 662 .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5 663 @*/ 664 PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 665 { 666 PetscErrorCode ierr; 667 668 PetscFunctionBegin; 669 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 670 ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "TSRKGetType" 676 /*@C 677 TSRKGetType - Get the type of RK scheme 678 679 Logically collective 680 681 Input Parameter: 682 . ts - timestepping context 683 684 Output Parameter: 685 . rktype - type of RK-scheme 686 687 Level: intermediate 688 689 .seealso: TSRKGetType() 690 @*/ 691 PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 692 { 693 PetscErrorCode ierr; 694 695 PetscFunctionBegin; 696 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 697 ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 #undef __FUNCT__ 702 #define __FUNCT__ "TSRKGetType_RK" 703 PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 704 { 705 TS_RK *rk = (TS_RK*)ts->data; 706 PetscErrorCode ierr; 707 708 PetscFunctionBegin; 709 if (!rk->tableau) { 710 ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 711 } 712 *rktype = rk->tableau->name; 713 PetscFunctionReturn(0); 714 } 715 #undef __FUNCT__ 716 #define __FUNCT__ "TSRKSetType_RK" 717 PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 718 { 719 TS_RK *rk = (TS_RK*)ts->data; 720 PetscErrorCode ierr; 721 PetscBool match; 722 RKTableauLink link; 723 724 PetscFunctionBegin; 725 if (rk->tableau) { 726 ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 727 if (match) PetscFunctionReturn(0); 728 } 729 for (link = RKTableauList; link; link=link->next) { 730 ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 731 if (match) { 732 ierr = TSReset_RK(ts);CHKERRQ(ierr); 733 rk->tableau = &link->tab; 734 PetscFunctionReturn(0); 735 } 736 } 737 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 738 PetscFunctionReturn(0); 739 } 740 741 /* ------------------------------------------------------------ */ 742 /*MC 743 TSRK - ODE and DAE solver using Runge-Kutta schemes 744 745 The user should provide the right hand side of the equation 746 using TSSetRHSFunction(). 747 748 Notes: 749 The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type 750 751 Level: beginner 752 753 .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 754 TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 755 756 M*/ 757 #undef __FUNCT__ 758 #define __FUNCT__ "TSCreate_RK" 759 PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 760 { 761 TS_RK *th; 762 PetscErrorCode ierr; 763 764 PetscFunctionBegin; 765 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 766 ierr = TSRKInitializePackage();CHKERRQ(ierr); 767 #endif 768 769 ts->ops->reset = TSReset_RK; 770 ts->ops->destroy = TSDestroy_RK; 771 ts->ops->view = TSView_RK; 772 ts->ops->load = TSLoad_RK; 773 ts->ops->setup = TSSetUp_RK; 774 ts->ops->step = TSStep_RK; 775 ts->ops->interpolate = TSInterpolate_RK; 776 ts->ops->evaluatestep = TSEvaluateStep_RK; 777 ts->ops->setfromoptions = TSSetFromOptions_RK; 778 779 ierr = PetscNewLog(ts,TS_RK,&th);CHKERRQ(ierr); 780 ts->data = (void*)th; 781 782 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 783 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 784 PetscFunctionReturn(0); 785 } 786