1 /* 2 Code for timestepping with additive Runge-Kutta IMEX method 3 4 Notes: 5 The general system is written as 6 7 F(t,X,Xdot) = G(t,X) 8 9 where F represents the stiff part of the physics and G represents the non-stiff part. 10 11 */ 12 #include <private/tsimpl.h> /*I "petscts.h" I*/ 13 14 static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2E; 15 static PetscBool TSARKIMEXRegisterAllCalled; 16 static PetscBool TSARKIMEXPackageInitialized; 17 18 typedef struct _ARKTableau *ARKTableau; 19 struct _ARKTableau { 20 char *name; 21 PetscInt order; /* Classical approximation order of the method */ 22 PetscInt s; /* Number of stages */ 23 PetscInt pinterp; /* Interpolation order */ 24 PetscReal *At,*bt,*ct; /* Stiff tableau */ 25 PetscReal *A,*b,*c; /* Non-stiff tableau */ 26 PetscReal *binterpt,*binterp; /* Dense output formula */ 27 }; 28 typedef struct _ARKTableauLink *ARKTableauLink; 29 struct _ARKTableauLink { 30 struct _ARKTableau tab; 31 ARKTableauLink next; 32 }; 33 static ARKTableauLink ARKTableauList; 34 35 typedef struct { 36 ARKTableau tableau; 37 Vec *Y; /* States computed during the step */ 38 Vec *YdotI; /* Time derivatives for the stiff part */ 39 Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 40 Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 41 Vec Work; /* Generic work vector */ 42 Vec Z; /* Ydot = shift(Y-Z) */ 43 PetscScalar *work; /* Scalar work */ 44 PetscReal shift; 45 PetscReal stage_time; 46 PetscBool imex; 47 } TS_ARKIMEX; 48 49 /*MC 50 TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 51 52 This method has one explicit stage and two implicit stages. 53 54 Level: advanced 55 56 .seealso: TSARKIMEX 57 M*/ 58 /*MC 59 TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 60 61 This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 62 63 Level: advanced 64 65 .seealso: TSARKIMEX 66 M*/ 67 /*MC 68 TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 69 70 This method has three implicit stages. 71 72 References: 73 L. Pareschi, G. Russo, Implicit-Explicit Runge-Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005, pp. 129-155 74 75 This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375 76 77 Level: advanced 78 79 .seealso: TSARKIMEX 80 M*/ 81 /*MC 82 TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 83 84 This method has one explicit stage and three implicit stages. 85 86 References: 87 Kennedy and Carpenter 2003. 88 89 Level: advanced 90 91 .seealso: TSARKIMEX 92 M*/ 93 /*MC 94 TSARKIMEXARS443 - Third order ARK IMEX scheme. 95 96 This method has one explicit stage and four implicit stages. 97 98 References: 99 U. Ascher, S. Ruuth, R. J. Spitheri, Implicit-explicit Runge-Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997), pp. 151–167. 100 101 This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375 102 103 Level: advanced 104 105 .seealso: TSARKIMEX 106 M*/ 107 /*MC 108 TSARKIMEXBPR3 - Third order ARK IMEX scheme. 109 110 This method has one explicit stage and four implicit stages. 111 112 References: 113 This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375 114 115 Level: advanced 116 117 .seealso: TSARKIMEX 118 M*/ 119 /*MC 120 TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 121 122 This method has one explicit stage and four implicit stages. 123 124 References: 125 Kennedy and Carpenter 2003. 126 127 Level: advanced 128 129 .seealso: TSARKIMEX 130 M*/ 131 /*MC 132 TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 133 134 This method has one explicit stage and five implicit stages. 135 136 References: 137 Kennedy and Carpenter 2003. 138 139 Level: advanced 140 141 .seealso: TSARKIMEX 142 M*/ 143 144 #undef __FUNCT__ 145 #define __FUNCT__ "TSARKIMEXRegisterAll" 146 /*@C 147 TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 148 149 Not Collective, but should be called by all processes which will need the schemes to be registered 150 151 Level: advanced 152 153 .keywords: TS, TSARKIMEX, register, all 154 155 .seealso: TSARKIMEXRegisterDestroy() 156 @*/ 157 PetscErrorCode TSARKIMEXRegisterAll(void) 158 { 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 163 TSARKIMEXRegisterAllCalled = PETSC_TRUE; 164 165 { 166 const PetscReal 167 A[3][3] = {{0,0,0}, 168 {0.41421356237309504880,0,0}, 169 {0.75,0.25,0}}, 170 At[3][3] = {{0,0,0}, 171 {0.12132034355964257320,0.29289321881345247560,0}, 172 {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, 173 binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 174 ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 175 } 176 { /* Optimal for linear implicit part */ 177 const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 178 A[3][3] = {{0,0,0}, 179 {2-s2,0,0}, 180 {(3-2*s2)/6,(3+2*s2)/6,0}}, 181 At[3][3] = {{0,0,0}, 182 {1-1/s2,1-1/s2,0}, 183 {1/(2*s2),1/(2*s2),1-1/s2}}, 184 binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 185 ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 186 } 187 { /* Optimal for linear implicit part */ 188 const PetscReal 189 A[3][3] = {{0,0,0}, 190 {0.5,0,0}, 191 {0.5,0.5,0}}, 192 At[3][3] = {{0.25,0,0}, 193 {0,0.25,0}, 194 {1./3,1./3,1./3}}; 195 ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 196 } 197 { 198 const PetscReal 199 A[4][4] = {{0,0,0,0}, 200 {1767732205903./2027836641118.,0,0,0}, 201 {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 202 {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 203 At[4][4] = {{0,0,0,0}, 204 {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 205 {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 206 {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 207 binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 208 {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 209 {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 210 {584795268549./6622622206610., 2508943948391./7218656332882.}}; 211 ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 212 } 213 { 214 const PetscReal 215 A[5][5] = {{0,0,0,0}, 216 {1./2,0,0,0,0}, 217 {11./18,1./18,0,0,0}, 218 {5./6,-5./6,.5,0,0}, 219 {1./4,7./4,3./4,-7./4,0}}, 220 At[5][5] = {{0,0,0,0,0}, 221 {0,1./2,0,0,0}, 222 {0,1./6,1./2,0,0}, 223 {0,-1./2,1./2,1./2,0}, 224 {0,3./2,-3./2,1./2,1./2}}; 225 ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 226 } 227 { 228 const PetscReal 229 A[5][5] = {{0,0,0,0}, 230 {1,0,0,0,0}, 231 {4./9,2./9,0,0,0}, 232 {1./4,0,3./4,0,0}, 233 {1./4,0,3./5,0,0}}, 234 At[5][5] = {{0,0,0,0}, 235 {.5,.5,0,0,0}, 236 {5./18,-1./9,.5,0,0}, 237 {.5,0,0,.5,0}, 238 {.25,0,.75,-.5,.5}}; 239 ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 240 } 241 { 242 const PetscReal 243 A[6][6] = {{0,0,0,0,0,0}, 244 {1./2,0,0,0,0,0}, 245 {13861./62500.,6889./62500.,0,0,0,0}, 246 {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 247 {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 248 {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 249 At[6][6] = {{0,0,0,0,0,0}, 250 {1./4,1./4,0,0,0,0}, 251 {8611./62500.,-1743./31250.,1./4,0,0,0}, 252 {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 253 {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 254 {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 255 binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 256 {0,0,0}, 257 {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 258 {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 259 {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 260 {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 261 ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 262 } 263 { 264 const PetscReal 265 A[8][8] = {{0,0,0,0,0,0,0,0}, 266 {41./100,0,0,0,0,0,0,0}, 267 {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 268 {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 269 {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 270 {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 271 {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 272 {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 273 At[8][8] = {{0,0,0,0,0,0,0,0}, 274 {41./200.,41./200.,0,0,0,0,0,0}, 275 {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 276 {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 277 {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 278 {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 279 {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 280 {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 281 binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.}, 282 {0 , 0 , 0 }, 283 {0 , 0 , 0 }, 284 {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.}, 285 {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.}, 286 {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.}, 287 {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.}, 288 {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }}; 289 ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 290 } 291 292 PetscFunctionReturn(0); 293 } 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "TSARKIMEXRegisterDestroy" 297 /*@C 298 TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 299 300 Not Collective 301 302 Level: advanced 303 304 .keywords: TSARKIMEX, register, destroy 305 .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic() 306 @*/ 307 PetscErrorCode TSARKIMEXRegisterDestroy(void) 308 { 309 PetscErrorCode ierr; 310 ARKTableauLink link; 311 312 PetscFunctionBegin; 313 while ((link = ARKTableauList)) { 314 ARKTableau t = &link->tab; 315 ARKTableauList = link->next; 316 ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 317 ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 318 ierr = PetscFree(t->name);CHKERRQ(ierr); 319 ierr = PetscFree(link);CHKERRQ(ierr); 320 } 321 TSARKIMEXRegisterAllCalled = PETSC_FALSE; 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "TSARKIMEXInitializePackage" 327 /*@C 328 TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 329 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 330 when using static libraries. 331 332 Input Parameter: 333 path - The dynamic library path, or PETSC_NULL 334 335 Level: developer 336 337 .keywords: TS, TSARKIMEX, initialize, package 338 .seealso: PetscInitialize() 339 @*/ 340 PetscErrorCode TSARKIMEXInitializePackage(const char path[]) 341 { 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 346 TSARKIMEXPackageInitialized = PETSC_TRUE; 347 ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 348 ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 349 PetscFunctionReturn(0); 350 } 351 352 #undef __FUNCT__ 353 #define __FUNCT__ "TSARKIMEXFinalizePackage" 354 /*@C 355 TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 356 called from PetscFinalize(). 357 358 Level: developer 359 360 .keywords: Petsc, destroy, package 361 .seealso: PetscFinalize() 362 @*/ 363 PetscErrorCode TSARKIMEXFinalizePackage(void) 364 { 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 TSARKIMEXPackageInitialized = PETSC_FALSE; 369 ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 370 PetscFunctionReturn(0); 371 } 372 373 #undef __FUNCT__ 374 #define __FUNCT__ "TSARKIMEXRegister" 375 /*@C 376 TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 377 378 Not Collective, but the same schemes should be registered on all processes on which they will be used 379 380 Input Parameters: 381 + name - identifier for method 382 . order - approximation order of method 383 . s - number of stages, this is the dimension of the matrices below 384 . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 385 . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At) 386 . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At) 387 . A - Non-stiff stage coefficients (dimension s*s, row-major) 388 . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At) 389 . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A) 390 . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 391 . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 392 - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt) 393 394 Notes: 395 Several ARK IMEX methods are provided, this function is only needed to create new methods. 396 397 Level: advanced 398 399 .keywords: TS, register 400 401 .seealso: TSARKIMEX 402 @*/ 403 PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s, 404 const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 405 const PetscReal A[],const PetscReal b[],const PetscReal c[], 406 PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 407 { 408 PetscErrorCode ierr; 409 ARKTableauLink link; 410 ARKTableau t; 411 PetscInt i,j; 412 413 PetscFunctionBegin; 414 ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 415 ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 416 t = &link->tab; 417 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 418 t->order = order; 419 t->s = s; 420 ierr = PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 421 ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 422 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 423 if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);} 424 else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 425 if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);} 426 else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i]; 427 if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);} 428 else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 429 if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);} 430 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 431 t->pinterp = pinterp; 432 ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr); 433 ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 434 ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 435 link->next = ARKTableauList; 436 ARKTableauList = link; 437 PetscFunctionReturn(0); 438 } 439 440 #undef __FUNCT__ 441 #define __FUNCT__ "TSStep_ARKIMEX" 442 static PetscErrorCode TSStep_ARKIMEX(TS ts) 443 { 444 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 445 ARKTableau tab = ark->tableau; 446 const PetscInt s = tab->s; 447 const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c; 448 PetscScalar *w = ark->work; 449 Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z; 450 SNES snes; 451 SNESConvergedReason snesreason; 452 PetscInt i,j,its,lits; 453 PetscReal next_time_step; 454 PetscReal h,t; 455 PetscErrorCode ierr; 456 457 PetscFunctionBegin; 458 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 459 next_time_step = ts->time_step; 460 h = ts->time_step; 461 t = ts->ptime; 462 463 for (i=0; i<s; i++) { 464 if (At[i*s+i] == 0) { /* This stage is explicit */ 465 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 466 for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 467 ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 468 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 469 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 470 } else { 471 ark->stage_time = t + h*ct[i]; 472 ark->shift = 1./(h*At[i*s+i]); 473 /* Affine part */ 474 ierr = VecZeroEntries(W);CHKERRQ(ierr); 475 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 476 ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 477 /* Ydot = shift*(Y-Z) */ 478 ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 479 for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 480 ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 481 /* Initial guess taken from last stage */ 482 ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr); 483 ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 484 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 485 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 486 ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 487 ts->nonlinear_its += its; ts->linear_its += lits; 488 if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 489 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 490 ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 } 494 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 495 ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr); 496 if (ark->imex) { 497 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 498 } else { 499 ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 500 } 501 } 502 for (j=0; j<s; j++) w[j] = -h*bt[j]; 503 ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 504 for (j=0; j<s; j++) w[j] = h*b[j]; 505 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 506 507 ts->ptime += ts->time_step; 508 ts->time_step = next_time_step; 509 ts->steps++; 510 PetscFunctionReturn(0); 511 } 512 513 #undef __FUNCT__ 514 #define __FUNCT__ "TSInterpolate_ARKIMEX" 515 static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 516 { 517 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 518 PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 519 PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */ 520 PetscScalar *bt,*b; 521 const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 522 PetscErrorCode ierr; 523 524 PetscFunctionBegin; 525 if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 526 ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr); 527 for (i=0; i<s; i++) bt[i] = b[i] = 0; 528 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 529 for (i=0; i<s; i++) { 530 bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0; 531 b[i] += ts->time_step * B[i*pinterp+j] * tt; 532 } 533 } 534 if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved"); 535 ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 536 ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 537 ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 538 ierr = PetscFree2(bt,b);CHKERRQ(ierr); 539 PetscFunctionReturn(0); 540 } 541 542 /*------------------------------------------------------------*/ 543 #undef __FUNCT__ 544 #define __FUNCT__ "TSReset_ARKIMEX" 545 static PetscErrorCode TSReset_ARKIMEX(TS ts) 546 { 547 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 548 PetscInt s; 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 if (!ark->tableau) PetscFunctionReturn(0); 553 s = ark->tableau->s; 554 ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 555 ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 556 ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 557 ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 558 ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 559 ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 560 ierr = PetscFree(ark->work);CHKERRQ(ierr); 561 PetscFunctionReturn(0); 562 } 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "TSDestroy_ARKIMEX" 566 static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 567 { 568 PetscErrorCode ierr; 569 570 PetscFunctionBegin; 571 ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 572 ierr = PetscFree(ts->data);CHKERRQ(ierr); 573 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr); 574 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr); 575 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",PETSC_NULL);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 /* 580 This defines the nonlinear equation that is to be solved with SNES 581 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 582 */ 583 #undef __FUNCT__ 584 #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 585 static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 586 { 587 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 592 ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 #undef __FUNCT__ 597 #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 598 static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 599 { 600 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 605 ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 606 PetscFunctionReturn(0); 607 } 608 609 #undef __FUNCT__ 610 #define __FUNCT__ "TSSetUp_ARKIMEX" 611 static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 612 { 613 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 614 ARKTableau tab = ark->tableau; 615 PetscInt s = tab->s; 616 PetscErrorCode ierr; 617 618 PetscFunctionBegin; 619 if (!ark->tableau) { 620 ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 621 } 622 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 623 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 624 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 625 ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 626 ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 627 ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 628 ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr); 629 PetscFunctionReturn(0); 630 } 631 /*------------------------------------------------------------*/ 632 633 #undef __FUNCT__ 634 #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 635 static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts) 636 { 637 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 638 PetscErrorCode ierr; 639 char arktype[256]; 640 641 PetscFunctionBegin; 642 ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr); 643 { 644 ARKTableauLink link; 645 PetscInt count,choice; 646 PetscBool flg; 647 const char **namelist; 648 ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr); 649 for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 650 ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 651 for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 652 ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 653 ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 654 ierr = PetscFree(namelist);CHKERRQ(ierr); 655 flg = (PetscBool)!ark->imex; 656 ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 657 ark->imex = (PetscBool)!flg; 658 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 659 } 660 ierr = PetscOptionsTail();CHKERRQ(ierr); 661 PetscFunctionReturn(0); 662 } 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "PetscFormatRealArray" 666 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 667 { 668 PetscErrorCode ierr; 669 PetscInt i; 670 size_t left,count; 671 char *p; 672 673 PetscFunctionBegin; 674 for (i=0,p=buf,left=len; i<n; i++) { 675 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 676 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 677 left -= count; 678 p += count; 679 *p++ = ' '; 680 } 681 p[i ? 0 : -1] = 0; 682 PetscFunctionReturn(0); 683 } 684 685 #undef __FUNCT__ 686 #define __FUNCT__ "TSView_ARKIMEX" 687 static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 688 { 689 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 690 ARKTableau tab = ark->tableau; 691 PetscBool iascii; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 696 if (iascii) { 697 const TSARKIMEXType arktype; 698 char buf[512]; 699 ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 700 ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 701 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 702 ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 703 ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 704 ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 705 } 706 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 707 PetscFunctionReturn(0); 708 } 709 710 #undef __FUNCT__ 711 #define __FUNCT__ "TSARKIMEXSetType" 712 /*@C 713 TSARKIMEXSetType - Set the type of ARK IMEX scheme 714 715 Logically collective 716 717 Input Parameter: 718 + ts - timestepping context 719 - arktype - type of ARK-IMEX scheme 720 721 Level: intermediate 722 723 .seealso: TSARKIMEXGetType() 724 @*/ 725 PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype) 726 { 727 PetscErrorCode ierr; 728 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 731 ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 732 PetscFunctionReturn(0); 733 } 734 735 #undef __FUNCT__ 736 #define __FUNCT__ "TSARKIMEXGetType" 737 /*@C 738 TSARKIMEXGetType - Get the type of ARK IMEX scheme 739 740 Logically collective 741 742 Input Parameter: 743 . ts - timestepping context 744 745 Output Parameter: 746 . arktype - type of ARK-IMEX scheme 747 748 Level: intermediate 749 750 .seealso: TSARKIMEXGetType() 751 @*/ 752 PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype) 753 { 754 PetscErrorCode ierr; 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 758 ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "TSARKIMEXSetFullyImplicit" 764 /*@C 765 TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 766 767 Logically collective 768 769 Input Parameter: 770 + ts - timestepping context 771 - flg - PETSC_TRUE for fully implicit 772 773 Level: intermediate 774 775 .seealso: TSARKIMEXGetType() 776 @*/ 777 PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 778 { 779 PetscErrorCode ierr; 780 781 PetscFunctionBegin; 782 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 783 ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 784 PetscFunctionReturn(0); 785 } 786 787 EXTERN_C_BEGIN 788 #undef __FUNCT__ 789 #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 790 PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype) 791 { 792 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 793 PetscErrorCode ierr; 794 795 PetscFunctionBegin; 796 if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);} 797 *arktype = ark->tableau->name; 798 PetscFunctionReturn(0); 799 } 800 #undef __FUNCT__ 801 #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 802 PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype) 803 { 804 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 805 PetscErrorCode ierr; 806 PetscBool match; 807 ARKTableauLink link; 808 809 PetscFunctionBegin; 810 if (ark->tableau) { 811 ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 812 if (match) PetscFunctionReturn(0); 813 } 814 for (link = ARKTableauList; link; link=link->next) { 815 ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 816 if (match) { 817 ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 818 ark->tableau = &link->tab; 819 PetscFunctionReturn(0); 820 } 821 } 822 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 823 PetscFunctionReturn(0); 824 } 825 #undef __FUNCT__ 826 #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX" 827 PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 828 { 829 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 830 831 PetscFunctionBegin; 832 ark->imex = (PetscBool)!flg; 833 PetscFunctionReturn(0); 834 } 835 EXTERN_C_END 836 837 /* ------------------------------------------------------------ */ 838 /*MC 839 TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes 840 841 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 842 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 843 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 844 845 Notes: 846 The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 847 848 This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 849 850 Level: beginner 851 852 .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3, 853 TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, TSARKIMEXRegister() 854 855 M*/ 856 EXTERN_C_BEGIN 857 #undef __FUNCT__ 858 #define __FUNCT__ "TSCreate_ARKIMEX" 859 PetscErrorCode TSCreate_ARKIMEX(TS ts) 860 { 861 TS_ARKIMEX *th; 862 PetscErrorCode ierr; 863 864 PetscFunctionBegin; 865 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 866 ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr); 867 #endif 868 869 ts->ops->reset = TSReset_ARKIMEX; 870 ts->ops->destroy = TSDestroy_ARKIMEX; 871 ts->ops->view = TSView_ARKIMEX; 872 ts->ops->setup = TSSetUp_ARKIMEX; 873 ts->ops->step = TSStep_ARKIMEX; 874 ts->ops->interpolate = TSInterpolate_ARKIMEX; 875 ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 876 ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 877 ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 878 879 ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr); 880 ts->data = (void*)th; 881 th->imex = PETSC_TRUE; 882 883 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 884 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 885 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 886 PetscFunctionReturn(0); 887 } 888 EXTERN_C_END 889