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