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