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