1 /* 2 Code for timestepping with Rosenbrock W methods 3 4 Notes: 5 The general system is written as 6 7 G(t,X,Xdot) = F(t,X) 8 9 where G represents the stiff part of the physics and F represents the non-stiff part. 10 This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian. 11 12 */ 13 #include <private/tsimpl.h> /*I "petscts.h" I*/ 14 15 #include <../src/mat/blockinvert.h> 16 17 static const TSRosWType TSRosWDefault = TSROSWRA34PW2; 18 static PetscBool TSRosWRegisterAllCalled; 19 static PetscBool TSRosWPackageInitialized; 20 21 typedef struct _RosWTableau *RosWTableau; 22 struct _RosWTableau { 23 char *name; 24 PetscInt order; /* Classical approximation order of the method */ 25 PetscInt s; /* Number of stages */ 26 PetscReal *A; /* Propagation table, strictly lower triangular */ 27 PetscReal *Gamma; /* Stage table, lower triangular with nonzero diagonal */ 28 PetscBool *GammaZeroDiag; /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */ 29 PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/ 30 PetscReal *b; /* Step completion table */ 31 PetscReal *bembed; /* Step completion table for embedded method of order one less */ 32 PetscReal *ASum; /* Row sum of A */ 33 PetscReal *GammaSum; /* Row sum of Gamma, only needed for non-autonomous systems */ 34 PetscReal *At; /* Propagation table in transformed variables */ 35 PetscReal *bt; /* Step completion table in transformed variables */ 36 PetscReal *bembedt; /* Step completion table of order one less in transformed variables */ 37 PetscReal *GammaInv; /* Inverse of Gamma, used for transformed variables */ 38 PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 39 }; 40 typedef struct _RosWTableauLink *RosWTableauLink; 41 struct _RosWTableauLink { 42 struct _RosWTableau tab; 43 RosWTableauLink next; 44 }; 45 static RosWTableauLink RosWTableauList; 46 47 typedef struct { 48 RosWTableau tableau; 49 Vec *Y; /* States computed during the step, used to complete the step */ 50 Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 51 Vec Ystage; /* Work vector for the state value at each stage */ 52 Vec Zdot; /* Ydot = Zdot + shift*Y */ 53 Vec Zstage; /* Y = Zstage + Y */ 54 PetscScalar *work; /* Scalar work space of length number of stages, used to prepare VecMAXPY() */ 55 PetscReal shift; 56 PetscReal stage_time; 57 PetscReal stage_explicit; /* Flag indicates that the current stage is explicit */ 58 PetscBool recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */ 59 TSStepStatus status; 60 } TS_RosW; 61 62 /*MC 63 TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method). 64 65 Only an approximate Jacobian is needed. 66 67 Level: intermediate 68 69 .seealso: TSROSW 70 M*/ 71 72 /*MC 73 TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method). 74 75 Only an approximate Jacobian is needed. 76 77 Level: intermediate 78 79 .seealso: TSROSW 80 M*/ 81 82 /*MC 83 TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme. 84 85 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P. 86 87 Level: intermediate 88 89 .seealso: TSROSW 90 M*/ 91 92 /*MC 93 TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme. 94 95 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M. 96 97 Level: intermediate 98 99 .seealso: TSROSW 100 M*/ 101 102 /*MC 103 TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1. 104 105 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 106 107 This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73. 108 109 References: 110 Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 111 112 Level: intermediate 113 114 .seealso: TSROSW 115 M*/ 116 117 /*MC 118 TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1. 119 120 Only an approximate Jacobian is needed. By default, it is only recomputed once per step. 121 122 This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48. 123 124 References: 125 Rang and Angermann, New Rosenbrock-W methods of order 3 for partial differential algebraic equations of index 1, 2005. 126 127 Level: intermediate 128 129 .seealso: TSROSW 130 M*/ 131 132 /*MC 133 TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme 134 135 By default, the Jacobian is only recomputed once per step. 136 137 Both the third order and embedded second order methods are stiffly accurate and L-stable. 138 139 References: 140 Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 141 142 Level: intermediate 143 144 .seealso: TSROSW, TSROSWSANDU3 145 M*/ 146 147 /*MC 148 TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme 149 150 By default, the Jacobian is only recomputed once per step. 151 152 The third order method is L-stable, but not stiffly accurate. 153 The second order embedded method is strongly A-stable with R(infty) = 0.5. 154 The internal stages are L-stable. 155 This method is called ROS3 in the paper. 156 157 References: 158 Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997. 159 160 Level: intermediate 161 162 .seealso: TSROSW, TSROSWRODAS3 163 M*/ 164 165 /*MC 166 TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages 167 168 By default, the Jacobian is only recomputed once per step. 169 170 A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3) 171 172 References: 173 Emil Constantinescu 174 175 Level: intermediate 176 177 .seealso: TSROSW, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, SSP 178 M*/ 179 180 /*MC 181 TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 182 183 By default, the Jacobian is only recomputed once per step. 184 185 L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 186 187 References: 188 Emil Constantinescu 189 190 Level: intermediate 191 192 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLLSSP3P4S2C, TSSSP 193 M*/ 194 195 /*MC 196 TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, three stages 197 198 By default, the Jacobian is only recomputed once per step. 199 200 L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2) 201 202 References: 203 Emil Constantinescu 204 205 Level: intermediate 206 207 .seealso: TSROSW, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSSSP 208 M*/ 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "TSRosWRegisterAll" 212 /*@C 213 TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW 214 215 Not Collective, but should be called by all processes which will need the schemes to be registered 216 217 Level: advanced 218 219 .keywords: TS, TSRosW, register, all 220 221 .seealso: TSRosWRegisterDestroy() 222 @*/ 223 PetscErrorCode TSRosWRegisterAll(void) 224 { 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 229 TSRosWRegisterAllCalled = PETSC_TRUE; 230 231 { 232 const PetscReal 233 A = 0, 234 Gamma = 1, 235 b = 1; 236 ierr = TSRosWRegister(TSROSWTHETA1,1,1,&A,&Gamma,&b,PETSC_NULL);CHKERRQ(ierr); 237 } 238 239 { 240 const PetscReal 241 A= 0, 242 Gamma = 0.5, 243 b = 1; 244 ierr = TSRosWRegister(TSROSWTHETA2,2,1,&A,&Gamma,&b,PETSC_NULL);CHKERRQ(ierr); 245 } 246 247 { 248 const PetscReal g = 1. + 1./PetscSqrtReal(2.0); 249 const PetscReal 250 A[2][2] = {{0,0}, {1.,0}}, 251 Gamma[2][2] = {{g,0}, {-2.*g,g}}, 252 b[2] = {0.5,0.5}, 253 b1[2] = {1.0,0.0}; 254 ierr = TSRosWRegister(TSROSW2P,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 255 } 256 { 257 const PetscReal g = 1. - 1./PetscSqrtReal(2.0); 258 const PetscReal 259 A[2][2] = {{0,0}, {1.,0}}, 260 Gamma[2][2] = {{g,0}, {-2.*g,g}}, 261 b[2] = {0.5,0.5}, 262 b1[2] = {1.0,0.0}; 263 ierr = TSRosWRegister(TSROSW2M,2,2,&A[0][0],&Gamma[0][0],b,b1);CHKERRQ(ierr); 264 } 265 { 266 const PetscReal g = 7.8867513459481287e-01; 267 const PetscReal 268 A[3][3] = {{0,0,0}, 269 {1.5773502691896257e+00,0,0}, 270 {0.5,0,0}}, 271 Gamma[3][3] = {{g,0,0}, 272 {-1.5773502691896257e+00,g,0}, 273 {-6.7075317547305480e-01,1.7075317547305482e-01,g}}, 274 b[3] = {1.0566243270259355e-01,4.9038105676657971e-02,8.4529946162074843e-01}, 275 b2[3] = {-1.7863279495408180e-01,1./3.,8.4529946162074843e-01}; 276 ierr = TSRosWRegister(TSROSWRA3PW,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 277 } 278 { 279 const PetscReal g = 4.3586652150845900e-01; 280 const PetscReal 281 A[4][4] = {{0,0,0,0}, 282 {8.7173304301691801e-01,0,0,0}, 283 {8.4457060015369423e-01,-1.1299064236484185e-01,0,0}, 284 {0,0,1.,0}}, 285 Gamma[4][4] = {{g,0,0,0}, 286 {-8.7173304301691801e-01,g,0,0}, 287 {-9.0338057013044082e-01,5.4180672388095326e-02,g,0}, 288 {2.4212380706095346e-01,-1.2232505839045147e+00,5.4526025533510214e-01,g}}, 289 b[4] = {2.4212380706095346e-01,-1.2232505839045147e+00,1.5452602553351020e+00,4.3586652150845900e-01}, 290 b2[4] = {3.7810903145819369e-01,-9.6042292212423178e-02,5.0000000000000000e-01,2.1793326075422950e-01}; 291 ierr = TSRosWRegister(TSROSWRA34PW2,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 292 } 293 { 294 const PetscReal g = 0.5; 295 const PetscReal 296 A[4][4] = {{0,0,0,0}, 297 {0,0,0,0}, 298 {1.,0,0,0}, 299 {0.75,-0.25,0.5,0}}, 300 Gamma[4][4] = {{g,0,0,0}, 301 {1.,g,0,0}, 302 {-0.25,-0.25,g,0}, 303 {1./12,1./12,-2./3,g}}, 304 b[4] = {5./6,-1./6,-1./6,0.5}, 305 b2[4] = {0.75,-0.25,0.5,0}; 306 ierr = TSRosWRegister(TSROSWRODAS3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 307 } 308 { 309 const PetscReal g = 0.43586652150845899941601945119356; 310 const PetscReal 311 A[3][3] = {{0,0,0}, 312 {g,0,0}, 313 {g,0,0}}, 314 Gamma[3][3] = {{g,0,0}, 315 {-0.19294655696029095575009695436041,g,0}, 316 {0,1.74927148125794685173529749738960,g}}, 317 b[3] = {-0.75457412385404315829818998646589,1.94100407061964420292840123379419,-0.18642994676560104463021124732829}, 318 b2[3] = {-1.53358745784149585370766523913002,2.81745131148625772213931745457622,-0.28386385364476186843165221544619}; 319 ierr = TSRosWRegister(TSROSWSANDU3,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 320 } 321 { 322 const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0; 323 const PetscReal 324 A[3][3] = {{0,0,0}, 325 {1,0,0}, 326 {0.25,0.25,0}}, 327 Gamma[3][3] = {{0,0,0}, 328 {(-3.0-s3)/6.0,g,0}, 329 {(-3.0-s3)/24.0,(-3.0-s3)/8.0,g}}, 330 b[3] = {1./6.,1./6.,2./3.}, 331 b2[3] = {1./4.,1./4.,1./2.}; 332 ierr = TSRosWRegister(TSROSWASSP3P3S1C,3,3,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 333 } 334 335 { 336 const PetscReal 337 A[4][4] = {{0,0,0,0}, 338 {1./2.,0,0,0}, 339 {1./2.,1./2.,0,0}, 340 {1./6.,1./6.,1./6.,0}}, 341 Gamma[4][4] = {{1./2.,0,0,0}, 342 {0.0,1./4.,0,0}, 343 {-2.,-2./3.,2./3.,0}, 344 {1./2.,5./36.,-2./9,0}}, 345 b[4] = {1./6.,1./6.,1./6.,1./2.}, 346 b2[4] = {1./8.,3./4.,1./8.,0}; 347 ierr = TSRosWRegister(TSROSWLASSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 348 } 349 350 { 351 const PetscReal 352 A[4][4] = {{0,0,0,0}, 353 {1./2.,0,0,0}, 354 {1./2.,1./2.,0,0}, 355 {1./6.,1./6.,1./6.,0}}, 356 Gamma[4][4] = {{1./2.,0,0,0}, 357 {0.0,3./4.,0,0}, 358 {-2./3.,-23./9.,2./9.,0}, 359 {1./18.,65./108.,-2./27,0}}, 360 b[4] = {1./6.,1./6.,1./6.,1./2.}, 361 b2[4] = {3./16.,10./16.,3./16.,0}; 362 ierr = TSRosWRegister(TSROSWLLSSP3P4S2C,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 363 } 364 365 { 366 PetscReal A[4][4],Gamma[4][4],b[4],b2[4]; 367 368 Gamma[0][0]=0.4358665215084589994160194475295062513822671686978816; 369 Gamma[0][1]=0; Gamma[0][2]=0; Gamma[0][3]=0; 370 Gamma[1][0]=-1.997527830934941248426324674704153457289527280554476; 371 Gamma[1][1]=0.4358665215084589994160194475295062513822671686978816; 372 Gamma[1][2]=0; Gamma[1][3]=0; 373 Gamma[2][0]=-1.007948511795029620852002345345404191008352770119903; 374 Gamma[2][1]=-0.004648958462629345562774289390054679806993396798458131; 375 Gamma[2][2]=0.4358665215084589994160194475295062513822671686978816; 376 Gamma[2][3]=0; 377 Gamma[3][0]=-0.6685429734233467180451604600279552604364311322650783; 378 Gamma[3][1]=0.6056625986449338476089525334450053439525178740492984; 379 Gamma[3][2]=-0.9717899277217721234705114616271378792182450260943198; 380 Gamma[3][3]=0; 381 382 A[0][0]=0; A[0][1]=0; A[0][2]=0; A[0][3]=0; 383 A[1][0]=0.8717330430169179988320388950590125027645343373957631; 384 A[1][1]=0; A[1][2]=0; A[1][3]=0; 385 A[2][0]=0.5275890119763004115618079766722914408876108660811028; 386 A[2][1]=0.07241098802369958843819203208518599088698057726988732; 387 A[2][2]=0; A[2][3]=0; 388 A[3][0]=0.3990960076760701320627260685975778145384666450351314; 389 A[3][1]=-0.4375576546135194437228463747348862825846903771419953; 390 A[3][2]=1.038461646937449311660120300601880176655352737312713; 391 A[3][3]=0; 392 393 b[0]=0.1876410243467238251612921333138006734899663569186926; 394 b[1]=-0.5952974735769549480478230473706443582188442040780541; 395 b[2]=0.9717899277217721234705114616271378792182450260943198; 396 b[3]=0.4358665215084589994160194475295062513822671686978816; 397 398 b2[0]=0.2147402862233891404862383521089097657790734483804460; 399 b2[1]=-0.4851622638849390928209050538171743017757490232519684; 400 b2[2]=0.8687250025203875511662123688667549217531982787600080; 401 b2[3]=0.4016969751411624011684543450940068201770721128357014; 402 403 ierr = TSRosWRegister(TSROSWARK3,3,4,&A[0][0],&Gamma[0][0],b,b2);CHKERRQ(ierr); 404 } 405 406 PetscFunctionReturn(0); 407 } 408 409 #undef __FUNCT__ 410 #define __FUNCT__ "TSRosWRegisterDestroy" 411 /*@C 412 TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 413 414 Not Collective 415 416 Level: advanced 417 418 .keywords: TSRosW, register, destroy 419 .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 420 @*/ 421 PetscErrorCode TSRosWRegisterDestroy(void) 422 { 423 PetscErrorCode ierr; 424 RosWTableauLink link; 425 426 PetscFunctionBegin; 427 while ((link = RosWTableauList)) { 428 RosWTableau t = &link->tab; 429 RosWTableauList = link->next; 430 ierr = PetscFree5(t->A,t->Gamma,t->b,t->ASum,t->GammaSum);CHKERRQ(ierr); 431 ierr = PetscFree5(t->At,t->bt,t->GammaInv,t->GammaZeroDiag,t->GammaExplicitCorr);CHKERRQ(ierr); 432 ierr = PetscFree2(t->bembed,t->bembedt);CHKERRQ(ierr); 433 ierr = PetscFree(t->name);CHKERRQ(ierr); 434 ierr = PetscFree(link);CHKERRQ(ierr); 435 } 436 TSRosWRegisterAllCalled = PETSC_FALSE; 437 PetscFunctionReturn(0); 438 } 439 440 #undef __FUNCT__ 441 #define __FUNCT__ "TSRosWInitializePackage" 442 /*@C 443 TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 444 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 445 when using static libraries. 446 447 Input Parameter: 448 path - The dynamic library path, or PETSC_NULL 449 450 Level: developer 451 452 .keywords: TS, TSRosW, initialize, package 453 .seealso: PetscInitialize() 454 @*/ 455 PetscErrorCode TSRosWInitializePackage(const char path[]) 456 { 457 PetscErrorCode ierr; 458 459 PetscFunctionBegin; 460 if (TSRosWPackageInitialized) PetscFunctionReturn(0); 461 TSRosWPackageInitialized = PETSC_TRUE; 462 ierr = TSRosWRegisterAll();CHKERRQ(ierr); 463 ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "TSRosWFinalizePackage" 469 /*@C 470 TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 471 called from PetscFinalize(). 472 473 Level: developer 474 475 .keywords: Petsc, destroy, package 476 .seealso: PetscFinalize() 477 @*/ 478 PetscErrorCode TSRosWFinalizePackage(void) 479 { 480 PetscErrorCode ierr; 481 482 PetscFunctionBegin; 483 TSRosWPackageInitialized = PETSC_FALSE; 484 ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "TSRosWRegister" 490 /*@C 491 TSRosWRegister - register a Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 492 493 Not Collective, but the same schemes should be registered on all processes on which they will be used 494 495 Input Parameters: 496 + name - identifier for method 497 . order - approximation order of method 498 . s - number of stages, this is the dimension of the matrices below 499 . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular 500 . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal 501 . b - Step completion table (dimension s) 502 - bembed - Step completion table for a scheme of order one less (dimension s, PETSC_NULL if no embedded scheme is available) 503 504 Notes: 505 Several Rosenbrock W methods are provided, this function is only needed to create new methods. 506 507 Level: advanced 508 509 .keywords: TS, register 510 511 .seealso: TSRosW 512 @*/ 513 PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 514 const PetscReal A[],const PetscReal Gamma[],const PetscReal b[],const PetscReal bembed[]) 515 { 516 PetscErrorCode ierr; 517 RosWTableauLink link; 518 RosWTableau t; 519 PetscInt i,j,k; 520 PetscScalar *GammaInv; 521 522 PetscFunctionBegin; 523 PetscValidCharPointer(name,1); 524 PetscValidPointer(A,4); 525 PetscValidPointer(Gamma,5); 526 PetscValidPointer(b,6); 527 if (bembed) PetscValidPointer(bembed,7); 528 529 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 530 ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 531 t = &link->tab; 532 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 533 t->order = order; 534 t->s = s; 535 ierr = PetscMalloc5(s*s,PetscReal,&t->A,s*s,PetscReal,&t->Gamma,s,PetscReal,&t->b,s,PetscReal,&t->ASum,s,PetscReal,&t->GammaSum);CHKERRQ(ierr); 536 ierr = PetscMalloc5(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s*s,PetscReal,&t->GammaInv,s,PetscBool,&t->GammaZeroDiag,s*s,PetscReal,&t->GammaExplicitCorr);CHKERRQ(ierr); 537 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 538 ierr = PetscMemcpy(t->Gamma,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 539 ierr = PetscMemcpy(t->GammaExplicitCorr,Gamma,s*s*sizeof(Gamma[0]));CHKERRQ(ierr); 540 ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 541 if (bembed) { 542 ierr = PetscMalloc2(s,PetscReal,&t->bembed,s,PetscReal,&t->bembedt);CHKERRQ(ierr); 543 ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 544 } 545 for (i=0; i<s; i++) { 546 t->ASum[i] = 0; 547 t->GammaSum[i] = 0; 548 for (j=0; j<s; j++) { 549 t->ASum[i] += A[i*s+j]; 550 t->GammaSum[i] += Gamma[i*s+j]; 551 } 552 } 553 ierr = PetscMalloc(s*s*sizeof(PetscScalar),&GammaInv);CHKERRQ(ierr); /* Need to use Scalar for inverse, then convert back to Real */ 554 for (i=0; i<s*s; i++) GammaInv[i] = Gamma[i]; 555 for (i=0; i<s; i++) { 556 if (Gamma[i*s+i] == 0.0) { 557 GammaInv[i*s+i] = 1.0; 558 t->GammaZeroDiag[i] = PETSC_TRUE; 559 } else { 560 t->GammaZeroDiag[i] = PETSC_FALSE; 561 } 562 } 563 564 switch (s) { 565 case 1: GammaInv[0] = 1./GammaInv[0]; break; 566 case 2: ierr = Kernel_A_gets_inverse_A_2(GammaInv,0);CHKERRQ(ierr); break; 567 case 3: ierr = Kernel_A_gets_inverse_A_3(GammaInv,0);CHKERRQ(ierr); break; 568 case 4: ierr = Kernel_A_gets_inverse_A_4(GammaInv,0);CHKERRQ(ierr); break; 569 case 5: { 570 PetscInt ipvt5[5]; 571 MatScalar work5[5*5]; 572 ierr = Kernel_A_gets_inverse_A_5(GammaInv,ipvt5,work5,0);CHKERRQ(ierr); break; 573 } 574 case 6: ierr = Kernel_A_gets_inverse_A_6(GammaInv,0);CHKERRQ(ierr); break; 575 case 7: ierr = Kernel_A_gets_inverse_A_7(GammaInv,0);CHKERRQ(ierr); break; 576 default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for %D stages",s); 577 } 578 for (i=0; i<s*s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]); 579 ierr = PetscFree(GammaInv);CHKERRQ(ierr); 580 581 for (i=0; i<s; i++) { 582 for (k=0; k<i+1; k++) { 583 t->GammaExplicitCorr[i*s+k]=(t->GammaExplicitCorr[i*s+k])*(t->GammaInv[k*s+k]); 584 for (j=k+1; j<i+1; j++) { 585 t->GammaExplicitCorr[i*s+k]+=(t->GammaExplicitCorr[i*s+j])*(t->GammaInv[j*s+k]); 586 } 587 } 588 } 589 590 for (i=0; i<s; i++) { 591 for (j=0; j<s; j++) { 592 t->At[i*s+j] = 0; 593 for (k=0; k<s; k++) { 594 t->At[i*s+j] += t->A[i*s+k] * t->GammaInv[k*s+j]; 595 } 596 } 597 t->bt[i] = 0; 598 for (j=0; j<s; j++) { 599 t->bt[i] += t->b[j] * t->GammaInv[j*s+i]; 600 } 601 if (bembed) { 602 t->bembedt[i] = 0; 603 for (j=0; j<s; j++) { 604 t->bembedt[i] += t->bembed[j] * t->GammaInv[j*s+i]; 605 } 606 } 607 } 608 t->ccfl = 1.0; /* Fix this */ 609 610 link->next = RosWTableauList; 611 RosWTableauList = link; 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "TSEvaluateStep_RosW" 617 /* 618 The step completion formula is 619 620 x1 = x0 + b^T Y 621 622 where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been 623 updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write 624 625 x1e = x0 + be^T Y 626 = x1 - b^T Y + be^T Y 627 = x1 + (be - b)^T Y 628 629 so we can evaluate the method of different order even after the step has been optimistically completed. 630 */ 631 static PetscErrorCode TSEvaluateStep_RosW(TS ts,PetscInt order,Vec X,PetscBool *done) 632 { 633 TS_RosW *ros = (TS_RosW*)ts->data; 634 RosWTableau tab = ros->tableau; 635 PetscScalar *w = ros->work; 636 PetscInt i; 637 PetscErrorCode ierr; 638 639 PetscFunctionBegin; 640 if (order == tab->order) { 641 if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */ 642 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 643 for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 644 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 645 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 646 if (done) *done = PETSC_TRUE; 647 PetscFunctionReturn(0); 648 } else if (order == tab->order-1) { 649 if (!tab->bembedt) goto unavailable; 650 if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */ 651 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 652 for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 653 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 654 } else { /* Use rollback-and-recomplete formula (bembedt - bt) */ 655 for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 656 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 657 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 658 } 659 if (done) *done = PETSC_TRUE; 660 PetscFunctionReturn(0); 661 } 662 unavailable: 663 if (done) *done = PETSC_FALSE; 664 else SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_SUP,"Rosenbrock-W '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 665 PetscFunctionReturn(0); 666 } 667 668 #undef __FUNCT__ 669 #define __FUNCT__ "TSStep_RosW" 670 static PetscErrorCode TSStep_RosW(TS ts) 671 { 672 TS_RosW *ros = (TS_RosW*)ts->data; 673 RosWTableau tab = ros->tableau; 674 const PetscInt s = tab->s; 675 const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 676 const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr; 677 const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 678 PetscScalar *w = ros->work; 679 Vec *Y = ros->Y,Ydot = ros->Ydot,Zdot = ros->Zdot,Zstage = ros->Zstage; 680 SNES snes; 681 TSAdapt adapt; 682 PetscInt i,j,its,lits,reject,next_scheme; 683 PetscReal next_time_step; 684 PetscBool accept; 685 PetscErrorCode ierr; 686 MatStructure str; 687 688 PetscFunctionBegin; 689 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 690 next_time_step = ts->time_step; 691 accept = PETSC_TRUE; 692 ros->status = TS_STEP_INCOMPLETE; 693 694 for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 695 const PetscReal h = ts->time_step; 696 for (i=0; i<s; i++) { 697 ros->stage_time = ts->ptime + h*ASum[i]; 698 if (GammaZeroDiag[i]) { 699 ros->stage_explicit = PETSC_TRUE; 700 ros->shift = 1./h; 701 } else { 702 ros->stage_explicit = PETSC_FALSE; 703 ros->shift = 1./(h*Gamma[i*s+i]); 704 } 705 706 ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 707 for (j=0; j<i; j++) w[j] = At[i*s+j]; 708 ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 709 710 for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 711 ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 712 ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 713 714 /* Initial guess taken from last stage */ 715 ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 716 717 if (!ros->stage_explicit) { 718 if (!ros->recompute_jacobian && !i) { 719 ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 720 } 721 ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 722 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 723 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 724 ts->nonlinear_its += its; ts->linear_its += lits; 725 } else { 726 Mat J,Jp; 727 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */ 728 ierr = TSComputeIFunction(ts,ros->stage_time,Zstage,Ydot,Y[i],PETSC_FALSE);CHKERRQ(ierr); 729 ierr = VecScale(Y[i],-1.0); 730 ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); /*Y[i]=F(Zstage)-Zdot[=GammaInv*Y]*/ 731 732 ierr = VecZeroEntries(Zstage);CHKERRQ(ierr); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */ 733 for (j=0; j<i; j++) w[j] = GammaExplicitCorr[i*s+j]; 734 ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 735 /*Y[i] += Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */ 736 str = SAME_NONZERO_PATTERN; 737 ierr = TSGetIJacobian(ts,&J,&Jp,PETSC_NULL,PETSC_NULL); 738 ierr = TSComputeIJacobian(ts,ros->stage_time,ts->vec_sol,Ydot,0,&J,&Jp,&str,PETSC_FALSE);CHKERRQ(ierr); 739 ierr = MatMult(J,Zstage,Zdot); 740 741 ierr = VecAXPY(Y[i],-1.0,Zdot);CHKERRQ(ierr); 742 ierr = VecScale(Y[i],h); 743 ts->linear_its += 1; 744 } 745 } 746 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 747 ros->status = TS_STEP_PENDING; 748 749 /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 750 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 751 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 752 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 753 ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 754 if (accept) { 755 /* ignore next_scheme for now */ 756 ts->ptime += ts->time_step; 757 ts->time_step = next_time_step; 758 ts->steps++; 759 ros->status = TS_STEP_COMPLETE; 760 break; 761 } else { /* Roll back the current step */ 762 for (i=0; i<s; i++) w[i] = -tab->bt[i]; 763 ierr = VecMAXPY(ts->vec_sol,s,w,Y);CHKERRQ(ierr); 764 ts->time_step = next_time_step; 765 ros->status = TS_STEP_INCOMPLETE; 766 } 767 } 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "TSInterpolate_RosW" 773 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 774 { 775 TS_RosW *ros = (TS_RosW*)ts->data; 776 777 PetscFunctionBegin; 778 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 779 PetscFunctionReturn(0); 780 } 781 782 /*------------------------------------------------------------*/ 783 #undef __FUNCT__ 784 #define __FUNCT__ "TSReset_RosW" 785 static PetscErrorCode TSReset_RosW(TS ts) 786 { 787 TS_RosW *ros = (TS_RosW*)ts->data; 788 PetscInt s; 789 PetscErrorCode ierr; 790 791 PetscFunctionBegin; 792 if (!ros->tableau) PetscFunctionReturn(0); 793 s = ros->tableau->s; 794 ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 795 ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 796 ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 797 ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 798 ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 799 ierr = PetscFree(ros->work);CHKERRQ(ierr); 800 PetscFunctionReturn(0); 801 } 802 803 #undef __FUNCT__ 804 #define __FUNCT__ "TSDestroy_RosW" 805 static PetscErrorCode TSDestroy_RosW(TS ts) 806 { 807 PetscErrorCode ierr; 808 809 PetscFunctionBegin; 810 ierr = TSReset_RosW(ts);CHKERRQ(ierr); 811 ierr = PetscFree(ts->data);CHKERRQ(ierr); 812 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 813 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 814 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 815 PetscFunctionReturn(0); 816 } 817 818 /* 819 This defines the nonlinear equation that is to be solved with SNES 820 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 821 */ 822 #undef __FUNCT__ 823 #define __FUNCT__ "SNESTSFormFunction_RosW" 824 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 825 { 826 TS_RosW *ros = (TS_RosW*)ts->data; 827 PetscErrorCode ierr; 828 829 PetscFunctionBegin; 830 ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 831 ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 832 ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 833 PetscFunctionReturn(0); 834 } 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "SNESTSFormJacobian_RosW" 838 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 839 { 840 TS_RosW *ros = (TS_RosW*)ts->data; 841 PetscErrorCode ierr; 842 843 PetscFunctionBegin; 844 /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 845 ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 846 PetscFunctionReturn(0); 847 } 848 849 #undef __FUNCT__ 850 #define __FUNCT__ "TSSetUp_RosW" 851 static PetscErrorCode TSSetUp_RosW(TS ts) 852 { 853 TS_RosW *ros = (TS_RosW*)ts->data; 854 RosWTableau tab = ros->tableau; 855 PetscInt s = tab->s; 856 PetscErrorCode ierr; 857 858 PetscFunctionBegin; 859 if (!ros->tableau) { 860 ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 861 } 862 ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 863 ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 864 ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 865 ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 866 ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 867 ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 868 PetscFunctionReturn(0); 869 } 870 /*------------------------------------------------------------*/ 871 872 #undef __FUNCT__ 873 #define __FUNCT__ "TSSetFromOptions_RosW" 874 static PetscErrorCode TSSetFromOptions_RosW(TS ts) 875 { 876 TS_RosW *ros = (TS_RosW*)ts->data; 877 PetscErrorCode ierr; 878 char rostype[256]; 879 880 PetscFunctionBegin; 881 ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 882 { 883 RosWTableauLink link; 884 PetscInt count,choice; 885 PetscBool flg; 886 const char **namelist; 887 SNES snes; 888 889 ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 890 for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 891 ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 892 for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 893 ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 894 ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 895 ierr = PetscFree(namelist);CHKERRQ(ierr); 896 897 ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 898 899 /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 900 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 901 if (!((PetscObject)snes)->type_name) { 902 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 903 } 904 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 905 } 906 ierr = PetscOptionsTail();CHKERRQ(ierr); 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "PetscFormatRealArray" 912 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 913 { 914 PetscErrorCode ierr; 915 PetscInt i; 916 size_t left,count; 917 char *p; 918 919 PetscFunctionBegin; 920 for (i=0,p=buf,left=len; i<n; i++) { 921 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 922 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 923 left -= count; 924 p += count; 925 *p++ = ' '; 926 } 927 p[i ? 0 : -1] = 0; 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "TSView_RosW" 933 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 934 { 935 TS_RosW *ros = (TS_RosW*)ts->data; 936 RosWTableau tab = ros->tableau; 937 PetscBool iascii; 938 PetscErrorCode ierr; 939 940 PetscFunctionBegin; 941 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 942 if (iascii) { 943 const TSRosWType rostype; 944 PetscInt i; 945 PetscReal abscissa[512]; 946 char buf[512]; 947 ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 948 ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 949 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 950 ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 951 for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 952 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 953 ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 954 } 955 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 956 PetscFunctionReturn(0); 957 } 958 959 #undef __FUNCT__ 960 #define __FUNCT__ "TSRosWSetType" 961 /*@C 962 TSRosWSetType - Set the type of Rosenbrock-W scheme 963 964 Logically collective 965 966 Input Parameter: 967 + ts - timestepping context 968 - rostype - type of Rosenbrock-W scheme 969 970 Level: beginner 971 972 .seealso: TSRosWGetType(), TSROSW, TSROSW2M, TSROSW2P, TSROSWRA3PW, TSROSWRA34PW2, TSROSWRODAS3, TSROSWSANDU3, TSROSWASSP3P3S1C, TSROSWLASSP3P4S2C, TSROSWLLSSP3P4S2C, TSROSWARK3 973 @*/ 974 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 975 { 976 PetscErrorCode ierr; 977 978 PetscFunctionBegin; 979 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 980 ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 981 PetscFunctionReturn(0); 982 } 983 984 #undef __FUNCT__ 985 #define __FUNCT__ "TSRosWGetType" 986 /*@C 987 TSRosWGetType - Get the type of Rosenbrock-W scheme 988 989 Logically collective 990 991 Input Parameter: 992 . ts - timestepping context 993 994 Output Parameter: 995 . rostype - type of Rosenbrock-W scheme 996 997 Level: intermediate 998 999 .seealso: TSRosWGetType() 1000 @*/ 1001 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 1002 { 1003 PetscErrorCode ierr; 1004 1005 PetscFunctionBegin; 1006 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1007 ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 1008 PetscFunctionReturn(0); 1009 } 1010 1011 #undef __FUNCT__ 1012 #define __FUNCT__ "TSRosWSetRecomputeJacobian" 1013 /*@C 1014 TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 1015 1016 Logically collective 1017 1018 Input Parameter: 1019 + ts - timestepping context 1020 - flg - PETSC_TRUE to recompute the Jacobian at each stage 1021 1022 Level: intermediate 1023 1024 .seealso: TSRosWGetType() 1025 @*/ 1026 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 1027 { 1028 PetscErrorCode ierr; 1029 1030 PetscFunctionBegin; 1031 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1032 ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1033 PetscFunctionReturn(0); 1034 } 1035 1036 EXTERN_C_BEGIN 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "TSRosWGetType_RosW" 1039 PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 1040 { 1041 TS_RosW *ros = (TS_RosW*)ts->data; 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 1046 *rostype = ros->tableau->name; 1047 PetscFunctionReturn(0); 1048 } 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "TSRosWSetType_RosW" 1051 PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 1052 { 1053 TS_RosW *ros = (TS_RosW*)ts->data; 1054 PetscErrorCode ierr; 1055 PetscBool match; 1056 RosWTableauLink link; 1057 1058 PetscFunctionBegin; 1059 if (ros->tableau) { 1060 ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 1061 if (match) PetscFunctionReturn(0); 1062 } 1063 for (link = RosWTableauList; link; link=link->next) { 1064 ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1065 if (match) { 1066 ierr = TSReset_RosW(ts);CHKERRQ(ierr); 1067 ros->tableau = &link->tab; 1068 PetscFunctionReturn(0); 1069 } 1070 } 1071 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1072 PetscFunctionReturn(0); 1073 } 1074 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 1077 PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1078 { 1079 TS_RosW *ros = (TS_RosW*)ts->data; 1080 1081 PetscFunctionBegin; 1082 ros->recompute_jacobian = flg; 1083 PetscFunctionReturn(0); 1084 } 1085 EXTERN_C_END 1086 1087 /* ------------------------------------------------------------ */ 1088 /*MC 1089 TSROSW - ODE solver using Rosenbrock-W schemes 1090 1091 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1092 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1093 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1094 1095 Notes: 1096 This method currently only works with autonomous ODE and DAE. 1097 1098 Developer notes: 1099 Rosenbrock-W methods are typically specified for autonomous ODE 1100 1101 $ xdot = f(x) 1102 1103 by the stage equations 1104 1105 $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 1106 1107 and step completion formula 1108 1109 $ x_1 = x_0 + sum_j b_j k_j 1110 1111 with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 1112 and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 1113 we define new variables for the stage equations 1114 1115 $ y_i = gamma_ij k_j 1116 1117 The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 1118 1119 $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 1120 1121 to rewrite the method as 1122 1123 $ [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j 1124 $ x_1 = x_0 + sum_j bt_j y_j 1125 1126 where we have introduced the mass matrix M. Continue by defining 1127 1128 $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 1129 1130 or, more compactly in tensor notation 1131 1132 $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 1133 1134 Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 1135 stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 1136 equation 1137 1138 $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 1139 1140 with initial guess y_i = 0. 1141 1142 Level: beginner 1143 1144 .seealso: TSCreate(), TS, TSSetType(), TSRosWSetType(), TSRosWRegister() 1145 1146 M*/ 1147 EXTERN_C_BEGIN 1148 #undef __FUNCT__ 1149 #define __FUNCT__ "TSCreate_RosW" 1150 PetscErrorCode TSCreate_RosW(TS ts) 1151 { 1152 TS_RosW *ros; 1153 PetscErrorCode ierr; 1154 1155 PetscFunctionBegin; 1156 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1157 ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1158 #endif 1159 1160 ts->ops->reset = TSReset_RosW; 1161 ts->ops->destroy = TSDestroy_RosW; 1162 ts->ops->view = TSView_RosW; 1163 ts->ops->setup = TSSetUp_RosW; 1164 ts->ops->step = TSStep_RosW; 1165 ts->ops->interpolate = TSInterpolate_RosW; 1166 ts->ops->evaluatestep = TSEvaluateStep_RosW; 1167 ts->ops->setfromoptions = TSSetFromOptions_RosW; 1168 ts->ops->snesfunction = SNESTSFormFunction_RosW; 1169 ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1170 1171 ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 1172 ts->data = (void*)ros; 1173 1174 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1175 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 1176 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1177 PetscFunctionReturn(0); 1178 } 1179 EXTERN_C_END 1180