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