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