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,U,Udot) = G(t,U) 8 9 where F represents the stiff part of the physics and G represents the non-stiff part. 10 11 */ 12 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 13 #include <petscdm.h> 14 15 static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3; 16 static PetscBool TSARKIMEXRegisterAllCalled; 17 static PetscBool TSARKIMEXPackageInitialized; 18 static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec); 19 20 typedef struct _ARKTableau *ARKTableau; 21 struct _ARKTableau { 22 char *name; 23 PetscInt order; /* Classical approximation order of the method */ 24 PetscInt s; /* Number of stages */ 25 PetscBool stiffly_accurate; /* The implicit part is stiffly accurate*/ 26 PetscBool FSAL_implicit; /* The implicit part is FSAL*/ 27 PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/ 28 PetscInt pinterp; /* Interpolation order */ 29 PetscReal *At,*bt,*ct; /* Stiff tableau */ 30 PetscReal *A,*b,*c; /* Non-stiff tableau */ 31 PetscReal *bembedt,*bembed; /* Embedded formula of order one less (order-1) */ 32 PetscReal *binterpt,*binterp; /* Dense output formula */ 33 PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 34 }; 35 typedef struct _ARKTableauLink *ARKTableauLink; 36 struct _ARKTableauLink { 37 struct _ARKTableau tab; 38 ARKTableauLink next; 39 }; 40 static ARKTableauLink ARKTableauList; 41 42 typedef struct { 43 ARKTableau tableau; 44 Vec *Y; /* States computed during the step */ 45 Vec *YdotI; /* Time derivatives for the stiff part */ 46 Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 47 Vec *Y_prev; /* States computed during the previous time step */ 48 Vec *YdotI_prev; /* Time derivatives for the stiff part for the previous time step*/ 49 Vec *YdotRHS_prev; /* Function evaluations for the non-stiff part for the previous time step*/ 50 Vec Ydot0; /* Holds the slope from the previous step in FSAL case */ 51 Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 52 Vec Z; /* Ydot = shift(Y-Z) */ 53 PetscScalar *work; /* Scalar work */ 54 PetscReal scoeff; /* shift = scoeff/dt */ 55 PetscReal stage_time; 56 PetscBool imex; 57 PetscBool extrapolate; /* Extrapolate initial guess from previous time-step stage values */ 58 TSStepStatus status; 59 } TS_ARKIMEX; 60 /*MC 61 TSARKIMEXARS122 - Second order ARK IMEX scheme. 62 63 This method has one explicit stage and one implicit stage. 64 65 Options Database: 66 . -ts_arkimex_type ars122 67 68 References: 69 . 1. - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 70 71 Level: advanced 72 73 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 74 M*/ 75 /*MC 76 TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 77 78 This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu. 79 80 Options Database: 81 . -ts_arkimex_type a2 82 83 Level: advanced 84 85 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 86 M*/ 87 /*MC 88 TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 89 90 This method has two implicit stages, and L-stable implicit scheme. 91 92 Options Database: 93 . -ts_arkimex_type l2 94 95 References: 96 . 1. - 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. 97 98 Level: advanced 99 100 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 101 M*/ 102 /*MC 103 TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method. 104 105 This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used. 106 107 Options Database: 108 . -ts_arkimex_type 1bee 109 110 Level: advanced 111 112 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 113 M*/ 114 /*MC 115 TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 116 117 This method has one explicit stage and two implicit stages. The implicit part is the same as in TSARKIMEX2D and TSARKIMEX2E, but the explicit part has a larger stability region on the negative real axis. This method was provided by Emil Constantinescu. 118 119 Options Database: 120 . -ts_arkimex_type 2c 121 122 Level: advanced 123 124 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 125 M*/ 126 /*MC 127 TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 128 129 This method has one explicit stage and two implicit stages. The stability function is independent of the explicit part in the infinity limit of the implict component. This method was provided by Emil Constantinescu. 130 131 Options Database: 132 . -ts_arkimex_type 2d 133 134 Level: advanced 135 136 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 137 M*/ 138 /*MC 139 TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 140 141 This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 142 143 Options Database: 144 . -ts_arkimex_type 2e 145 146 Level: advanced 147 148 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 149 M*/ 150 /*MC 151 TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 152 153 This method has three implicit stages. 154 155 References: 156 . 1. - 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. 157 158 This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375 159 160 Options Database: 161 . -ts_arkimex_type prssp2 162 163 Level: advanced 164 165 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 166 M*/ 167 /*MC 168 TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 169 170 This method has one explicit stage and three implicit stages. 171 172 Options Database: 173 . -ts_arkimex_type 3 174 175 References: 176 . 1. - Kennedy and Carpenter 2003. 177 178 Level: advanced 179 180 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 181 M*/ 182 /*MC 183 TSARKIMEXARS443 - Third order ARK IMEX scheme. 184 185 This method has one explicit stage and four implicit stages. 186 187 Options Database: 188 . -ts_arkimex_type ars443 189 190 References: 191 + 1. - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 192 - 2. - This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375 193 194 Level: advanced 195 196 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 197 M*/ 198 /*MC 199 TSARKIMEXBPR3 - Third order ARK IMEX scheme. 200 201 This method has one explicit stage and four implicit stages. 202 203 Options Database: 204 . -ts_arkimex_type bpr3 205 206 References: 207 . This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375 208 209 Level: advanced 210 211 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 212 M*/ 213 /*MC 214 TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 215 216 This method has one explicit stage and four implicit stages. 217 218 Options Database: 219 . -ts_arkimex_type 4 220 221 References: 222 . Kennedy and Carpenter 2003. 223 224 Level: advanced 225 226 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 227 M*/ 228 /*MC 229 TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 230 231 This method has one explicit stage and five implicit stages. 232 233 Options Database: 234 . -ts_arkimex_type 5 235 236 References: 237 . Kennedy and Carpenter 2003. 238 239 Level: advanced 240 241 .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType() 242 M*/ 243 244 /*@C 245 TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 246 247 Not Collective, but should be called by all processes which will need the schemes to be registered 248 249 Level: advanced 250 251 .keywords: TS, TSARKIMEX, register, all 252 253 .seealso: TSARKIMEXRegisterDestroy() 254 @*/ 255 PetscErrorCode TSARKIMEXRegisterAll(void) 256 { 257 PetscErrorCode ierr; 258 259 PetscFunctionBegin; 260 if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 261 TSARKIMEXRegisterAllCalled = PETSC_TRUE; 262 263 { 264 const PetscReal 265 A[3][3] = {{0.0,0.0,0.0}, 266 {0.0,0.0,0.0}, 267 {0.0,0.5,0.0}}, 268 At[3][3] = {{1.0,0.0,0.0}, 269 {0.0,0.5,0.0}, 270 {0.0,0.5,0.5}}, 271 b[3] = {0.0,0.5,0.5}, 272 bembedt[3] = {1.0,0.0,0.0}; 273 ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 274 } 275 { 276 const PetscReal 277 A[2][2] = {{0.0,0.0}, 278 {0.5,0.0}}, 279 At[2][2] = {{0.0,0.0}, 280 {0.0,0.5}}, 281 b[2] = {0.0,1.0}, 282 bembedt[2] = {0.5,0.5}; 283 /* binterpt[2][2] = {{1.0,-1.0},{0.0,1.0}}; second order dense output has poor stability properties and hence it is not currently in use*/ 284 ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 285 } 286 { 287 const PetscReal 288 A[2][2] = {{0.0,0.0}, 289 {1.0,0.0}}, 290 At[2][2] = {{0.0,0.0}, 291 {0.5,0.5}}, 292 b[2] = {0.5,0.5}, 293 bembedt[2] = {0.0,1.0}; 294 /* binterpt[2][2] = {{1.0,-0.5},{0.0,0.5}} second order dense output has poor stability properties and hence it is not currently in use*/ 295 ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 296 } 297 { 298 /* const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0); Direct evaluation: 0.2928932188134524755992. Used below to ensure all values are available at compile time */ 299 const PetscReal 300 A[2][2] = {{0.0,0.0}, 301 {1.0,0.0}}, 302 At[2][2] = {{0.2928932188134524755992,0.0}, 303 {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}}, 304 b[2] = {0.5,0.5}, 305 bembedt[2] = {0.0,1.0}, 306 binterpt[2][2] = {{ (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}, 307 {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}}, 308 binterp[2][2] = {{1.0,-0.5},{0.0,0.5}}; 309 ierr = TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);CHKERRQ(ierr); 310 } 311 { 312 /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 313 const PetscReal 314 A[3][3] = {{0,0,0}, 315 {2-1.414213562373095048802,0,0}, 316 {0.5,0.5,0}}, 317 At[3][3] = {{0,0,0}, 318 {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 319 {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 320 bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 321 binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 322 {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 323 {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 324 ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 325 } 326 { 327 /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 328 const PetscReal 329 A[3][3] = {{0,0,0}, 330 {2-1.414213562373095048802,0,0}, 331 {0.75,0.25,0}}, 332 At[3][3] = {{0,0,0}, 333 {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 334 {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 335 bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 336 binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 337 {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 338 {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 339 ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 340 } 341 { /* Optimal for linear implicit part */ 342 /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 343 const PetscReal 344 A[3][3] = {{0,0,0}, 345 {2-1.414213562373095048802,0,0}, 346 {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}}, 347 At[3][3] = {{0,0,0}, 348 {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 349 {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 350 bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 351 binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 352 {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 353 {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 354 ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 355 } 356 { /* Optimal for linear implicit part */ 357 const PetscReal 358 A[3][3] = {{0,0,0}, 359 {0.5,0,0}, 360 {0.5,0.5,0}}, 361 At[3][3] = {{0.25,0,0}, 362 {0,0.25,0}, 363 {1./3,1./3,1./3}}; 364 ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr); 365 } 366 { 367 const PetscReal 368 A[4][4] = {{0,0,0,0}, 369 {1767732205903./2027836641118.,0,0,0}, 370 {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 371 {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 372 At[4][4] = {{0,0,0,0}, 373 {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 374 {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 375 {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 376 bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.}, 377 binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 378 {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 379 {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 380 {584795268549./6622622206610., 2508943948391./7218656332882.}}; 381 ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 382 } 383 { 384 const PetscReal 385 A[5][5] = {{0,0,0,0,0}, 386 {1./2,0,0,0,0}, 387 {11./18,1./18,0,0,0}, 388 {5./6,-5./6,.5,0,0}, 389 {1./4,7./4,3./4,-7./4,0}}, 390 At[5][5] = {{0,0,0,0,0}, 391 {0,1./2,0,0,0}, 392 {0,1./6,1./2,0,0}, 393 {0,-1./2,1./2,1./2,0}, 394 {0,3./2,-3./2,1./2,1./2}}, 395 *bembedt = NULL; 396 ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 397 } 398 { 399 const PetscReal 400 A[5][5] = {{0,0,0,0,0}, 401 {1,0,0,0,0}, 402 {4./9,2./9,0,0,0}, 403 {1./4,0,3./4,0,0}, 404 {1./4,0,3./5,0,0}}, 405 At[5][5] = {{0,0,0,0,0}, 406 {.5,.5,0,0,0}, 407 {5./18,-1./9,.5,0,0}, 408 {.5,0,0,.5,0}, 409 {.25,0,.75,-.5,.5}}, 410 *bembedt = NULL; 411 ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 412 } 413 { 414 const PetscReal 415 A[6][6] = {{0,0,0,0,0,0}, 416 {1./2,0,0,0,0,0}, 417 {13861./62500.,6889./62500.,0,0,0,0}, 418 {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 419 {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 420 {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 421 At[6][6] = {{0,0,0,0,0,0}, 422 {1./4,1./4,0,0,0,0}, 423 {8611./62500.,-1743./31250.,1./4,0,0,0}, 424 {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 425 {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 426 {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 427 bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.}, 428 binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 429 {0,0,0}, 430 {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 431 {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 432 {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 433 {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 434 ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 435 } 436 { 437 const PetscReal 438 A[8][8] = {{0,0,0,0,0,0,0,0}, 439 {41./100,0,0,0,0,0,0,0}, 440 {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 441 {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 442 {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 443 {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 444 {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 445 {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 446 At[8][8] = {{0,0,0,0,0,0,0,0}, 447 {41./200.,41./200.,0,0,0,0,0,0}, 448 {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 449 {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 450 {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 451 {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 452 {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 453 {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 454 bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.}, 455 binterpt[8][3] = {{-17674230611817./10670229744614., 43486358583215./12773830924787., -9257016797708./5021505065439.}, 456 {0, 0, 0 }, 457 {0, 0, 0 }, 458 {65168852399939./7868540260826., -91478233927265./11067650958493., 26096422576131./11239449250142.}, 459 {15494834004392./5936557850923., -79368583304911./10890268929626., 92396832856987./20362823103730.}, 460 {-99329723586156./26959484932159., -12239297817655./9152339842473., 30029262896817./10175596800299.}, 461 {-19024464361622./5461577185407., 115839755401235./10719374521269., -26136350496073./3983972220547.}, 462 {-6511271360970./6095937251113., 5843115559534./2180450260947., -5289405421727./3760307252460. }}; 463 ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 464 } 465 PetscFunctionReturn(0); 466 } 467 468 /*@C 469 TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 470 471 Not Collective 472 473 Level: advanced 474 475 .keywords: TSARKIMEX, register, destroy 476 .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll() 477 @*/ 478 PetscErrorCode TSARKIMEXRegisterDestroy(void) 479 { 480 PetscErrorCode ierr; 481 ARKTableauLink link; 482 483 PetscFunctionBegin; 484 while ((link = ARKTableauList)) { 485 ARKTableau t = &link->tab; 486 ARKTableauList = link->next; 487 ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 488 ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr); 489 ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 490 ierr = PetscFree(t->name);CHKERRQ(ierr); 491 ierr = PetscFree(link);CHKERRQ(ierr); 492 } 493 TSARKIMEXRegisterAllCalled = PETSC_FALSE; 494 PetscFunctionReturn(0); 495 } 496 497 /*@C 498 TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 499 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 500 when using static libraries. 501 502 Level: developer 503 504 .keywords: TS, TSARKIMEX, initialize, package 505 .seealso: PetscInitialize() 506 @*/ 507 PetscErrorCode TSARKIMEXInitializePackage(void) 508 { 509 PetscErrorCode ierr; 510 511 PetscFunctionBegin; 512 if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 513 TSARKIMEXPackageInitialized = PETSC_TRUE; 514 ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 515 ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 516 PetscFunctionReturn(0); 517 } 518 519 /*@C 520 TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 521 called from PetscFinalize(). 522 523 Level: developer 524 525 .keywords: Petsc, destroy, package 526 .seealso: PetscFinalize() 527 @*/ 528 PetscErrorCode TSARKIMEXFinalizePackage(void) 529 { 530 PetscErrorCode ierr; 531 532 PetscFunctionBegin; 533 TSARKIMEXPackageInitialized = PETSC_FALSE; 534 ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 535 PetscFunctionReturn(0); 536 } 537 538 /*@C 539 TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 540 541 Not Collective, but the same schemes should be registered on all processes on which they will be used 542 543 Input Parameters: 544 + name - identifier for method 545 . order - approximation order of method 546 . s - number of stages, this is the dimension of the matrices below 547 . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 548 . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 549 . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 550 . A - Non-stiff stage coefficients (dimension s*s, row-major) 551 . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 552 . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 553 . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 554 . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 555 . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 556 . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 557 - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 558 559 Notes: 560 Several ARK IMEX methods are provided, this function is only needed to create new methods. 561 562 Level: advanced 563 564 .keywords: TS, register 565 566 .seealso: TSARKIMEX 567 @*/ 568 PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s, 569 const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 570 const PetscReal A[],const PetscReal b[],const PetscReal c[], 571 const PetscReal bembedt[],const PetscReal bembed[], 572 PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 573 { 574 PetscErrorCode ierr; 575 ARKTableauLink link; 576 ARKTableau t; 577 PetscInt i,j; 578 579 PetscFunctionBegin; 580 ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 581 ierr = PetscCalloc1(1,&link);CHKERRQ(ierr); 582 t = &link->tab; 583 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 584 t->order = order; 585 t->s = s; 586 ierr = PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 587 ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 588 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 589 if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); } 590 else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 591 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 592 else for (i=0; i<s; i++) t->b[i] = t->bt[i]; 593 if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); } 594 else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 595 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 596 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 597 t->stiffly_accurate = PETSC_TRUE; 598 for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 599 t->explicit_first_stage = PETSC_TRUE; 600 for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 601 /*def of FSAL can be made more precise*/ 602 t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 603 if (bembedt) { 604 ierr = PetscMalloc2(s,&t->bembedt,s,&t->bembed);CHKERRQ(ierr); 605 ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr); 606 ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr); 607 } 608 609 t->pinterp = pinterp; 610 ierr = PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);CHKERRQ(ierr); 611 ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 612 ierr = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 613 link->next = ARKTableauList; 614 ARKTableauList = link; 615 PetscFunctionReturn(0); 616 } 617 618 /* 619 The step completion formula is 620 621 x1 = x0 - h bt^T YdotI + h b^T YdotRHS 622 623 This function can be called before or after ts->vec_sol has been updated. 624 Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 625 We can write 626 627 x1e = x0 - h bet^T YdotI + h be^T YdotRHS 628 = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 629 = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 630 631 so we can evaluate the method with different order even after the step has been optimistically completed. 632 */ 633 static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 634 { 635 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 636 ARKTableau tab = ark->tableau; 637 PetscScalar *w = ark->work; 638 PetscReal h; 639 PetscInt s = tab->s,j; 640 PetscErrorCode ierr; 641 642 PetscFunctionBegin; 643 switch (ark->status) { 644 case TS_STEP_INCOMPLETE: 645 case TS_STEP_PENDING: 646 h = ts->time_step; break; 647 case TS_STEP_COMPLETE: 648 h = ts->ptime - ts->ptime_prev; break; 649 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 650 } 651 if (order == tab->order) { 652 if (ark->status == TS_STEP_INCOMPLETE) { 653 if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 654 ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr); 655 } else { /* Use the standard completion formula (bt,b) */ 656 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 657 for (j=0; j<s; j++) w[j] = h*tab->bt[j]; 658 ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 659 if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 660 for (j=0; j<s; j++) w[j] = h*tab->b[j]; 661 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 662 } 663 } 664 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 665 if (done) *done = PETSC_TRUE; 666 PetscFunctionReturn(0); 667 } else if (order == tab->order-1) { 668 if (!tab->bembedt) goto unavailable; 669 if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 670 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 671 for (j=0; j<s; j++) w[j] = h*tab->bembedt[j]; 672 ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 673 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 674 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 675 } else { /* Rollback and re-complete using (bet-be,be-b) */ 676 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 677 for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]); 678 ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr); 679 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 680 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 681 } 682 if (done) *done = PETSC_TRUE; 683 PetscFunctionReturn(0); 684 } 685 unavailable: 686 if (done) *done = PETSC_FALSE; 687 else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 688 PetscFunctionReturn(0); 689 } 690 691 static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 692 { 693 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 694 ARKTableau tab = ark->tableau; 695 const PetscInt s = tab->s; 696 const PetscReal *bt = tab->bt,*b = tab->b; 697 PetscScalar *w = ark->work; 698 Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 699 PetscInt j; 700 PetscReal h; 701 PetscErrorCode ierr; 702 703 PetscFunctionBegin; 704 switch (ark->status) { 705 case TS_STEP_INCOMPLETE: 706 case TS_STEP_PENDING: 707 h = ts->time_step; break; 708 case TS_STEP_COMPLETE: 709 h = ts->ptime - ts->ptime_prev; break; 710 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 711 } 712 for (j=0; j<s; j++) w[j] = -h*bt[j]; 713 ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 714 for (j=0; j<s; j++) w[j] = -h*b[j]; 715 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 static PetscErrorCode TSStep_ARKIMEX(TS ts) 720 { 721 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 722 ARKTableau tab = ark->tableau; 723 const PetscInt s = tab->s; 724 const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 725 PetscScalar *w = ark->work; 726 Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z; 727 PetscBool extrapolate = ark->extrapolate; 728 TSAdapt adapt; 729 SNES snes; 730 PetscInt i,j,its,lits; 731 PetscInt rejections = 0; 732 PetscBool stageok,accept = PETSC_TRUE; 733 PetscReal next_time_step = ts->time_step; 734 PetscErrorCode ierr; 735 736 PetscFunctionBegin; 737 if (ark->extrapolate && !ark->Y_prev) { 738 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 739 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 740 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 741 } 742 743 if (!ts->steprollback) { 744 if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 745 ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr); 746 } 747 if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 748 for (i = 0; i<s; i++) { 749 ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr); 750 ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr); 751 ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr); 752 } 753 } 754 } 755 756 if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 757 TS ts_start; 758 ierr = TSClone(ts,&ts_start);CHKERRQ(ierr); 759 ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr); 760 ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr); 761 ierr = TSSetMaxSteps(ts_start,ts->steps+1);CHKERRQ(ierr); 762 ierr = TSSetMaxTime(ts_start,ts->ptime+ts->time_step);CHKERRQ(ierr); 763 ierr = TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 764 ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr); 765 ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr); 766 ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr); 767 ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr); 768 769 ierr = TSRestartStep(ts_start);CHKERRQ(ierr); 770 ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr); 771 ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr); 772 ierr = TSGetTimeStep(ts_start,&ts->time_step);CHKERRQ(ierr); 773 774 { /* Save the initial slope for the next step */ 775 TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data; 776 ierr = VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);CHKERRQ(ierr); 777 } 778 ts->steps++; 779 780 /* Set the correct TS in SNES */ 781 /* We'll try to bypass this by changing the method on the fly */ 782 { 783 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 784 ierr = TSSetSNES(ts,snes);CHKERRQ(ierr); 785 } 786 ierr = TSDestroy(&ts_start);CHKERRQ(ierr); 787 } 788 789 ark->status = TS_STEP_INCOMPLETE; 790 while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 791 PetscReal t = ts->ptime; 792 PetscReal h = ts->time_step; 793 for (i=0; i<s; i++) { 794 ark->stage_time = t + h*ct[i]; 795 ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr); 796 if (At[i*s+i] == 0) { /* This stage is explicit */ 797 if (i!=0 && ts->equation_type >= TS_EQ_IMPLICIT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Explicit stages other than the first one are not supported for implicit problems"); 798 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 799 for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 800 ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 801 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 802 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 803 } else { 804 ark->scoeff = 1./At[i*s+i]; 805 /* Ydot = shift*(Y-Z) */ 806 ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 807 for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 808 ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 809 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 810 ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr); 811 if (extrapolate && !ts->steprestart) { 812 /* Initial guess extrapolated from previous time step stage values */ 813 ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr); 814 } else { 815 /* Initial guess taken from last stage */ 816 ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr); 817 } 818 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 819 ierr = SNESSolve(snes,NULL,Y[i]);CHKERRQ(ierr); 820 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 821 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 822 ts->snes_its += its; ts->ksp_its += lits; 823 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 824 ierr = TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);CHKERRQ(ierr); 825 if (!stageok) { 826 /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 827 * use extrapolation to initialize the solves on the next attempt. */ 828 extrapolate = PETSC_FALSE; 829 goto reject_step; 830 } 831 } 832 if (ts->equation_type >= TS_EQ_IMPLICIT) { 833 if (i==0 && tab->explicit_first_stage) { 834 if (!tab->stiffly_accurate ) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated",ark->tableau->name); 835 ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr); /* YdotI = YdotI(tn-1) */ 836 } else { 837 ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 838 } 839 } else { 840 if (i==0 && tab->explicit_first_stage) { 841 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 842 ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);/* YdotI = -G(t,Y,0) */ 843 ierr = VecScale(YdotI[i],-1.0);CHKERRQ(ierr); 844 } else { 845 ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 846 } 847 if (ark->imex) { 848 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 849 } else { 850 ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 851 } 852 } 853 ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr); 854 } 855 856 ark->status = TS_STEP_INCOMPLETE; 857 ierr = TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 858 ark->status = TS_STEP_PENDING; 859 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 860 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 861 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 862 ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 863 ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 864 if (!accept) { /* Roll back the current step */ 865 ierr = TSRollBack_ARKIMEX(ts);CHKERRQ(ierr); 866 ts->time_step = next_time_step; 867 goto reject_step; 868 } 869 870 ts->ptime += ts->time_step; 871 ts->time_step = next_time_step; 872 break; 873 874 reject_step: 875 ts->reject++; accept = PETSC_FALSE; 876 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 877 ts->reason = TS_DIVERGED_STEP_REJECTED; 878 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 879 } 880 } 881 PetscFunctionReturn(0); 882 } 883 884 static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 885 { 886 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 887 PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 888 PetscReal h; 889 PetscReal tt,t; 890 PetscScalar *bt,*b; 891 const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 892 PetscErrorCode ierr; 893 894 PetscFunctionBegin; 895 if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 896 switch (ark->status) { 897 case TS_STEP_INCOMPLETE: 898 case TS_STEP_PENDING: 899 h = ts->time_step; 900 t = (itime - ts->ptime)/h; 901 break; 902 case TS_STEP_COMPLETE: 903 h = ts->ptime - ts->ptime_prev; 904 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 905 break; 906 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 907 } 908 ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 909 for (i=0; i<s; i++) bt[i] = b[i] = 0; 910 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 911 for (i=0; i<s; i++) { 912 bt[i] += h * Bt[i*pinterp+j] * tt; 913 b[i] += h * B[i*pinterp+j] * tt; 914 } 915 } 916 ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 917 ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 918 ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 919 ierr = PetscFree2(bt,b);CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 923 static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 924 { 925 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 926 PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 927 PetscReal h,h_prev,t,tt; 928 PetscScalar *bt,*b; 929 const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 930 PetscErrorCode ierr; 931 932 PetscFunctionBegin; 933 if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 934 ierr = PetscCalloc2(s,&bt,s,&b);CHKERRQ(ierr); 935 h = ts->time_step; 936 h_prev = ts->ptime - ts->ptime_prev; 937 t = 1 + h/h_prev*c; 938 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 939 for (i=0; i<s; i++) { 940 bt[i] += h * Bt[i*pinterp+j] * tt; 941 b[i] += h * B[i*pinterp+j] * tt; 942 } 943 } 944 if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 945 ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr); 946 ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr); 947 ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr); 948 ierr = PetscFree2(bt,b);CHKERRQ(ierr); 949 PetscFunctionReturn(0); 950 } 951 952 /*------------------------------------------------------------*/ 953 954 static PetscErrorCode TSARKIMEXTableauReset(TS ts) 955 { 956 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 957 ARKTableau tab = ark->tableau; 958 PetscErrorCode ierr; 959 960 PetscFunctionBegin; 961 if (!tab) PetscFunctionReturn(0); 962 ierr = PetscFree(ark->work);CHKERRQ(ierr); 963 ierr = VecDestroyVecs(tab->s,&ark->Y);CHKERRQ(ierr); 964 ierr = VecDestroyVecs(tab->s,&ark->YdotI);CHKERRQ(ierr); 965 ierr = VecDestroyVecs(tab->s,&ark->YdotRHS);CHKERRQ(ierr); 966 ierr = VecDestroyVecs(tab->s,&ark->Y_prev);CHKERRQ(ierr); 967 ierr = VecDestroyVecs(tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 968 ierr = VecDestroyVecs(tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 969 PetscFunctionReturn(0); 970 } 971 972 static PetscErrorCode TSReset_ARKIMEX(TS ts) 973 { 974 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 975 PetscErrorCode ierr; 976 977 PetscFunctionBegin; 978 ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr); 979 ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 980 ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr); 981 ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 982 PetscFunctionReturn(0); 983 } 984 985 static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 986 { 987 TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 988 PetscErrorCode ierr; 989 990 PetscFunctionBegin; 991 if (Z) { 992 if (dm && dm != ts->dm) { 993 ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 994 } else *Z = ax->Z; 995 } 996 if (Ydot) { 997 if (dm && dm != ts->dm) { 998 ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 999 } else *Ydot = ax->Ydot; 1000 } 1001 PetscFunctionReturn(0); 1002 } 1003 1004 1005 static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1006 { 1007 PetscErrorCode ierr; 1008 1009 PetscFunctionBegin; 1010 if (Z) { 1011 if (dm && dm != ts->dm) { 1012 ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1013 } 1014 } 1015 if (Ydot) { 1016 if (dm && dm != ts->dm) { 1017 ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1018 } 1019 } 1020 PetscFunctionReturn(0); 1021 } 1022 1023 /* 1024 This defines the nonlinear equation that is to be solved with SNES 1025 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 1026 */ 1027 static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 1028 { 1029 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1030 DM dm,dmsave; 1031 Vec Z,Ydot; 1032 PetscReal shift = ark->scoeff / ts->time_step; 1033 PetscErrorCode ierr; 1034 1035 PetscFunctionBegin; 1036 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1037 ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1038 ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 1039 dmsave = ts->dm; 1040 ts->dm = dm; 1041 1042 ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr); 1043 1044 ts->dm = dmsave; 1045 ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 1050 { 1051 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1052 DM dm,dmsave; 1053 Vec Ydot; 1054 PetscReal shift = ark->scoeff / ts->time_step; 1055 PetscErrorCode ierr; 1056 1057 PetscFunctionBegin; 1058 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1059 ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1060 /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1061 dmsave = ts->dm; 1062 ts->dm = dm; 1063 1064 ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr); 1065 1066 ts->dm = dmsave; 1067 ierr = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1068 PetscFunctionReturn(0); 1069 } 1070 1071 static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1072 { 1073 PetscFunctionBegin; 1074 PetscFunctionReturn(0); 1075 } 1076 1077 static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1078 { 1079 TS ts = (TS)ctx; 1080 PetscErrorCode ierr; 1081 Vec Z,Z_c; 1082 1083 PetscFunctionBegin; 1084 ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 1085 ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1086 ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 1087 ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 1088 ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 1089 ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1090 PetscFunctionReturn(0); 1091 } 1092 1093 1094 static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1095 { 1096 PetscFunctionBegin; 1097 PetscFunctionReturn(0); 1098 } 1099 1100 static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1101 { 1102 TS ts = (TS)ctx; 1103 PetscErrorCode ierr; 1104 Vec Z,Z_c; 1105 1106 PetscFunctionBegin; 1107 ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 1108 ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1109 1110 ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1111 ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1112 1113 ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 1114 ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1115 PetscFunctionReturn(0); 1116 } 1117 1118 static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 1119 { 1120 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1121 ARKTableau tab = ark->tableau; 1122 PetscErrorCode ierr; 1123 1124 PetscFunctionBegin; 1125 ierr = PetscMalloc1(tab->s,&ark->work);CHKERRQ(ierr); 1126 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);CHKERRQ(ierr); 1127 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);CHKERRQ(ierr); 1128 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);CHKERRQ(ierr); 1129 if (ark->extrapolate) { 1130 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 1131 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 1132 ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 1133 } 1134 PetscFunctionReturn(0); 1135 } 1136 1137 static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 1138 { 1139 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1140 PetscErrorCode ierr; 1141 DM dm; 1142 SNES snes; 1143 1144 PetscFunctionBegin; 1145 ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr); 1146 ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 1147 ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr); 1148 ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 1149 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1150 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1151 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1152 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1153 PetscFunctionReturn(0); 1154 } 1155 /*------------------------------------------------------------*/ 1156 1157 static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 1158 { 1159 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1160 PetscErrorCode ierr; 1161 1162 PetscFunctionBegin; 1163 ierr = PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");CHKERRQ(ierr); 1164 { 1165 ARKTableauLink link; 1166 PetscInt count,choice; 1167 PetscBool flg; 1168 const char **namelist; 1169 for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 1170 ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 1171 for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 1172 ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);CHKERRQ(ierr); 1173 if (flg) {ierr = TSARKIMEXSetType(ts,namelist[choice]);CHKERRQ(ierr);} 1174 ierr = PetscFree(namelist);CHKERRQ(ierr); 1175 1176 flg = (PetscBool) !ark->imex; 1177 ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr); 1178 ark->imex = (PetscBool) !flg; 1179 ierr = PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->extrapolate,&ark->extrapolate,NULL);CHKERRQ(ierr); 1180 } 1181 ierr = PetscOptionsTail();CHKERRQ(ierr); 1182 PetscFunctionReturn(0); 1183 } 1184 1185 static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 1186 { 1187 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1188 PetscBool iascii; 1189 PetscErrorCode ierr; 1190 1191 PetscFunctionBegin; 1192 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1193 if (iascii) { 1194 ARKTableau tab = ark->tableau; 1195 TSARKIMEXType arktype; 1196 char buf[512]; 1197 ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 1198 ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 1199 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 1200 ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 1201 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1202 ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 1203 ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr); 1204 ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr); 1205 ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 1206 } 1207 PetscFunctionReturn(0); 1208 } 1209 1210 static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1211 { 1212 PetscErrorCode ierr; 1213 SNES snes; 1214 TSAdapt adapt; 1215 1216 PetscFunctionBegin; 1217 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 1218 ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1219 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1220 ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 1221 /* function and Jacobian context for SNES when used with TS is always ts object */ 1222 ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 1223 ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 1224 PetscFunctionReturn(0); 1225 } 1226 1227 /*@C 1228 TSARKIMEXSetType - Set the type of ARK IMEX scheme 1229 1230 Logically collective 1231 1232 Input Parameter: 1233 + ts - timestepping context 1234 - arktype - type of ARK-IMEX scheme 1235 1236 Options Database: 1237 . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> 1238 1239 Level: intermediate 1240 1241 .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, 1242 TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 1243 @*/ 1244 PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 1245 { 1246 PetscErrorCode ierr; 1247 1248 PetscFunctionBegin; 1249 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1250 PetscValidCharPointer(arktype,2); 1251 ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 1252 PetscFunctionReturn(0); 1253 } 1254 1255 /*@C 1256 TSARKIMEXGetType - Get the type of ARK IMEX scheme 1257 1258 Logically collective 1259 1260 Input Parameter: 1261 . ts - timestepping context 1262 1263 Output Parameter: 1264 . arktype - type of ARK-IMEX scheme 1265 1266 Level: intermediate 1267 1268 .seealso: TSARKIMEXGetType() 1269 @*/ 1270 PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 1271 { 1272 PetscErrorCode ierr; 1273 1274 PetscFunctionBegin; 1275 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1276 ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 1277 PetscFunctionReturn(0); 1278 } 1279 1280 /*@ 1281 TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 1282 1283 Logically collective 1284 1285 Input Parameter: 1286 + ts - timestepping context 1287 - flg - PETSC_TRUE for fully implicit 1288 1289 Level: intermediate 1290 1291 .seealso: TSARKIMEXGetType() 1292 @*/ 1293 PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 1294 { 1295 PetscErrorCode ierr; 1296 1297 PetscFunctionBegin; 1298 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1299 ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1300 PetscFunctionReturn(0); 1301 } 1302 1303 static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 1304 { 1305 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1306 1307 PetscFunctionBegin; 1308 *arktype = ark->tableau->name; 1309 PetscFunctionReturn(0); 1310 } 1311 static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype) 1312 { 1313 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1314 PetscErrorCode ierr; 1315 PetscBool match; 1316 ARKTableauLink link; 1317 1318 PetscFunctionBegin; 1319 if (ark->tableau) { 1320 ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 1321 if (match) PetscFunctionReturn(0); 1322 } 1323 for (link = ARKTableauList; link; link=link->next) { 1324 ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 1325 if (match) { 1326 if (ts->setupcalled) {ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr);} 1327 ark->tableau = &link->tab; 1328 if (ts->setupcalled) {ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr);} 1329 ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1330 PetscFunctionReturn(0); 1331 } 1332 } 1333 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 1334 PetscFunctionReturn(0); 1335 } 1336 1337 static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 1338 { 1339 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1340 1341 PetscFunctionBegin; 1342 ark->imex = (PetscBool)!flg; 1343 PetscFunctionReturn(0); 1344 } 1345 1346 static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1347 { 1348 PetscErrorCode ierr; 1349 1350 PetscFunctionBegin; 1351 ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 1352 if (ts->dm) { 1353 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1354 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1355 } 1356 ierr = PetscFree(ts->data);CHKERRQ(ierr); 1357 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr); 1358 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr); 1359 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 1360 PetscFunctionReturn(0); 1361 } 1362 1363 /* ------------------------------------------------------------ */ 1364 /*MC 1365 TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes 1366 1367 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1368 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1369 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1370 1371 Notes: 1372 The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1373 1374 If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 1375 1376 Methods with an explicit stage can only be used with ODE in which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 1377 1378 Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1379 1380 Level: beginner 1381 1382 .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX1BEE, 1383 TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122, 1384 TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister() 1385 1386 M*/ 1387 PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 1388 { 1389 TS_ARKIMEX *th; 1390 PetscErrorCode ierr; 1391 1392 PetscFunctionBegin; 1393 ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 1394 1395 ts->ops->reset = TSReset_ARKIMEX; 1396 ts->ops->destroy = TSDestroy_ARKIMEX; 1397 ts->ops->view = TSView_ARKIMEX; 1398 ts->ops->load = TSLoad_ARKIMEX; 1399 ts->ops->setup = TSSetUp_ARKIMEX; 1400 ts->ops->step = TSStep_ARKIMEX; 1401 ts->ops->interpolate = TSInterpolate_ARKIMEX; 1402 ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 1403 ts->ops->rollback = TSRollBack_ARKIMEX; 1404 ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 1405 ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 1406 ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 1407 1408 ts->usessnes = PETSC_TRUE; 1409 1410 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1411 ts->data = (void*)th; 1412 th->imex = PETSC_TRUE; 1413 1414 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 1415 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 1416 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 1417 1418 ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 1419 PetscFunctionReturn(0); 1420 } 1421