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