18a381b04SJed Brown /* 28a381b04SJed Brown Code for timestepping with additive Runge-Kutta IMEX method 38a381b04SJed Brown 48a381b04SJed Brown Notes: 58a381b04SJed Brown The general system is written as 68a381b04SJed Brown 7f9c1d6abSBarry Smith F(t,U,Udot) = G(t,U) 88a381b04SJed Brown 98a381b04SJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 108a381b04SJed Brown 118a381b04SJed Brown */ 12af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 131e25c274SJed Brown #include <petscdm.h> 148a381b04SJed Brown 1519fd82e9SBarry Smith static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3; 168a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled; 178a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized; 1856dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS, PetscReal, Vec); 198a381b04SJed Brown 208a381b04SJed Brown typedef struct _ARKTableau *ARKTableau; 218a381b04SJed Brown struct _ARKTableau { 228a381b04SJed Brown char *name; 234f385281SJed Brown PetscInt order; /* Classical approximation order of the method */ 244f385281SJed Brown PetscInt s; /* Number of stages */ 25e817cc15SEmil Constantinescu PetscBool stiffly_accurate; /* The implicit part is stiffly accurate*/ 26e817cc15SEmil Constantinescu PetscBool FSAL_implicit; /* The implicit part is FSAL*/ 27e817cc15SEmil Constantinescu PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/ 284f385281SJed Brown PetscInt pinterp; /* Interpolation order */ 294f385281SJed Brown PetscReal *At, *bt, *ct; /* Stiff tableau */ 308a381b04SJed Brown PetscReal *A, *b, *c; /* Non-stiff tableau */ 31108c343cSJed Brown PetscReal *bembedt, *bembed; /* Embedded formula of order one less (order-1) */ 32cd652676SJed Brown PetscReal *binterpt, *binterp; /* Dense output formula */ 33108c343cSJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 348a381b04SJed Brown }; 358a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 368a381b04SJed Brown struct _ARKTableauLink { 378a381b04SJed Brown struct _ARKTableau tab; 388a381b04SJed Brown ARKTableauLink next; 398a381b04SJed Brown }; 408a381b04SJed Brown static ARKTableauLink ARKTableauList; 418a381b04SJed Brown 428a381b04SJed Brown typedef struct { 438a381b04SJed Brown ARKTableau tableau; 448a381b04SJed Brown Vec *Y; /* States computed during the step */ 458a381b04SJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 468a381b04SJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 4756dcabbaSDebojyoti Ghosh Vec *Y_prev; /* States computed during the previous time step */ 4856dcabbaSDebojyoti Ghosh Vec *YdotI_prev; /* Time derivatives for the stiff part for the previous time step*/ 4956dcabbaSDebojyoti Ghosh Vec *YdotRHS_prev; /* Function evaluations for the non-stiff part for the previous time step*/ 50e817cc15SEmil Constantinescu Vec Ydot0; /* Holds the slope from the previous step in FSAL case */ 518a381b04SJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 528a381b04SJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 538a381b04SJed Brown PetscScalar *work; /* Scalar work */ 54b296d7d5SJed Brown PetscReal scoeff; /* shift = scoeff/dt */ 558a381b04SJed Brown PetscReal stage_time; 564cc180ffSJed Brown PetscBool imex; 5796400bd6SLisandro Dalcin PetscBool extrapolate; /* Extrapolate initial guess from previous time-step stage values */ 58108c343cSJed Brown TSStepStatus status; 598a381b04SJed Brown } TS_ARKIMEX; 601f80e275SEmil Constantinescu /*MC 611f80e275SEmil Constantinescu TSARKIMEXARS122 - Second order ARK IMEX scheme. 628a381b04SJed Brown 631f80e275SEmil Constantinescu This method has one explicit stage and one implicit stage. 641f80e275SEmil Constantinescu 65*bcf0153eSBarry Smith Options Database Key: 6667b8a455SSatish Balay . -ts_arkimex_type ars122 - set arkimex type to ars122 679bd3e852SBarry Smith 68*bcf0153eSBarry Smith Level: advanced 69*bcf0153eSBarry Smith 701f80e275SEmil Constantinescu References: 71606c0280SSatish Balay . * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 721f80e275SEmil Constantinescu 73*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 741f80e275SEmil Constantinescu M*/ 751f80e275SEmil Constantinescu /*MC 761f80e275SEmil Constantinescu TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 771f80e275SEmil Constantinescu 781f80e275SEmil Constantinescu This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu. 791f80e275SEmil Constantinescu 80*bcf0153eSBarry Smith Options Database Key: 8167b8a455SSatish Balay . -ts_arkimex_type a2 - set arkimex type to a2 829bd3e852SBarry Smith 831f80e275SEmil Constantinescu Level: advanced 841f80e275SEmil Constantinescu 85*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 861f80e275SEmil Constantinescu M*/ 871f80e275SEmil Constantinescu /*MC 881f80e275SEmil Constantinescu TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 891f80e275SEmil Constantinescu 901f80e275SEmil Constantinescu This method has two implicit stages, and L-stable implicit scheme. 911f80e275SEmil Constantinescu 92*bcf0153eSBarry Smith Options Database Key: 9367b8a455SSatish Balay . -ts_arkimex_type l2 - set arkimex type to l2 949bd3e852SBarry Smith 95*bcf0153eSBarry Smith Level: advanced 96*bcf0153eSBarry Smith 971f80e275SEmil Constantinescu References: 98606c0280SSatish Balay . * - 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. 991f80e275SEmil Constantinescu 100*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1011f80e275SEmil Constantinescu M*/ 1021f80e275SEmil Constantinescu /*MC 103c79dcfbdSBarry Smith TSARKIMEX1BEE - First order backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method. 104e817cc15SEmil Constantinescu 105e817cc15SEmil Constantinescu This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used. 106e817cc15SEmil Constantinescu 107*bcf0153eSBarry Smith Options Database Key: 10867b8a455SSatish Balay . -ts_arkimex_type 1bee - set arkimex type to 1bee 1099bd3e852SBarry Smith 110e817cc15SEmil Constantinescu Level: advanced 111e817cc15SEmil Constantinescu 112*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 113e817cc15SEmil Constantinescu M*/ 114e817cc15SEmil Constantinescu /*MC 1151f80e275SEmil Constantinescu TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 1161f80e275SEmil Constantinescu 1171f80e275SEmil Constantinescu 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. 1181f80e275SEmil Constantinescu 119*bcf0153eSBarry Smith Options Database Key: 12067b8a455SSatish Balay . -ts_arkimex_type 2c - set arkimex type to 2c 1219bd3e852SBarry Smith 1221f80e275SEmil Constantinescu Level: advanced 1231f80e275SEmil Constantinescu 124*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1251f80e275SEmil Constantinescu M*/ 12664f491ddSJed Brown /*MC 12764f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 12864f491ddSJed Brown 129617a39beSEmil Constantinescu 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. 13064f491ddSJed Brown 131*bcf0153eSBarry Smith Options Database Key: 13267b8a455SSatish Balay . -ts_arkimex_type 2d - set arkimex type to 2d 1339bd3e852SBarry Smith 134b330ce4dSSatish Balay Level: advanced 135b330ce4dSSatish Balay 136*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 13764f491ddSJed Brown M*/ 13864f491ddSJed Brown /*MC 13964f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 14064f491ddSJed Brown 14164f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 14264f491ddSJed Brown 143*bcf0153eSBarry Smith Options Database Key: 14467b8a455SSatish Balay . -ts_arkimex_type 2e - set arkimex type to 2e 1459bd3e852SBarry Smith 146b330ce4dSSatish Balay Level: advanced 147b330ce4dSSatish Balay 148*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 14964f491ddSJed Brown M*/ 15064f491ddSJed Brown /*MC 1516cf0794eSJed Brown TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 1526cf0794eSJed Brown 1536cf0794eSJed Brown This method has three implicit stages. 1546cf0794eSJed Brown 1556cf0794eSJed Brown References: 156606c0280SSatish Balay . * - 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. 1576cf0794eSJed Brown 158a8d69d7bSBarry Smith This method is referred to as SSP2-(3,3,2) in https://arxiv.org/abs/1110.4375 1596cf0794eSJed Brown 160*bcf0153eSBarry Smith Options Database Key: 16167b8a455SSatish Balay . -ts_arkimex_type prssp2 - set arkimex type to prssp2 1629bd3e852SBarry Smith 1636cf0794eSJed Brown Level: advanced 1646cf0794eSJed Brown 165*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1666cf0794eSJed Brown M*/ 1676cf0794eSJed Brown /*MC 16864f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 16964f491ddSJed Brown 17064f491ddSJed Brown This method has one explicit stage and three implicit stages. 17164f491ddSJed Brown 172*bcf0153eSBarry Smith Options Database Key: 17367b8a455SSatish Balay . -ts_arkimex_type 3 - set arkimex type to 3 1749bd3e852SBarry Smith 175*bcf0153eSBarry Smith Level: advanced 176*bcf0153eSBarry Smith 17764f491ddSJed Brown References: 178606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 17964f491ddSJed Brown 180*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 18164f491ddSJed Brown M*/ 18264f491ddSJed Brown /*MC 1836cf0794eSJed Brown TSARKIMEXARS443 - Third order ARK IMEX scheme. 1846cf0794eSJed Brown 1856cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1866cf0794eSJed Brown 187*bcf0153eSBarry Smith Options Database Key: 18867b8a455SSatish Balay . -ts_arkimex_type ars443 - set arkimex type to ars443 1899bd3e852SBarry Smith 190*bcf0153eSBarry Smith Level: advanced 191*bcf0153eSBarry Smith 1926cf0794eSJed Brown References: 193606c0280SSatish Balay + * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 194606c0280SSatish Balay - * - This method is referred to as ARS(4,4,3) in https://arxiv.org/abs/1110.4375 1956cf0794eSJed Brown 196*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1976cf0794eSJed Brown M*/ 1986cf0794eSJed Brown /*MC 1996cf0794eSJed Brown TSARKIMEXBPR3 - Third order ARK IMEX scheme. 2006cf0794eSJed Brown 2016cf0794eSJed Brown This method has one explicit stage and four implicit stages. 2026cf0794eSJed Brown 203*bcf0153eSBarry Smith Options Database Key: 20467b8a455SSatish Balay . -ts_arkimex_type bpr3 - set arkimex type to bpr3 2059bd3e852SBarry Smith 206*bcf0153eSBarry Smith Level: advanced 207*bcf0153eSBarry Smith 2086cf0794eSJed Brown References: 209606c0280SSatish Balay . * - This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375 2106cf0794eSJed Brown 211*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 2126cf0794eSJed Brown M*/ 2136cf0794eSJed Brown /*MC 21464f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 21564f491ddSJed Brown 21664f491ddSJed Brown This method has one explicit stage and four implicit stages. 21764f491ddSJed Brown 218*bcf0153eSBarry Smith Options Database Key: 21967b8a455SSatish Balay . -ts_arkimex_type 4 - set arkimex type to4 2209bd3e852SBarry Smith 221*bcf0153eSBarry Smith Level: advanced 222*bcf0153eSBarry Smith 22364f491ddSJed Brown References: 224606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 22564f491ddSJed Brown 226*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 22764f491ddSJed Brown M*/ 22864f491ddSJed Brown /*MC 22964f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 23064f491ddSJed Brown 23164f491ddSJed Brown This method has one explicit stage and five implicit stages. 23264f491ddSJed Brown 233*bcf0153eSBarry Smith Options Database Key: 23467b8a455SSatish Balay . -ts_arkimex_type 5 - set arkimex type to 5 2359bd3e852SBarry Smith 236*bcf0153eSBarry Smith Level: advanced 237*bcf0153eSBarry Smith 23864f491ddSJed Brown References: 239606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 24064f491ddSJed Brown 241*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 24264f491ddSJed Brown M*/ 24364f491ddSJed Brown 2448a381b04SJed Brown /*@C 245*bcf0153eSBarry Smith TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in `TSARKIMEX` 2468a381b04SJed Brown 247fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 2488a381b04SJed Brown 2498a381b04SJed Brown Level: advanced 2508a381b04SJed Brown 251*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSARKIMEX`, `TSARKIMEXRegisterDestroy()` 2528a381b04SJed Brown @*/ 253d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXRegisterAll(void) 254d71ae5a4SJacob Faibussowitsch { 2558a381b04SJed Brown PetscFunctionBegin; 2568a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2578a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 258e817cc15SEmil Constantinescu 259e817cc15SEmil Constantinescu { 2609371c9d4SSatish Balay const PetscReal A[3][3] = 2619371c9d4SSatish Balay { 262e817cc15SEmil Constantinescu {0.0, 0.0, 0.0}, 2639371c9d4SSatish Balay {0.0, 0.0, 0.0}, 2649371c9d4SSatish Balay {0.0, 0.5, 0.0} 2659371c9d4SSatish Balay }, 2669371c9d4SSatish Balay At[3][3] = {{1.0, 0.0, 0.0}, {0.0, 0.5, 0.0}, {0.0, 0.5, 0.5}}, b[3] = {0.0, 0.5, 0.5}, bembedt[3] = {1.0, 0.0, 0.0}; 2679566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX1BEE, 2, 3, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 268e817cc15SEmil Constantinescu } 2698a381b04SJed Brown { 2709371c9d4SSatish Balay const PetscReal A[2][2] = 2719371c9d4SSatish Balay { 2729371c9d4SSatish Balay {0.0, 0.0}, 2739371c9d4SSatish Balay {0.5, 0.0} 2749371c9d4SSatish Balay }, 2759371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.0, 0.5}}, b[2] = {0.0, 1.0}, bembedt[2] = {0.5, 0.5}; 2761f80e275SEmil Constantinescu /* 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*/ 2779566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS122, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2781f80e275SEmil Constantinescu } 2791f80e275SEmil Constantinescu { 2809371c9d4SSatish Balay const PetscReal A[2][2] = 2819371c9d4SSatish Balay { 2829371c9d4SSatish Balay {0.0, 0.0}, 2839371c9d4SSatish Balay {1.0, 0.0} 2849371c9d4SSatish Balay }, 2859371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.5, 0.5}}, b[2] = {0.5, 0.5}, bembedt[2] = {0.0, 1.0}; 2861f80e275SEmil Constantinescu /* 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*/ 2879566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXA2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2881f80e275SEmil Constantinescu } 2891f80e275SEmil Constantinescu { 290da80777bSKarl Rupp /* const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0); Direct evaluation: 0.2928932188134524755992. Used below to ensure all values are available at compile time */ 2919371c9d4SSatish Balay const PetscReal A[2][2] = 2929371c9d4SSatish Balay { 2939371c9d4SSatish Balay {0.0, 0.0}, 2949371c9d4SSatish Balay {1.0, 0.0} 2959371c9d4SSatish Balay }, 2969371c9d4SSatish Balay At[2][2] = {{0.2928932188134524755992, 0.0}, {1.0 - 2.0 * 0.2928932188134524755992, 0.2928932188134524755992}}, b[2] = {0.5, 0.5}, bembedt[2] = {0.0, 1.0}, binterpt[2][2] = {{(0.2928932188134524755992 - 1.0) / (2.0 * 0.2928932188134524755992 - 1.0), -1 / (2.0 * (1.0 - 2.0 * 0.2928932188134524755992))}, {1 - (0.2928932188134524755992 - 1.0) / (2.0 * 0.2928932188134524755992 - 1.0), -1 / (2.0 * (1.0 - 2.0 * 0.2928932188134524755992))}}, binterp[2][2] = {{1.0, -0.5}, {0.0, 0.5}}; 2979566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXL2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 2, binterpt[0], binterp[0])); 2981f80e275SEmil Constantinescu } 2991f80e275SEmil Constantinescu { 300da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3019371c9d4SSatish Balay const PetscReal A[3][3] = 3029371c9d4SSatish Balay { 3039371c9d4SSatish Balay {0, 0, 0}, 304da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3059371c9d4SSatish Balay {0.5, 0.5, 0} 3069371c9d4SSatish Balay }, 3079371c9d4SSatish Balay At[3][3] = {{0, 0, 0}, {1 - 1 / 1.414213562373095048802, 1 - 1 / 1.414213562373095048802, 0}, {1 / (2 * 1.414213562373095048802), 1 / (2 * 1.414213562373095048802), 1 - 1 / 1.414213562373095048802}}, bembedt[3] = {(4. - 1.414213562373095048802) / 8., (4. - 1.414213562373095048802) / 8., 1 / (2. * 1.414213562373095048802)}, binterpt[3][2] = {{1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 - 1.414213562373095048802, 1.0 / 1.414213562373095048802}}; 3089566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2C, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3091f80e275SEmil Constantinescu } 3101f80e275SEmil Constantinescu { 311da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3129371c9d4SSatish Balay const PetscReal A[3][3] = 3139371c9d4SSatish Balay { 3149371c9d4SSatish Balay {0, 0, 0}, 315da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3169371c9d4SSatish Balay {0.75, 0.25, 0} 3179371c9d4SSatish Balay }, 3189371c9d4SSatish Balay At[3][3] = {{0, 0, 0}, {1 - 1 / 1.414213562373095048802, 1 - 1 / 1.414213562373095048802, 0}, {1 / (2 * 1.414213562373095048802), 1 / (2 * 1.414213562373095048802), 1 - 1 / 1.414213562373095048802}}, bembedt[3] = {(4. - 1.414213562373095048802) / 8., (4. - 1.414213562373095048802) / 8., 1 / (2. * 1.414213562373095048802)}, binterpt[3][2] = {{1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 - 1.414213562373095048802, 1.0 / 1.414213562373095048802}}; 3199566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2D, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3208a381b04SJed Brown } 32106db7b1cSJed Brown { /* Optimal for linear implicit part */ 322da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3239371c9d4SSatish Balay const PetscReal A[3][3] = 3249371c9d4SSatish Balay { 3259371c9d4SSatish Balay {0, 0, 0}, 326da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3279371c9d4SSatish Balay {(3 - 2 * 1.414213562373095048802) / 6, (3 + 2 * 1.414213562373095048802) / 6, 0} 3289371c9d4SSatish Balay }, 3299371c9d4SSatish Balay At[3][3] = {{0, 0, 0}, {1 - 1 / 1.414213562373095048802, 1 - 1 / 1.414213562373095048802, 0}, {1 / (2 * 1.414213562373095048802), 1 / (2 * 1.414213562373095048802), 1 - 1 / 1.414213562373095048802}}, bembedt[3] = {(4. - 1.414213562373095048802) / 8., (4. - 1.414213562373095048802) / 8., 1 / (2. * 1.414213562373095048802)}, binterpt[3][2] = {{1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 / 1.414213562373095048802, -1.0 / (2.0 * 1.414213562373095048802)}, {1.0 - 1.414213562373095048802, 1.0 / 1.414213562373095048802}}; 3309566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2E, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 331a3a57f36SJed Brown } 3326cf0794eSJed Brown { /* Optimal for linear implicit part */ 3339371c9d4SSatish Balay const PetscReal A[3][3] = 3349371c9d4SSatish Balay { 3359371c9d4SSatish Balay {0, 0, 0}, 3366cf0794eSJed Brown {0.5, 0, 0}, 3379371c9d4SSatish Balay {0.5, 0.5, 0} 3389371c9d4SSatish Balay }, 3399371c9d4SSatish Balay At[3][3] = {{0.25, 0, 0}, {0, 0.25, 0}, {1. / 3, 1. / 3, 1. / 3}}; 3409566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXPRSSP2, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, NULL, NULL, 0, NULL, NULL)); 3416cf0794eSJed Brown } 342a3a57f36SJed Brown { 3439371c9d4SSatish Balay const PetscReal A[4][4] = 3449371c9d4SSatish Balay { 3459371c9d4SSatish Balay {0, 0, 0, 0}, 3464040e9f2SJed Brown {1767732205903. / 2027836641118., 0, 0, 0}, 3474040e9f2SJed Brown {5535828885825. / 10492691773637., 788022342437. / 10882634858940., 0, 0}, 3489371c9d4SSatish Balay {6485989280629. / 16251701735622., -4246266847089. / 9704473918619., 10755448449292. / 10357097424841., 0} 3499371c9d4SSatish Balay }, 3509371c9d4SSatish Balay At[4][4] = {{0, 0, 0, 0}, {1767732205903. / 4055673282236., 1767732205903. / 4055673282236., 0, 0}, {2746238789719. / 10658868560708., -640167445237. / 6845629431997., 1767732205903. / 4055673282236., 0}, {1471266399579. / 7840856788654., -4482444167858. / 7529755066697., 11266239266428. / 11593286722821., 1767732205903. / 4055673282236.}}, bembedt[4] = {2756255671327. / 12835298489170., -10771552573575. / 22201958757719., 9247589265047. / 10645013368117., 2193209047091. / 5459859503100.}, binterpt[4][2] = {{4655552711362. / 22874653954995., -215264564351. / 13552729205753.}, {-18682724506714. / 9892148508045., 17870216137069. / 13817060693119.}, {34259539580243. / 13192909600954., -28141676662227. / 17317692491321.}, {584795268549. / 6622622206610., 2508943948391. / 7218656332882.}}; 3519566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX3, 3, 4, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 352a3a57f36SJed Brown } 353a3a57f36SJed Brown { 3549371c9d4SSatish Balay const PetscReal A[5][5] = 3559371c9d4SSatish Balay { 3569371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3576cf0794eSJed Brown {1. / 2, 0, 0, 0, 0}, 3586cf0794eSJed Brown {11. / 18, 1. / 18, 0, 0, 0}, 3596cf0794eSJed Brown {5. / 6, -5. / 6, .5, 0, 0}, 3609371c9d4SSatish Balay {1. / 4, 7. / 4, 3. / 4, -7. / 4, 0} 3619371c9d4SSatish Balay }, 3629371c9d4SSatish Balay At[5][5] = {{0, 0, 0, 0, 0}, {0, 1. / 2, 0, 0, 0}, {0, 1. / 6, 1. / 2, 0, 0}, {0, -1. / 2, 1. / 2, 1. / 2, 0}, {0, 3. / 2, -3. / 2, 1. / 2, 1. / 2}}, *bembedt = NULL; 3639566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS443, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3646cf0794eSJed Brown } 3656cf0794eSJed Brown { 3669371c9d4SSatish Balay const PetscReal A[5][5] = 3679371c9d4SSatish Balay { 3689371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3696cf0794eSJed Brown {1, 0, 0, 0, 0}, 3706cf0794eSJed Brown {4. / 9, 2. / 9, 0, 0, 0}, 3716cf0794eSJed Brown {1. / 4, 0, 3. / 4, 0, 0}, 3729371c9d4SSatish Balay {1. / 4, 0, 3. / 5, 0, 0} 3739371c9d4SSatish Balay }, 3749371c9d4SSatish Balay At[5][5] = {{0, 0, 0, 0, 0}, {.5, .5, 0, 0, 0}, {5. / 18, -1. / 9, .5, 0, 0}, {.5, 0, 0, .5, 0}, {.25, 0, .75, -.5, .5}}, *bembedt = NULL; 3759566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXBPR3, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3766cf0794eSJed Brown } 3776cf0794eSJed Brown { 3789371c9d4SSatish Balay const PetscReal A[6][6] = 3799371c9d4SSatish Balay { 3809371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 381a3a57f36SJed Brown {1. / 2, 0, 0, 0, 0, 0}, 3824040e9f2SJed Brown {13861. / 62500., 6889. / 62500., 0, 0, 0, 0}, 3834040e9f2SJed Brown {-116923316275. / 2393684061468., -2731218467317. / 15368042101831., 9408046702089. / 11113171139209., 0, 0, 0}, 3844040e9f2SJed Brown {-451086348788. / 2902428689909., -2682348792572. / 7519795681897., 12662868775082. / 11960479115383., 3355817975965. / 11060851509271., 0, 0}, 3859371c9d4SSatish Balay {647845179188. / 3216320057751., 73281519250. / 8382639484533., 552539513391. / 3454668386233., 3354512671639. / 8306763924573., 4040. / 17871., 0} 3869371c9d4SSatish Balay }, 3879371c9d4SSatish Balay At[6][6] = {{0, 0, 0, 0, 0, 0}, {1. / 4, 1. / 4, 0, 0, 0, 0}, {8611. / 62500., -1743. / 31250., 1. / 4, 0, 0, 0}, {5012029. / 34652500., -654441. / 2922500., 174375. / 388108., 1. / 4, 0, 0}, {15267082809. / 155376265600., -71443401. / 120774400., 730878875. / 902184768., 2285395. / 8070912., 1. / 4, 0}, {82889. / 524892., 0, 15625. / 83664., 69875. / 102672., -2260. / 8211, 1. / 4}}, bembedt[6] = {4586570599. / 29645900160., 0, 178811875. / 945068544., 814220225. / 1159782912., -3700637. / 11593932., 61727. / 225920.}, binterpt[6][3] = {{6943876665148. / 7220017795957., -54480133. / 30881146., 6818779379841. / 7100303317025.}, {0, 0, 0}, 3889371c9d4SSatish Balay {7640104374378. / 9702883013639., -11436875. / 14766696., 2173542590792. / 12501825683035.}, {-20649996744609. / 7521556579894., 174696575. / 18121608., -31592104683404. / 5083833661969.}, 3899371c9d4SSatish Balay {8854892464581. / 2390941311638., -12120380. / 966161., 61146701046299. / 7138195549469.}, {-11397109935349. / 6675773540249., 3843. / 706., -17219254887155. / 4939391667607.}}; 3909566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX4, 4, 6, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 391a3a57f36SJed Brown } 392a3a57f36SJed Brown { 3939371c9d4SSatish Balay const PetscReal A[8][8] = 3949371c9d4SSatish Balay { 3959371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 396a3a57f36SJed Brown {41. / 100, 0, 0, 0, 0, 0, 0, 0}, 3974040e9f2SJed Brown {367902744464. / 2072280473677., 677623207551. / 8224143866563., 0, 0, 0, 0, 0, 0}, 3984040e9f2SJed Brown {1268023523408. / 10340822734521., 0, 1029933939417. / 13636558850479., 0, 0, 0, 0, 0}, 3994040e9f2SJed Brown {14463281900351. / 6315353703477., 0, 66114435211212. / 5879490589093., -54053170152839. / 4284798021562., 0, 0, 0, 0}, 4004040e9f2SJed Brown {14090043504691. / 34967701212078., 0, 15191511035443. / 11219624916014., -18461159152457. / 12425892160975., -281667163811. / 9011619295870., 0, 0, 0}, 4014040e9f2SJed Brown {19230459214898. / 13134317526959., 0, 21275331358303. / 2942455364971., -38145345988419. / 4862620318723., -1. / 8, -1. / 8, 0, 0}, 4029371c9d4SSatish Balay {-19977161125411. / 11928030595625., 0, -40795976796054. / 6384907823539., 177454434618887. / 12078138498510., 782672205425. / 8267701900261., -69563011059811. / 9646580694205., 7356628210526. / 4942186776405., 0} 4039371c9d4SSatish Balay }, 4049371c9d4SSatish Balay At[8][8] = {{0, 0, 0, 0, 0, 0, 0, 0}, {41. / 200., 41. / 200., 0, 0, 0, 0, 0, 0}, {41. / 400., -567603406766. / 11931857230679., 41. / 200., 0, 0, 0, 0, 0}, {683785636431. / 9252920307686., 0, -110385047103. / 1367015193373., 41. / 200., 0, 0, 0, 0}, {3016520224154. / 10081342136671., 0, 30586259806659. / 12414158314087., -22760509404356. / 11113319521817., 41. / 200., 0, 0, 0}, {218866479029. / 1489978393911., 0, 638256894668. / 5436446318841., -1179710474555. / 5321154724896., -60928119172. / 8023461067671., 41. / 200., 0, 0}, {1020004230633. / 5715676835656., 0, 25762820946817. / 25263940353407., -2161375909145. / 9755907335909., -211217309593. / 5846859502534., -4269925059573. / 7827059040749., 41. / 200, 0}, {-872700587467. / 9133579230613., 0, 0, 22348218063261. / 9555858737531., -1143369518992. / 8141816002931., -39379526789629. / 19018526304540., 32727382324388. / 42900044865799., 41. / 200.}}, bembedt[8] = {-975461918565. / 9796059967033., 0, 0, 78070527104295. / 32432590147079., -548382580838. / 3424219808633., -33438840321285. / 15594753105479., 3629800801594. / 4656183773603., 4035322873751. / 18575991585200.}, binterpt[8][3] = {{-17674230611817. / 10670229744614., 43486358583215. / 12773830924787., -9257016797708. / 5021505065439.}, {0, 0, 0}, {0, 0, 0}, {65168852399939. / 7868540260826., -91478233927265. / 11067650958493., 26096422576131. / 11239449250142.}, {15494834004392. / 5936557850923., -79368583304911. / 10890268929626., 92396832856987. / 20362823103730.}, {-99329723586156. / 26959484932159., -12239297817655. / 9152339842473., 30029262896817. / 10175596800299.}, {-19024464361622. / 5461577185407., 115839755401235. / 10719374521269., -26136350496073. / 3983972220547.}, {-6511271360970. / 6095937251113., 5843115559534. / 2180450260947., -5289405421727. / 3760307252460.}}; 4059566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX5, 5, 8, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 406a3a57f36SJed Brown } 4078a381b04SJed Brown PetscFunctionReturn(0); 4088a381b04SJed Brown } 4098a381b04SJed Brown 4108a381b04SJed Brown /*@C 411*bcf0153eSBarry Smith TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by `TSARKIMEXRegister()`. 4128a381b04SJed Brown 4138a381b04SJed Brown Not Collective 4148a381b04SJed Brown 4158a381b04SJed Brown Level: advanced 4168a381b04SJed Brown 417*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXRegister()`, `TSARKIMEXRegisterAll()` 4188a381b04SJed Brown @*/ 419d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXRegisterDestroy(void) 420d71ae5a4SJacob Faibussowitsch { 4218a381b04SJed Brown ARKTableauLink link; 4228a381b04SJed Brown 4238a381b04SJed Brown PetscFunctionBegin; 4248a381b04SJed Brown while ((link = ARKTableauList)) { 4258a381b04SJed Brown ARKTableau t = &link->tab; 4268a381b04SJed Brown ARKTableauList = link->next; 4279566063dSJacob Faibussowitsch PetscCall(PetscFree6(t->At, t->bt, t->ct, t->A, t->b, t->c)); 4289566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembedt, t->bembed)); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->binterpt, t->binterp)); 4309566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 4319566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 4328a381b04SJed Brown } 4338a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4348a381b04SJed Brown PetscFunctionReturn(0); 4358a381b04SJed Brown } 4368a381b04SJed Brown 4378a381b04SJed Brown /*@C 438*bcf0153eSBarry Smith TSARKIMEXInitializePackage - This function initializes everything in the `TSARKIMEX` package. It is called 439*bcf0153eSBarry Smith from `TSInitializePackage()`. 4408a381b04SJed Brown 4418a381b04SJed Brown Level: developer 4428a381b04SJed Brown 443*bcf0153eSBarry Smith .seealso: [](chapter_ts), `PetscInitialize()`, `TSARKIMEXFinalizePackage()` 4448a381b04SJed Brown @*/ 445d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXInitializePackage(void) 446d71ae5a4SJacob Faibussowitsch { 4478a381b04SJed Brown PetscFunctionBegin; 4488a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 4498a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 4509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterAll()); 4519566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSARKIMEXFinalizePackage)); 4528a381b04SJed Brown PetscFunctionReturn(0); 4538a381b04SJed Brown } 4548a381b04SJed Brown 4558a381b04SJed Brown /*@C 456*bcf0153eSBarry Smith TSARKIMEXFinalizePackage - This function destroys everything in the `TSARKIMEX` package. It is 457*bcf0153eSBarry Smith called from `PetscFinalize()`. 4588a381b04SJed Brown 4598a381b04SJed Brown Level: developer 4608a381b04SJed Brown 461*bcf0153eSBarry Smith .seealso: [](chapter_ts), `PetscFinalize()`, `TSARKIMEXInitializePackage()` 4628a381b04SJed Brown @*/ 463d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXFinalizePackage(void) 464d71ae5a4SJacob Faibussowitsch { 4658a381b04SJed Brown PetscFunctionBegin; 4668a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 4679566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterDestroy()); 4688a381b04SJed Brown PetscFunctionReturn(0); 4698a381b04SJed Brown } 4708a381b04SJed Brown 471cd652676SJed Brown /*@C 472*bcf0153eSBarry Smith TSARKIMEXRegister - register a `TSARKIMEX` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 473cd652676SJed Brown 474cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 475cd652676SJed Brown 476cd652676SJed Brown Input Parameters: 477cd652676SJed Brown + name - identifier for method 478cd652676SJed Brown . order - approximation order of method 479cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 480cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 4810298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 4820298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 483cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 4840298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 4850298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 4860298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 4870298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 488cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 489cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 4900298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 491cd652676SJed Brown 492cd652676SJed Brown Level: advanced 493cd652676SJed Brown 494*bcf0153eSBarry Smith Note: 495*bcf0153eSBarry Smith Several `TSARKIMEX` methods are provided, this function is only needed to create new methods. 496*bcf0153eSBarry Smith 497*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSType`, `TS` 498cd652676SJed Brown @*/ 499d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name, PetscInt order, PetscInt s, const PetscReal At[], const PetscReal bt[], const PetscReal ct[], const PetscReal A[], const PetscReal b[], const PetscReal c[], const PetscReal bembedt[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[], const PetscReal binterp[]) 500d71ae5a4SJacob Faibussowitsch { 5018a381b04SJed Brown ARKTableauLink link; 5028a381b04SJed Brown ARKTableau t; 5038a381b04SJed Brown PetscInt i, j; 5048a381b04SJed Brown 5058a381b04SJed Brown PetscFunctionBegin; 5069566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 5079566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 5088a381b04SJed Brown t = &link->tab; 5099566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 5108a381b04SJed Brown t->order = order; 5118a381b04SJed Brown t->s = s; 5129566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(s * s, &t->At, s, &t->bt, s, &t->ct, s * s, &t->A, s, &t->b, s, &t->c)); 5139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->At, At, s * s)); 5149566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 5159566063dSJacob Faibussowitsch if (bt) PetscCall(PetscArraycpy(t->bt, bt, s)); 5169371c9d4SSatish Balay else 5179371c9d4SSatish Balay for (i = 0; i < s; i++) t->bt[i] = At[(s - 1) * s + i]; 5189566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 5199371c9d4SSatish Balay else 5209371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = t->bt[i]; 5219566063dSJacob Faibussowitsch if (ct) PetscCall(PetscArraycpy(t->ct, ct, s)); 5229371c9d4SSatish Balay else 5239371c9d4SSatish Balay for (i = 0; i < s; i++) 5249371c9d4SSatish Balay for (j = 0, t->ct[i] = 0; j < s; j++) t->ct[i] += At[i * s + j]; 5259566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 5269371c9d4SSatish Balay else 5279371c9d4SSatish Balay for (i = 0; i < s; i++) 5289371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 529e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 5309371c9d4SSatish Balay for (i = 0; i < s; i++) 5319371c9d4SSatish Balay if (t->At[(s - 1) * s + i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 532e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 5339371c9d4SSatish Balay for (i = 0; i < s; i++) 5349371c9d4SSatish Balay if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 535e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5364e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 537108c343cSJed Brown if (bembedt) { 5389566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &t->bembedt, s, &t->bembed)); 5399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembedt, bembedt, s)); 5409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed ? bembed : bembedt, s)); 541108c343cSJed Brown } 542108c343cSJed Brown 5434f385281SJed Brown t->pinterp = pinterp; 5449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s * pinterp, &t->binterpt, s * pinterp, &t->binterp)); 5459566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 5469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp ? binterp : binterpt, s * pinterp)); 5478a381b04SJed Brown link->next = ARKTableauList; 5488a381b04SJed Brown ARKTableauList = link; 5498a381b04SJed Brown PetscFunctionReturn(0); 5508a381b04SJed Brown } 5518a381b04SJed Brown 552108c343cSJed Brown /* 553108c343cSJed Brown The step completion formula is 554108c343cSJed Brown 555108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 556108c343cSJed Brown 557108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 558108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 559108c343cSJed Brown We can write 560108c343cSJed Brown 561108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 562108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 563108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 564108c343cSJed Brown 565108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 566108c343cSJed Brown */ 567d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts, PetscInt order, Vec X, PetscBool *done) 568d71ae5a4SJacob Faibussowitsch { 569108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 570108c343cSJed Brown ARKTableau tab = ark->tableau; 571108c343cSJed Brown PetscScalar *w = ark->work; 572108c343cSJed Brown PetscReal h; 573108c343cSJed Brown PetscInt s = tab->s, j; 574108c343cSJed Brown 575108c343cSJed Brown PetscFunctionBegin; 576108c343cSJed Brown switch (ark->status) { 577108c343cSJed Brown case TS_STEP_INCOMPLETE: 578d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 579d71ae5a4SJacob Faibussowitsch h = ts->time_step; 580d71ae5a4SJacob Faibussowitsch break; 581d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 582d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 583d71ae5a4SJacob Faibussowitsch break; 584d71ae5a4SJacob Faibussowitsch default: 585d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 586108c343cSJed Brown } 587108c343cSJed Brown if (order == tab->order) { 588e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 589740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 5909566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[s - 1], X)); 591e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 5929566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 593e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bt[j]; 5949566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 595e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 596108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->b[j]; 5979566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 598e817cc15SEmil Constantinescu } 599e817cc15SEmil Constantinescu } 6009566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 601108c343cSJed Brown if (done) *done = PETSC_TRUE; 602108c343cSJed Brown PetscFunctionReturn(0); 603108c343cSJed Brown } else if (order == tab->order - 1) { 604108c343cSJed Brown if (!tab->bembedt) goto unavailable; 605108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 6069566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 607e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembedt[j]; 6089566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 609108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 6109566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 611108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 6129566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 613e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembedt[j] - tab->bt[j]); 6149566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, tab->s, w, ark->YdotI)); 615108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 6169566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 617108c343cSJed Brown } 618108c343cSJed Brown if (done) *done = PETSC_TRUE; 619108c343cSJed Brown PetscFunctionReturn(0); 620108c343cSJed Brown } 621108c343cSJed Brown unavailable: 622108c343cSJed Brown if (done) *done = PETSC_FALSE; 6239371c9d4SSatish Balay else 6249371c9d4SSatish Balay SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "ARKIMEX '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name, 6259371c9d4SSatish Balay tab->order, order); 626108c343cSJed Brown PetscFunctionReturn(0); 627108c343cSJed Brown } 628108c343cSJed Brown 629d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts, PetscBool *id) 630d71ae5a4SJacob Faibussowitsch { 631c79dcfbdSBarry Smith Vec Udot, Y1, Y2; 632c79dcfbdSBarry Smith TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 633c79dcfbdSBarry Smith PetscReal norm; 634c79dcfbdSBarry Smith 635c79dcfbdSBarry Smith PetscFunctionBegin; 6369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Udot)); 6379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y1)); 6389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y2)); 6399566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y1, ark->imex)); 6409566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Udot, NULL)); 6419566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y2, ark->imex)); 6429566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Y1)); 6439566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Udot)); 6449566063dSJacob Faibussowitsch PetscCall(VecNorm(Y2, NORM_2, &norm)); 645c79dcfbdSBarry Smith if (norm < 100.0 * PETSC_MACHINE_EPSILON) { 646c79dcfbdSBarry Smith *id = PETSC_TRUE; 647c79dcfbdSBarry Smith } else { 648c79dcfbdSBarry Smith *id = PETSC_FALSE; 6499566063dSJacob Faibussowitsch PetscCall(PetscInfo((PetscObject)ts, "IFunction(Udot = random) - IFunction(Udot = 0) is not near Udot, %g, suspect mass matrix implied in IFunction() is not the identity as required\n", (double)norm)); 650c79dcfbdSBarry Smith } 6519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Udot)); 6529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 6539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y2)); 654c79dcfbdSBarry Smith PetscFunctionReturn(0); 655c79dcfbdSBarry Smith } 656c79dcfbdSBarry Smith 657d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 658d71ae5a4SJacob Faibussowitsch { 65924655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 66024655328SShri ARKTableau tab = ark->tableau; 66124655328SShri const PetscInt s = tab->s; 66224655328SShri const PetscReal *bt = tab->bt, *b = tab->b; 66324655328SShri PetscScalar *w = ark->work; 66424655328SShri Vec *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS; 66524655328SShri PetscInt j; 666be5899b3SLisandro Dalcin PetscReal h; 66724655328SShri 66824655328SShri PetscFunctionBegin; 669be5899b3SLisandro Dalcin switch (ark->status) { 670be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 671d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 672d71ae5a4SJacob Faibussowitsch h = ts->time_step; 673d71ae5a4SJacob Faibussowitsch break; 674d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 675d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 676d71ae5a4SJacob Faibussowitsch break; 677d71ae5a4SJacob Faibussowitsch default: 678d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 679be5899b3SLisandro Dalcin } 68024655328SShri for (j = 0; j < s; j++) w[j] = -h * bt[j]; 6819566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotI)); 68224655328SShri for (j = 0; j < s; j++) w[j] = -h * b[j]; 6839566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 68424655328SShri PetscFunctionReturn(0); 68524655328SShri } 68624655328SShri 687d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_ARKIMEX(TS ts) 688d71ae5a4SJacob Faibussowitsch { 6898a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 6908a381b04SJed Brown ARKTableau tab = ark->tableau; 6918a381b04SJed Brown const PetscInt s = tab->s; 69224655328SShri const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct, *c = tab->c; 693406d0ec2SJed Brown PetscScalar *w = ark->work; 6941297b224SEmil Constantinescu Vec *Y = ark->Y, *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS, Ydot = ark->Ydot, Ydot0 = ark->Ydot0, Z = ark->Z; 69596400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 696108c343cSJed Brown TSAdapt adapt; 6978a381b04SJed Brown SNES snes; 698fecfb714SLisandro Dalcin PetscInt i, j, its, lits; 699be5899b3SLisandro Dalcin PetscInt rejections = 0; 70096400bd6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 70196400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7028a381b04SJed Brown 7038a381b04SJed Brown PetscFunctionBegin; 70496400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 7059566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 7069566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 7079566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 70896400bd6SLisandro Dalcin } 70996400bd6SLisandro Dalcin 71096400bd6SLisandro Dalcin if (!ts->steprollback) { 71196400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 7129566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[s - 1], Ydot0)); 71396400bd6SLisandro Dalcin } 714fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 71596400bd6SLisandro Dalcin for (i = 0; i < s; i++) { 7169566063dSJacob Faibussowitsch PetscCall(VecCopy(Y[i], ark->Y_prev[i])); 7179566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotRHS[i], ark->YdotRHS_prev[i])); 7189566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[i], ark->YdotI_prev[i])); 71996400bd6SLisandro Dalcin } 72096400bd6SLisandro Dalcin } 72196400bd6SLisandro Dalcin } 72296400bd6SLisandro Dalcin 723fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 72496400bd6SLisandro Dalcin TS ts_start; 725c79dcfbdSBarry Smith if (PetscDefined(USE_DEBUG)) { 726c79dcfbdSBarry Smith PetscBool id = PETSC_FALSE; 7279566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTestMassIdentity(ts, &id)); 7283c633725SBarry Smith PetscCheck(id, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "This scheme requires an identity mass matrix, however the TSIFunction you provide does not utilize an identity mass matrix"); 729c79dcfbdSBarry Smith } 7309566063dSJacob Faibussowitsch PetscCall(TSClone(ts, &ts_start)); 7319566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts_start, ts->vec_sol)); 7329566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts_start, ts->ptime)); 7339566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts_start, ts->steps + 1)); 7349566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts_start, ts->ptime + ts->time_step)); 7359566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts_start, TS_EXACTFINALTIME_STEPOVER)); 7369566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts_start, ts->time_step)); 7379566063dSJacob Faibussowitsch PetscCall(TSSetType(ts_start, TSARKIMEX)); 7389566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts_start, PETSC_TRUE)); 7399566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts_start, TSARKIMEX1BEE)); 74034561852SEmil Constantinescu 7419566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts_start)); 7429566063dSJacob Faibussowitsch PetscCall(TSSolve(ts_start, ts->vec_sol)); 7439566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts_start, &ts->ptime)); 7449566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts_start, &ts->time_step)); 745bbd56ea5SKarl Rupp 74685fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 74785fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX *)ts_start->data; 7489566063dSJacob Faibussowitsch PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s - 1], Ydot0)); 74985fc7851SLisandro Dalcin } 75096400bd6SLisandro Dalcin ts->steps++; 75134561852SEmil Constantinescu 752d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 753d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 75496400bd6SLisandro Dalcin { 7559566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7569566063dSJacob Faibussowitsch PetscCall(TSSetSNES(ts, snes)); 75796400bd6SLisandro Dalcin } 7589566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts_start)); 759e817cc15SEmil Constantinescu } 760e817cc15SEmil Constantinescu 761108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 76296400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 76396400bd6SLisandro Dalcin PetscReal t = ts->ptime; 764108c343cSJed Brown PetscReal h = ts->time_step; 7658a381b04SJed Brown for (i = 0; i < s; i++) { 7669be3e283SDebojyoti Ghosh ark->stage_time = t + h * ct[i]; 7679566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ark->stage_time)); 7688a381b04SJed Brown if (At[i * s + i] == 0) { /* This stage is explicit */ 7693c633725SBarry Smith PetscCheck(i == 0 || ts->equation_type < TS_EQ_IMPLICIT, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Explicit stages other than the first one are not supported for implicit problems"); 7709566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 771e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotI)); 7738a381b04SJed Brown for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7749566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 7758a381b04SJed Brown } else { 776b296d7d5SJed Brown ark->scoeff = 1. / At[i * s + i]; 7778a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 7789566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Z)); 779e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7809566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotI)); 781c58d1302SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7829566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotRHS)); 783fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 78456dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 7859566063dSJacob Faibussowitsch PetscCall(TSExtrapolate_ARKIMEX(ts, c[i], Y[i])); 78656dcabbaSDebojyoti Ghosh } else { 7878a381b04SJed Brown /* Initial guess taken from last stage */ 7889566063dSJacob Faibussowitsch PetscCall(VecCopy(i > 0 ? Y[i - 1] : ts->vec_sol, Y[i])); 78956dcabbaSDebojyoti Ghosh } 7909566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7919566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, Y[i])); 7929566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 7939566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 7949371c9d4SSatish Balay ts->snes_its += its; 7959371c9d4SSatish Balay ts->ksp_its += lits; 7969566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 7979566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, ark->stage_time, Y[i], &stageok)); 79896400bd6SLisandro Dalcin if (!stageok) { 7991be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8001be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 80196400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8021be93e3eSJed Brown goto reject_step; 8031be93e3eSJed Brown } 8048a381b04SJed Brown } 805e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 806e817cc15SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 8079371c9d4SSatish Balay PetscCheck(tab->stiffly_accurate, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated", 8089371c9d4SSatish Balay ark->tableau->name); 8099566063dSJacob Faibussowitsch PetscCall(VecCopy(Ydot0, YdotI[0])); /* YdotI = YdotI(tn-1) */ 810e817cc15SEmil Constantinescu } else { 8119566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 812e817cc15SEmil Constantinescu } 813e817cc15SEmil Constantinescu } else { 8145eca1a21SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 8159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); 8169566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, t + h * ct[i], Y[i], Ydot, YdotI[i], ark->imex)); /* YdotI = -G(t,Y,0) */ 8179566063dSJacob Faibussowitsch PetscCall(VecScale(YdotI[i], -1.0)); 8185eca1a21SEmil Constantinescu } else { 8199566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 8205eca1a21SEmil Constantinescu } 8214cc180ffSJed Brown if (ark->imex) { 8229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 8234cc180ffSJed Brown } else { 8249566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(YdotRHS[i])); 8254cc180ffSJed Brown } 8268a381b04SJed Brown } 8279566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ark->stage_time, i, Y)); 828e817cc15SEmil Constantinescu } 82996400bd6SLisandro Dalcin 830be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 8319566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_ARKIMEX(ts, tab->order, ts->vec_sol, NULL)); 832108c343cSJed Brown ark->status = TS_STEP_PENDING; 8339566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8349566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8359566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8369566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 83796400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 83896400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8399566063dSJacob Faibussowitsch PetscCall(TSRollBack_ARKIMEX(ts)); 840be5899b3SLisandro Dalcin ts->time_step = next_time_step; 841be5899b3SLisandro Dalcin goto reject_step; 84296400bd6SLisandro Dalcin } 84396400bd6SLisandro Dalcin 8448a381b04SJed Brown ts->ptime += ts->time_step; 845cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 846108c343cSJed Brown break; 84796400bd6SLisandro Dalcin 84896400bd6SLisandro Dalcin reject_step: 8499371c9d4SSatish Balay ts->reject++; 8509371c9d4SSatish Balay accept = PETSC_FALSE; 851be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 85296400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 85363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 854108c343cSJed Brown } 855f85781f1SEmil Constantinescu } 8568a381b04SJed Brown PetscFunctionReturn(0); 8578a381b04SJed Brown } 8588a381b04SJed Brown 859d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_ARKIMEX(TS ts, PetscReal itime, Vec X) 860d71ae5a4SJacob Faibussowitsch { 861cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 8624f385281SJed Brown PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 863108c343cSJed Brown PetscReal h; 864108c343cSJed Brown PetscReal tt, t; 865cd652676SJed Brown PetscScalar *bt, *b; 866cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 867cd652676SJed Brown 868cd652676SJed Brown PetscFunctionBegin; 8693c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 870108c343cSJed Brown switch (ark->status) { 871108c343cSJed Brown case TS_STEP_INCOMPLETE: 872108c343cSJed Brown case TS_STEP_PENDING: 873108c343cSJed Brown h = ts->time_step; 874108c343cSJed Brown t = (itime - ts->ptime) / h; 875108c343cSJed Brown break; 876108c343cSJed Brown case TS_STEP_COMPLETE: 877be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 878108c343cSJed Brown t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 879108c343cSJed Brown break; 880d71ae5a4SJacob Faibussowitsch default: 881d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 882108c343cSJed Brown } 8839566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &bt, s, &b)); 884cd652676SJed Brown for (i = 0; i < s; i++) bt[i] = b[i] = 0; 8854f385281SJed Brown for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 886cd652676SJed Brown for (i = 0; i < s; i++) { 887c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 888108c343cSJed Brown b[i] += h * B[i * pinterp + j] * tt; 889cd652676SJed Brown } 890cd652676SJed Brown } 8919566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[0], X)); 8929566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI)); 8939566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS)); 8949566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 895cd652676SJed Brown PetscFunctionReturn(0); 896cd652676SJed Brown } 897cd652676SJed Brown 898d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts, PetscReal c, Vec X) 899d71ae5a4SJacob Faibussowitsch { 90056dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 90156dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 902be5899b3SLisandro Dalcin PetscReal h, h_prev, t, tt; 90356dcabbaSDebojyoti Ghosh PetscScalar *bt, *b; 90456dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 90556dcabbaSDebojyoti Ghosh 90656dcabbaSDebojyoti Ghosh PetscFunctionBegin; 9073c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 9089566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(s, &bt, s, &b)); 90981d12688SDebojyoti Ghosh h = ts->time_step; 910be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 911be5899b3SLisandro Dalcin t = 1 + h / h_prev * c; 91256dcabbaSDebojyoti Ghosh for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 91356dcabbaSDebojyoti Ghosh for (i = 0; i < s; i++) { 91481d12688SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 91556dcabbaSDebojyoti Ghosh b[i] += h * B[i * pinterp + j] * tt; 91656dcabbaSDebojyoti Ghosh } 91756dcabbaSDebojyoti Ghosh } 9183c633725SBarry Smith PetscCheck(ark->Y_prev, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Stages from previous step have not been stored"); 9199566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y_prev[0], X)); 9209566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI_prev)); 9219566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS_prev)); 9229566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 92356dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 92456dcabbaSDebojyoti Ghosh } 92556dcabbaSDebojyoti Ghosh 9268a381b04SJed Brown /*------------------------------------------------------------*/ 92796400bd6SLisandro Dalcin 928d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTableauReset(TS ts) 929d71ae5a4SJacob Faibussowitsch { 93096400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 93196400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 93296400bd6SLisandro Dalcin 93396400bd6SLisandro Dalcin PetscFunctionBegin; 93496400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 9359566063dSJacob Faibussowitsch PetscCall(PetscFree(ark->work)); 9369566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y)); 9379566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI)); 9389566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS)); 9399566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev)); 9409566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev)); 9419566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev)); 94296400bd6SLisandro Dalcin PetscFunctionReturn(0); 94396400bd6SLisandro Dalcin } 94496400bd6SLisandro Dalcin 945d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_ARKIMEX(TS ts) 946d71ae5a4SJacob Faibussowitsch { 9478a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 9488a381b04SJed Brown 9498a381b04SJed Brown PetscFunctionBegin; 9509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauReset(ts)); 9519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot)); 9529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot0)); 9539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Z)); 9548a381b04SJed Brown PetscFunctionReturn(0); 9558a381b04SJed Brown } 9568a381b04SJed Brown 957d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) 958d71ae5a4SJacob Faibussowitsch { 959d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX *)ts->data; 960d5e6173cSPeter Brune 961d5e6173cSPeter Brune PetscFunctionBegin; 962d5e6173cSPeter Brune if (Z) { 963d5e6173cSPeter Brune if (dm && dm != ts->dm) { 9649566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 965d5e6173cSPeter Brune } else *Z = ax->Z; 966d5e6173cSPeter Brune } 967d5e6173cSPeter Brune if (Ydot) { 968d5e6173cSPeter Brune if (dm && dm != ts->dm) { 9699566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 970d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 971d5e6173cSPeter Brune } 972d5e6173cSPeter Brune PetscFunctionReturn(0); 973d5e6173cSPeter Brune } 974d5e6173cSPeter Brune 975d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXRestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) 976d71ae5a4SJacob Faibussowitsch { 977d5e6173cSPeter Brune PetscFunctionBegin; 978d5e6173cSPeter Brune if (Z) { 97948a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 980d5e6173cSPeter Brune } 981d5e6173cSPeter Brune if (Ydot) { 98248a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 983d5e6173cSPeter Brune } 984d5e6173cSPeter Brune PetscFunctionReturn(0); 985d5e6173cSPeter Brune } 986d5e6173cSPeter Brune 9878a381b04SJed Brown /* 9888a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 9898a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 9908a381b04SJed Brown */ 991d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes, Vec X, Vec F, TS ts) 992d71ae5a4SJacob Faibussowitsch { 9938a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 994d5e6173cSPeter Brune DM dm, dmsave; 995d5e6173cSPeter Brune Vec Z, Ydot; 996b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 9978a381b04SJed Brown 9988a381b04SJed Brown PetscFunctionBegin; 9999566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 10009566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, &Ydot)); 10019566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Ydot, -shift, shift, 0, Z, X)); /* Ydot = shift*(X-Z) */ 1002d5e6173cSPeter Brune dmsave = ts->dm; 1003d5e6173cSPeter Brune ts->dm = dm; 1004740132f1SEmil Constantinescu 10059566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ark->stage_time, X, Ydot, F, ark->imex)); 1006e817cc15SEmil Constantinescu 1007d5e6173cSPeter Brune ts->dm = dmsave; 10089566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, &Ydot)); 10098a381b04SJed Brown PetscFunctionReturn(0); 10108a381b04SJed Brown } 10118a381b04SJed Brown 1012d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes, Vec X, Mat A, Mat B, TS ts) 1013d71ae5a4SJacob Faibussowitsch { 10148a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1015d5e6173cSPeter Brune DM dm, dmsave; 1016d5e6173cSPeter Brune Vec Ydot; 1017b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10188a381b04SJed Brown 10198a381b04SJed Brown PetscFunctionBegin; 10209566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 10219566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, NULL, &Ydot)); 10228a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1023d5e6173cSPeter Brune dmsave = ts->dm; 1024d5e6173cSPeter Brune ts->dm = dm; 1025740132f1SEmil Constantinescu 10269566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ark->stage_time, X, Ydot, shift, A, B, ark->imex)); 1027740132f1SEmil Constantinescu 1028d5e6173cSPeter Brune ts->dm = dmsave; 10299566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, NULL, &Ydot)); 1030d5e6173cSPeter Brune PetscFunctionReturn(0); 1031d5e6173cSPeter Brune } 1032d5e6173cSPeter Brune 1033d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine, DM coarse, void *ctx) 1034d71ae5a4SJacob Faibussowitsch { 1035d5e6173cSPeter Brune PetscFunctionBegin; 1036d5e6173cSPeter Brune PetscFunctionReturn(0); 1037d5e6173cSPeter Brune } 1038d5e6173cSPeter Brune 1039d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1040d71ae5a4SJacob Faibussowitsch { 1041d5e6173cSPeter Brune TS ts = (TS)ctx; 1042d5e6173cSPeter Brune Vec Z, Z_c; 1043d5e6173cSPeter Brune 1044d5e6173cSPeter Brune PetscFunctionBegin; 10459566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, fine, &Z, NULL)); 10469566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, coarse, &Z_c, NULL)); 10479566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Z, Z_c)); 10489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Z_c, rscale, Z_c)); 10499566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, fine, &Z, NULL)); 10509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, coarse, &Z_c, NULL)); 10518a381b04SJed Brown PetscFunctionReturn(0); 10528a381b04SJed Brown } 10538a381b04SJed Brown 1054d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm, DM subdm, void *ctx) 1055d71ae5a4SJacob Faibussowitsch { 1056cdb298fcSPeter Brune PetscFunctionBegin; 1057cdb298fcSPeter Brune PetscFunctionReturn(0); 1058cdb298fcSPeter Brune } 1059cdb298fcSPeter Brune 1060d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1061d71ae5a4SJacob Faibussowitsch { 1062cdb298fcSPeter Brune TS ts = (TS)ctx; 1063cdb298fcSPeter Brune Vec Z, Z_c; 1064cdb298fcSPeter Brune 1065cdb298fcSPeter Brune PetscFunctionBegin; 10669566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, NULL)); 10679566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, subdm, &Z_c, NULL)); 1068cdb298fcSPeter Brune 10699566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 10709566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 1071cdb298fcSPeter Brune 10729566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, NULL)); 10739566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, subdm, &Z_c, NULL)); 1074cdb298fcSPeter Brune PetscFunctionReturn(0); 1075cdb298fcSPeter Brune } 1076cdb298fcSPeter Brune 1077d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 1078d71ae5a4SJacob Faibussowitsch { 107996400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 108096400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 108196400bd6SLisandro Dalcin 108296400bd6SLisandro Dalcin PetscFunctionBegin; 10839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &ark->work)); 10849566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y)); 10859566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI)); 10869566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS)); 108796400bd6SLisandro Dalcin if (ark->extrapolate) { 10889566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 10899566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 10909566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 109196400bd6SLisandro Dalcin } 109296400bd6SLisandro Dalcin PetscFunctionReturn(0); 109396400bd6SLisandro Dalcin } 109496400bd6SLisandro Dalcin 1095d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 1096d71ae5a4SJacob Faibussowitsch { 10978a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1098d5e6173cSPeter Brune DM dm; 109996400bd6SLisandro Dalcin SNES snes; 1100f9c1d6abSBarry Smith 11018a381b04SJed Brown PetscFunctionBegin; 11029566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauSetUp(ts)); 11039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot)); 11049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot0)); 11059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Z)); 11069566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11079566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 11089566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 11099566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 11108a381b04SJed Brown PetscFunctionReturn(0); 11118a381b04SJed Brown } 11128a381b04SJed Brown /*------------------------------------------------------------*/ 11138a381b04SJed Brown 1114d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts, PetscOptionItems *PetscOptionsObject) 1115d71ae5a4SJacob Faibussowitsch { 11164cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 11178a381b04SJed Brown 11188a381b04SJed Brown PetscFunctionBegin; 1119d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ARKIMEX ODE solver options"); 11208a381b04SJed Brown { 11218a381b04SJed Brown ARKTableauLink link; 11228a381b04SJed Brown PetscInt count, choice; 11238a381b04SJed Brown PetscBool flg; 11248a381b04SJed Brown const char **namelist; 11259371c9d4SSatish Balay for (link = ARKTableauList, count = 0; link; link = link->next, count++) 11269371c9d4SSatish Balay ; 11279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 11288a381b04SJed Brown for (link = ARKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 11299566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_arkimex_type", "Family of ARK IMEX method", "TSARKIMEXSetType", (const char *const *)namelist, count, ark->tableau->name, &choice, &flg)); 11309566063dSJacob Faibussowitsch if (flg) PetscCall(TSARKIMEXSetType(ts, namelist[choice])); 11319566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 113296400bd6SLisandro Dalcin 11334cc180ffSJed Brown flg = (PetscBool)!ark->imex; 11349566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_arkimex_fully_implicit", "Solve the problem fully implicitly", "TSARKIMEXSetFullyImplicit", flg, &flg, NULL)); 11354cc180ffSJed Brown ark->imex = (PetscBool)!flg; 11369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate", "Extrapolate the initial guess for the stage solution from stage values of the previous time step", "", ark->extrapolate, &ark->extrapolate, NULL)); 11378a381b04SJed Brown } 1138d0609cedSBarry Smith PetscOptionsHeadEnd(); 11398a381b04SJed Brown PetscFunctionReturn(0); 11408a381b04SJed Brown } 11418a381b04SJed Brown 1142d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_ARKIMEX(TS ts, PetscViewer viewer) 1143d71ae5a4SJacob Faibussowitsch { 11448a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 11458a381b04SJed Brown PetscBool iascii; 11468a381b04SJed Brown 11478a381b04SJed Brown PetscFunctionBegin; 11489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 11498a381b04SJed Brown if (iascii) { 11509c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 115119fd82e9SBarry Smith TSARKIMEXType arktype; 11528a381b04SJed Brown char buf[512]; 11533a28c0e5SStefano Zampini PetscBool flg; 11543a28c0e5SStefano Zampini 11559566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetType(ts, &arktype)); 11569566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetFullyImplicit(ts, &flg)); 11579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ARK IMEX %s\n", arktype)); 11589566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ct)); 11599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Stiff abscissa ct = %s\n", buf)); 11609566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c)); 11619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Fully implicit: %s\n", flg ? "yes" : "no")); 11629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s\n", tab->stiffly_accurate ? "yes" : "no")); 11639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Explicit first stage: %s\n", tab->explicit_first_stage ? "yes" : "no")); 11649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "FSAL property: %s\n", tab->FSAL_implicit ? "yes" : "no")); 11659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Nonstiff abscissa c = %s\n", buf)); 11668a381b04SJed Brown } 11678a381b04SJed Brown PetscFunctionReturn(0); 11688a381b04SJed Brown } 11698a381b04SJed Brown 1170d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_ARKIMEX(TS ts, PetscViewer viewer) 1171d71ae5a4SJacob Faibussowitsch { 1172f2c2a1b9SBarry Smith SNES snes; 11739c334d8fSLisandro Dalcin TSAdapt adapt; 1174f2c2a1b9SBarry Smith 1175f2c2a1b9SBarry Smith PetscFunctionBegin; 11769566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 11779566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 11789566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 11799566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes, viewer)); 1180ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 11819566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 11829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 1183f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1184f2c2a1b9SBarry Smith } 1185f2c2a1b9SBarry Smith 11868a381b04SJed Brown /*@C 1187*bcf0153eSBarry Smith TSARKIMEXSetType - Set the type of `TSARKIMEX` scheme 11888a381b04SJed Brown 11898a381b04SJed Brown Logically collective 11908a381b04SJed Brown 1191d8d19677SJose E. Roman Input Parameters: 11928a381b04SJed Brown + ts - timestepping context 1193*bcf0153eSBarry Smith - arktype - type of `TSARKIMEX` scheme 11948a381b04SJed Brown 1195*bcf0153eSBarry Smith Options Database Key: 1196*bcf0153eSBarry Smith . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> - set `TSARKIMEX` scheme type 11979bd3e852SBarry Smith 11988a381b04SJed Brown Level: intermediate 11998a381b04SJed Brown 1200*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEXGetType()`, `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEX1BEE`, `TSARKIMEXA2`, `TSARKIMEXL2`, `TSARKIMEXARS122`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEXPRSSP2`, 1201db781477SPatrick Sanan `TSARKIMEX3`, `TSARKIMEXBPR3`, `TSARKIMEXARS443`, `TSARKIMEX4`, `TSARKIMEX5` 12028a381b04SJed Brown @*/ 1203d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType arktype) 1204d71ae5a4SJacob Faibussowitsch { 12058a381b04SJed Brown PetscFunctionBegin; 12068a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1207b92453a8SLisandro Dalcin PetscValidCharPointer(arktype, 2); 1208cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetType_C", (TS, TSARKIMEXType), (ts, arktype)); 12098a381b04SJed Brown PetscFunctionReturn(0); 12108a381b04SJed Brown } 12118a381b04SJed Brown 12128a381b04SJed Brown /*@C 1213*bcf0153eSBarry Smith TSARKIMEXGetType - Get the type of `TSARKIMEX` scheme 12148a381b04SJed Brown 12158a381b04SJed Brown Logically collective 12168a381b04SJed Brown 12178a381b04SJed Brown Input Parameter: 12188a381b04SJed Brown . ts - timestepping context 12198a381b04SJed Brown 12208a381b04SJed Brown Output Parameter: 1221*bcf0153eSBarry Smith . arktype - type of `TSARKIMEX` scheme 12228a381b04SJed Brown 12238a381b04SJed Brown Level: intermediate 12248a381b04SJed Brown 1225*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`c, `TSARKIMEXGetType()` 12268a381b04SJed Brown @*/ 1227d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *arktype) 1228d71ae5a4SJacob Faibussowitsch { 12298a381b04SJed Brown PetscFunctionBegin; 12308a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1231cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetType_C", (TS, TSARKIMEXType *), (ts, arktype)); 12328a381b04SJed Brown PetscFunctionReturn(0); 12338a381b04SJed Brown } 12348a381b04SJed Brown 123516353aafSBarry Smith /*@ 1236*bcf0153eSBarry Smith TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly, including the part that is normally solved explicitly 12374cc180ffSJed Brown 12384cc180ffSJed Brown Logically collective 12394cc180ffSJed Brown 1240d8d19677SJose E. Roman Input Parameters: 12414cc180ffSJed Brown + ts - timestepping context 1242*bcf0153eSBarry Smith - flg - `PETSC_TRUE` for fully implicit 12434cc180ffSJed Brown 12444cc180ffSJed Brown Level: intermediate 12454cc180ffSJed Brown 1246*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEX`, `TSARKIMEXGetType()`, `TSARKIMEXGetFullyImplicit()` 12474cc180ffSJed Brown @*/ 1248d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts, PetscBool flg) 1249d71ae5a4SJacob Faibussowitsch { 12504cc180ffSJed Brown PetscFunctionBegin; 12514cc180ffSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12523a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts, flg, 2); 1253cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetFullyImplicit_C", (TS, PetscBool), (ts, flg)); 12544cc180ffSJed Brown PetscFunctionReturn(0); 12554cc180ffSJed Brown } 12564cc180ffSJed Brown 12573a28c0e5SStefano Zampini /*@ 12583a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 12593a28c0e5SStefano Zampini 12603a28c0e5SStefano Zampini Logically collective 12613a28c0e5SStefano Zampini 12623a28c0e5SStefano Zampini Input Parameter: 12633a28c0e5SStefano Zampini . ts - timestepping context 12643a28c0e5SStefano Zampini 12657a7aea1fSJed Brown Output Parameter: 1266*bcf0153eSBarry Smith . flg - `PETSC_TRUE` for fully implicit 12673a28c0e5SStefano Zampini 12683a28c0e5SStefano Zampini Level: intermediate 12693a28c0e5SStefano Zampini 1270*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()` 12713a28c0e5SStefano Zampini @*/ 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts, PetscBool *flg) 1273d71ae5a4SJacob Faibussowitsch { 12743a28c0e5SStefano Zampini PetscFunctionBegin; 12753a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1276dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 2); 1277cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetFullyImplicit_C", (TS, PetscBool *), (ts, flg)); 12783a28c0e5SStefano Zampini PetscFunctionReturn(0); 12793a28c0e5SStefano Zampini } 12803a28c0e5SStefano Zampini 1281d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts, TSARKIMEXType *arktype) 1282d71ae5a4SJacob Faibussowitsch { 12838a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12848a381b04SJed Brown 12858a381b04SJed Brown PetscFunctionBegin; 12868a381b04SJed Brown *arktype = ark->tableau->name; 12878a381b04SJed Brown PetscFunctionReturn(0); 12888a381b04SJed Brown } 1289d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts, TSARKIMEXType arktype) 1290d71ae5a4SJacob Faibussowitsch { 12918a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12928a381b04SJed Brown PetscBool match; 12938a381b04SJed Brown ARKTableauLink link; 12948a381b04SJed Brown 12958a381b04SJed Brown PetscFunctionBegin; 12968a381b04SJed Brown if (ark->tableau) { 12979566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ark->tableau->name, arktype, &match)); 12988a381b04SJed Brown if (match) PetscFunctionReturn(0); 12998a381b04SJed Brown } 13008a381b04SJed Brown for (link = ARKTableauList; link; link = link->next) { 13019566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, arktype, &match)); 13028a381b04SJed Brown if (match) { 13039566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauReset(ts)); 13048a381b04SJed Brown ark->tableau = &link->tab; 13059566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauSetUp(ts)); 13062ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13078a381b04SJed Brown PetscFunctionReturn(0); 13088a381b04SJed Brown } 13098a381b04SJed Brown } 131098921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", arktype); 13118a381b04SJed Brown } 1312e0877f53SBarry Smith 1313d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts, PetscBool flg) 1314d71ae5a4SJacob Faibussowitsch { 13154cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 13164cc180ffSJed Brown 13174cc180ffSJed Brown PetscFunctionBegin; 13184cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13194cc180ffSJed Brown PetscFunctionReturn(0); 13204cc180ffSJed Brown } 13218a381b04SJed Brown 1322d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts, PetscBool *flg) 1323d71ae5a4SJacob Faibussowitsch { 13243a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 13253a28c0e5SStefano Zampini 13263a28c0e5SStefano Zampini PetscFunctionBegin; 13273a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 13283a28c0e5SStefano Zampini PetscFunctionReturn(0); 13293a28c0e5SStefano Zampini } 13303a28c0e5SStefano Zampini 1331d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1332d71ae5a4SJacob Faibussowitsch { 1333b3a6b972SJed Brown PetscFunctionBegin; 13349566063dSJacob Faibussowitsch PetscCall(TSReset_ARKIMEX(ts)); 1335b3a6b972SJed Brown if (ts->dm) { 13369566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 13379566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 1338b3a6b972SJed Brown } 13399566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 13409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", NULL)); 13419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", NULL)); 13429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", NULL)); 13432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", NULL)); 1344b3a6b972SJed Brown PetscFunctionReturn(0); 1345b3a6b972SJed Brown } 1346b3a6b972SJed Brown 13478a381b04SJed Brown /* ------------------------------------------------------------ */ 13488a381b04SJed Brown /*MC 1349c79dcfbdSBarry Smith TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes 13508a381b04SJed Brown 1351fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1352fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1353*bcf0153eSBarry Smith of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`. 1354d0685a90SJed Brown 13558a381b04SJed Brown Level: beginner 13568a381b04SJed Brown 1357*bcf0153eSBarry Smith Notes: 1358*bcf0153eSBarry Smith The default is `TSARKIMEX3`, it can be changed with `TSARKIMEXSetType()` or -ts_arkimex_type 13598a381b04SJed Brown 1360*bcf0153eSBarry Smith If the equation is implicit or a DAE, then `TSSetEquationType()` needs to be set accordingly. Refer to the manual for further information. 1361*bcf0153eSBarry Smith 1362*bcf0153eSBarry Smith 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). 1363*bcf0153eSBarry Smith 1364*bcf0153eSBarry Smith Consider trying `TSROSW` if the stiff part is linear or weakly nonlinear. 1365*bcf0153eSBarry Smith 1366*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSARKIMEXSetType()`, `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()`, `TSARKIMEXGetFullyImplicit()`, 1367*bcf0153eSBarry Smith `TSARKIMEX1BEE`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEX3`, `TSARKIMEXL2`, `TSARKIMEXA2`, `TSARKIMEXARS122`, 1368*bcf0153eSBarry Smith `TSARKIMEX4`, `TSARKIMEX5`, `TSARKIMEXPRSSP2`, `TSARKIMEXARS443`, `TSARKIMEXBPR3`, `TSARKIMEXType`, `TSARKIMEXRegister()`, `TSType` 13698a381b04SJed Brown M*/ 1370d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 1371d71ae5a4SJacob Faibussowitsch { 13728a381b04SJed Brown TS_ARKIMEX *th; 13738a381b04SJed Brown 13748a381b04SJed Brown PetscFunctionBegin; 13759566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 13768a381b04SJed Brown 13778a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 13788a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 13798a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1380f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 13818a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 13828a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1383cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1384108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 138524655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 13868a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 13878a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 13888a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 13898a381b04SJed Brown 1390efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1391efd4aadfSBarry Smith 13924dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 13938a381b04SJed Brown ts->data = (void *)th; 13944cc180ffSJed Brown th->imex = PETSC_TRUE; 13958a381b04SJed Brown 13969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", TSARKIMEXGetType_ARKIMEX)); 13979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", TSARKIMEXSetType_ARKIMEX)); 13989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", TSARKIMEXSetFullyImplicit_ARKIMEX)); 13999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", TSARKIMEXGetFullyImplicit_ARKIMEX)); 140096400bd6SLisandro Dalcin 14019566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts, TSARKIMEXDefault)); 14028a381b04SJed Brown PetscFunctionReturn(0); 14038a381b04SJed Brown } 1404