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