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