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 = PetscCalloc1(1,&link);CHKERRQ(ierr); 557 t = &link->tab; 558 ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 559 t->order = order; 560 t->s = s; 561 ierr = PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 562 ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 563 ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 564 if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); } 565 else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 566 if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 567 else for (i=0; i<s; i++) t->b[i] = t->bt[i]; 568 if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); } 569 else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 570 if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 571 else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 572 t->stiffly_accurate = PETSC_TRUE; 573 for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 574 t->explicit_first_stage = PETSC_TRUE; 575 for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 576 /*def of FSAL can be made more precise*/ 577 t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 578 if (bembedt) { 579 ierr = PetscMalloc2(s,&t->bembedt,s,&t->bembed);CHKERRQ(ierr); 580 ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr); 581 ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr); 582 } 583 584 t->pinterp = pinterp; 585 ierr = PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);CHKERRQ(ierr); 586 ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 587 ierr = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 588 link->next = ARKTableauList; 589 ARKTableauList = link; 590 PetscFunctionReturn(0); 591 } 592 593 #undef __FUNCT__ 594 #define __FUNCT__ "TSEvaluateStep_ARKIMEX" 595 /* 596 The step completion formula is 597 598 x1 = x0 - h bt^T YdotI + h b^T YdotRHS 599 600 This function can be called before or after ts->vec_sol has been updated. 601 Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 602 We can write 603 604 x1e = x0 - h bet^T YdotI + h be^T YdotRHS 605 = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 606 = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 607 608 so we can evaluate the method with different order even after the step has been optimistically completed. 609 */ 610 static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 611 { 612 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 613 ARKTableau tab = ark->tableau; 614 PetscScalar *w = ark->work; 615 PetscReal h; 616 PetscInt s = tab->s,j; 617 PetscErrorCode ierr; 618 619 PetscFunctionBegin; 620 switch (ark->status) { 621 case TS_STEP_INCOMPLETE: 622 case TS_STEP_PENDING: 623 h = ts->time_step; break; 624 case TS_STEP_COMPLETE: 625 h = ts->time_step_prev; break; 626 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 627 } 628 if (order == tab->order) { 629 if (ark->status == TS_STEP_INCOMPLETE) { 630 if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 631 ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr); 632 } else { /* Use the standard completion formula (bt,b) */ 633 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 634 for (j=0; j<s; j++) w[j] = h*tab->bt[j]; 635 ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 636 if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 637 for (j=0; j<s; j++) w[j] = h*tab->b[j]; 638 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 639 } 640 } 641 } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 642 if (done) *done = PETSC_TRUE; 643 PetscFunctionReturn(0); 644 } else if (order == tab->order-1) { 645 if (!tab->bembedt) goto unavailable; 646 if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 647 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 648 for (j=0; j<s; j++) w[j] = h*tab->bembedt[j]; 649 ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 650 for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 651 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 652 } else { /* Rollback and re-complete using (bet-be,be-b) */ 653 ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 654 for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]); 655 ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr); 656 for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 657 ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 658 } 659 if (done) *done = PETSC_TRUE; 660 PetscFunctionReturn(0); 661 } 662 unavailable: 663 if (done) *done = PETSC_FALSE; 664 else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 665 PetscFunctionReturn(0); 666 } 667 668 #undef __FUNCT__ 669 #define __FUNCT__ "TSRollBack_ARKIMEX" 670 static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 671 { 672 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 673 ARKTableau tab = ark->tableau; 674 const PetscInt s = tab->s; 675 const PetscReal *bt = tab->bt,*b = tab->b; 676 PetscScalar *w = ark->work; 677 Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 678 PetscInt j; 679 PetscReal h=ts->time_step; 680 PetscErrorCode ierr; 681 682 PetscFunctionBegin; 683 for (j=0; j<s; j++) w[j] = -h*bt[j]; 684 ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 685 for (j=0; j<s; j++) w[j] = -h*b[j]; 686 ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 687 ark->status = TS_STEP_INCOMPLETE; 688 PetscFunctionReturn(0); 689 } 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "TSStep_ARKIMEX" 693 static PetscErrorCode TSStep_ARKIMEX(TS ts) 694 { 695 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 696 ARKTableau tab = ark->tableau; 697 const PetscInt s = tab->s; 698 const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 699 PetscScalar *w = ark->work; 700 Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,W = ark->Work,Z = ark->Z; 701 PetscBool init_guess_extrp = ark->init_guess_extrp; 702 TSAdapt adapt; 703 SNES snes; 704 PetscInt i,j,its,lits,reject,next_scheme; 705 PetscReal t; 706 PetscReal next_time_step; 707 PetscBool accept; 708 PetscErrorCode ierr; 709 710 PetscFunctionBegin; 711 if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage) { 712 PetscReal valid_time; 713 PetscBool isvalid; 714 ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,valid_time,isvalid);CHKERRQ(ierr); 715 if (!isvalid || valid_time != ts->ptime) { 716 TS ts_start; 717 SNES snes_start; 718 DM dm; 719 PetscReal atol; 720 Vec vatol; 721 PetscReal rtol; 722 Vec vrtol; 723 724 ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts_start);CHKERRQ(ierr); 725 ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr); 726 ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr); 727 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 728 ierr = TSSetDM(ts_start,dm);CHKERRQ(ierr); 729 730 ts_start->adapt=ts->adapt; 731 PetscObjectReference((PetscObject)ts_start->adapt); 732 733 ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr); 734 ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr); 735 ierr = TSSetDuration(ts_start,1,ts->ptime+ts->time_step);CHKERRQ(ierr); 736 ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr); 737 ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr); 738 ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr); 739 ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr); 740 ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr); 741 ierr = TSGetTolerances(ts,&atol,&vatol,&rtol,&vrtol);CHKERRQ(ierr); 742 ierr = TSSetTolerances(ts_start,atol,vatol,rtol,vrtol);CHKERRQ(ierr); 743 ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr); 744 ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr); 745 746 ts->time_step = ts_start->time_step; 747 ts->steps++; 748 ierr = VecCopy(((TS_ARKIMEX*)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr); 749 ts_start->snes=NULL; 750 ierr = TSSetSNES(ts,snes_start);CHKERRQ(ierr); 751 ierr = SNESDestroy(&snes_start);CHKERRQ(ierr); 752 ierr = TSDestroy(&ts_start);CHKERRQ(ierr); 753 } 754 } 755 756 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 757 t = ts->ptime; 758 next_time_step = ts->time_step; 759 accept = PETSC_TRUE; 760 ark->status = TS_STEP_INCOMPLETE; 761 762 763 for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 764 PetscReal h = ts->time_step; 765 ierr = TSPreStep(ts);CHKERRQ(ierr); 766 for (i=0; i<s; i++) { 767 ark->stage_time = t + h*ct[i]; 768 if (At[i*s+i] == 0) { /* This stage is explicit */ 769 ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 770 for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 771 ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 772 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 773 ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 774 } else { 775 ark->scoeff = 1./At[i*s+i]; 776 ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr); 777 /* Affine part */ 778 ierr = VecZeroEntries(W);CHKERRQ(ierr); 779 for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 780 ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 781 ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr); 782 783 /* Ydot = shift*(Y-Z) */ 784 ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 785 for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 786 ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 787 788 if (init_guess_extrp && ark->prev_step_valid) { 789 /* Initial guess extrapolated from previous time step stage values */ 790 ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr); 791 } else { 792 /* Initial guess taken from last stage */ 793 ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr); 794 } 795 ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 796 ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 797 ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 798 ts->snes_its += its; ts->ksp_its += lits; 799 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 800 ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 801 if (!accept) { 802 /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 803 * use extrapolation to initialize the solves on the next attempt. */ 804 ark->prev_step_valid = PETSC_FALSE; 805 goto reject_step; 806 } 807 } 808 ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr); 809 if (ts->equation_type>=TS_EQ_IMPLICIT) { 810 if (i==0 && tab->explicit_first_stage) { 811 ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr); 812 } else { 813 ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 814 } 815 } else { 816 ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 817 ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr); 818 ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr); 819 if (ark->imex) { 820 ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 821 } else { 822 ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 823 } 824 } 825 } 826 ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 827 ark->status = TS_STEP_PENDING; 828 829 /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 830 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 831 ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 832 ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 833 ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 834 if (accept) { 835 /* ignore next_scheme for now */ 836 ts->ptime += ts->time_step; 837 ts->time_step = next_time_step; 838 ts->steps++; 839 if (ts->equation_type>=TS_EQ_IMPLICIT) { /* save the initial slope for the next step*/ 840 ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr); 841 } 842 ark->status = TS_STEP_COMPLETE; 843 if (tab->explicit_first_stage) { 844 ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 845 } 846 /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 847 if (ark->init_guess_extrp) { 848 for (i = 0; i<s; i++) { 849 ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr); 850 ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr); 851 ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr); 852 } 853 ark->prev_step_valid = PETSC_TRUE; 854 } 855 break; 856 } else { /* Roll back the current step */ 857 ts->ptime += next_time_step; /* This will be undone in rollback */ 858 ark->status = TS_STEP_INCOMPLETE; 859 ierr = TSRollBack(ts);CHKERRQ(ierr); 860 } 861 reject_step: continue; 862 } 863 if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 864 PetscFunctionReturn(0); 865 } 866 867 #undef __FUNCT__ 868 #define __FUNCT__ "TSInterpolate_ARKIMEX" 869 static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 870 { 871 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 872 PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 873 PetscReal h; 874 PetscReal tt,t; 875 PetscScalar *bt,*b; 876 const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 877 PetscErrorCode ierr; 878 879 PetscFunctionBegin; 880 if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 881 switch (ark->status) { 882 case TS_STEP_INCOMPLETE: 883 case TS_STEP_PENDING: 884 h = ts->time_step; 885 t = (itime - ts->ptime)/h; 886 break; 887 case TS_STEP_COMPLETE: 888 h = ts->time_step_prev; 889 t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 890 break; 891 default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 892 } 893 ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 894 for (i=0; i<s; i++) bt[i] = b[i] = 0; 895 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 896 for (i=0; i<s; i++) { 897 bt[i] += h * Bt[i*pinterp+j] * tt; 898 b[i] += h * B[i*pinterp+j] * tt; 899 } 900 } 901 ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 902 ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 903 ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 904 ierr = PetscFree2(bt,b);CHKERRQ(ierr); 905 PetscFunctionReturn(0); 906 } 907 908 #undef __FUNCT__ 909 #define __FUNCT__ "TSExtrapolate_ARKIMEX" 910 static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 911 { 912 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 913 PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 914 PetscReal h; 915 PetscReal tt,t; 916 PetscScalar *bt,*b; 917 const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 918 PetscErrorCode ierr; 919 920 PetscFunctionBegin; 921 if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 922 t = 1.0 + (ts->time_step/ts->time_step_prev)*c; 923 h = ts->time_step; 924 ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 925 for (i=0; i<s; i++) bt[i] = b[i] = 0; 926 for (j=0,tt=t; j<pinterp; j++,tt*=t) { 927 for (i=0; i<s; i++) { 928 bt[i] += h * Bt[i*pinterp+j] * tt; 929 b[i] += h * B[i*pinterp+j] * tt; 930 } 931 } 932 if (!ark->prev_step_valid) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 933 ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr); 934 ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr); 935 ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr); 936 ierr = PetscFree2(bt,b);CHKERRQ(ierr); 937 PetscFunctionReturn(0); 938 } 939 940 /*------------------------------------------------------------*/ 941 #undef __FUNCT__ 942 #define __FUNCT__ "TSReset_ARKIMEX" 943 static PetscErrorCode TSReset_ARKIMEX(TS ts) 944 { 945 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 946 PetscInt s; 947 PetscErrorCode ierr; 948 949 PetscFunctionBegin; 950 if (!ark->tableau) PetscFunctionReturn(0); 951 s = ark->tableau->s; 952 ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 953 ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 954 ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 955 if (&ark->init_guess_extrp) { 956 ierr = VecDestroyVecs(s,&ark->Y_prev);CHKERRQ(ierr); 957 ierr = VecDestroyVecs(s,&ark->YdotI_prev);CHKERRQ(ierr); 958 ierr = VecDestroyVecs(s,&ark->YdotRHS_prev);CHKERRQ(ierr); 959 } 960 ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 961 ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 962 ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr); 963 ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 964 ierr = PetscFree(ark->work);CHKERRQ(ierr); 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "TSDestroy_ARKIMEX" 970 static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 971 { 972 PetscErrorCode ierr; 973 974 PetscFunctionBegin; 975 ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 976 ierr = PetscFree(ts->data);CHKERRQ(ierr); 977 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr); 978 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr); 979 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 980 PetscFunctionReturn(0); 981 } 982 983 984 #undef __FUNCT__ 985 #define __FUNCT__ "TSARKIMEXGetVecs" 986 static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 987 { 988 TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 989 PetscErrorCode ierr; 990 991 PetscFunctionBegin; 992 if (Z) { 993 if (dm && dm != ts->dm) { 994 ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 995 } else *Z = ax->Z; 996 } 997 if (Ydot) { 998 if (dm && dm != ts->dm) { 999 ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1000 } else *Ydot = ax->Ydot; 1001 } 1002 PetscFunctionReturn(0); 1003 } 1004 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "TSARKIMEXRestoreVecs" 1008 static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1009 { 1010 PetscErrorCode ierr; 1011 1012 PetscFunctionBegin; 1013 if (Z) { 1014 if (dm && dm != ts->dm) { 1015 ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1016 } 1017 } 1018 if (Ydot) { 1019 if (dm && dm != ts->dm) { 1020 ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1021 } 1022 } 1023 PetscFunctionReturn(0); 1024 } 1025 1026 /* 1027 This defines the nonlinear equation that is to be solved with SNES 1028 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 1029 */ 1030 #undef __FUNCT__ 1031 #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 1032 static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 1033 { 1034 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1035 DM dm,dmsave; 1036 Vec Z,Ydot; 1037 PetscReal shift = ark->scoeff / ts->time_step; 1038 PetscErrorCode ierr; 1039 1040 PetscFunctionBegin; 1041 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1042 ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1043 ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 1044 dmsave = ts->dm; 1045 ts->dm = dm; 1046 1047 ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr); 1048 1049 ts->dm = dmsave; 1050 ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1051 PetscFunctionReturn(0); 1052 } 1053 1054 #undef __FUNCT__ 1055 #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 1056 static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 1057 { 1058 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1059 DM dm,dmsave; 1060 Vec Ydot; 1061 PetscReal shift = ark->scoeff / ts->time_step; 1062 PetscErrorCode ierr; 1063 1064 PetscFunctionBegin; 1065 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1066 ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1067 /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1068 dmsave = ts->dm; 1069 ts->dm = dm; 1070 1071 ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr); 1072 1073 ts->dm = dmsave; 1074 ierr = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1075 PetscFunctionReturn(0); 1076 } 1077 1078 #undef __FUNCT__ 1079 #define __FUNCT__ "DMCoarsenHook_TSARKIMEX" 1080 static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1081 { 1082 PetscFunctionBegin; 1083 PetscFunctionReturn(0); 1084 } 1085 1086 #undef __FUNCT__ 1087 #define __FUNCT__ "DMRestrictHook_TSARKIMEX" 1088 static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1089 { 1090 TS ts = (TS)ctx; 1091 PetscErrorCode ierr; 1092 Vec Z,Z_c; 1093 1094 PetscFunctionBegin; 1095 ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 1096 ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1097 ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 1098 ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 1099 ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 1100 ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1101 PetscFunctionReturn(0); 1102 } 1103 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "DMSubDomainHook_TSARKIMEX" 1107 static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1108 { 1109 PetscFunctionBegin; 1110 PetscFunctionReturn(0); 1111 } 1112 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX" 1115 static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1116 { 1117 TS ts = (TS)ctx; 1118 PetscErrorCode ierr; 1119 Vec Z,Z_c; 1120 1121 PetscFunctionBegin; 1122 ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 1123 ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1124 1125 ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1126 ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1127 1128 ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 1129 ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1130 PetscFunctionReturn(0); 1131 } 1132 1133 #undef __FUNCT__ 1134 #define __FUNCT__ "TSSetUp_ARKIMEX" 1135 static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 1136 { 1137 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1138 ARKTableau tab; 1139 PetscInt s; 1140 PetscErrorCode ierr; 1141 DM dm; 1142 1143 PetscFunctionBegin; 1144 if (!ark->tableau) { 1145 ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 1146 } 1147 tab = ark->tableau; 1148 s = tab->s; 1149 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 1150 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 1151 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 1152 if (ark->init_guess_extrp) { 1153 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y_prev);CHKERRQ(ierr); 1154 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI_prev);CHKERRQ(ierr); 1155 ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS_prev);CHKERRQ(ierr); 1156 } 1157 ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 1158 ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 1159 ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr); 1160 ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 1161 ierr = PetscMalloc1(s,&ark->work);CHKERRQ(ierr); 1162 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1163 if (dm) { 1164 ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1165 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1166 } 1167 PetscFunctionReturn(0); 1168 } 1169 /*------------------------------------------------------------*/ 1170 1171 #undef __FUNCT__ 1172 #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 1173 static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts) 1174 { 1175 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1176 PetscErrorCode ierr; 1177 char arktype[256]; 1178 1179 PetscFunctionBegin; 1180 ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr); 1181 { 1182 ARKTableauLink link; 1183 PetscInt count,choice; 1184 PetscBool flg; 1185 const char **namelist; 1186 ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr); 1187 for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 1188 ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 1189 for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 1190 ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 1191 ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 1192 ierr = PetscFree(namelist);CHKERRQ(ierr); 1193 flg = (PetscBool) !ark->imex; 1194 ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr); 1195 ark->imex = (PetscBool) !flg; 1196 ark->init_guess_extrp = PETSC_FALSE; 1197 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); 1198 ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 1199 } 1200 ierr = PetscOptionsTail();CHKERRQ(ierr); 1201 PetscFunctionReturn(0); 1202 } 1203 1204 #undef __FUNCT__ 1205 #define __FUNCT__ "PetscFormatRealArray" 1206 static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 1207 { 1208 PetscErrorCode ierr; 1209 PetscInt i; 1210 size_t left,count; 1211 char *p; 1212 1213 PetscFunctionBegin; 1214 for (i=0,p=buf,left=len; i<n; i++) { 1215 ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 1216 if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 1217 left -= count; 1218 p += count; 1219 *p++ = ' '; 1220 } 1221 p[i ? 0 : -1] = 0; 1222 PetscFunctionReturn(0); 1223 } 1224 1225 #undef __FUNCT__ 1226 #define __FUNCT__ "TSView_ARKIMEX" 1227 static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 1228 { 1229 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1230 ARKTableau tab = ark->tableau; 1231 PetscBool iascii; 1232 PetscErrorCode ierr; 1233 TSAdapt adapt; 1234 1235 PetscFunctionBegin; 1236 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1237 if (iascii) { 1238 TSARKIMEXType arktype; 1239 char buf[512]; 1240 ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 1241 ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 1242 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 1243 ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 1244 ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1245 ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 1246 ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr); 1247 ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr); 1248 ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 1249 } 1250 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 1251 ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr); 1252 ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 1253 PetscFunctionReturn(0); 1254 } 1255 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "TSLoad_ARKIMEX" 1258 static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1259 { 1260 PetscErrorCode ierr; 1261 SNES snes; 1262 TSAdapt tsadapt; 1263 1264 PetscFunctionBegin; 1265 ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr); 1266 ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr); 1267 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1268 ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 1269 /* function and Jacobian context for SNES when used with TS is always ts object */ 1270 ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 1271 ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 1272 PetscFunctionReturn(0); 1273 } 1274 1275 #undef __FUNCT__ 1276 #define __FUNCT__ "TSARKIMEXSetType" 1277 /*@C 1278 TSARKIMEXSetType - Set the type of ARK IMEX scheme 1279 1280 Logically collective 1281 1282 Input Parameter: 1283 + ts - timestepping context 1284 - arktype - type of ARK-IMEX scheme 1285 1286 Level: intermediate 1287 1288 .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 1289 @*/ 1290 PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 1291 { 1292 PetscErrorCode ierr; 1293 1294 PetscFunctionBegin; 1295 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1296 ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 1297 PetscFunctionReturn(0); 1298 } 1299 1300 #undef __FUNCT__ 1301 #define __FUNCT__ "TSARKIMEXGetType" 1302 /*@C 1303 TSARKIMEXGetType - Get the type of ARK IMEX scheme 1304 1305 Logically collective 1306 1307 Input Parameter: 1308 . ts - timestepping context 1309 1310 Output Parameter: 1311 . arktype - type of ARK-IMEX scheme 1312 1313 Level: intermediate 1314 1315 .seealso: TSARKIMEXGetType() 1316 @*/ 1317 PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 1318 { 1319 PetscErrorCode ierr; 1320 1321 PetscFunctionBegin; 1322 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1323 ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 1324 PetscFunctionReturn(0); 1325 } 1326 1327 #undef __FUNCT__ 1328 #define __FUNCT__ "TSARKIMEXSetFullyImplicit" 1329 /*@C 1330 TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 1331 1332 Logically collective 1333 1334 Input Parameter: 1335 + ts - timestepping context 1336 - flg - PETSC_TRUE for fully implicit 1337 1338 Level: intermediate 1339 1340 .seealso: TSARKIMEXGetType() 1341 @*/ 1342 PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 1343 { 1344 PetscErrorCode ierr; 1345 1346 PetscFunctionBegin; 1347 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1348 ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 1354 PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 1355 { 1356 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1357 PetscErrorCode ierr; 1358 1359 PetscFunctionBegin; 1360 if (!ark->tableau) { 1361 ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 1362 } 1363 *arktype = ark->tableau->name; 1364 PetscFunctionReturn(0); 1365 } 1366 #undef __FUNCT__ 1367 #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 1368 PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype) 1369 { 1370 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1371 PetscErrorCode ierr; 1372 PetscBool match; 1373 ARKTableauLink link; 1374 1375 PetscFunctionBegin; 1376 if (ark->tableau) { 1377 ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 1378 if (match) PetscFunctionReturn(0); 1379 } 1380 for (link = ARKTableauList; link; link=link->next) { 1381 ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 1382 if (match) { 1383 ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 1384 ark->tableau = &link->tab; 1385 PetscFunctionReturn(0); 1386 } 1387 } 1388 SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 1389 PetscFunctionReturn(0); 1390 } 1391 #undef __FUNCT__ 1392 #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX" 1393 PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 1394 { 1395 TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1396 1397 PetscFunctionBegin; 1398 ark->imex = (PetscBool)!flg; 1399 PetscFunctionReturn(0); 1400 } 1401 1402 /* ------------------------------------------------------------ */ 1403 /*MC 1404 TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes 1405 1406 These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1407 nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1408 of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1409 1410 Notes: 1411 The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1412 1413 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). 1414 1415 Level: beginner 1416 1417 .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3, 1418 TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister() 1419 1420 M*/ 1421 #undef __FUNCT__ 1422 #define __FUNCT__ "TSCreate_ARKIMEX" 1423 PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 1424 { 1425 TS_ARKIMEX *th; 1426 PetscErrorCode ierr; 1427 1428 PetscFunctionBegin; 1429 ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 1430 1431 ts->ops->reset = TSReset_ARKIMEX; 1432 ts->ops->destroy = TSDestroy_ARKIMEX; 1433 ts->ops->view = TSView_ARKIMEX; 1434 ts->ops->load = TSLoad_ARKIMEX; 1435 ts->ops->setup = TSSetUp_ARKIMEX; 1436 ts->ops->step = TSStep_ARKIMEX; 1437 ts->ops->interpolate = TSInterpolate_ARKIMEX; 1438 ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 1439 ts->ops->rollback = TSRollBack_ARKIMEX; 1440 ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 1441 ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 1442 ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 1443 1444 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1445 ts->data = (void*)th; 1446 th->imex = PETSC_TRUE; 1447 1448 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 1449 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 1450 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 1451 PetscFunctionReturn(0); 1452 } 1453