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 PetscBool step_taken; /* ts->vec_sol has been advanced to the end of the current time step */ 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->step_taken) {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 594 else { 595 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 596 for (i=0; i<tab->s; i++) w[i] = tab->bt[i]; 597 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 598 } 599 if (done) *done = PETSC_TRUE; 600 PetscFunctionReturn(0); 601 } else if (order == tab->order-1) { 602 if (!tab->bembedt) goto unavailable; 603 if (ros->step_taken) { 604 for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i]; 605 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 606 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 607 } else { 608 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 609 for (i=0; i<tab->s; i++) w[i] = tab->bembedt[i]; 610 ierr = VecMAXPY(X,tab->s,w,ros->Y);CHKERRQ(ierr); 611 } 612 if (done) *done = PETSC_TRUE; 613 PetscFunctionReturn(0); 614 } 615 unavailable: 616 if (done) *done = PETSC_FALSE; 617 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); 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "TSStep_RosW" 623 static PetscErrorCode TSStep_RosW(TS ts) 624 { 625 TS_RosW *ros = (TS_RosW*)ts->data; 626 RosWTableau tab = ros->tableau; 627 const PetscInt s = tab->s; 628 const PetscReal *At = tab->At,*Gamma = tab->Gamma,*ASum = tab->ASum,*GammaInv = tab->GammaInv; 629 const PetscBool *GammaZeroDiag = tab->GammaZeroDiag; 630 PetscScalar *w = ros->work; 631 Vec *Y = ros->Y,Zdot = ros->Zdot,Zstage = ros->Zstage; 632 SNES snes; 633 TSAdapt adapt; 634 PetscInt i,j,its,lits,reject,next_scheme; 635 PetscReal next_time_step; 636 PetscBool accept; 637 PetscErrorCode ierr; 638 639 PetscFunctionBegin; 640 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 641 next_time_step = ts->time_step; 642 accept = PETSC_TRUE; 643 ros->step_taken = PETSC_FALSE; 644 645 for (reject=0; reject<ts->max_reject; reject++,ts->reject++) { 646 const PetscReal h = ts->time_step; 647 for (i=0; i<s; i++) { 648 ros->stage_time = ts->ptime + h*ASum[i]; 649 if (GammaZeroDiag[i]) { 650 ros->stage_explicit = PETSC_TRUE; 651 ros->shift = 1./h; 652 } else { 653 ros->stage_explicit = PETSC_FALSE; 654 ros->shift = 1./(h*Gamma[i*s+i]); 655 } 656 657 ierr = VecCopy(ts->vec_sol,Zstage);CHKERRQ(ierr); 658 for (j=0; j<i; j++) w[j] = At[i*s+j]; 659 ierr = VecMAXPY(Zstage,i,w,Y);CHKERRQ(ierr); 660 661 for (j=0; j<i; j++) w[j] = 1./h * GammaInv[i*s+j]; 662 ierr = VecZeroEntries(Zdot);CHKERRQ(ierr); 663 ierr = VecMAXPY(Zdot,i,w,Y);CHKERRQ(ierr); 664 665 /* Initial guess taken from last stage */ 666 ierr = VecZeroEntries(Y[i]);CHKERRQ(ierr); 667 668 if (!ros->recompute_jacobian && !i) { 669 ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 670 } 671 672 ierr = SNESSolve(snes,PETSC_NULL,Y[i]);CHKERRQ(ierr); 673 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 674 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 675 ts->nonlinear_its += its; ts->linear_its += lits; 676 } 677 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 678 ros->step_taken = PETSC_TRUE; 679 680 /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 681 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 682 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 683 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 684 ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 685 if (accept) { 686 /* ignore next_scheme for now */ 687 ts->ptime += ts->time_step; 688 ts->time_step = next_time_step; 689 ts->steps++; 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->step_taken = PETSC_FALSE; 696 } 697 } 698 699 PetscFunctionReturn(0); 700 } 701 702 #undef __FUNCT__ 703 #define __FUNCT__ "TSInterpolate_RosW" 704 static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 705 { 706 TS_RosW *ros = (TS_RosW*)ts->data; 707 708 PetscFunctionBegin; 709 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ros->tableau->name); 710 PetscFunctionReturn(0); 711 } 712 713 /*------------------------------------------------------------*/ 714 #undef __FUNCT__ 715 #define __FUNCT__ "TSReset_RosW" 716 static PetscErrorCode TSReset_RosW(TS ts) 717 { 718 TS_RosW *ros = (TS_RosW*)ts->data; 719 PetscInt s; 720 PetscErrorCode ierr; 721 722 PetscFunctionBegin; 723 if (!ros->tableau) PetscFunctionReturn(0); 724 s = ros->tableau->s; 725 ierr = VecDestroyVecs(s,&ros->Y);CHKERRQ(ierr); 726 ierr = VecDestroy(&ros->Ydot);CHKERRQ(ierr); 727 ierr = VecDestroy(&ros->Ystage);CHKERRQ(ierr); 728 ierr = VecDestroy(&ros->Zdot);CHKERRQ(ierr); 729 ierr = VecDestroy(&ros->Zstage);CHKERRQ(ierr); 730 ierr = PetscFree(ros->work);CHKERRQ(ierr); 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNCT__ 735 #define __FUNCT__ "TSDestroy_RosW" 736 static PetscErrorCode TSDestroy_RosW(TS ts) 737 { 738 PetscErrorCode ierr; 739 740 PetscFunctionBegin; 741 ierr = TSReset_RosW(ts);CHKERRQ(ierr); 742 ierr = PetscFree(ts->data);CHKERRQ(ierr); 743 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 744 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 745 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","",PETSC_NULL);CHKERRQ(ierr); 746 PetscFunctionReturn(0); 747 } 748 749 /* 750 This defines the nonlinear equation that is to be solved with SNES 751 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 752 */ 753 #undef __FUNCT__ 754 #define __FUNCT__ "SNESTSFormFunction_RosW" 755 static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 756 { 757 TS_RosW *ros = (TS_RosW*)ts->data; 758 PetscErrorCode ierr; 759 760 PetscFunctionBegin; 761 if (ros->stage_explicit) { 762 ierr = VecAXPBY(ros->Ydot,ros->shift,0.0,X);CHKERRQ(ierr); /* Ydot = shift*X*/ 763 } else { 764 ierr = VecWAXPY(ros->Ydot,ros->shift,X,ros->Zdot);CHKERRQ(ierr); /* Ydot = shift*X + Zdot */ 765 } 766 ierr = VecWAXPY(ros->Ystage,1.0,X,ros->Zstage);CHKERRQ(ierr); /* Ystage = X + Zstage */ 767 ierr = TSComputeIFunction(ts,ros->stage_time,ros->Ystage,ros->Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "SNESTSFormJacobian_RosW" 773 static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 774 { 775 TS_RosW *ros = (TS_RosW*)ts->data; 776 PetscErrorCode ierr; 777 778 PetscFunctionBegin; 779 /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 780 ierr = TSComputeIJacobian(ts,ros->stage_time,ros->Ystage,ros->Ydot,ros->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 781 PetscFunctionReturn(0); 782 } 783 784 #undef __FUNCT__ 785 #define __FUNCT__ "TSSetUp_RosW" 786 static PetscErrorCode TSSetUp_RosW(TS ts) 787 { 788 TS_RosW *ros = (TS_RosW*)ts->data; 789 RosWTableau tab = ros->tableau; 790 PetscInt s = tab->s; 791 PetscErrorCode ierr; 792 793 PetscFunctionBegin; 794 if (!ros->tableau) { 795 ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 796 } 797 ierr = VecDuplicateVecs(ts->vec_sol,s,&ros->Y);CHKERRQ(ierr); 798 ierr = VecDuplicate(ts->vec_sol,&ros->Ydot);CHKERRQ(ierr); 799 ierr = VecDuplicate(ts->vec_sol,&ros->Ystage);CHKERRQ(ierr); 800 ierr = VecDuplicate(ts->vec_sol,&ros->Zdot);CHKERRQ(ierr); 801 ierr = VecDuplicate(ts->vec_sol,&ros->Zstage);CHKERRQ(ierr); 802 ierr = PetscMalloc(s*sizeof(ros->work[0]),&ros->work);CHKERRQ(ierr); 803 PetscFunctionReturn(0); 804 } 805 /*------------------------------------------------------------*/ 806 807 #undef __FUNCT__ 808 #define __FUNCT__ "TSSetFromOptions_RosW" 809 static PetscErrorCode TSSetFromOptions_RosW(TS ts) 810 { 811 TS_RosW *ros = (TS_RosW*)ts->data; 812 PetscErrorCode ierr; 813 char rostype[256]; 814 815 PetscFunctionBegin; 816 ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 817 { 818 RosWTableauLink link; 819 PetscInt count,choice; 820 PetscBool flg; 821 const char **namelist; 822 SNES snes; 823 824 ierr = PetscStrncpy(rostype,TSRosWDefault,sizeof rostype);CHKERRQ(ierr); 825 for (link=RosWTableauList,count=0; link; link=link->next,count++) ; 826 ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 827 for (link=RosWTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 828 ierr = PetscOptionsEList("-ts_rosw_type","Family of Rosenbrock-W method","TSRosWSetType",(const char*const*)namelist,count,rostype,&choice,&flg);CHKERRQ(ierr); 829 ierr = TSRosWSetType(ts,flg ? namelist[choice] : rostype);CHKERRQ(ierr); 830 ierr = PetscFree(namelist);CHKERRQ(ierr); 831 832 ierr = PetscOptionsBool("-ts_rosw_recompute_jacobian","Recompute the Jacobian at each stage","TSRosWSetRecomputeJacobian",ros->recompute_jacobian,&ros->recompute_jacobian,PETSC_NULL);CHKERRQ(ierr); 833 834 /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */ 835 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 836 if (!((PetscObject)snes)->type_name) { 837 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 838 } 839 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 840 } 841 ierr = PetscOptionsTail();CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845 #undef __FUNCT__ 846 #define __FUNCT__ "PetscFormatRealArray" 847 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 848 { 849 PetscErrorCode ierr; 850 PetscInt i; 851 size_t left,count; 852 char *p; 853 854 PetscFunctionBegin; 855 for (i=0,p=buf,left=len; i<n; i++) { 856 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 857 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 858 left -= count; 859 p += count; 860 *p++ = ' '; 861 } 862 p[i ? 0 : -1] = 0; 863 PetscFunctionReturn(0); 864 } 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "TSView_RosW" 868 static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 869 { 870 TS_RosW *ros = (TS_RosW*)ts->data; 871 RosWTableau tab = ros->tableau; 872 PetscBool iascii; 873 PetscErrorCode ierr; 874 875 PetscFunctionBegin; 876 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 877 if (iascii) { 878 const TSRosWType rostype; 879 PetscInt i; 880 PetscReal abscissa[512]; 881 char buf[512]; 882 ierr = TSRosWGetType(ts,&rostype);CHKERRQ(ierr); 883 ierr = PetscViewerASCIIPrintf(viewer," Rosenbrock-W %s\n",rostype);CHKERRQ(ierr); 884 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ASum);CHKERRQ(ierr); 885 ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A = %s\n",buf);CHKERRQ(ierr); 886 for (i=0; i<tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i]; 887 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,abscissa);CHKERRQ(ierr); 888 ierr = PetscViewerASCIIPrintf(viewer," Abscissa of A+Gamma = %s\n",buf);CHKERRQ(ierr); 889 } 890 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 891 PetscFunctionReturn(0); 892 } 893 894 #undef __FUNCT__ 895 #define __FUNCT__ "TSRosWSetType" 896 /*@C 897 TSRosWSetType - Set the type of Rosenbrock-W scheme 898 899 Logically collective 900 901 Input Parameter: 902 + ts - timestepping context 903 - rostype - type of Rosenbrock-W scheme 904 905 Level: intermediate 906 907 .seealso: TSRosWGetType() 908 @*/ 909 PetscErrorCode TSRosWSetType(TS ts,const TSRosWType rostype) 910 { 911 PetscErrorCode ierr; 912 913 PetscFunctionBegin; 914 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 915 ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,rostype));CHKERRQ(ierr); 916 PetscFunctionReturn(0); 917 } 918 919 #undef __FUNCT__ 920 #define __FUNCT__ "TSRosWGetType" 921 /*@C 922 TSRosWGetType - Get the type of Rosenbrock-W scheme 923 924 Logically collective 925 926 Input Parameter: 927 . ts - timestepping context 928 929 Output Parameter: 930 . rostype - type of Rosenbrock-W scheme 931 932 Level: intermediate 933 934 .seealso: TSRosWGetType() 935 @*/ 936 PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *rostype) 937 { 938 PetscErrorCode ierr; 939 940 PetscFunctionBegin; 941 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 942 ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,rostype));CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "TSRosWSetRecomputeJacobian" 948 /*@C 949 TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step. 950 951 Logically collective 952 953 Input Parameter: 954 + ts - timestepping context 955 - flg - PETSC_TRUE to recompute the Jacobian at each stage 956 957 Level: intermediate 958 959 .seealso: TSRosWGetType() 960 @*/ 961 PetscErrorCode TSRosWSetRecomputeJacobian(TS ts,PetscBool flg) 962 { 963 PetscErrorCode ierr; 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 967 ierr = PetscTryMethod(ts,"TSRosWSetRecomputeJacobian_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 EXTERN_C_BEGIN 972 #undef __FUNCT__ 973 #define __FUNCT__ "TSRosWGetType_RosW" 974 PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *rostype) 975 { 976 TS_RosW *ros = (TS_RosW*)ts->data; 977 PetscErrorCode ierr; 978 979 PetscFunctionBegin; 980 if (!ros->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 981 *rostype = ros->tableau->name; 982 PetscFunctionReturn(0); 983 } 984 #undef __FUNCT__ 985 #define __FUNCT__ "TSRosWSetType_RosW" 986 PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType rostype) 987 { 988 TS_RosW *ros = (TS_RosW*)ts->data; 989 PetscErrorCode ierr; 990 PetscBool match; 991 RosWTableauLink link; 992 993 PetscFunctionBegin; 994 if (ros->tableau) { 995 ierr = PetscStrcmp(ros->tableau->name,rostype,&match);CHKERRQ(ierr); 996 if (match) PetscFunctionReturn(0); 997 } 998 for (link = RosWTableauList; link; link=link->next) { 999 ierr = PetscStrcmp(link->tab.name,rostype,&match);CHKERRQ(ierr); 1000 if (match) { 1001 ierr = TSReset_RosW(ts);CHKERRQ(ierr); 1002 ros->tableau = &link->tab; 1003 PetscFunctionReturn(0); 1004 } 1005 } 1006 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rostype); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 #undef __FUNCT__ 1011 #define __FUNCT__ "TSRosWSetRecomputeJacobian_RosW" 1012 PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts,PetscBool flg) 1013 { 1014 TS_RosW *ros = (TS_RosW*)ts->data; 1015 1016 PetscFunctionBegin; 1017 ros->recompute_jacobian = flg; 1018 PetscFunctionReturn(0); 1019 } 1020 EXTERN_C_END 1021 1022 /* ------------------------------------------------------------ */ 1023 /*MC 1024 TSRosW - ODE solver using Rosenbrock-W schemes 1025 1026 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1027 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1028 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1029 1030 Notes: 1031 This method currently only works with autonomous ODE and DAE. 1032 1033 Developer notes: 1034 Rosenbrock-W methods are typically specified for autonomous ODE 1035 1036 $ xdot = f(x) 1037 1038 by the stage equations 1039 1040 $ k_i = h f(x_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j 1041 1042 and step completion formula 1043 1044 $ x_1 = x_0 + sum_j b_j k_j 1045 1046 with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(x) 1047 and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner, 1048 we define new variables for the stage equations 1049 1050 $ y_i = gamma_ij k_j 1051 1052 The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define 1053 1054 $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-i} 1055 1056 to rewrite the method as 1057 1058 $ [M/(h gamma_ii) - J] y_i = f(x_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j 1059 $ x_1 = x_0 + sum_j bt_j y_j 1060 1061 where we have introduced the mass matrix M. Continue by defining 1062 1063 $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j 1064 1065 or, more compactly in tensor notation 1066 1067 $ Ydot = 1/h (Gamma^{-1} \otimes I) Y . 1068 1069 Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current 1070 stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the 1071 equation 1072 1073 $ g(x_0 + sum_j a_ij y_j + y_i, ydot_i) = 0 1074 1075 with initial guess y_i = 0. 1076 1077 Level: beginner 1078 1079 .seealso: TSCreate(), TS, TSSetType(), TSRosWRegister() 1080 1081 M*/ 1082 EXTERN_C_BEGIN 1083 #undef __FUNCT__ 1084 #define __FUNCT__ "TSCreate_RosW" 1085 PetscErrorCode TSCreate_RosW(TS ts) 1086 { 1087 TS_RosW *ros; 1088 PetscErrorCode ierr; 1089 1090 PetscFunctionBegin; 1091 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 1092 ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 1093 #endif 1094 1095 ts->ops->reset = TSReset_RosW; 1096 ts->ops->destroy = TSDestroy_RosW; 1097 ts->ops->view = TSView_RosW; 1098 ts->ops->setup = TSSetUp_RosW; 1099 ts->ops->step = TSStep_RosW; 1100 ts->ops->interpolate = TSInterpolate_RosW; 1101 ts->ops->evaluatestep = TSEvaluateStep_RosW; 1102 ts->ops->setfromoptions = TSSetFromOptions_RosW; 1103 ts->ops->snesfunction = SNESTSFormFunction_RosW; 1104 ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 1105 1106 ierr = PetscNewLog(ts,TS_RosW,&ros);CHKERRQ(ierr); 1107 ts->data = (void*)ros; 1108 1109 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 1110 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 1111 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetRecomputeJacobian_C","TSRosWSetRecomputeJacobian_RosW",TSRosWSetRecomputeJacobian_RosW);CHKERRQ(ierr); 1112 PetscFunctionReturn(0); 1113 } 1114 EXTERN_C_END 1115