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 659bd3e852SBarry Smith Options Database: 6667b8a455SSatish Balay . -ts_arkimex_type ars122 - set arkimex type to ars122 679bd3e852SBarry Smith 681f80e275SEmil Constantinescu References: 69606c0280SSatish Balay . * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 701f80e275SEmil Constantinescu 711f80e275SEmil Constantinescu Level: advanced 721f80e275SEmil Constantinescu 73db781477SPatrick Sanan .seealso: `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 809bd3e852SBarry Smith Options Database: 8167b8a455SSatish Balay . -ts_arkimex_type a2 - set arkimex type to a2 829bd3e852SBarry Smith 831f80e275SEmil Constantinescu Level: advanced 841f80e275SEmil Constantinescu 85db781477SPatrick Sanan .seealso: `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 929bd3e852SBarry Smith Options Database: 9367b8a455SSatish Balay . -ts_arkimex_type l2 - set arkimex type to l2 949bd3e852SBarry Smith 951f80e275SEmil Constantinescu References: 96606c0280SSatish 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. 971f80e275SEmil Constantinescu 981f80e275SEmil Constantinescu Level: advanced 991f80e275SEmil Constantinescu 100db781477SPatrick Sanan .seealso: `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 1079bd3e852SBarry Smith Options Database: 10867b8a455SSatish Balay . -ts_arkimex_type 1bee - set arkimex type to 1bee 1099bd3e852SBarry Smith 110e817cc15SEmil Constantinescu Level: advanced 111e817cc15SEmil Constantinescu 112db781477SPatrick Sanan .seealso: `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 1199bd3e852SBarry Smith Options Database: 12067b8a455SSatish Balay . -ts_arkimex_type 2c - set arkimex type to 2c 1219bd3e852SBarry Smith 1221f80e275SEmil Constantinescu Level: advanced 1231f80e275SEmil Constantinescu 124db781477SPatrick Sanan .seealso: `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 1319bd3e852SBarry Smith Options Database: 13267b8a455SSatish Balay . -ts_arkimex_type 2d - set arkimex type to 2d 1339bd3e852SBarry Smith 134b330ce4dSSatish Balay Level: advanced 135b330ce4dSSatish Balay 136db781477SPatrick Sanan .seealso: `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 1439bd3e852SBarry Smith Options Database: 14467b8a455SSatish Balay . -ts_arkimex_type 2e - set arkimex type to 2e 1459bd3e852SBarry Smith 146b330ce4dSSatish Balay Level: advanced 147b330ce4dSSatish Balay 148db781477SPatrick Sanan .seealso: `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 1609bd3e852SBarry Smith Options Database: 16167b8a455SSatish Balay . -ts_arkimex_type prssp2 - set arkimex type to prssp2 1629bd3e852SBarry Smith 1636cf0794eSJed Brown Level: advanced 1646cf0794eSJed Brown 165db781477SPatrick Sanan .seealso: `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 1729bd3e852SBarry Smith Options Database: 17367b8a455SSatish Balay . -ts_arkimex_type 3 - set arkimex type to 3 1749bd3e852SBarry Smith 17564f491ddSJed Brown References: 176606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 17764f491ddSJed Brown 178b330ce4dSSatish Balay Level: advanced 179b330ce4dSSatish Balay 180db781477SPatrick Sanan .seealso: `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 1879bd3e852SBarry Smith Options Database: 18867b8a455SSatish Balay . -ts_arkimex_type ars443 - set arkimex type to ars443 1899bd3e852SBarry Smith 1906cf0794eSJed Brown References: 191606c0280SSatish Balay + * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 192606c0280SSatish Balay - * - This method is referred to as ARS(4,4,3) in https://arxiv.org/abs/1110.4375 1936cf0794eSJed Brown 1946cf0794eSJed Brown Level: advanced 1956cf0794eSJed Brown 196db781477SPatrick Sanan .seealso: `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 2039bd3e852SBarry Smith Options Database: 20467b8a455SSatish Balay . -ts_arkimex_type bpr3 - set arkimex type to bpr3 2059bd3e852SBarry Smith 2066cf0794eSJed Brown References: 207606c0280SSatish Balay . * - This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375 2086cf0794eSJed Brown 2096cf0794eSJed Brown Level: advanced 2106cf0794eSJed Brown 211db781477SPatrick Sanan .seealso: `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 2189bd3e852SBarry Smith Options Database: 21967b8a455SSatish Balay . -ts_arkimex_type 4 - set arkimex type to4 2209bd3e852SBarry Smith 22164f491ddSJed Brown References: 222606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 22364f491ddSJed Brown 224b330ce4dSSatish Balay Level: advanced 225b330ce4dSSatish Balay 226db781477SPatrick Sanan .seealso: `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 2339bd3e852SBarry Smith Options Database: 23467b8a455SSatish Balay . -ts_arkimex_type 5 - set arkimex type to 5 2359bd3e852SBarry Smith 23664f491ddSJed Brown References: 237606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 23864f491ddSJed Brown 239b330ce4dSSatish Balay Level: advanced 240b330ce4dSSatish Balay 241db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 24264f491ddSJed Brown M*/ 24364f491ddSJed Brown 2448a381b04SJed Brown /*@C 2458a381b04SJed Brown 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 251db781477SPatrick Sanan .seealso: `TSARKIMEXRegisterDestroy()` 2528a381b04SJed Brown @*/ 2539371c9d4SSatish Balay PetscErrorCode TSARKIMEXRegisterAll(void) { 2548a381b04SJed Brown PetscFunctionBegin; 2558a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2568a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 257e817cc15SEmil Constantinescu 258e817cc15SEmil Constantinescu { 2599371c9d4SSatish Balay const PetscReal A[3][3] = 2609371c9d4SSatish Balay { 261e817cc15SEmil Constantinescu {0.0, 0.0, 0.0}, 2629371c9d4SSatish Balay {0.0, 0.0, 0.0}, 2639371c9d4SSatish Balay {0.0, 0.5, 0.0} 2649371c9d4SSatish Balay }, 2659371c9d4SSatish 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}; 2669566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX1BEE, 2, 3, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 267e817cc15SEmil Constantinescu } 2688a381b04SJed Brown { 2699371c9d4SSatish Balay const PetscReal A[2][2] = 2709371c9d4SSatish Balay { 2719371c9d4SSatish Balay {0.0, 0.0}, 2729371c9d4SSatish Balay {0.5, 0.0} 2739371c9d4SSatish Balay }, 2749371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.0, 0.5}}, b[2] = {0.0, 1.0}, bembedt[2] = {0.5, 0.5}; 2751f80e275SEmil 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*/ 2769566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS122, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2771f80e275SEmil Constantinescu } 2781f80e275SEmil Constantinescu { 2799371c9d4SSatish Balay const PetscReal A[2][2] = 2809371c9d4SSatish Balay { 2819371c9d4SSatish Balay {0.0, 0.0}, 2829371c9d4SSatish Balay {1.0, 0.0} 2839371c9d4SSatish Balay }, 2849371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.5, 0.5}}, b[2] = {0.5, 0.5}, bembedt[2] = {0.0, 1.0}; 2851f80e275SEmil 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*/ 2869566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXA2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2871f80e275SEmil Constantinescu } 2881f80e275SEmil Constantinescu { 289da80777bSKarl 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 */ 2909371c9d4SSatish Balay const PetscReal A[2][2] = 2919371c9d4SSatish Balay { 2929371c9d4SSatish Balay {0.0, 0.0}, 2939371c9d4SSatish Balay {1.0, 0.0} 2949371c9d4SSatish Balay }, 2959371c9d4SSatish 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}}; 2969566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXL2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 2, binterpt[0], binterp[0])); 2971f80e275SEmil Constantinescu } 2981f80e275SEmil Constantinescu { 299da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3009371c9d4SSatish Balay const PetscReal A[3][3] = 3019371c9d4SSatish Balay { 3029371c9d4SSatish Balay {0, 0, 0}, 303da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3049371c9d4SSatish Balay {0.5, 0.5, 0} 3059371c9d4SSatish Balay }, 3069371c9d4SSatish 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}}; 3079566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2C, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3081f80e275SEmil Constantinescu } 3091f80e275SEmil Constantinescu { 310da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3119371c9d4SSatish Balay const PetscReal A[3][3] = 3129371c9d4SSatish Balay { 3139371c9d4SSatish Balay {0, 0, 0}, 314da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3159371c9d4SSatish Balay {0.75, 0.25, 0} 3169371c9d4SSatish Balay }, 3179371c9d4SSatish 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}}; 3189566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2D, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3198a381b04SJed Brown } 32006db7b1cSJed Brown { /* Optimal for linear implicit part */ 321da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3229371c9d4SSatish Balay const PetscReal A[3][3] = 3239371c9d4SSatish Balay { 3249371c9d4SSatish Balay {0, 0, 0}, 325da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3269371c9d4SSatish Balay {(3 - 2 * 1.414213562373095048802) / 6, (3 + 2 * 1.414213562373095048802) / 6, 0} 3279371c9d4SSatish Balay }, 3289371c9d4SSatish 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}}; 3299566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2E, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 330a3a57f36SJed Brown } 3316cf0794eSJed Brown { /* Optimal for linear implicit part */ 3329371c9d4SSatish Balay const PetscReal A[3][3] = 3339371c9d4SSatish Balay { 3349371c9d4SSatish Balay {0, 0, 0}, 3356cf0794eSJed Brown {0.5, 0, 0}, 3369371c9d4SSatish Balay {0.5, 0.5, 0} 3379371c9d4SSatish Balay }, 3389371c9d4SSatish Balay At[3][3] = {{0.25, 0, 0}, {0, 0.25, 0}, {1. / 3, 1. / 3, 1. / 3}}; 3399566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXPRSSP2, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, NULL, NULL, 0, NULL, NULL)); 3406cf0794eSJed Brown } 341a3a57f36SJed Brown { 3429371c9d4SSatish Balay const PetscReal A[4][4] = 3439371c9d4SSatish Balay { 3449371c9d4SSatish Balay {0, 0, 0, 0}, 3454040e9f2SJed Brown {1767732205903. / 2027836641118., 0, 0, 0}, 3464040e9f2SJed Brown {5535828885825. / 10492691773637., 788022342437. / 10882634858940., 0, 0}, 3479371c9d4SSatish Balay {6485989280629. / 16251701735622., -4246266847089. / 9704473918619., 10755448449292. / 10357097424841., 0} 3489371c9d4SSatish Balay }, 3499371c9d4SSatish 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.}}; 3509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX3, 3, 4, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 351a3a57f36SJed Brown } 352a3a57f36SJed Brown { 3539371c9d4SSatish Balay const PetscReal A[5][5] = 3549371c9d4SSatish Balay { 3559371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3566cf0794eSJed Brown {1. / 2, 0, 0, 0, 0}, 3576cf0794eSJed Brown {11. / 18, 1. / 18, 0, 0, 0}, 3586cf0794eSJed Brown {5. / 6, -5. / 6, .5, 0, 0}, 3599371c9d4SSatish Balay {1. / 4, 7. / 4, 3. / 4, -7. / 4, 0} 3609371c9d4SSatish Balay }, 3619371c9d4SSatish 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; 3629566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS443, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3636cf0794eSJed Brown } 3646cf0794eSJed Brown { 3659371c9d4SSatish Balay const PetscReal A[5][5] = 3669371c9d4SSatish Balay { 3679371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3686cf0794eSJed Brown {1, 0, 0, 0, 0}, 3696cf0794eSJed Brown {4. / 9, 2. / 9, 0, 0, 0}, 3706cf0794eSJed Brown {1. / 4, 0, 3. / 4, 0, 0}, 3719371c9d4SSatish Balay {1. / 4, 0, 3. / 5, 0, 0} 3729371c9d4SSatish Balay }, 3739371c9d4SSatish 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; 3749566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXBPR3, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3756cf0794eSJed Brown } 3766cf0794eSJed Brown { 3779371c9d4SSatish Balay const PetscReal A[6][6] = 3789371c9d4SSatish Balay { 3799371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 380a3a57f36SJed Brown {1. / 2, 0, 0, 0, 0, 0}, 3814040e9f2SJed Brown {13861. / 62500., 6889. / 62500., 0, 0, 0, 0}, 3824040e9f2SJed Brown {-116923316275. / 2393684061468., -2731218467317. / 15368042101831., 9408046702089. / 11113171139209., 0, 0, 0}, 3834040e9f2SJed Brown {-451086348788. / 2902428689909., -2682348792572. / 7519795681897., 12662868775082. / 11960479115383., 3355817975965. / 11060851509271., 0, 0}, 3849371c9d4SSatish Balay {647845179188. / 3216320057751., 73281519250. / 8382639484533., 552539513391. / 3454668386233., 3354512671639. / 8306763924573., 4040. / 17871., 0} 3859371c9d4SSatish Balay }, 3869371c9d4SSatish 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}, 3879371c9d4SSatish Balay {7640104374378. / 9702883013639., -11436875. / 14766696., 2173542590792. / 12501825683035.}, {-20649996744609. / 7521556579894., 174696575. / 18121608., -31592104683404. / 5083833661969.}, 3889371c9d4SSatish Balay {8854892464581. / 2390941311638., -12120380. / 966161., 61146701046299. / 7138195549469.}, {-11397109935349. / 6675773540249., 3843. / 706., -17219254887155. / 4939391667607.}}; 3899566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX4, 4, 6, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 390a3a57f36SJed Brown } 391a3a57f36SJed Brown { 3929371c9d4SSatish Balay const PetscReal A[8][8] = 3939371c9d4SSatish Balay { 3949371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 395a3a57f36SJed Brown {41. / 100, 0, 0, 0, 0, 0, 0, 0}, 3964040e9f2SJed Brown {367902744464. / 2072280473677., 677623207551. / 8224143866563., 0, 0, 0, 0, 0, 0}, 3974040e9f2SJed Brown {1268023523408. / 10340822734521., 0, 1029933939417. / 13636558850479., 0, 0, 0, 0, 0}, 3984040e9f2SJed Brown {14463281900351. / 6315353703477., 0, 66114435211212. / 5879490589093., -54053170152839. / 4284798021562., 0, 0, 0, 0}, 3994040e9f2SJed Brown {14090043504691. / 34967701212078., 0, 15191511035443. / 11219624916014., -18461159152457. / 12425892160975., -281667163811. / 9011619295870., 0, 0, 0}, 4004040e9f2SJed Brown {19230459214898. / 13134317526959., 0, 21275331358303. / 2942455364971., -38145345988419. / 4862620318723., -1. / 8, -1. / 8, 0, 0}, 4019371c9d4SSatish Balay {-19977161125411. / 11928030595625., 0, -40795976796054. / 6384907823539., 177454434618887. / 12078138498510., 782672205425. / 8267701900261., -69563011059811. / 9646580694205., 7356628210526. / 4942186776405., 0} 4029371c9d4SSatish Balay }, 4039371c9d4SSatish 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.}}; 4049566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX5, 5, 8, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 405a3a57f36SJed Brown } 4068a381b04SJed Brown PetscFunctionReturn(0); 4078a381b04SJed Brown } 4088a381b04SJed Brown 4098a381b04SJed Brown /*@C 4108a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 4118a381b04SJed Brown 4128a381b04SJed Brown Not Collective 4138a381b04SJed Brown 4148a381b04SJed Brown Level: advanced 4158a381b04SJed Brown 416db781477SPatrick Sanan .seealso: `TSARKIMEXRegister()`, `TSARKIMEXRegisterAll()` 4178a381b04SJed Brown @*/ 4189371c9d4SSatish Balay PetscErrorCode TSARKIMEXRegisterDestroy(void) { 4198a381b04SJed Brown ARKTableauLink link; 4208a381b04SJed Brown 4218a381b04SJed Brown PetscFunctionBegin; 4228a381b04SJed Brown while ((link = ARKTableauList)) { 4238a381b04SJed Brown ARKTableau t = &link->tab; 4248a381b04SJed Brown ARKTableauList = link->next; 4259566063dSJacob Faibussowitsch PetscCall(PetscFree6(t->At, t->bt, t->ct, t->A, t->b, t->c)); 4269566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembedt, t->bembed)); 4279566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->binterpt, t->binterp)); 4289566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 4299566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 4308a381b04SJed Brown } 4318a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4328a381b04SJed Brown PetscFunctionReturn(0); 4338a381b04SJed Brown } 4348a381b04SJed Brown 4358a381b04SJed Brown /*@C 4368a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 4378a690491SBarry Smith from TSInitializePackage(). 4388a381b04SJed Brown 4398a381b04SJed Brown Level: developer 4408a381b04SJed Brown 441db781477SPatrick Sanan .seealso: `PetscInitialize()` 4428a381b04SJed Brown @*/ 4439371c9d4SSatish Balay PetscErrorCode TSARKIMEXInitializePackage(void) { 4448a381b04SJed Brown PetscFunctionBegin; 4458a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 4468a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 4479566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterAll()); 4489566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSARKIMEXFinalizePackage)); 4498a381b04SJed Brown PetscFunctionReturn(0); 4508a381b04SJed Brown } 4518a381b04SJed Brown 4528a381b04SJed Brown /*@C 4538a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 4548a381b04SJed Brown called from PetscFinalize(). 4558a381b04SJed Brown 4568a381b04SJed Brown Level: developer 4578a381b04SJed Brown 458db781477SPatrick Sanan .seealso: `PetscFinalize()` 4598a381b04SJed Brown @*/ 4609371c9d4SSatish Balay PetscErrorCode TSARKIMEXFinalizePackage(void) { 4618a381b04SJed Brown PetscFunctionBegin; 4628a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 4639566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterDestroy()); 4648a381b04SJed Brown PetscFunctionReturn(0); 4658a381b04SJed Brown } 4668a381b04SJed Brown 467cd652676SJed Brown /*@C 468cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 469cd652676SJed Brown 470cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 471cd652676SJed Brown 472cd652676SJed Brown Input Parameters: 473cd652676SJed Brown + name - identifier for method 474cd652676SJed Brown . order - approximation order of method 475cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 476cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 4770298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 4780298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 479cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 4800298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 4810298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 4820298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 4830298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 484cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 485cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 4860298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 487cd652676SJed Brown 488cd652676SJed Brown Notes: 489cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 490cd652676SJed Brown 491cd652676SJed Brown Level: advanced 492cd652676SJed Brown 493db781477SPatrick Sanan .seealso: `TSARKIMEX` 494cd652676SJed Brown @*/ 4959371c9d4SSatish Balay 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[]) { 4968a381b04SJed Brown ARKTableauLink link; 4978a381b04SJed Brown ARKTableau t; 4988a381b04SJed Brown PetscInt i, j; 4998a381b04SJed Brown 5008a381b04SJed Brown PetscFunctionBegin; 5019566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 5029566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 5038a381b04SJed Brown t = &link->tab; 5049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 5058a381b04SJed Brown t->order = order; 5068a381b04SJed Brown t->s = s; 5079566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(s * s, &t->At, s, &t->bt, s, &t->ct, s * s, &t->A, s, &t->b, s, &t->c)); 5089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->At, At, s * s)); 5099566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 5109566063dSJacob Faibussowitsch if (bt) PetscCall(PetscArraycpy(t->bt, bt, s)); 5119371c9d4SSatish Balay else 5129371c9d4SSatish Balay for (i = 0; i < s; i++) t->bt[i] = At[(s - 1) * s + i]; 5139566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 5149371c9d4SSatish Balay else 5159371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = t->bt[i]; 5169566063dSJacob Faibussowitsch if (ct) PetscCall(PetscArraycpy(t->ct, ct, s)); 5179371c9d4SSatish Balay else 5189371c9d4SSatish Balay for (i = 0; i < s; i++) 5199371c9d4SSatish Balay for (j = 0, t->ct[i] = 0; j < s; j++) t->ct[i] += At[i * s + j]; 5209566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 5219371c9d4SSatish Balay else 5229371c9d4SSatish Balay for (i = 0; i < s; i++) 5239371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 524e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 5259371c9d4SSatish Balay for (i = 0; i < s; i++) 5269371c9d4SSatish Balay if (t->At[(s - 1) * s + i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 527e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 5289371c9d4SSatish Balay for (i = 0; i < s; i++) 5299371c9d4SSatish Balay if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 530e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5314e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 532108c343cSJed Brown if (bembedt) { 5339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &t->bembedt, s, &t->bembed)); 5349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembedt, bembedt, s)); 5359566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed ? bembed : bembedt, s)); 536108c343cSJed Brown } 537108c343cSJed Brown 5384f385281SJed Brown t->pinterp = pinterp; 5399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s * pinterp, &t->binterpt, s * pinterp, &t->binterp)); 5409566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 5419566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp ? binterp : binterpt, s * pinterp)); 5428a381b04SJed Brown link->next = ARKTableauList; 5438a381b04SJed Brown ARKTableauList = link; 5448a381b04SJed Brown PetscFunctionReturn(0); 5458a381b04SJed Brown } 5468a381b04SJed Brown 547108c343cSJed Brown /* 548108c343cSJed Brown The step completion formula is 549108c343cSJed Brown 550108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 551108c343cSJed Brown 552108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 553108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 554108c343cSJed Brown We can write 555108c343cSJed Brown 556108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 557108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 558108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 559108c343cSJed Brown 560108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 561108c343cSJed Brown */ 5629371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts, PetscInt order, Vec X, PetscBool *done) { 563108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 564108c343cSJed Brown ARKTableau tab = ark->tableau; 565108c343cSJed Brown PetscScalar *w = ark->work; 566108c343cSJed Brown PetscReal h; 567108c343cSJed Brown PetscInt s = tab->s, j; 568108c343cSJed Brown 569108c343cSJed Brown PetscFunctionBegin; 570108c343cSJed Brown switch (ark->status) { 571108c343cSJed Brown case TS_STEP_INCOMPLETE: 5729371c9d4SSatish Balay case TS_STEP_PENDING: h = ts->time_step; break; 5739371c9d4SSatish Balay case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break; 574ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 575108c343cSJed Brown } 576108c343cSJed Brown if (order == tab->order) { 577e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 578740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 5799566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[s - 1], X)); 580e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 5819566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 582e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bt[j]; 5839566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 584e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 585108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->b[j]; 5869566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 587e817cc15SEmil Constantinescu } 588e817cc15SEmil Constantinescu } 5899566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 590108c343cSJed Brown if (done) *done = PETSC_TRUE; 591108c343cSJed Brown PetscFunctionReturn(0); 592108c343cSJed Brown } else if (order == tab->order - 1) { 593108c343cSJed Brown if (!tab->bembedt) goto unavailable; 594108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 5959566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 596e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembedt[j]; 5979566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 598108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 5999566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 600108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 6019566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 602e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembedt[j] - tab->bt[j]); 6039566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, tab->s, w, ark->YdotI)); 604108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 6059566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 606108c343cSJed Brown } 607108c343cSJed Brown if (done) *done = PETSC_TRUE; 608108c343cSJed Brown PetscFunctionReturn(0); 609108c343cSJed Brown } 610108c343cSJed Brown unavailable: 611108c343cSJed Brown if (done) *done = PETSC_FALSE; 6129371c9d4SSatish Balay else 6139371c9d4SSatish 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, 6149371c9d4SSatish Balay tab->order, order); 615108c343cSJed Brown PetscFunctionReturn(0); 616108c343cSJed Brown } 617108c343cSJed Brown 6189371c9d4SSatish Balay static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts, PetscBool *id) { 619c79dcfbdSBarry Smith Vec Udot, Y1, Y2; 620c79dcfbdSBarry Smith TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 621c79dcfbdSBarry Smith PetscReal norm; 622c79dcfbdSBarry Smith 623c79dcfbdSBarry Smith PetscFunctionBegin; 6249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Udot)); 6259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y1)); 6269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y2)); 6279566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y1, ark->imex)); 6289566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Udot, NULL)); 6299566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y2, ark->imex)); 6309566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Y1)); 6319566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Udot)); 6329566063dSJacob Faibussowitsch PetscCall(VecNorm(Y2, NORM_2, &norm)); 633c79dcfbdSBarry Smith if (norm < 100.0 * PETSC_MACHINE_EPSILON) { 634c79dcfbdSBarry Smith *id = PETSC_TRUE; 635c79dcfbdSBarry Smith } else { 636c79dcfbdSBarry Smith *id = PETSC_FALSE; 6379566063dSJacob 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)); 638c79dcfbdSBarry Smith } 6399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Udot)); 6409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 6419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y2)); 642c79dcfbdSBarry Smith PetscFunctionReturn(0); 643c79dcfbdSBarry Smith } 644c79dcfbdSBarry Smith 6459371c9d4SSatish Balay static PetscErrorCode TSRollBack_ARKIMEX(TS ts) { 64624655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 64724655328SShri ARKTableau tab = ark->tableau; 64824655328SShri const PetscInt s = tab->s; 64924655328SShri const PetscReal *bt = tab->bt, *b = tab->b; 65024655328SShri PetscScalar *w = ark->work; 65124655328SShri Vec *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS; 65224655328SShri PetscInt j; 653be5899b3SLisandro Dalcin PetscReal h; 65424655328SShri 65524655328SShri PetscFunctionBegin; 656be5899b3SLisandro Dalcin switch (ark->status) { 657be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 6589371c9d4SSatish Balay case TS_STEP_PENDING: h = ts->time_step; break; 6599371c9d4SSatish Balay case TS_STEP_COMPLETE: h = ts->ptime - ts->ptime_prev; break; 660be5899b3SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 661be5899b3SLisandro Dalcin } 66224655328SShri for (j = 0; j < s; j++) w[j] = -h * bt[j]; 6639566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotI)); 66424655328SShri for (j = 0; j < s; j++) w[j] = -h * b[j]; 6659566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 66624655328SShri PetscFunctionReturn(0); 66724655328SShri } 66824655328SShri 6699371c9d4SSatish Balay static PetscErrorCode TSStep_ARKIMEX(TS ts) { 6708a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 6718a381b04SJed Brown ARKTableau tab = ark->tableau; 6728a381b04SJed Brown const PetscInt s = tab->s; 67324655328SShri const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct, *c = tab->c; 674406d0ec2SJed Brown PetscScalar *w = ark->work; 6751297b224SEmil Constantinescu Vec *Y = ark->Y, *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS, Ydot = ark->Ydot, Ydot0 = ark->Ydot0, Z = ark->Z; 67696400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 677108c343cSJed Brown TSAdapt adapt; 6788a381b04SJed Brown SNES snes; 679fecfb714SLisandro Dalcin PetscInt i, j, its, lits; 680be5899b3SLisandro Dalcin PetscInt rejections = 0; 68196400bd6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 68296400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 6838a381b04SJed Brown 6848a381b04SJed Brown PetscFunctionBegin; 68596400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 6869566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 6879566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 6889566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 68996400bd6SLisandro Dalcin } 69096400bd6SLisandro Dalcin 69196400bd6SLisandro Dalcin if (!ts->steprollback) { 69296400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 6939566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[s - 1], Ydot0)); 69496400bd6SLisandro Dalcin } 695fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 69696400bd6SLisandro Dalcin for (i = 0; i < s; i++) { 6979566063dSJacob Faibussowitsch PetscCall(VecCopy(Y[i], ark->Y_prev[i])); 6989566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotRHS[i], ark->YdotRHS_prev[i])); 6999566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[i], ark->YdotI_prev[i])); 70096400bd6SLisandro Dalcin } 70196400bd6SLisandro Dalcin } 70296400bd6SLisandro Dalcin } 70396400bd6SLisandro Dalcin 704fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 70596400bd6SLisandro Dalcin TS ts_start; 706c79dcfbdSBarry Smith if (PetscDefined(USE_DEBUG)) { 707c79dcfbdSBarry Smith PetscBool id = PETSC_FALSE; 7089566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTestMassIdentity(ts, &id)); 7093c633725SBarry 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"); 710c79dcfbdSBarry Smith } 7119566063dSJacob Faibussowitsch PetscCall(TSClone(ts, &ts_start)); 7129566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts_start, ts->vec_sol)); 7139566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts_start, ts->ptime)); 7149566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts_start, ts->steps + 1)); 7159566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts_start, ts->ptime + ts->time_step)); 7169566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts_start, TS_EXACTFINALTIME_STEPOVER)); 7179566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts_start, ts->time_step)); 7189566063dSJacob Faibussowitsch PetscCall(TSSetType(ts_start, TSARKIMEX)); 7199566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts_start, PETSC_TRUE)); 7209566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts_start, TSARKIMEX1BEE)); 72134561852SEmil Constantinescu 7229566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts_start)); 7239566063dSJacob Faibussowitsch PetscCall(TSSolve(ts_start, ts->vec_sol)); 7249566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts_start, &ts->ptime)); 7259566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts_start, &ts->time_step)); 726bbd56ea5SKarl Rupp 72785fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 72885fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX *)ts_start->data; 7299566063dSJacob Faibussowitsch PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s - 1], Ydot0)); 73085fc7851SLisandro Dalcin } 73196400bd6SLisandro Dalcin ts->steps++; 73234561852SEmil Constantinescu 733d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 734d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 73596400bd6SLisandro Dalcin { 7369566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7379566063dSJacob Faibussowitsch PetscCall(TSSetSNES(ts, snes)); 73896400bd6SLisandro Dalcin } 7399566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts_start)); 740e817cc15SEmil Constantinescu } 741e817cc15SEmil Constantinescu 742108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 74396400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 74496400bd6SLisandro Dalcin PetscReal t = ts->ptime; 745108c343cSJed Brown PetscReal h = ts->time_step; 7468a381b04SJed Brown for (i = 0; i < s; i++) { 7479be3e283SDebojyoti Ghosh ark->stage_time = t + h * ct[i]; 7489566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ark->stage_time)); 7498a381b04SJed Brown if (At[i * s + i] == 0) { /* This stage is explicit */ 7503c633725SBarry 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"); 7519566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 752e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7539566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotI)); 7548a381b04SJed Brown for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 7568a381b04SJed Brown } else { 757b296d7d5SJed Brown ark->scoeff = 1. / At[i * s + i]; 7588a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 7599566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Z)); 760e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7619566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotI)); 762c58d1302SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7639566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotRHS)); 764fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 76556dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 7669566063dSJacob Faibussowitsch PetscCall(TSExtrapolate_ARKIMEX(ts, c[i], Y[i])); 76756dcabbaSDebojyoti Ghosh } else { 7688a381b04SJed Brown /* Initial guess taken from last stage */ 7699566063dSJacob Faibussowitsch PetscCall(VecCopy(i > 0 ? Y[i - 1] : ts->vec_sol, Y[i])); 77056dcabbaSDebojyoti Ghosh } 7719566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7729566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, Y[i])); 7739566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 7749566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 7759371c9d4SSatish Balay ts->snes_its += its; 7769371c9d4SSatish Balay ts->ksp_its += lits; 7779566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 7789566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, ark->stage_time, Y[i], &stageok)); 77996400bd6SLisandro Dalcin if (!stageok) { 7801be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 7811be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 78296400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 7831be93e3eSJed Brown goto reject_step; 7841be93e3eSJed Brown } 7858a381b04SJed Brown } 786e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 787e817cc15SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 7889371c9d4SSatish 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", 7899371c9d4SSatish Balay ark->tableau->name); 7909566063dSJacob Faibussowitsch PetscCall(VecCopy(Ydot0, YdotI[0])); /* YdotI = YdotI(tn-1) */ 791e817cc15SEmil Constantinescu } else { 7929566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 793e817cc15SEmil Constantinescu } 794e817cc15SEmil Constantinescu } else { 7955eca1a21SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 7969566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); 7979566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, t + h * ct[i], Y[i], Ydot, YdotI[i], ark->imex)); /* YdotI = -G(t,Y,0) */ 7989566063dSJacob Faibussowitsch PetscCall(VecScale(YdotI[i], -1.0)); 7995eca1a21SEmil Constantinescu } else { 8009566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 8015eca1a21SEmil Constantinescu } 8024cc180ffSJed Brown if (ark->imex) { 8039566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 8044cc180ffSJed Brown } else { 8059566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(YdotRHS[i])); 8064cc180ffSJed Brown } 8078a381b04SJed Brown } 8089566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ark->stage_time, i, Y)); 809e817cc15SEmil Constantinescu } 81096400bd6SLisandro Dalcin 811be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 8129566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_ARKIMEX(ts, tab->order, ts->vec_sol, NULL)); 813108c343cSJed Brown ark->status = TS_STEP_PENDING; 8149566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8159566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8169566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8179566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 81896400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 81996400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8209566063dSJacob Faibussowitsch PetscCall(TSRollBack_ARKIMEX(ts)); 821be5899b3SLisandro Dalcin ts->time_step = next_time_step; 822be5899b3SLisandro Dalcin goto reject_step; 82396400bd6SLisandro Dalcin } 82496400bd6SLisandro Dalcin 8258a381b04SJed Brown ts->ptime += ts->time_step; 826cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 827108c343cSJed Brown break; 82896400bd6SLisandro Dalcin 82996400bd6SLisandro Dalcin reject_step: 8309371c9d4SSatish Balay ts->reject++; 8319371c9d4SSatish Balay accept = PETSC_FALSE; 832be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 83396400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 83463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 835108c343cSJed Brown } 836f85781f1SEmil Constantinescu } 8378a381b04SJed Brown PetscFunctionReturn(0); 8388a381b04SJed Brown } 8398a381b04SJed Brown 8409371c9d4SSatish Balay static PetscErrorCode TSInterpolate_ARKIMEX(TS ts, PetscReal itime, Vec X) { 841cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 8424f385281SJed Brown PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 843108c343cSJed Brown PetscReal h; 844108c343cSJed Brown PetscReal tt, t; 845cd652676SJed Brown PetscScalar *bt, *b; 846cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 847cd652676SJed Brown 848cd652676SJed Brown PetscFunctionBegin; 8493c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 850108c343cSJed Brown switch (ark->status) { 851108c343cSJed Brown case TS_STEP_INCOMPLETE: 852108c343cSJed Brown case TS_STEP_PENDING: 853108c343cSJed Brown h = ts->time_step; 854108c343cSJed Brown t = (itime - ts->ptime) / h; 855108c343cSJed Brown break; 856108c343cSJed Brown case TS_STEP_COMPLETE: 857be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 858108c343cSJed Brown t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 859108c343cSJed Brown break; 860ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 861108c343cSJed Brown } 8629566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &bt, s, &b)); 863cd652676SJed Brown for (i = 0; i < s; i++) bt[i] = b[i] = 0; 8644f385281SJed Brown for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 865cd652676SJed Brown for (i = 0; i < s; i++) { 866c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 867108c343cSJed Brown b[i] += h * B[i * pinterp + j] * tt; 868cd652676SJed Brown } 869cd652676SJed Brown } 8709566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[0], X)); 8719566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI)); 8729566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS)); 8739566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 874cd652676SJed Brown PetscFunctionReturn(0); 875cd652676SJed Brown } 876cd652676SJed Brown 8779371c9d4SSatish Balay static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts, PetscReal c, Vec X) { 87856dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 87956dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 880be5899b3SLisandro Dalcin PetscReal h, h_prev, t, tt; 88156dcabbaSDebojyoti Ghosh PetscScalar *bt, *b; 88256dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 88356dcabbaSDebojyoti Ghosh 88456dcabbaSDebojyoti Ghosh PetscFunctionBegin; 8853c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 8869566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(s, &bt, s, &b)); 88781d12688SDebojyoti Ghosh h = ts->time_step; 888be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 889be5899b3SLisandro Dalcin t = 1 + h / h_prev * c; 89056dcabbaSDebojyoti Ghosh for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 89156dcabbaSDebojyoti Ghosh for (i = 0; i < s; i++) { 89281d12688SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 89356dcabbaSDebojyoti Ghosh b[i] += h * B[i * pinterp + j] * tt; 89456dcabbaSDebojyoti Ghosh } 89556dcabbaSDebojyoti Ghosh } 8963c633725SBarry Smith PetscCheck(ark->Y_prev, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Stages from previous step have not been stored"); 8979566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y_prev[0], X)); 8989566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI_prev)); 8999566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS_prev)); 9009566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 90156dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 90256dcabbaSDebojyoti Ghosh } 90356dcabbaSDebojyoti Ghosh 9048a381b04SJed Brown /*------------------------------------------------------------*/ 90596400bd6SLisandro Dalcin 9069371c9d4SSatish Balay static PetscErrorCode TSARKIMEXTableauReset(TS ts) { 90796400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 90896400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 90996400bd6SLisandro Dalcin 91096400bd6SLisandro Dalcin PetscFunctionBegin; 91196400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 9129566063dSJacob Faibussowitsch PetscCall(PetscFree(ark->work)); 9139566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y)); 9149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI)); 9159566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS)); 9169566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev)); 9179566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev)); 9189566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev)); 91996400bd6SLisandro Dalcin PetscFunctionReturn(0); 92096400bd6SLisandro Dalcin } 92196400bd6SLisandro Dalcin 9229371c9d4SSatish Balay static PetscErrorCode TSReset_ARKIMEX(TS ts) { 9238a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 9248a381b04SJed Brown 9258a381b04SJed Brown PetscFunctionBegin; 9269566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauReset(ts)); 9279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot)); 9289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot0)); 9299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Z)); 9308a381b04SJed Brown PetscFunctionReturn(0); 9318a381b04SJed Brown } 9328a381b04SJed Brown 9339371c9d4SSatish Balay static PetscErrorCode TSARKIMEXGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) { 934d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX *)ts->data; 935d5e6173cSPeter Brune 936d5e6173cSPeter Brune PetscFunctionBegin; 937d5e6173cSPeter Brune if (Z) { 938d5e6173cSPeter Brune if (dm && dm != ts->dm) { 9399566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 940d5e6173cSPeter Brune } else *Z = ax->Z; 941d5e6173cSPeter Brune } 942d5e6173cSPeter Brune if (Ydot) { 943d5e6173cSPeter Brune if (dm && dm != ts->dm) { 9449566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 945d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 946d5e6173cSPeter Brune } 947d5e6173cSPeter Brune PetscFunctionReturn(0); 948d5e6173cSPeter Brune } 949d5e6173cSPeter Brune 9509371c9d4SSatish Balay static PetscErrorCode TSARKIMEXRestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) { 951d5e6173cSPeter Brune PetscFunctionBegin; 952d5e6173cSPeter Brune if (Z) { 95348a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 954d5e6173cSPeter Brune } 955d5e6173cSPeter Brune if (Ydot) { 95648a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 957d5e6173cSPeter Brune } 958d5e6173cSPeter Brune PetscFunctionReturn(0); 959d5e6173cSPeter Brune } 960d5e6173cSPeter Brune 9618a381b04SJed Brown /* 9628a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 9638a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 9648a381b04SJed Brown */ 9659371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes, Vec X, Vec F, TS ts) { 9668a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 967d5e6173cSPeter Brune DM dm, dmsave; 968d5e6173cSPeter Brune Vec Z, Ydot; 969b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 9708a381b04SJed Brown 9718a381b04SJed Brown PetscFunctionBegin; 9729566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 9739566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, &Ydot)); 9749566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Ydot, -shift, shift, 0, Z, X)); /* Ydot = shift*(X-Z) */ 975d5e6173cSPeter Brune dmsave = ts->dm; 976d5e6173cSPeter Brune ts->dm = dm; 977740132f1SEmil Constantinescu 9789566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ark->stage_time, X, Ydot, F, ark->imex)); 979e817cc15SEmil Constantinescu 980d5e6173cSPeter Brune ts->dm = dmsave; 9819566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, &Ydot)); 9828a381b04SJed Brown PetscFunctionReturn(0); 9838a381b04SJed Brown } 9848a381b04SJed Brown 9859371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes, Vec X, Mat A, Mat B, TS ts) { 9868a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 987d5e6173cSPeter Brune DM dm, dmsave; 988d5e6173cSPeter Brune Vec Ydot; 989b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 9908a381b04SJed Brown 9918a381b04SJed Brown PetscFunctionBegin; 9929566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 9939566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, NULL, &Ydot)); 9948a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 995d5e6173cSPeter Brune dmsave = ts->dm; 996d5e6173cSPeter Brune ts->dm = dm; 997740132f1SEmil Constantinescu 9989566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ark->stage_time, X, Ydot, shift, A, B, ark->imex)); 999740132f1SEmil Constantinescu 1000d5e6173cSPeter Brune ts->dm = dmsave; 10019566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, NULL, &Ydot)); 1002d5e6173cSPeter Brune PetscFunctionReturn(0); 1003d5e6173cSPeter Brune } 1004d5e6173cSPeter Brune 10059371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine, DM coarse, void *ctx) { 1006d5e6173cSPeter Brune PetscFunctionBegin; 1007d5e6173cSPeter Brune PetscFunctionReturn(0); 1008d5e6173cSPeter Brune } 1009d5e6173cSPeter Brune 10109371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 1011d5e6173cSPeter Brune TS ts = (TS)ctx; 1012d5e6173cSPeter Brune Vec Z, Z_c; 1013d5e6173cSPeter Brune 1014d5e6173cSPeter Brune PetscFunctionBegin; 10159566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, fine, &Z, NULL)); 10169566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, coarse, &Z_c, NULL)); 10179566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Z, Z_c)); 10189566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Z_c, rscale, Z_c)); 10199566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, fine, &Z, NULL)); 10209566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, coarse, &Z_c, NULL)); 10218a381b04SJed Brown PetscFunctionReturn(0); 10228a381b04SJed Brown } 10238a381b04SJed Brown 10249371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm, DM subdm, void *ctx) { 1025cdb298fcSPeter Brune PetscFunctionBegin; 1026cdb298fcSPeter Brune PetscFunctionReturn(0); 1027cdb298fcSPeter Brune } 1028cdb298fcSPeter Brune 10299371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 1030cdb298fcSPeter Brune TS ts = (TS)ctx; 1031cdb298fcSPeter Brune Vec Z, Z_c; 1032cdb298fcSPeter Brune 1033cdb298fcSPeter Brune PetscFunctionBegin; 10349566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, NULL)); 10359566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, subdm, &Z_c, NULL)); 1036cdb298fcSPeter Brune 10379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 10389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 1039cdb298fcSPeter Brune 10409566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, NULL)); 10419566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, subdm, &Z_c, NULL)); 1042cdb298fcSPeter Brune PetscFunctionReturn(0); 1043cdb298fcSPeter Brune } 1044cdb298fcSPeter Brune 10459371c9d4SSatish Balay static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) { 104696400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 104796400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 104896400bd6SLisandro Dalcin 104996400bd6SLisandro Dalcin PetscFunctionBegin; 10509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &ark->work)); 10519566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y)); 10529566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI)); 10539566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS)); 105496400bd6SLisandro Dalcin if (ark->extrapolate) { 10559566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 10569566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 10579566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 105896400bd6SLisandro Dalcin } 105996400bd6SLisandro Dalcin PetscFunctionReturn(0); 106096400bd6SLisandro Dalcin } 106196400bd6SLisandro Dalcin 10629371c9d4SSatish Balay static PetscErrorCode TSSetUp_ARKIMEX(TS ts) { 10638a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1064d5e6173cSPeter Brune DM dm; 106596400bd6SLisandro Dalcin SNES snes; 1066f9c1d6abSBarry Smith 10678a381b04SJed Brown PetscFunctionBegin; 10689566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauSetUp(ts)); 10699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot)); 10709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot0)); 10719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Z)); 10729566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 10739566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 10749566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 10759566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 10768a381b04SJed Brown PetscFunctionReturn(0); 10778a381b04SJed Brown } 10788a381b04SJed Brown /*------------------------------------------------------------*/ 10798a381b04SJed Brown 10809371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts, PetscOptionItems *PetscOptionsObject) { 10814cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 10828a381b04SJed Brown 10838a381b04SJed Brown PetscFunctionBegin; 1084d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ARKIMEX ODE solver options"); 10858a381b04SJed Brown { 10868a381b04SJed Brown ARKTableauLink link; 10878a381b04SJed Brown PetscInt count, choice; 10888a381b04SJed Brown PetscBool flg; 10898a381b04SJed Brown const char **namelist; 10909371c9d4SSatish Balay for (link = ARKTableauList, count = 0; link; link = link->next, count++) 10919371c9d4SSatish Balay ; 10929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 10938a381b04SJed Brown for (link = ARKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 10949566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_arkimex_type", "Family of ARK IMEX method", "TSARKIMEXSetType", (const char *const *)namelist, count, ark->tableau->name, &choice, &flg)); 10959566063dSJacob Faibussowitsch if (flg) PetscCall(TSARKIMEXSetType(ts, namelist[choice])); 10969566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 109796400bd6SLisandro Dalcin 10984cc180ffSJed Brown flg = (PetscBool)!ark->imex; 10999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_arkimex_fully_implicit", "Solve the problem fully implicitly", "TSARKIMEXSetFullyImplicit", flg, &flg, NULL)); 11004cc180ffSJed Brown ark->imex = (PetscBool)!flg; 11019566063dSJacob 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)); 11028a381b04SJed Brown } 1103d0609cedSBarry Smith PetscOptionsHeadEnd(); 11048a381b04SJed Brown PetscFunctionReturn(0); 11058a381b04SJed Brown } 11068a381b04SJed Brown 11079371c9d4SSatish Balay static PetscErrorCode TSView_ARKIMEX(TS ts, PetscViewer viewer) { 11088a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 11098a381b04SJed Brown PetscBool iascii; 11108a381b04SJed Brown 11118a381b04SJed Brown PetscFunctionBegin; 11129566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 11138a381b04SJed Brown if (iascii) { 11149c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 111519fd82e9SBarry Smith TSARKIMEXType arktype; 11168a381b04SJed Brown char buf[512]; 11173a28c0e5SStefano Zampini PetscBool flg; 11183a28c0e5SStefano Zampini 11199566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetType(ts, &arktype)); 11209566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetFullyImplicit(ts, &flg)); 11219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ARK IMEX %s\n", arktype)); 11229566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ct)); 11239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Stiff abscissa ct = %s\n", buf)); 11249566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c)); 11259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Fully implicit: %s\n", flg ? "yes" : "no")); 11269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s\n", tab->stiffly_accurate ? "yes" : "no")); 11279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Explicit first stage: %s\n", tab->explicit_first_stage ? "yes" : "no")); 11289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "FSAL property: %s\n", tab->FSAL_implicit ? "yes" : "no")); 11299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Nonstiff abscissa c = %s\n", buf)); 11308a381b04SJed Brown } 11318a381b04SJed Brown PetscFunctionReturn(0); 11328a381b04SJed Brown } 11338a381b04SJed Brown 11349371c9d4SSatish Balay static PetscErrorCode TSLoad_ARKIMEX(TS ts, PetscViewer viewer) { 1135f2c2a1b9SBarry Smith SNES snes; 11369c334d8fSLisandro Dalcin TSAdapt adapt; 1137f2c2a1b9SBarry Smith 1138f2c2a1b9SBarry Smith PetscFunctionBegin; 11399566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 11409566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 11419566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 11429566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes, viewer)); 1143ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 11449566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 11459566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 1146f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1147f2c2a1b9SBarry Smith } 1148f2c2a1b9SBarry Smith 11498a381b04SJed Brown /*@C 11508a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 11518a381b04SJed Brown 11528a381b04SJed Brown Logically collective 11538a381b04SJed Brown 1154d8d19677SJose E. Roman Input Parameters: 11558a381b04SJed Brown + ts - timestepping context 11568a381b04SJed Brown - arktype - type of ARK-IMEX scheme 11578a381b04SJed Brown 11589bd3e852SBarry Smith Options Database: 115967b8a455SSatish Balay . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> - set ARK IMEX scheme type 11609bd3e852SBarry Smith 11618a381b04SJed Brown Level: intermediate 11628a381b04SJed Brown 1163db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEX1BEE`, `TSARKIMEXA2`, `TSARKIMEXL2`, `TSARKIMEXARS122`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEXPRSSP2`, 1164db781477SPatrick Sanan `TSARKIMEX3`, `TSARKIMEXBPR3`, `TSARKIMEXARS443`, `TSARKIMEX4`, `TSARKIMEX5` 11658a381b04SJed Brown @*/ 11669371c9d4SSatish Balay PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType arktype) { 11678a381b04SJed Brown PetscFunctionBegin; 11688a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1169b92453a8SLisandro Dalcin PetscValidCharPointer(arktype, 2); 1170cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetType_C", (TS, TSARKIMEXType), (ts, arktype)); 11718a381b04SJed Brown PetscFunctionReturn(0); 11728a381b04SJed Brown } 11738a381b04SJed Brown 11748a381b04SJed Brown /*@C 11758a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 11768a381b04SJed Brown 11778a381b04SJed Brown Logically collective 11788a381b04SJed Brown 11798a381b04SJed Brown Input Parameter: 11808a381b04SJed Brown . ts - timestepping context 11818a381b04SJed Brown 11828a381b04SJed Brown Output Parameter: 11838a381b04SJed Brown . arktype - type of ARK-IMEX scheme 11848a381b04SJed Brown 11858a381b04SJed Brown Level: intermediate 11868a381b04SJed Brown 1187db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()` 11888a381b04SJed Brown @*/ 11899371c9d4SSatish Balay PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *arktype) { 11908a381b04SJed Brown PetscFunctionBegin; 11918a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1192cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetType_C", (TS, TSARKIMEXType *), (ts, arktype)); 11938a381b04SJed Brown PetscFunctionReturn(0); 11948a381b04SJed Brown } 11958a381b04SJed Brown 119616353aafSBarry Smith /*@ 11974cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 11984cc180ffSJed Brown 11994cc180ffSJed Brown Logically collective 12004cc180ffSJed Brown 1201d8d19677SJose E. Roman Input Parameters: 12024cc180ffSJed Brown + ts - timestepping context 12034cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 12044cc180ffSJed Brown 12054cc180ffSJed Brown Level: intermediate 12064cc180ffSJed Brown 1207db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEXGetFullyImplicit()` 12084cc180ffSJed Brown @*/ 12099371c9d4SSatish Balay PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts, PetscBool flg) { 12104cc180ffSJed Brown PetscFunctionBegin; 12114cc180ffSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 12123a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts, flg, 2); 1213cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetFullyImplicit_C", (TS, PetscBool), (ts, flg)); 12144cc180ffSJed Brown PetscFunctionReturn(0); 12154cc180ffSJed Brown } 12164cc180ffSJed Brown 12173a28c0e5SStefano Zampini /*@ 12183a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 12193a28c0e5SStefano Zampini 12203a28c0e5SStefano Zampini Logically collective 12213a28c0e5SStefano Zampini 12223a28c0e5SStefano Zampini Input Parameter: 12233a28c0e5SStefano Zampini . ts - timestepping context 12243a28c0e5SStefano Zampini 12257a7aea1fSJed Brown Output Parameter: 12263a28c0e5SStefano Zampini . flg - PETSC_TRUE for fully implicit 12273a28c0e5SStefano Zampini 12283a28c0e5SStefano Zampini Level: intermediate 12293a28c0e5SStefano Zampini 1230db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()` 12313a28c0e5SStefano Zampini @*/ 12329371c9d4SSatish Balay PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts, PetscBool *flg) { 12333a28c0e5SStefano Zampini PetscFunctionBegin; 12343a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1235dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 2); 1236cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetFullyImplicit_C", (TS, PetscBool *), (ts, flg)); 12373a28c0e5SStefano Zampini PetscFunctionReturn(0); 12383a28c0e5SStefano Zampini } 12393a28c0e5SStefano Zampini 12409371c9d4SSatish Balay static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts, TSARKIMEXType *arktype) { 12418a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12428a381b04SJed Brown 12438a381b04SJed Brown PetscFunctionBegin; 12448a381b04SJed Brown *arktype = ark->tableau->name; 12458a381b04SJed Brown PetscFunctionReturn(0); 12468a381b04SJed Brown } 12479371c9d4SSatish Balay static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts, TSARKIMEXType arktype) { 12488a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12498a381b04SJed Brown PetscBool match; 12508a381b04SJed Brown ARKTableauLink link; 12518a381b04SJed Brown 12528a381b04SJed Brown PetscFunctionBegin; 12538a381b04SJed Brown if (ark->tableau) { 12549566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ark->tableau->name, arktype, &match)); 12558a381b04SJed Brown if (match) PetscFunctionReturn(0); 12568a381b04SJed Brown } 12578a381b04SJed Brown for (link = ARKTableauList; link; link = link->next) { 12589566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, arktype, &match)); 12598a381b04SJed Brown if (match) { 12609566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauReset(ts)); 12618a381b04SJed Brown ark->tableau = &link->tab; 12629566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauSetUp(ts)); 12632ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 12648a381b04SJed Brown PetscFunctionReturn(0); 12658a381b04SJed Brown } 12668a381b04SJed Brown } 126798921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", arktype); 12688a381b04SJed Brown } 1269e0877f53SBarry Smith 12709371c9d4SSatish Balay static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts, PetscBool flg) { 12714cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12724cc180ffSJed Brown 12734cc180ffSJed Brown PetscFunctionBegin; 12744cc180ffSJed Brown ark->imex = (PetscBool)!flg; 12754cc180ffSJed Brown PetscFunctionReturn(0); 12764cc180ffSJed Brown } 12778a381b04SJed Brown 12789371c9d4SSatish Balay static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts, PetscBool *flg) { 12793a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 12803a28c0e5SStefano Zampini 12813a28c0e5SStefano Zampini PetscFunctionBegin; 12823a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 12833a28c0e5SStefano Zampini PetscFunctionReturn(0); 12843a28c0e5SStefano Zampini } 12853a28c0e5SStefano Zampini 12869371c9d4SSatish Balay static PetscErrorCode TSDestroy_ARKIMEX(TS ts) { 1287b3a6b972SJed Brown PetscFunctionBegin; 12889566063dSJacob Faibussowitsch PetscCall(TSReset_ARKIMEX(ts)); 1289b3a6b972SJed Brown if (ts->dm) { 12909566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 12919566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 1292b3a6b972SJed Brown } 12939566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 12949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", NULL)); 12959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", NULL)); 12969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", NULL)); 12972e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", NULL)); 1298b3a6b972SJed Brown PetscFunctionReturn(0); 1299b3a6b972SJed Brown } 1300b3a6b972SJed Brown 13018a381b04SJed Brown /* ------------------------------------------------------------ */ 13028a381b04SJed Brown /*MC 1303c79dcfbdSBarry Smith TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes 13048a381b04SJed Brown 1305fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1306fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1307fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1308fca742c7SJed Brown 1309fca742c7SJed Brown Notes: 1310a4386c9eSJed Brown The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1311c8058688SBarry Smith 13125eca1a21SEmil Constantinescu If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 13135eca1a21SEmil Constantinescu 1314a4386c9eSJed Brown 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). 1315fca742c7SJed Brown 1316d0685a90SJed Brown Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1317d0685a90SJed Brown 13188a381b04SJed Brown Level: beginner 13198a381b04SJed Brown 1320db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSARKIMEXSetType()`, `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()`, `TSARKIMEXGetFullyImplicit()`, 1321db781477SPatrick Sanan `TSARKIMEX1BEE`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEX3`, `TSARKIMEXL2`, `TSARKIMEXA2`, `TSARKIMEXARS122`, 1322db781477SPatrick Sanan `TSARKIMEX4`, `TSARKIMEX5`, `TSARKIMEXPRSSP2`, `TSARKIMEXARS443`, `TSARKIMEXBPR3`, `TSARKIMEXType`, `TSARKIMEXRegister()` 13238a381b04SJed Brown 13248a381b04SJed Brown M*/ 13259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) { 13268a381b04SJed Brown TS_ARKIMEX *th; 13278a381b04SJed Brown 13288a381b04SJed Brown PetscFunctionBegin; 13299566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 13308a381b04SJed Brown 13318a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 13328a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 13338a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1334f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 13358a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 13368a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1337cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1338108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 133924655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 13408a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 13418a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 13428a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 13438a381b04SJed Brown 1344efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1345efd4aadfSBarry Smith 1346*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 13478a381b04SJed Brown ts->data = (void *)th; 13484cc180ffSJed Brown th->imex = PETSC_TRUE; 13498a381b04SJed Brown 13509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", TSARKIMEXGetType_ARKIMEX)); 13519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", TSARKIMEXSetType_ARKIMEX)); 13529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", TSARKIMEXSetFullyImplicit_ARKIMEX)); 13539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", TSARKIMEXGetFullyImplicit_ARKIMEX)); 135496400bd6SLisandro Dalcin 13559566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts, TSARKIMEXDefault)); 13568a381b04SJed Brown PetscFunctionReturn(0); 13578a381b04SJed Brown } 1358