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