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