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