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; 5980ab5e31SHong Zhang 6080ab5e31SHong Zhang /* context for sensitivity analysis */ 6180ab5e31SHong Zhang Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 6280ab5e31SHong Zhang Vec *VecsSensiTemp; /* Vectors to be multiplied with Jacobian transpose */ 6380ab5e31SHong Zhang Vec *VecsSensiPTemp; /* Temporary Vectors to store JacobianP-transpose-vector product */ 648a381b04SJed Brown } TS_ARKIMEX; 651f80e275SEmil Constantinescu /*MC 661f80e275SEmil Constantinescu TSARKIMEXARS122 - Second order ARK IMEX scheme. 678a381b04SJed Brown 681f80e275SEmil Constantinescu This method has one explicit stage and one implicit stage. 691f80e275SEmil Constantinescu 70bcf0153eSBarry Smith Options Database Key: 7167b8a455SSatish Balay . -ts_arkimex_type ars122 - set arkimex type to ars122 729bd3e852SBarry Smith 73bcf0153eSBarry Smith Level: advanced 74bcf0153eSBarry Smith 751f80e275SEmil Constantinescu References: 76606c0280SSatish Balay . * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 771f80e275SEmil Constantinescu 78*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 791f80e275SEmil Constantinescu M*/ 801f80e275SEmil Constantinescu /*MC 811f80e275SEmil Constantinescu TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 821f80e275SEmil Constantinescu 831f80e275SEmil 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. 841f80e275SEmil Constantinescu 85bcf0153eSBarry Smith Options Database Key: 8667b8a455SSatish Balay . -ts_arkimex_type a2 - set arkimex type to a2 879bd3e852SBarry Smith 881f80e275SEmil Constantinescu Level: advanced 891f80e275SEmil Constantinescu 90*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 911f80e275SEmil Constantinescu M*/ 921f80e275SEmil Constantinescu /*MC 931f80e275SEmil Constantinescu TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 941f80e275SEmil Constantinescu 951f80e275SEmil Constantinescu This method has two implicit stages, and L-stable implicit scheme. 961f80e275SEmil Constantinescu 97bcf0153eSBarry Smith Options Database Key: 9867b8a455SSatish Balay . -ts_arkimex_type l2 - set arkimex type to l2 999bd3e852SBarry Smith 100bcf0153eSBarry Smith Level: advanced 101bcf0153eSBarry Smith 1021f80e275SEmil Constantinescu References: 103606c0280SSatish 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. 1041f80e275SEmil Constantinescu 105*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1061f80e275SEmil Constantinescu M*/ 1071f80e275SEmil Constantinescu /*MC 108c79dcfbdSBarry Smith TSARKIMEX1BEE - First order backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method. 109e817cc15SEmil Constantinescu 110e817cc15SEmil Constantinescu This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used. 111e817cc15SEmil Constantinescu 112bcf0153eSBarry Smith Options Database Key: 11367b8a455SSatish Balay . -ts_arkimex_type 1bee - set arkimex type to 1bee 1149bd3e852SBarry Smith 115e817cc15SEmil Constantinescu Level: advanced 116e817cc15SEmil Constantinescu 117*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 118e817cc15SEmil Constantinescu M*/ 119e817cc15SEmil Constantinescu /*MC 1201f80e275SEmil Constantinescu TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 1211f80e275SEmil Constantinescu 1221f80e275SEmil 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. 1231f80e275SEmil Constantinescu 124bcf0153eSBarry Smith Options Database Key: 12567b8a455SSatish Balay . -ts_arkimex_type 2c - set arkimex type to 2c 1269bd3e852SBarry Smith 1271f80e275SEmil Constantinescu Level: advanced 1281f80e275SEmil Constantinescu 129*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1301f80e275SEmil Constantinescu M*/ 13164f491ddSJed Brown /*MC 13264f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 13364f491ddSJed Brown 134da81f932SPierre Jolivet 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 implicit component. This method was provided by Emil Constantinescu. 13564f491ddSJed Brown 136bcf0153eSBarry Smith Options Database Key: 13767b8a455SSatish Balay . -ts_arkimex_type 2d - set arkimex type to 2d 1389bd3e852SBarry Smith 139b330ce4dSSatish Balay Level: advanced 140b330ce4dSSatish Balay 141*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 14264f491ddSJed Brown M*/ 14364f491ddSJed Brown /*MC 14464f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 14564f491ddSJed Brown 14664f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 14764f491ddSJed Brown 148bcf0153eSBarry Smith Options Database Key: 14967b8a455SSatish Balay . -ts_arkimex_type 2e - set arkimex type to 2e 1509bd3e852SBarry Smith 151b330ce4dSSatish Balay Level: advanced 152b330ce4dSSatish Balay 153*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 15464f491ddSJed Brown M*/ 15564f491ddSJed Brown /*MC 1566cf0794eSJed Brown TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 1576cf0794eSJed Brown 1586cf0794eSJed Brown This method has three implicit stages. 1596cf0794eSJed Brown 1606cf0794eSJed Brown References: 161606c0280SSatish 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. 1626cf0794eSJed Brown 163a8d69d7bSBarry Smith This method is referred to as SSP2-(3,3,2) in https://arxiv.org/abs/1110.4375 1646cf0794eSJed Brown 165bcf0153eSBarry Smith Options Database Key: 16667b8a455SSatish Balay . -ts_arkimex_type prssp2 - set arkimex type to prssp2 1679bd3e852SBarry Smith 1686cf0794eSJed Brown Level: advanced 1696cf0794eSJed Brown 170*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1716cf0794eSJed Brown M*/ 1726cf0794eSJed Brown /*MC 17364f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 17464f491ddSJed Brown 17564f491ddSJed Brown This method has one explicit stage and three implicit stages. 17664f491ddSJed Brown 177bcf0153eSBarry Smith Options Database Key: 17867b8a455SSatish Balay . -ts_arkimex_type 3 - set arkimex type to 3 1799bd3e852SBarry Smith 180bcf0153eSBarry Smith Level: advanced 181bcf0153eSBarry Smith 18264f491ddSJed Brown References: 183606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 18464f491ddSJed Brown 185*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 18664f491ddSJed Brown M*/ 18764f491ddSJed Brown /*MC 1886cf0794eSJed Brown TSARKIMEXARS443 - Third order ARK IMEX scheme. 1896cf0794eSJed Brown 1906cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1916cf0794eSJed Brown 192bcf0153eSBarry Smith Options Database Key: 19367b8a455SSatish Balay . -ts_arkimex_type ars443 - set arkimex type to ars443 1949bd3e852SBarry Smith 195bcf0153eSBarry Smith Level: advanced 196bcf0153eSBarry Smith 1976cf0794eSJed Brown References: 198606c0280SSatish Balay + * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 199606c0280SSatish Balay - * - This method is referred to as ARS(4,4,3) in https://arxiv.org/abs/1110.4375 2006cf0794eSJed Brown 201*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 2026cf0794eSJed Brown M*/ 2036cf0794eSJed Brown /*MC 2046cf0794eSJed Brown TSARKIMEXBPR3 - Third order ARK IMEX scheme. 2056cf0794eSJed Brown 2066cf0794eSJed Brown This method has one explicit stage and four implicit stages. 2076cf0794eSJed Brown 208bcf0153eSBarry Smith Options Database Key: 20967b8a455SSatish Balay . -ts_arkimex_type bpr3 - set arkimex type to bpr3 2109bd3e852SBarry Smith 211bcf0153eSBarry Smith Level: advanced 212bcf0153eSBarry Smith 2136cf0794eSJed Brown References: 214606c0280SSatish Balay . * - This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375 2156cf0794eSJed Brown 216*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 2176cf0794eSJed Brown M*/ 2186cf0794eSJed Brown /*MC 21964f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 22064f491ddSJed Brown 22164f491ddSJed Brown This method has one explicit stage and four implicit stages. 22264f491ddSJed Brown 223bcf0153eSBarry Smith Options Database Key: 22467b8a455SSatish Balay . -ts_arkimex_type 4 - set arkimex type to4 2259bd3e852SBarry Smith 226bcf0153eSBarry Smith Level: advanced 227bcf0153eSBarry Smith 22864f491ddSJed Brown References: 229606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 23064f491ddSJed Brown 231*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 23264f491ddSJed Brown M*/ 23364f491ddSJed Brown /*MC 23464f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 23564f491ddSJed Brown 23664f491ddSJed Brown This method has one explicit stage and five implicit stages. 23764f491ddSJed Brown 238bcf0153eSBarry Smith Options Database Key: 23967b8a455SSatish Balay . -ts_arkimex_type 5 - set arkimex type to 5 2409bd3e852SBarry Smith 241bcf0153eSBarry Smith Level: advanced 242bcf0153eSBarry Smith 24364f491ddSJed Brown References: 244606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 24564f491ddSJed Brown 246*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 24764f491ddSJed Brown M*/ 24864f491ddSJed Brown 2498a381b04SJed Brown /*@C 250bcf0153eSBarry Smith TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in `TSARKIMEX` 2518a381b04SJed Brown 252fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 2538a381b04SJed Brown 2548a381b04SJed Brown Level: advanced 2558a381b04SJed Brown 256*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSARKIMEX`, `TSARKIMEXRegisterDestroy()` 2578a381b04SJed Brown @*/ 258d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXRegisterAll(void) 259d71ae5a4SJacob Faibussowitsch { 2608a381b04SJed Brown PetscFunctionBegin; 2613ba16761SJacob Faibussowitsch if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 2628a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 263e817cc15SEmil Constantinescu 264e817cc15SEmil Constantinescu { 2659371c9d4SSatish Balay const PetscReal A[3][3] = 2669371c9d4SSatish Balay { 267e817cc15SEmil Constantinescu {0.0, 0.0, 0.0}, 2689371c9d4SSatish Balay {0.0, 0.0, 0.0}, 2699371c9d4SSatish Balay {0.0, 0.5, 0.0} 2709371c9d4SSatish Balay }, 2719371c9d4SSatish 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}; 2729566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX1BEE, 2, 3, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 273e817cc15SEmil Constantinescu } 2748a381b04SJed Brown { 2759371c9d4SSatish Balay const PetscReal A[2][2] = 2769371c9d4SSatish Balay { 2779371c9d4SSatish Balay {0.0, 0.0}, 2789371c9d4SSatish Balay {0.5, 0.0} 2799371c9d4SSatish Balay }, 2809371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.0, 0.5}}, b[2] = {0.0, 1.0}, bembedt[2] = {0.5, 0.5}; 2811f80e275SEmil 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*/ 2829566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS122, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2831f80e275SEmil Constantinescu } 2841f80e275SEmil Constantinescu { 2859371c9d4SSatish Balay const PetscReal A[2][2] = 2869371c9d4SSatish Balay { 2879371c9d4SSatish Balay {0.0, 0.0}, 2889371c9d4SSatish Balay {1.0, 0.0} 2899371c9d4SSatish Balay }, 2909371c9d4SSatish Balay At[2][2] = {{0.0, 0.0}, {0.5, 0.5}}, b[2] = {0.5, 0.5}, bembedt[2] = {0.0, 1.0}; 2911f80e275SEmil 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*/ 2929566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXA2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 1, b, NULL)); 2931f80e275SEmil Constantinescu } 2941f80e275SEmil Constantinescu { 295da80777bSKarl 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 */ 2969371c9d4SSatish Balay const PetscReal A[2][2] = 2979371c9d4SSatish Balay { 2989371c9d4SSatish Balay {0.0, 0.0}, 2999371c9d4SSatish Balay {1.0, 0.0} 3009371c9d4SSatish Balay }, 3019371c9d4SSatish 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}}; 3029566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXL2, 2, 2, &At[0][0], b, NULL, &A[0][0], b, NULL, bembedt, bembedt, 2, binterpt[0], binterp[0])); 3031f80e275SEmil Constantinescu } 3041f80e275SEmil Constantinescu { 305da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3069371c9d4SSatish Balay const PetscReal A[3][3] = 3079371c9d4SSatish Balay { 3089371c9d4SSatish Balay {0, 0, 0}, 309da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3109371c9d4SSatish Balay {0.5, 0.5, 0} 3119371c9d4SSatish Balay }, 3129371c9d4SSatish 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}}; 3139566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2C, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3141f80e275SEmil Constantinescu } 3151f80e275SEmil Constantinescu { 316da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3179371c9d4SSatish Balay const PetscReal A[3][3] = 3189371c9d4SSatish Balay { 3199371c9d4SSatish Balay {0, 0, 0}, 320da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3219371c9d4SSatish Balay {0.75, 0.25, 0} 3229371c9d4SSatish Balay }, 3239371c9d4SSatish 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}}; 3249566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2D, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 3258a381b04SJed Brown } 32606db7b1cSJed Brown { /* Optimal for linear implicit part */ 327da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 3289371c9d4SSatish Balay const PetscReal A[3][3] = 3299371c9d4SSatish Balay { 3309371c9d4SSatish Balay {0, 0, 0}, 331da80777bSKarl Rupp {2 - 1.414213562373095048802, 0, 0}, 3329371c9d4SSatish Balay {(3 - 2 * 1.414213562373095048802) / 6, (3 + 2 * 1.414213562373095048802) / 6, 0} 3339371c9d4SSatish Balay }, 3349371c9d4SSatish 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}}; 3359566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2E, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 336a3a57f36SJed Brown } 3376cf0794eSJed Brown { /* Optimal for linear implicit part */ 3389371c9d4SSatish Balay const PetscReal A[3][3] = 3399371c9d4SSatish Balay { 3409371c9d4SSatish Balay {0, 0, 0}, 3416cf0794eSJed Brown {0.5, 0, 0}, 3429371c9d4SSatish Balay {0.5, 0.5, 0} 3439371c9d4SSatish Balay }, 3449371c9d4SSatish Balay At[3][3] = {{0.25, 0, 0}, {0, 0.25, 0}, {1. / 3, 1. / 3, 1. / 3}}; 3459566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXPRSSP2, 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, NULL, NULL, 0, NULL, NULL)); 3466cf0794eSJed Brown } 347a3a57f36SJed Brown { 3489371c9d4SSatish Balay const PetscReal A[4][4] = 3499371c9d4SSatish Balay { 3509371c9d4SSatish Balay {0, 0, 0, 0}, 3514040e9f2SJed Brown {1767732205903. / 2027836641118., 0, 0, 0}, 3524040e9f2SJed Brown {5535828885825. / 10492691773637., 788022342437. / 10882634858940., 0, 0}, 3539371c9d4SSatish Balay {6485989280629. / 16251701735622., -4246266847089. / 9704473918619., 10755448449292. / 10357097424841., 0} 3549371c9d4SSatish Balay }, 3559371c9d4SSatish 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.}}; 3569566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX3, 3, 4, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 2, binterpt[0], NULL)); 357a3a57f36SJed Brown } 358a3a57f36SJed Brown { 3599371c9d4SSatish Balay const PetscReal A[5][5] = 3609371c9d4SSatish Balay { 3619371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3626cf0794eSJed Brown {1. / 2, 0, 0, 0, 0}, 3636cf0794eSJed Brown {11. / 18, 1. / 18, 0, 0, 0}, 3646cf0794eSJed Brown {5. / 6, -5. / 6, .5, 0, 0}, 3659371c9d4SSatish Balay {1. / 4, 7. / 4, 3. / 4, -7. / 4, 0} 3669371c9d4SSatish Balay }, 3679371c9d4SSatish 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; 3689566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS443, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3696cf0794eSJed Brown } 3706cf0794eSJed Brown { 3719371c9d4SSatish Balay const PetscReal A[5][5] = 3729371c9d4SSatish Balay { 3739371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3746cf0794eSJed Brown {1, 0, 0, 0, 0}, 3756cf0794eSJed Brown {4. / 9, 2. / 9, 0, 0, 0}, 3766cf0794eSJed Brown {1. / 4, 0, 3. / 4, 0, 0}, 3779371c9d4SSatish Balay {1. / 4, 0, 3. / 5, 0, 0} 3789371c9d4SSatish Balay }, 3799371c9d4SSatish 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; 3809566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXBPR3, 3, 5, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 0, NULL, NULL)); 3816cf0794eSJed Brown } 3826cf0794eSJed Brown { 3839371c9d4SSatish Balay const PetscReal A[6][6] = 3849371c9d4SSatish Balay { 3859371c9d4SSatish Balay {0, 0, 0, 0, 0, 0}, 386a3a57f36SJed Brown {1. / 2, 0, 0, 0, 0, 0}, 3874040e9f2SJed Brown {13861. / 62500., 6889. / 62500., 0, 0, 0, 0}, 3884040e9f2SJed Brown {-116923316275. / 2393684061468., -2731218467317. / 15368042101831., 9408046702089. / 11113171139209., 0, 0, 0}, 3894040e9f2SJed Brown {-451086348788. / 2902428689909., -2682348792572. / 7519795681897., 12662868775082. / 11960479115383., 3355817975965. / 11060851509271., 0, 0}, 3909371c9d4SSatish Balay {647845179188. / 3216320057751., 73281519250. / 8382639484533., 552539513391. / 3454668386233., 3354512671639. / 8306763924573., 4040. / 17871., 0} 3919371c9d4SSatish Balay }, 3929371c9d4SSatish 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}, 3939371c9d4SSatish Balay {7640104374378. / 9702883013639., -11436875. / 14766696., 2173542590792. / 12501825683035.}, {-20649996744609. / 7521556579894., 174696575. / 18121608., -31592104683404. / 5083833661969.}, 3949371c9d4SSatish Balay {8854892464581. / 2390941311638., -12120380. / 966161., 61146701046299. / 7138195549469.}, {-11397109935349. / 6675773540249., 3843. / 706., -17219254887155. / 4939391667607.}}; 3959566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX4, 4, 6, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 396a3a57f36SJed Brown } 397a3a57f36SJed Brown { 3989371c9d4SSatish Balay const PetscReal A[8][8] = 3999371c9d4SSatish Balay { 4009371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0}, 401a3a57f36SJed Brown {41. / 100, 0, 0, 0, 0, 0, 0, 0}, 4024040e9f2SJed Brown {367902744464. / 2072280473677., 677623207551. / 8224143866563., 0, 0, 0, 0, 0, 0}, 4034040e9f2SJed Brown {1268023523408. / 10340822734521., 0, 1029933939417. / 13636558850479., 0, 0, 0, 0, 0}, 4044040e9f2SJed Brown {14463281900351. / 6315353703477., 0, 66114435211212. / 5879490589093., -54053170152839. / 4284798021562., 0, 0, 0, 0}, 4054040e9f2SJed Brown {14090043504691. / 34967701212078., 0, 15191511035443. / 11219624916014., -18461159152457. / 12425892160975., -281667163811. / 9011619295870., 0, 0, 0}, 4064040e9f2SJed Brown {19230459214898. / 13134317526959., 0, 21275331358303. / 2942455364971., -38145345988419. / 4862620318723., -1. / 8, -1. / 8, 0, 0}, 4079371c9d4SSatish Balay {-19977161125411. / 11928030595625., 0, -40795976796054. / 6384907823539., 177454434618887. / 12078138498510., 782672205425. / 8267701900261., -69563011059811. / 9646580694205., 7356628210526. / 4942186776405., 0} 4089371c9d4SSatish Balay }, 4099371c9d4SSatish 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.}}; 4109566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX5, 5, 8, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembedt, 3, binterpt[0], NULL)); 411a3a57f36SJed Brown } 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4138a381b04SJed Brown } 4148a381b04SJed Brown 4158a381b04SJed Brown /*@C 416bcf0153eSBarry Smith TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by `TSARKIMEXRegister()`. 4178a381b04SJed Brown 4188a381b04SJed Brown Not Collective 4198a381b04SJed Brown 4208a381b04SJed Brown Level: advanced 4218a381b04SJed Brown 422*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXRegister()`, `TSARKIMEXRegisterAll()` 4238a381b04SJed Brown @*/ 424d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXRegisterDestroy(void) 425d71ae5a4SJacob Faibussowitsch { 4268a381b04SJed Brown ARKTableauLink link; 4278a381b04SJed Brown 4288a381b04SJed Brown PetscFunctionBegin; 4298a381b04SJed Brown while ((link = ARKTableauList)) { 4308a381b04SJed Brown ARKTableau t = &link->tab; 4318a381b04SJed Brown ARKTableauList = link->next; 4329566063dSJacob Faibussowitsch PetscCall(PetscFree6(t->At, t->bt, t->ct, t->A, t->b, t->c)); 4339566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembedt, t->bembed)); 4349566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->binterpt, t->binterp)); 4359566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 4369566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 4378a381b04SJed Brown } 4388a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4408a381b04SJed Brown } 4418a381b04SJed Brown 4428a381b04SJed Brown /*@C 443bcf0153eSBarry Smith TSARKIMEXInitializePackage - This function initializes everything in the `TSARKIMEX` package. It is called 444bcf0153eSBarry Smith from `TSInitializePackage()`. 4458a381b04SJed Brown 4468a381b04SJed Brown Level: developer 4478a381b04SJed Brown 448*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`, `TSARKIMEXFinalizePackage()` 4498a381b04SJed Brown @*/ 450d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXInitializePackage(void) 451d71ae5a4SJacob Faibussowitsch { 4528a381b04SJed Brown PetscFunctionBegin; 4533ba16761SJacob Faibussowitsch if (TSARKIMEXPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 4548a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 4559566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterAll()); 4569566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSARKIMEXFinalizePackage)); 4573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4588a381b04SJed Brown } 4598a381b04SJed Brown 4608a381b04SJed Brown /*@C 461bcf0153eSBarry Smith TSARKIMEXFinalizePackage - This function destroys everything in the `TSARKIMEX` package. It is 462bcf0153eSBarry Smith called from `PetscFinalize()`. 4638a381b04SJed Brown 4648a381b04SJed Brown Level: developer 4658a381b04SJed Brown 466*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`, `TSARKIMEXInitializePackage()` 4678a381b04SJed Brown @*/ 468d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXFinalizePackage(void) 469d71ae5a4SJacob Faibussowitsch { 4708a381b04SJed Brown PetscFunctionBegin; 4718a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 4729566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterDestroy()); 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4748a381b04SJed Brown } 4758a381b04SJed Brown 476cd652676SJed Brown /*@C 477bcf0153eSBarry Smith TSARKIMEXRegister - register a `TSARKIMEX` scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 478cd652676SJed Brown 479cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 480cd652676SJed Brown 481cd652676SJed Brown Input Parameters: 482cd652676SJed Brown + name - identifier for method 483cd652676SJed Brown . order - approximation order of method 484cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 485cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 4860298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 4870298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 488cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 4890298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 4900298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 4910298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 4920298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 493cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 494cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 4950298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 496cd652676SJed Brown 497cd652676SJed Brown Level: advanced 498cd652676SJed Brown 499bcf0153eSBarry Smith Note: 500bcf0153eSBarry Smith Several `TSARKIMEX` methods are provided, this function is only needed to create new methods. 501bcf0153eSBarry Smith 502*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSType`, `TS` 503cd652676SJed Brown @*/ 504d71ae5a4SJacob 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[]) 505d71ae5a4SJacob Faibussowitsch { 5068a381b04SJed Brown ARKTableauLink link; 5078a381b04SJed Brown ARKTableau t; 5088a381b04SJed Brown PetscInt i, j; 5098a381b04SJed Brown 5108a381b04SJed Brown PetscFunctionBegin; 5119566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 5129566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 5138a381b04SJed Brown t = &link->tab; 5149566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 5158a381b04SJed Brown t->order = order; 5168a381b04SJed Brown t->s = s; 5179566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(s * s, &t->At, s, &t->bt, s, &t->ct, s * s, &t->A, s, &t->b, s, &t->c)); 5189566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->At, At, s * s)); 5199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A, A, s * s)); 5209566063dSJacob Faibussowitsch if (bt) PetscCall(PetscArraycpy(t->bt, bt, s)); 5219371c9d4SSatish Balay else 5229371c9d4SSatish Balay for (i = 0; i < s; i++) t->bt[i] = At[(s - 1) * s + i]; 5239566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b, b, s)); 5249371c9d4SSatish Balay else 5259371c9d4SSatish Balay for (i = 0; i < s; i++) t->b[i] = t->bt[i]; 5269566063dSJacob Faibussowitsch if (ct) PetscCall(PetscArraycpy(t->ct, ct, s)); 5279371c9d4SSatish Balay else 5289371c9d4SSatish Balay for (i = 0; i < s; i++) 5299371c9d4SSatish Balay for (j = 0, t->ct[i] = 0; j < s; j++) t->ct[i] += At[i * s + j]; 5309566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c, c, s)); 5319371c9d4SSatish Balay else 5329371c9d4SSatish Balay for (i = 0; i < s; i++) 5339371c9d4SSatish Balay for (j = 0, t->c[i] = 0; j < s; j++) t->c[i] += A[i * s + j]; 534e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 5359371c9d4SSatish Balay for (i = 0; i < s; i++) 5369371c9d4SSatish Balay if (t->At[(s - 1) * s + i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 537e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 5389371c9d4SSatish Balay for (i = 0; i < s; i++) 5399371c9d4SSatish Balay if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 540e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5414e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 542108c343cSJed Brown if (bembedt) { 5439566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &t->bembedt, s, &t->bembed)); 5449566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembedt, bembedt, s)); 5459566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed, bembed ? bembed : bembedt, s)); 546108c343cSJed Brown } 547108c343cSJed Brown 5484f385281SJed Brown t->pinterp = pinterp; 5499566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s * pinterp, &t->binterpt, s * pinterp, &t->binterp)); 5509566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp)); 5519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp, binterp ? binterp : binterpt, s * pinterp)); 5528a381b04SJed Brown link->next = ARKTableauList; 5538a381b04SJed Brown ARKTableauList = link; 5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5558a381b04SJed Brown } 5568a381b04SJed Brown 557108c343cSJed Brown /* 558108c343cSJed Brown The step completion formula is 559108c343cSJed Brown 560108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 561108c343cSJed Brown 562108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 563108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 564108c343cSJed Brown We can write 565108c343cSJed Brown 566108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 567108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 568108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 569108c343cSJed Brown 570108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 571108c343cSJed Brown */ 572d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts, PetscInt order, Vec X, PetscBool *done) 573d71ae5a4SJacob Faibussowitsch { 574108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 575108c343cSJed Brown ARKTableau tab = ark->tableau; 576108c343cSJed Brown PetscScalar *w = ark->work; 577108c343cSJed Brown PetscReal h; 578108c343cSJed Brown PetscInt s = tab->s, j; 579108c343cSJed Brown 580108c343cSJed Brown PetscFunctionBegin; 581108c343cSJed Brown switch (ark->status) { 582108c343cSJed Brown case TS_STEP_INCOMPLETE: 583d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 584d71ae5a4SJacob Faibussowitsch h = ts->time_step; 585d71ae5a4SJacob Faibussowitsch break; 586d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 587d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 588d71ae5a4SJacob Faibussowitsch break; 589d71ae5a4SJacob Faibussowitsch default: 590d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 591108c343cSJed Brown } 592108c343cSJed Brown if (order == tab->order) { 593e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 594740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 5959566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[s - 1], X)); 596e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 5979566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 598e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bt[j]; 5999566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 600e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 601108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->b[j]; 6029566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 603e817cc15SEmil Constantinescu } 604e817cc15SEmil Constantinescu } 6059566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol, X)); 606108c343cSJed Brown if (done) *done = PETSC_TRUE; 6073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 608108c343cSJed Brown } else if (order == tab->order - 1) { 609108c343cSJed Brown if (!tab->bembedt) goto unavailable; 610108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 6119566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 612e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * tab->bembedt[j]; 6139566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotI)); 614108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * tab->bembed[j]; 6159566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 616108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 6179566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 618e817cc15SEmil Constantinescu for (j = 0; j < s; j++) w[j] = h * (tab->bembedt[j] - tab->bt[j]); 6199566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, tab->s, w, ark->YdotI)); 620108c343cSJed Brown for (j = 0; j < s; j++) w[j] = h * (tab->bembed[j] - tab->b[j]); 6219566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, w, ark->YdotRHS)); 622108c343cSJed Brown } 623108c343cSJed Brown if (done) *done = PETSC_TRUE; 6243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 625108c343cSJed Brown } 626108c343cSJed Brown unavailable: 627108c343cSJed Brown if (done) *done = PETSC_FALSE; 6289371c9d4SSatish Balay else 6299371c9d4SSatish 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, 6309371c9d4SSatish Balay tab->order, order); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 632108c343cSJed Brown } 633108c343cSJed Brown 634d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts, PetscBool *id) 635d71ae5a4SJacob Faibussowitsch { 636c79dcfbdSBarry Smith Vec Udot, Y1, Y2; 637c79dcfbdSBarry Smith TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 638c79dcfbdSBarry Smith PetscReal norm; 639c79dcfbdSBarry Smith 640c79dcfbdSBarry Smith PetscFunctionBegin; 6419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Udot)); 6429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y1)); 6439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &Y2)); 6449566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y1, ark->imex)); 6459566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Udot, NULL)); 6469566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ts->ptime, ts->vec_sol, Udot, Y2, ark->imex)); 6479566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Y1)); 6489566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2, -1.0, Udot)); 6499566063dSJacob Faibussowitsch PetscCall(VecNorm(Y2, NORM_2, &norm)); 650c79dcfbdSBarry Smith if (norm < 100.0 * PETSC_MACHINE_EPSILON) { 651c79dcfbdSBarry Smith *id = PETSC_TRUE; 652c79dcfbdSBarry Smith } else { 653c79dcfbdSBarry Smith *id = PETSC_FALSE; 6549566063dSJacob 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)); 655c79dcfbdSBarry Smith } 6569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Udot)); 6579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 6589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y2)); 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 660c79dcfbdSBarry Smith } 661c79dcfbdSBarry Smith 662d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 663d71ae5a4SJacob Faibussowitsch { 66424655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 66524655328SShri ARKTableau tab = ark->tableau; 66624655328SShri const PetscInt s = tab->s; 66724655328SShri const PetscReal *bt = tab->bt, *b = tab->b; 66824655328SShri PetscScalar *w = ark->work; 66924655328SShri Vec *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS; 67024655328SShri PetscInt j; 671be5899b3SLisandro Dalcin PetscReal h; 67224655328SShri 67324655328SShri PetscFunctionBegin; 674be5899b3SLisandro Dalcin switch (ark->status) { 675be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 676d71ae5a4SJacob Faibussowitsch case TS_STEP_PENDING: 677d71ae5a4SJacob Faibussowitsch h = ts->time_step; 678d71ae5a4SJacob Faibussowitsch break; 679d71ae5a4SJacob Faibussowitsch case TS_STEP_COMPLETE: 680d71ae5a4SJacob Faibussowitsch h = ts->ptime - ts->ptime_prev; 681d71ae5a4SJacob Faibussowitsch break; 682d71ae5a4SJacob Faibussowitsch default: 683d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 684be5899b3SLisandro Dalcin } 68524655328SShri for (j = 0; j < s; j++) w[j] = -h * bt[j]; 6869566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotI)); 68724655328SShri for (j = 0; j < s; j++) w[j] = -h * b[j]; 6889566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol, s, w, YdotRHS)); 6893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69024655328SShri } 69124655328SShri 692d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_ARKIMEX(TS ts) 693d71ae5a4SJacob Faibussowitsch { 6948a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 6958a381b04SJed Brown ARKTableau tab = ark->tableau; 6968a381b04SJed Brown const PetscInt s = tab->s; 69724655328SShri const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct, *c = tab->c; 698406d0ec2SJed Brown PetscScalar *w = ark->work; 6991297b224SEmil Constantinescu Vec *Y = ark->Y, *YdotI = ark->YdotI, *YdotRHS = ark->YdotRHS, Ydot = ark->Ydot, Ydot0 = ark->Ydot0, Z = ark->Z; 70096400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 701108c343cSJed Brown TSAdapt adapt; 7028a381b04SJed Brown SNES snes; 703fecfb714SLisandro Dalcin PetscInt i, j, its, lits; 704be5899b3SLisandro Dalcin PetscInt rejections = 0; 70596400bd6SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 70696400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7078a381b04SJed Brown 7088a381b04SJed Brown PetscFunctionBegin; 70996400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 7109566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 7119566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 7129566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 71396400bd6SLisandro Dalcin } 71496400bd6SLisandro Dalcin 71596400bd6SLisandro Dalcin if (!ts->steprollback) { 71696400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 7179566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[s - 1], Ydot0)); 71896400bd6SLisandro Dalcin } 719fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 72096400bd6SLisandro Dalcin for (i = 0; i < s; i++) { 7219566063dSJacob Faibussowitsch PetscCall(VecCopy(Y[i], ark->Y_prev[i])); 7229566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotRHS[i], ark->YdotRHS_prev[i])); 7239566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[i], ark->YdotI_prev[i])); 72496400bd6SLisandro Dalcin } 72596400bd6SLisandro Dalcin } 72696400bd6SLisandro Dalcin } 72796400bd6SLisandro Dalcin 728fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 72996400bd6SLisandro Dalcin TS ts_start; 730c79dcfbdSBarry Smith if (PetscDefined(USE_DEBUG)) { 731c79dcfbdSBarry Smith PetscBool id = PETSC_FALSE; 7329566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTestMassIdentity(ts, &id)); 7333c633725SBarry 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"); 734c79dcfbdSBarry Smith } 7359566063dSJacob Faibussowitsch PetscCall(TSClone(ts, &ts_start)); 7369566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts_start, ts->vec_sol)); 7379566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts_start, ts->ptime)); 7389566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts_start, ts->steps + 1)); 7399566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts_start, ts->ptime + ts->time_step)); 7409566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts_start, TS_EXACTFINALTIME_STEPOVER)); 7419566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts_start, ts->time_step)); 7429566063dSJacob Faibussowitsch PetscCall(TSSetType(ts_start, TSARKIMEX)); 7439566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts_start, PETSC_TRUE)); 7449566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts_start, TSARKIMEX1BEE)); 74534561852SEmil Constantinescu 7469566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts_start)); 7479566063dSJacob Faibussowitsch PetscCall(TSSolve(ts_start, ts->vec_sol)); 7489566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts_start, &ts->ptime)); 7499566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts_start, &ts->time_step)); 750bbd56ea5SKarl Rupp 75185fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 75285fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX *)ts_start->data; 7539566063dSJacob Faibussowitsch PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s - 1], Ydot0)); 75485fc7851SLisandro Dalcin } 75596400bd6SLisandro Dalcin ts->steps++; 75634561852SEmil Constantinescu 757d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 758d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 75996400bd6SLisandro Dalcin { 7609566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7619566063dSJacob Faibussowitsch PetscCall(TSSetSNES(ts, snes)); 76296400bd6SLisandro Dalcin } 7639566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts_start)); 764e817cc15SEmil Constantinescu } 765e817cc15SEmil Constantinescu 766108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 76796400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 76896400bd6SLisandro Dalcin PetscReal t = ts->ptime; 769108c343cSJed Brown PetscReal h = ts->time_step; 7708a381b04SJed Brown for (i = 0; i < s; i++) { 7719be3e283SDebojyoti Ghosh ark->stage_time = t + h * ct[i]; 7729566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, ark->stage_time)); 7738a381b04SJed Brown if (At[i * s + i] == 0) { /* This stage is explicit */ 7743c633725SBarry 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"); 7759566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 776e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7779566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotI)); 7788a381b04SJed Brown for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7799566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i], i, w, YdotRHS)); 7808a381b04SJed Brown } else { 781b296d7d5SJed Brown ark->scoeff = 1. / At[i * s + i]; 7828a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 7839566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Z)); 784e817cc15SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * At[i * s + j]; 7859566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotI)); 786c58d1302SEmil Constantinescu for (j = 0; j < i; j++) w[j] = h * A[i * s + j]; 7879566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, i, w, YdotRHS)); 788fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 78956dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 7909566063dSJacob Faibussowitsch PetscCall(TSExtrapolate_ARKIMEX(ts, c[i], Y[i])); 79156dcabbaSDebojyoti Ghosh } else { 7928a381b04SJed Brown /* Initial guess taken from last stage */ 7939566063dSJacob Faibussowitsch PetscCall(VecCopy(i > 0 ? Y[i - 1] : ts->vec_sol, Y[i])); 79456dcabbaSDebojyoti Ghosh } 7959566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 7969566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, Y[i])); 7979566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 7989566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes, &lits)); 7999371c9d4SSatish Balay ts->snes_its += its; 8009371c9d4SSatish Balay ts->ksp_its += lits; 8019566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8029566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt, ts, ark->stage_time, Y[i], &stageok)); 80396400bd6SLisandro Dalcin if (!stageok) { 8041be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8051be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 80696400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8071be93e3eSJed Brown goto reject_step; 8081be93e3eSJed Brown } 8098a381b04SJed Brown } 810e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 811e817cc15SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 8129371c9d4SSatish 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", 8139371c9d4SSatish Balay ark->tableau->name); 8149566063dSJacob Faibussowitsch PetscCall(VecCopy(Ydot0, YdotI[0])); /* YdotI = YdotI(tn-1) */ 815e817cc15SEmil Constantinescu } else { 8169566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 817e817cc15SEmil Constantinescu } 818e817cc15SEmil Constantinescu } else { 8195eca1a21SEmil Constantinescu if (i == 0 && tab->explicit_first_stage) { 8209566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); 8219566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, t + h * ct[i], Y[i], Ydot, YdotI[i], ark->imex)); /* YdotI = -G(t,Y,0) */ 8229566063dSJacob Faibussowitsch PetscCall(VecScale(YdotI[i], -1.0)); 8235eca1a21SEmil Constantinescu } else { 8249566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i], -ark->scoeff / h, ark->scoeff / h, 0, Z, Y[i])); /* YdotI = shift*(X-Z) */ 8255eca1a21SEmil Constantinescu } 8264cc180ffSJed Brown if (ark->imex) { 8279566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * c[i], Y[i], YdotRHS[i])); 8284cc180ffSJed Brown } else { 8299566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(YdotRHS[i])); 8304cc180ffSJed Brown } 8318a381b04SJed Brown } 8329566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, ark->stage_time, i, Y)); 833e817cc15SEmil Constantinescu } 83496400bd6SLisandro Dalcin 835be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 8369566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_ARKIMEX(ts, tab->order, ts->vec_sol, NULL)); 837108c343cSJed Brown ark->status = TS_STEP_PENDING; 8389566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 8399566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8409566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE)); 8419566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 84296400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 84396400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8449566063dSJacob Faibussowitsch PetscCall(TSRollBack_ARKIMEX(ts)); 845be5899b3SLisandro Dalcin ts->time_step = next_time_step; 846be5899b3SLisandro Dalcin goto reject_step; 84796400bd6SLisandro Dalcin } 84896400bd6SLisandro Dalcin 8498a381b04SJed Brown ts->ptime += ts->time_step; 850cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 851108c343cSJed Brown break; 85296400bd6SLisandro Dalcin 85396400bd6SLisandro Dalcin reject_step: 8549371c9d4SSatish Balay ts->reject++; 8559371c9d4SSatish Balay accept = PETSC_FALSE; 856be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 85796400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 85863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 859108c343cSJed Brown } 860f85781f1SEmil Constantinescu } 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8628a381b04SJed Brown } 86380ab5e31SHong Zhang /* 86480ab5e31SHong Zhang This adjoint step function assumes the partitioned ODE system has an identity mass matrix and thus can be represented in the form 86580ab5e31SHong Zhang Udot = H(t,U) + G(t,U) 86680ab5e31SHong Zhang This corresponds to F(t,U,Udot) = Udot-H(t,U). 86780ab5e31SHong Zhang 86880ab5e31SHong Zhang The complete adjoint equations are 86980ab5e31SHong Zhang (shift*I - dHdu) lambda_s[i] = 1/at[i][i] ( 87080ab5e31SHong Zhang (b_i dGdU + bt[i] dHdU) lambda_{n+1} + dGdU \sum_{j=i+1}^s a[j][i] lambda_s[j] 87180ab5e31SHong Zhang + dHdU \sum_{j=i+1}^s at[j][i] lambda_s[j]), i = s-1,...,0 87280ab5e31SHong Zhang lambda_n = lambda_{n+1} + \sum_{j=1}^s lambda_s[j] 87380ab5e31SHong Zhang mu_{n+1}[i] = h (at[i][i] dHdP lambda_s[i] 87480ab5e31SHong Zhang + (b_i dGdP + bt[i] dHdP) lambda_{n+1} + dGdP \sum_{j=i+1}^s a[j][i] lambda_s[j] 87580ab5e31SHong Zhang + dHdP \sum_{j=i+1}^s at[j][i] lambda_s[j]), i = s-1,...,0 87680ab5e31SHong Zhang mu_n = mu_{n+1} + \sum_{j=1}^s mu_{n+1}[j] 87780ab5e31SHong Zhang where shift = 1/(at[i][i]*h) 87880ab5e31SHong Zhang 87980ab5e31SHong Zhang If at[i][i] is 0, the first equation falls back to 88080ab5e31SHong Zhang lambda_s[i] = h ( 88180ab5e31SHong Zhang (b_i dGdU + bt[i] dHdU) lambda_{n+1} + dGdU \sum_{j=i+1}^s a[j][i] lambda_s[j] 88280ab5e31SHong Zhang + dHdU \sum_{j=i+1}^s at[j][i] lambda_s[j]) 88380ab5e31SHong Zhang 88480ab5e31SHong Zhang */ 88580ab5e31SHong Zhang static PetscErrorCode TSAdjointStep_ARKIMEX(TS ts) 88680ab5e31SHong Zhang { 88780ab5e31SHong Zhang TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 88880ab5e31SHong Zhang ARKTableau tab = ark->tableau; 88980ab5e31SHong Zhang const PetscInt s = tab->s; 89080ab5e31SHong Zhang const PetscReal *At = tab->At, *A = tab->A, *ct = tab->ct, *c = tab->c, *b = tab->b, *bt = tab->bt; 89180ab5e31SHong Zhang PetscScalar *w = ark->work; 89280ab5e31SHong Zhang Vec *Y = ark->Y, Ydot = ark->Ydot, *VecsDeltaLam = ark->VecsDeltaLam, *VecsSensiTemp = ark->VecsSensiTemp, *VecsSensiPTemp = ark->VecsSensiPTemp; 89380ab5e31SHong Zhang Mat Jex, Jim, Jimpre; 89480ab5e31SHong Zhang PetscInt i, j, nadj; 89580ab5e31SHong Zhang PetscReal t = ts->ptime, stage_time_ex; 89680ab5e31SHong Zhang PetscReal adjoint_time_step = -ts->time_step; /* always positive since ts->time_step is negative */ 89780ab5e31SHong Zhang KSP ksp; 89880ab5e31SHong Zhang 89980ab5e31SHong Zhang PetscFunctionBegin; 90080ab5e31SHong Zhang ark->status = TS_STEP_INCOMPLETE; 90180ab5e31SHong Zhang PetscCall(SNESGetKSP(ts->snes, &ksp)); 90280ab5e31SHong Zhang PetscCall(TSGetRHSMats_Private(ts, &Jex, NULL)); 90380ab5e31SHong Zhang PetscCall(TSGetIJacobian(ts, &Jim, &Jimpre, NULL, NULL)); 90480ab5e31SHong Zhang 90580ab5e31SHong Zhang for (i = s - 1; i >= 0; i--) { 90680ab5e31SHong Zhang ark->stage_time = t - adjoint_time_step * (1.0 - ct[i]); 90780ab5e31SHong Zhang stage_time_ex = t - adjoint_time_step * (1.0 - c[i]); 90880ab5e31SHong Zhang if (At[i * s + i] == 0) { // This stage is explicit 90980ab5e31SHong Zhang ark->scoeff = 0.; 91080ab5e31SHong Zhang } else { 91180ab5e31SHong Zhang ark->scoeff = -1. / At[i * s + i]; // this makes shift=ark->scoeff/ts->time_step positive since ts->time_step is negative 91280ab5e31SHong Zhang } 91380ab5e31SHong Zhang PetscCall(TSComputeSNESJacobian(ts, Y[i], Jim, Jimpre)); 91480ab5e31SHong Zhang PetscCall(TSComputeRHSJacobian(ts, stage_time_ex, Y[i], Jex, Jex)); 91580ab5e31SHong Zhang if (ts->vecs_sensip) { 91680ab5e31SHong Zhang PetscCall(TSComputeIJacobianP(ts, ark->stage_time, Y[i], Ydot, ark->scoeff / adjoint_time_step, ts->Jacp, PETSC_TRUE)); // get dFdP (-dHdP), Ydot not really used since mass matrix is identity 91780ab5e31SHong Zhang PetscCall(TSComputeRHSJacobianP(ts, stage_time_ex, Y[i], ts->Jacprhs)); // get dGdP 91880ab5e31SHong Zhang } 91980ab5e31SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) { 92080ab5e31SHong Zhang /* Build RHS (stored in VecsDeltaLam) for first-order adjoint */ 92180ab5e31SHong Zhang if (s - i - 1 > 0) { 92280ab5e31SHong Zhang /* Temp = -\sum_{j=i+1}^s at[j][i] lambda_s[j] */ 92380ab5e31SHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = -At[j * s + i]; 92480ab5e31SHong Zhang PetscCall(VecSet(VecsSensiTemp[nadj], 0)); 92580ab5e31SHong Zhang PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 92680ab5e31SHong Zhang /* (shift I - dHdU) Temp */ 92780ab5e31SHong Zhang PetscCall(MatMultTranspose(Jim, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i])); 92880ab5e31SHong Zhang /* cancel out shift Temp where shift=-scoeff/h */ 92980ab5e31SHong Zhang PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], ark->scoeff / adjoint_time_step, VecsSensiTemp[nadj])); 93080ab5e31SHong Zhang if (ts->vecs_sensip) { 93180ab5e31SHong Zhang /* - dHdP Temp */ 93280ab5e31SHong Zhang PetscCall(MatMultTranspose(ts->Jacp, VecsSensiTemp[nadj], VecsSensiPTemp[nadj])); 93380ab5e31SHong Zhang /* mu_n += h dHdP Temp */ 93480ab5e31SHong Zhang PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsSensiPTemp[nadj])); 93580ab5e31SHong Zhang } 93680ab5e31SHong Zhang /* Temp = \sum_{j=i+1}^s a[j][i] lambda_s[j] */ 93780ab5e31SHong Zhang for (j = i + 1; j < s; j++) w[j - i - 1] = A[j * s + i]; 93880ab5e31SHong Zhang PetscCall(VecSet(VecsSensiTemp[nadj], 0)); 93980ab5e31SHong Zhang PetscCall(VecMAXPY(VecsSensiTemp[nadj], s - i - 1, w, &VecsDeltaLam[nadj * s + i + 1])); 94080ab5e31SHong Zhang /* dGdU Temp */ 94180ab5e31SHong Zhang PetscCall(MatMultTransposeAdd(Jex, VecsSensiTemp[nadj], VecsDeltaLam[nadj * s + i], VecsDeltaLam[nadj * s + i])); 94280ab5e31SHong Zhang if (ts->vecs_sensip) { 94380ab5e31SHong Zhang /* dGdP Temp */ 94480ab5e31SHong Zhang PetscCall(MatMultTranspose(ts->Jacprhs, VecsSensiTemp[nadj], VecsSensiPTemp[nadj])); 94580ab5e31SHong Zhang /* mu_n += h dGdP Temp */ 94680ab5e31SHong Zhang PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, VecsSensiPTemp[nadj])); 94780ab5e31SHong Zhang } 94880ab5e31SHong Zhang } else { 94980ab5e31SHong Zhang PetscCall(VecSet(VecsDeltaLam[nadj * s + i], 0)); // make sure it is initialized 95080ab5e31SHong Zhang } 95180ab5e31SHong Zhang if (bt[i]) { // bt[i] dHdU lambda_{n+1} 95280ab5e31SHong Zhang /* (shift I - dHdU)^T lambda_{n+1} */ 95380ab5e31SHong Zhang PetscCall(MatMultTranspose(Jim, ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 95480ab5e31SHong Zhang /* -bt[i] (shift I - dHdU)^T lambda_{n+1} */ 95580ab5e31SHong Zhang PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -bt[i], VecsSensiTemp[nadj])); 95680ab5e31SHong Zhang /* cancel out -bt[i] shift lambda_{n+1} */ 95780ab5e31SHong Zhang PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], -bt[i] * ark->scoeff / adjoint_time_step, ts->vecs_sensi[nadj])); 95880ab5e31SHong Zhang if (ts->vecs_sensip) { 95980ab5e31SHong Zhang /* -dHdP lambda_{n+1} */ 96080ab5e31SHong Zhang PetscCall(MatMultTranspose(ts->Jacp, ts->vecs_sensi[nadj], VecsSensiPTemp[nadj])); 96180ab5e31SHong Zhang /* mu_n += h bt[i] dHdP lambda_{n+1} */ 96280ab5e31SHong Zhang PetscCall(VecAXPY(ts->vecs_sensip[nadj], -bt[i] * adjoint_time_step, VecsSensiPTemp[nadj])); 96380ab5e31SHong Zhang } 96480ab5e31SHong Zhang } 96580ab5e31SHong Zhang if (b[i]) { // b[i] dGdU lambda_{n+1} 96680ab5e31SHong Zhang /* dGdU lambda_{n+1} */ 96780ab5e31SHong Zhang PetscCall(MatMultTranspose(Jex, ts->vecs_sensi[nadj], VecsSensiTemp[nadj])); 96880ab5e31SHong Zhang /* b[i] dGdU lambda_{n+1} */ 96980ab5e31SHong Zhang PetscCall(VecAXPY(VecsDeltaLam[nadj * s + i], b[i], VecsSensiTemp[nadj])); 97080ab5e31SHong Zhang if (ts->vecs_sensip) { 97180ab5e31SHong Zhang /* dGdP lambda_{n+1} */ 97280ab5e31SHong Zhang PetscCall(MatMultTranspose(ts->Jacprhs, ts->vecs_sensi[nadj], VecsSensiPTemp[nadj])); 97380ab5e31SHong Zhang /* mu_n += h b[i] dGdP lambda_{n+1} */ 97480ab5e31SHong Zhang PetscCall(VecAXPY(ts->vecs_sensip[nadj], b[i] * adjoint_time_step, VecsSensiPTemp[nadj])); 97580ab5e31SHong Zhang } 97680ab5e31SHong Zhang } 97780ab5e31SHong Zhang /* Build LHS for first-order adjoint */ 97880ab5e31SHong Zhang if (At[i * s + i] == 0) { // This stage is explicit 97980ab5e31SHong Zhang PetscCall(VecScale(VecsDeltaLam[nadj * s + i], adjoint_time_step)); 98080ab5e31SHong Zhang } else { 98180ab5e31SHong Zhang KSP ksp; 98280ab5e31SHong Zhang KSPConvergedReason kspreason; 98380ab5e31SHong Zhang PetscCall(SNESGetKSP(ts->snes, &ksp)); 98480ab5e31SHong Zhang PetscCall(KSPSetOperators(ksp, Jim, Jimpre)); 98580ab5e31SHong Zhang PetscCall(VecScale(VecsDeltaLam[nadj * s + i], 1. / At[i * s + i])); 98680ab5e31SHong Zhang PetscCall(KSPSolveTranspose(ksp, VecsDeltaLam[nadj * s + i], VecsDeltaLam[nadj * s + i])); 98780ab5e31SHong Zhang PetscCall(KSPGetConvergedReason(ksp, &kspreason)); 98880ab5e31SHong Zhang if (kspreason < 0) { 98980ab5e31SHong Zhang ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE; 99080ab5e31SHong Zhang PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj)); 99180ab5e31SHong Zhang } 99280ab5e31SHong Zhang if (ts->vecs_sensip) { 99380ab5e31SHong Zhang /* -dHdP lambda_s[i] */ 99480ab5e31SHong Zhang PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj * s + i], VecsSensiPTemp[nadj])); 99580ab5e31SHong Zhang /* mu_n += h at[i][i] dHdP lambda_s[i] */ 99680ab5e31SHong Zhang PetscCall(VecAXPY(ts->vecs_sensip[nadj], -At[i * s + i] * adjoint_time_step, VecsSensiPTemp[nadj])); 99780ab5e31SHong Zhang } 99880ab5e31SHong Zhang } 99980ab5e31SHong Zhang } 100080ab5e31SHong Zhang } 100180ab5e31SHong Zhang for (j = 0; j < s; j++) w[j] = 1.0; 100280ab5e31SHong Zhang for (nadj = 0; nadj < ts->numcost; nadj++) // no need to do this for mu's 100380ab5e31SHong Zhang PetscCall(VecMAXPY(ts->vecs_sensi[nadj], s, w, &VecsDeltaLam[nadj * s])); 100480ab5e31SHong Zhang ark->status = TS_STEP_COMPLETE; 100580ab5e31SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 100680ab5e31SHong Zhang } 10078a381b04SJed Brown 1008d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_ARKIMEX(TS ts, PetscReal itime, Vec X) 1009d71ae5a4SJacob Faibussowitsch { 1010cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 10114f385281SJed Brown PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 1012108c343cSJed Brown PetscReal h; 1013108c343cSJed Brown PetscReal tt, t; 1014cd652676SJed Brown PetscScalar *bt, *b; 1015cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 1016cd652676SJed Brown 1017cd652676SJed Brown PetscFunctionBegin; 10183c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 1019108c343cSJed Brown switch (ark->status) { 1020108c343cSJed Brown case TS_STEP_INCOMPLETE: 1021108c343cSJed Brown case TS_STEP_PENDING: 1022108c343cSJed Brown h = ts->time_step; 1023108c343cSJed Brown t = (itime - ts->ptime) / h; 1024108c343cSJed Brown break; 1025108c343cSJed Brown case TS_STEP_COMPLETE: 1026be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 1027108c343cSJed Brown t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */ 1028108c343cSJed Brown break; 1029d71ae5a4SJacob Faibussowitsch default: 1030d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus"); 1031108c343cSJed Brown } 10329566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s, &bt, s, &b)); 1033cd652676SJed Brown for (i = 0; i < s; i++) bt[i] = b[i] = 0; 10344f385281SJed Brown for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 1035cd652676SJed Brown for (i = 0; i < s; i++) { 1036c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 1037108c343cSJed Brown b[i] += h * B[i * pinterp + j] * tt; 1038cd652676SJed Brown } 1039cd652676SJed Brown } 10409566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[0], X)); 10419566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI)); 10429566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS)); 10439566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 10443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1045cd652676SJed Brown } 1046cd652676SJed Brown 1047d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts, PetscReal c, Vec X) 1048d71ae5a4SJacob Faibussowitsch { 104956dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 105056dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s, pinterp = ark->tableau->pinterp, i, j; 1051be5899b3SLisandro Dalcin PetscReal h, h_prev, t, tt; 105256dcabbaSDebojyoti Ghosh PetscScalar *bt, *b; 105356dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt, *B = ark->tableau->binterp; 105456dcabbaSDebojyoti Ghosh 105556dcabbaSDebojyoti Ghosh PetscFunctionBegin; 10563c633725SBarry Smith PetscCheck(Bt && B, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSARKIMEX %s does not have an interpolation formula", ark->tableau->name); 10579566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(s, &bt, s, &b)); 105881d12688SDebojyoti Ghosh h = ts->time_step; 1059be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 1060be5899b3SLisandro Dalcin t = 1 + h / h_prev * c; 106156dcabbaSDebojyoti Ghosh for (j = 0, tt = t; j < pinterp; j++, tt *= t) { 106256dcabbaSDebojyoti Ghosh for (i = 0; i < s; i++) { 106381d12688SDebojyoti Ghosh bt[i] += h * Bt[i * pinterp + j] * tt; 106456dcabbaSDebojyoti Ghosh b[i] += h * B[i * pinterp + j] * tt; 106556dcabbaSDebojyoti Ghosh } 106656dcabbaSDebojyoti Ghosh } 10673c633725SBarry Smith PetscCheck(ark->Y_prev, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Stages from previous step have not been stored"); 10689566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y_prev[0], X)); 10699566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, bt, ark->YdotI_prev)); 10709566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, b, ark->YdotRHS_prev)); 10719566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt, b)); 10723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107356dcabbaSDebojyoti Ghosh } 107456dcabbaSDebojyoti Ghosh 10758a381b04SJed Brown /*------------------------------------------------------------*/ 107696400bd6SLisandro Dalcin 1077d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTableauReset(TS ts) 1078d71ae5a4SJacob Faibussowitsch { 107996400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 108096400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 108196400bd6SLisandro Dalcin 108296400bd6SLisandro Dalcin PetscFunctionBegin; 10833ba16761SJacob Faibussowitsch if (!tab) PetscFunctionReturn(PETSC_SUCCESS); 10849566063dSJacob Faibussowitsch PetscCall(PetscFree(ark->work)); 10859566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y)); 10869566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI)); 10879566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS)); 10889566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->Y_prev)); 10899566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotI_prev)); 10909566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &ark->YdotRHS_prev)); 10913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109296400bd6SLisandro Dalcin } 109396400bd6SLisandro Dalcin 1094d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_ARKIMEX(TS ts) 1095d71ae5a4SJacob Faibussowitsch { 10968a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 10978a381b04SJed Brown 10988a381b04SJed Brown PetscFunctionBegin; 10999566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauReset(ts)); 11009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot)); 11019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot0)); 11029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Z)); 11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11048a381b04SJed Brown } 11058a381b04SJed Brown 110680ab5e31SHong Zhang static PetscErrorCode TSAdjointReset_ARKIMEX(TS ts) 110780ab5e31SHong Zhang { 110880ab5e31SHong Zhang TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 110980ab5e31SHong Zhang ARKTableau tab = ark->tableau; 111080ab5e31SHong Zhang 111180ab5e31SHong Zhang PetscFunctionBegin; 111280ab5e31SHong Zhang PetscCall(VecDestroyVecs(tab->s * ts->numcost, &ark->VecsDeltaLam)); 111380ab5e31SHong Zhang PetscCall(VecDestroyVecs(ts->numcost, &ark->VecsSensiTemp)); 111480ab5e31SHong Zhang PetscCall(VecDestroyVecs(ts->numcost, &ark->VecsSensiPTemp)); 111580ab5e31SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 111680ab5e31SHong Zhang } 111780ab5e31SHong Zhang 1118d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) 1119d71ae5a4SJacob Faibussowitsch { 1120d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX *)ts->data; 1121d5e6173cSPeter Brune 1122d5e6173cSPeter Brune PetscFunctionBegin; 1123d5e6173cSPeter Brune if (Z) { 1124d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11259566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 1126d5e6173cSPeter Brune } else *Z = ax->Z; 1127d5e6173cSPeter Brune } 1128d5e6173cSPeter Brune if (Ydot) { 1129d5e6173cSPeter Brune if (dm && dm != ts->dm) { 11309566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 1131d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 1132d5e6173cSPeter Brune } 11333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1134d5e6173cSPeter Brune } 1135d5e6173cSPeter Brune 1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXRestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydot) 1137d71ae5a4SJacob Faibussowitsch { 1138d5e6173cSPeter Brune PetscFunctionBegin; 1139d5e6173cSPeter Brune if (Z) { 114048a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Z", Z)); 1141d5e6173cSPeter Brune } 1142d5e6173cSPeter Brune if (Ydot) { 114348a46eb9SPierre Jolivet if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSARKIMEX_Ydot", Ydot)); 1144d5e6173cSPeter Brune } 11453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1146d5e6173cSPeter Brune } 1147d5e6173cSPeter Brune 11488a381b04SJed Brown /* 11498a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 11508a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 11518a381b04SJed Brown */ 1152d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes, Vec X, Vec F, TS ts) 1153d71ae5a4SJacob Faibussowitsch { 11548a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1155d5e6173cSPeter Brune DM dm, dmsave; 1156d5e6173cSPeter Brune Vec Z, Ydot; 1157b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 11588a381b04SJed Brown 11598a381b04SJed Brown PetscFunctionBegin; 11609566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 11619566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, &Ydot)); 11629566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Ydot, -shift, shift, 0, Z, X)); /* Ydot = shift*(X-Z) */ 1163d5e6173cSPeter Brune dmsave = ts->dm; 1164d5e6173cSPeter Brune ts->dm = dm; 1165740132f1SEmil Constantinescu 11669566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ark->stage_time, X, Ydot, F, ark->imex)); 1167e817cc15SEmil Constantinescu 1168d5e6173cSPeter Brune ts->dm = dmsave; 11699566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, &Ydot)); 11703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11718a381b04SJed Brown } 11728a381b04SJed Brown 1173d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes, Vec X, Mat A, Mat B, TS ts) 1174d71ae5a4SJacob Faibussowitsch { 11758a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1176d5e6173cSPeter Brune DM dm, dmsave; 1177d5e6173cSPeter Brune Vec Ydot; 1178b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 11798a381b04SJed Brown 11808a381b04SJed Brown PetscFunctionBegin; 11819566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 11829566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, NULL, &Ydot)); 11838a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1184d5e6173cSPeter Brune dmsave = ts->dm; 1185d5e6173cSPeter Brune ts->dm = dm; 1186740132f1SEmil Constantinescu 11879566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ark->stage_time, X, Ydot, shift, A, B, ark->imex)); 1188740132f1SEmil Constantinescu 1189d5e6173cSPeter Brune ts->dm = dmsave; 11909566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, NULL, &Ydot)); 11913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1192d5e6173cSPeter Brune } 1193d5e6173cSPeter Brune 119480ab5e31SHong Zhang static PetscErrorCode TSGetStages_ARKIMEX(TS ts, PetscInt *ns, Vec *Y[]) 119580ab5e31SHong Zhang { 119680ab5e31SHong Zhang TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 119780ab5e31SHong Zhang 119880ab5e31SHong Zhang PetscFunctionBegin; 119980ab5e31SHong Zhang if (ns) *ns = ark->tableau->s; 120080ab5e31SHong Zhang if (Y) *Y = ark->Y; 120180ab5e31SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 120280ab5e31SHong Zhang } 120380ab5e31SHong Zhang 1204d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine, DM coarse, void *ctx) 1205d71ae5a4SJacob Faibussowitsch { 1206d5e6173cSPeter Brune PetscFunctionBegin; 12073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1208d5e6173cSPeter Brune } 1209d5e6173cSPeter Brune 1210d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) 1211d71ae5a4SJacob Faibussowitsch { 1212d5e6173cSPeter Brune TS ts = (TS)ctx; 1213d5e6173cSPeter Brune Vec Z, Z_c; 1214d5e6173cSPeter Brune 1215d5e6173cSPeter Brune PetscFunctionBegin; 12169566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, fine, &Z, NULL)); 12179566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, coarse, &Z_c, NULL)); 12189566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, Z, Z_c)); 12199566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Z_c, rscale, Z_c)); 12209566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, fine, &Z, NULL)); 12219566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, coarse, &Z_c, NULL)); 12223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12238a381b04SJed Brown } 12248a381b04SJed Brown 1225d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm, DM subdm, void *ctx) 1226d71ae5a4SJacob Faibussowitsch { 1227cdb298fcSPeter Brune PetscFunctionBegin; 12283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1229cdb298fcSPeter Brune } 1230cdb298fcSPeter Brune 1231d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) 1232d71ae5a4SJacob Faibussowitsch { 1233cdb298fcSPeter Brune TS ts = (TS)ctx; 1234cdb298fcSPeter Brune Vec Z, Z_c; 1235cdb298fcSPeter Brune 1236cdb298fcSPeter Brune PetscFunctionBegin; 12379566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, dm, &Z, NULL)); 12389566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts, subdm, &Z_c, NULL)); 1239cdb298fcSPeter Brune 12409566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 12419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, Z, Z_c, INSERT_VALUES, SCATTER_FORWARD)); 1242cdb298fcSPeter Brune 12439566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, dm, &Z, NULL)); 12449566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts, subdm, &Z_c, NULL)); 12453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1246cdb298fcSPeter Brune } 1247cdb298fcSPeter Brune 1248d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 1249d71ae5a4SJacob Faibussowitsch { 125096400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 125196400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 125296400bd6SLisandro Dalcin 125396400bd6SLisandro Dalcin PetscFunctionBegin; 12549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &ark->work)); 12559566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y)); 12569566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI)); 12579566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS)); 125896400bd6SLisandro Dalcin if (ark->extrapolate) { 12599566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->Y_prev)); 12609566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotI_prev)); 12619566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ark->YdotRHS_prev)); 126296400bd6SLisandro Dalcin } 12633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126496400bd6SLisandro Dalcin } 126596400bd6SLisandro Dalcin 1266d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 1267d71ae5a4SJacob Faibussowitsch { 12688a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 1269d5e6173cSPeter Brune DM dm; 127096400bd6SLisandro Dalcin SNES snes; 1271f9c1d6abSBarry Smith 12728a381b04SJed Brown PetscFunctionBegin; 12739566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauSetUp(ts)); 12749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot)); 12759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Ydot0)); 12769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ark->Z)); 12779566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 12789566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 12799566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 12809566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 12813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12828a381b04SJed Brown } 12838a381b04SJed Brown /*------------------------------------------------------------*/ 12848a381b04SJed Brown 128580ab5e31SHong Zhang static PetscErrorCode TSAdjointSetUp_ARKIMEX(TS ts) 128680ab5e31SHong Zhang { 128780ab5e31SHong Zhang TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 128880ab5e31SHong Zhang ARKTableau tab = ark->tableau; 128980ab5e31SHong Zhang 129080ab5e31SHong Zhang PetscFunctionBegin; 129180ab5e31SHong Zhang PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], tab->s * ts->numcost, &ark->VecsDeltaLam)); 129280ab5e31SHong Zhang PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &ark->VecsSensiTemp)); 129380ab5e31SHong Zhang if (ts->vecs_sensip) { PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &ark->VecsSensiPTemp)); } 129480ab5e31SHong Zhang if (PetscDefined(USE_DEBUG)) { 129580ab5e31SHong Zhang PetscBool id = PETSC_FALSE; 129680ab5e31SHong Zhang PetscCall(TSARKIMEXTestMassIdentity(ts, &id)); 129780ab5e31SHong Zhang PetscCheck(id, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Adjoint ARKIMEX requires an identity mass matrix, however the TSIFunction you provide does not utilize an identity mass matrix"); 129880ab5e31SHong Zhang } 129980ab5e31SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 130080ab5e31SHong Zhang } 130180ab5e31SHong Zhang 1302d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts, PetscOptionItems *PetscOptionsObject) 1303d71ae5a4SJacob Faibussowitsch { 13044cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 13058a381b04SJed Brown 13068a381b04SJed Brown PetscFunctionBegin; 1307d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ARKIMEX ODE solver options"); 13088a381b04SJed Brown { 13098a381b04SJed Brown ARKTableauLink link; 13108a381b04SJed Brown PetscInt count, choice; 13118a381b04SJed Brown PetscBool flg; 13128a381b04SJed Brown const char **namelist; 13139371c9d4SSatish Balay for (link = ARKTableauList, count = 0; link; link = link->next, count++) 13149371c9d4SSatish Balay ; 13159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 13168a381b04SJed Brown for (link = ARKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 13179566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_arkimex_type", "Family of ARK IMEX method", "TSARKIMEXSetType", (const char *const *)namelist, count, ark->tableau->name, &choice, &flg)); 13189566063dSJacob Faibussowitsch if (flg) PetscCall(TSARKIMEXSetType(ts, namelist[choice])); 13199566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 132096400bd6SLisandro Dalcin 13214cc180ffSJed Brown flg = (PetscBool)!ark->imex; 13229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_arkimex_fully_implicit", "Solve the problem fully implicitly", "TSARKIMEXSetFullyImplicit", flg, &flg, NULL)); 13234cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13249566063dSJacob 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)); 13258a381b04SJed Brown } 1326d0609cedSBarry Smith PetscOptionsHeadEnd(); 13273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13288a381b04SJed Brown } 13298a381b04SJed Brown 1330d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_ARKIMEX(TS ts, PetscViewer viewer) 1331d71ae5a4SJacob Faibussowitsch { 13328a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 13338a381b04SJed Brown PetscBool iascii; 13348a381b04SJed Brown 13358a381b04SJed Brown PetscFunctionBegin; 13369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 13378a381b04SJed Brown if (iascii) { 13389c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 133919fd82e9SBarry Smith TSARKIMEXType arktype; 13408a381b04SJed Brown char buf[512]; 13413a28c0e5SStefano Zampini PetscBool flg; 13423a28c0e5SStefano Zampini 13439566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetType(ts, &arktype)); 13449566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetFullyImplicit(ts, &flg)); 13459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " ARK IMEX %s\n", arktype)); 13469566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ct)); 13479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Stiff abscissa ct = %s\n", buf)); 13489566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->c)); 13499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Fully implicit: %s\n", flg ? "yes" : "no")); 13509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s\n", tab->stiffly_accurate ? "yes" : "no")); 13519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Explicit first stage: %s\n", tab->explicit_first_stage ? "yes" : "no")); 13529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "FSAL property: %s\n", tab->FSAL_implicit ? "yes" : "no")); 13539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Nonstiff abscissa c = %s\n", buf)); 13548a381b04SJed Brown } 13553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13568a381b04SJed Brown } 13578a381b04SJed Brown 1358d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSLoad_ARKIMEX(TS ts, PetscViewer viewer) 1359d71ae5a4SJacob Faibussowitsch { 1360f2c2a1b9SBarry Smith SNES snes; 13619c334d8fSLisandro Dalcin TSAdapt adapt; 1362f2c2a1b9SBarry Smith 1363f2c2a1b9SBarry Smith PetscFunctionBegin; 13649566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 13659566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 13669566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 13679566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes, viewer)); 1368ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 13699566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, ts)); 13709566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts)); 13713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1372f2c2a1b9SBarry Smith } 1373f2c2a1b9SBarry Smith 13748a381b04SJed Brown /*@C 1375bcf0153eSBarry Smith TSARKIMEXSetType - Set the type of `TSARKIMEX` scheme 13768a381b04SJed Brown 137720f4b53cSBarry Smith Logically Collective 13788a381b04SJed Brown 1379d8d19677SJose E. Roman Input Parameters: 13808a381b04SJed Brown + ts - timestepping context 1381bcf0153eSBarry Smith - arktype - type of `TSARKIMEX` scheme 13828a381b04SJed Brown 1383bcf0153eSBarry Smith Options Database Key: 1384bcf0153eSBarry Smith . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> - set `TSARKIMEX` scheme type 13859bd3e852SBarry Smith 13868a381b04SJed Brown Level: intermediate 13878a381b04SJed Brown 1388*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXGetType()`, `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEX1BEE`, `TSARKIMEXA2`, `TSARKIMEXL2`, `TSARKIMEXARS122`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEXPRSSP2`, 1389db781477SPatrick Sanan `TSARKIMEX3`, `TSARKIMEXBPR3`, `TSARKIMEXARS443`, `TSARKIMEX4`, `TSARKIMEX5` 13908a381b04SJed Brown @*/ 1391d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType arktype) 1392d71ae5a4SJacob Faibussowitsch { 13938a381b04SJed Brown PetscFunctionBegin; 13948a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1395b92453a8SLisandro Dalcin PetscValidCharPointer(arktype, 2); 1396cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetType_C", (TS, TSARKIMEXType), (ts, arktype)); 13973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13988a381b04SJed Brown } 13998a381b04SJed Brown 14008a381b04SJed Brown /*@C 1401bcf0153eSBarry Smith TSARKIMEXGetType - Get the type of `TSARKIMEX` scheme 14028a381b04SJed Brown 140320f4b53cSBarry Smith Logically Collective 14048a381b04SJed Brown 14058a381b04SJed Brown Input Parameter: 14068a381b04SJed Brown . ts - timestepping context 14078a381b04SJed Brown 14088a381b04SJed Brown Output Parameter: 1409bcf0153eSBarry Smith . arktype - type of `TSARKIMEX` scheme 14108a381b04SJed Brown 14118a381b04SJed Brown Level: intermediate 14128a381b04SJed Brown 1413*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`c, `TSARKIMEXGetType()` 14148a381b04SJed Brown @*/ 1415d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *arktype) 1416d71ae5a4SJacob Faibussowitsch { 14178a381b04SJed Brown PetscFunctionBegin; 14188a381b04SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1419cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetType_C", (TS, TSARKIMEXType *), (ts, arktype)); 14203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14218a381b04SJed Brown } 14228a381b04SJed Brown 142316353aafSBarry Smith /*@ 1424bcf0153eSBarry Smith TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly, including the part that is normally solved explicitly 14254cc180ffSJed Brown 142620f4b53cSBarry Smith Logically Collective 14274cc180ffSJed Brown 1428d8d19677SJose E. Roman Input Parameters: 14294cc180ffSJed Brown + ts - timestepping context 1430bcf0153eSBarry Smith - flg - `PETSC_TRUE` for fully implicit 14314cc180ffSJed Brown 14324cc180ffSJed Brown Level: intermediate 14334cc180ffSJed Brown 1434*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEX`, `TSARKIMEXGetType()`, `TSARKIMEXGetFullyImplicit()` 14354cc180ffSJed Brown @*/ 1436d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts, PetscBool flg) 1437d71ae5a4SJacob Faibussowitsch { 14384cc180ffSJed Brown PetscFunctionBegin; 14394cc180ffSJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 14403a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts, flg, 2); 1441cac4c232SBarry Smith PetscTryMethod(ts, "TSARKIMEXSetFullyImplicit_C", (TS, PetscBool), (ts, flg)); 14423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14434cc180ffSJed Brown } 14444cc180ffSJed Brown 14453a28c0e5SStefano Zampini /*@ 14463a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 14473a28c0e5SStefano Zampini 144820f4b53cSBarry Smith Logically Collective 14493a28c0e5SStefano Zampini 14503a28c0e5SStefano Zampini Input Parameter: 14513a28c0e5SStefano Zampini . ts - timestepping context 14523a28c0e5SStefano Zampini 14537a7aea1fSJed Brown Output Parameter: 1454bcf0153eSBarry Smith . flg - `PETSC_TRUE` for fully implicit 14553a28c0e5SStefano Zampini 14563a28c0e5SStefano Zampini Level: intermediate 14573a28c0e5SStefano Zampini 1458*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()` 14593a28c0e5SStefano Zampini @*/ 1460d71ae5a4SJacob Faibussowitsch PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts, PetscBool *flg) 1461d71ae5a4SJacob Faibussowitsch { 14623a28c0e5SStefano Zampini PetscFunctionBegin; 14633a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1464dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg, 2); 1465cac4c232SBarry Smith PetscUseMethod(ts, "TSARKIMEXGetFullyImplicit_C", (TS, PetscBool *), (ts, flg)); 14663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14673a28c0e5SStefano Zampini } 14683a28c0e5SStefano Zampini 1469d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts, TSARKIMEXType *arktype) 1470d71ae5a4SJacob Faibussowitsch { 14718a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 14728a381b04SJed Brown 14738a381b04SJed Brown PetscFunctionBegin; 14748a381b04SJed Brown *arktype = ark->tableau->name; 14753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14768a381b04SJed Brown } 1477d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts, TSARKIMEXType arktype) 1478d71ae5a4SJacob Faibussowitsch { 14798a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 14808a381b04SJed Brown PetscBool match; 14818a381b04SJed Brown ARKTableauLink link; 14828a381b04SJed Brown 14838a381b04SJed Brown PetscFunctionBegin; 14848a381b04SJed Brown if (ark->tableau) { 14859566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ark->tableau->name, arktype, &match)); 14863ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS); 14878a381b04SJed Brown } 14888a381b04SJed Brown for (link = ARKTableauList; link; link = link->next) { 14899566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, arktype, &match)); 14908a381b04SJed Brown if (match) { 14919566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauReset(ts)); 14928a381b04SJed Brown ark->tableau = &link->tab; 14939566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauSetUp(ts)); 14942ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 14953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14968a381b04SJed Brown } 14978a381b04SJed Brown } 149898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", arktype); 14998a381b04SJed Brown } 1500e0877f53SBarry Smith 1501d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts, PetscBool flg) 1502d71ae5a4SJacob Faibussowitsch { 15034cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 15044cc180ffSJed Brown 15054cc180ffSJed Brown PetscFunctionBegin; 15064cc180ffSJed Brown ark->imex = (PetscBool)!flg; 15073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15084cc180ffSJed Brown } 15098a381b04SJed Brown 1510d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts, PetscBool *flg) 1511d71ae5a4SJacob Faibussowitsch { 15123a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX *)ts->data; 15133a28c0e5SStefano Zampini 15143a28c0e5SStefano Zampini PetscFunctionBegin; 15153a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 15163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15173a28c0e5SStefano Zampini } 15183a28c0e5SStefano Zampini 1519d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1520d71ae5a4SJacob Faibussowitsch { 1521b3a6b972SJed Brown PetscFunctionBegin; 15229566063dSJacob Faibussowitsch PetscCall(TSReset_ARKIMEX(ts)); 1523b3a6b972SJed Brown if (ts->dm) { 15249566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSARKIMEX, DMRestrictHook_TSARKIMEX, ts)); 15259566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSARKIMEX, DMSubDomainRestrictHook_TSARKIMEX, ts)); 1526b3a6b972SJed Brown } 15279566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 15289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", NULL)); 15299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", NULL)); 15309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", NULL)); 15312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", NULL)); 15323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1533b3a6b972SJed Brown } 1534b3a6b972SJed Brown 15358a381b04SJed Brown /* ------------------------------------------------------------ */ 15368a381b04SJed Brown /*MC 1537c79dcfbdSBarry Smith TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes 15388a381b04SJed Brown 1539fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1540fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1541bcf0153eSBarry Smith of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`. 1542d0685a90SJed Brown 15438a381b04SJed Brown Level: beginner 15448a381b04SJed Brown 1545bcf0153eSBarry Smith Notes: 1546bcf0153eSBarry Smith The default is `TSARKIMEX3`, it can be changed with `TSARKIMEXSetType()` or -ts_arkimex_type 15478a381b04SJed Brown 1548bcf0153eSBarry Smith If the equation is implicit or a DAE, then `TSSetEquationType()` needs to be set accordingly. Refer to the manual for further information. 1549bcf0153eSBarry Smith 1550bcf0153eSBarry 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). 1551bcf0153eSBarry Smith 1552bcf0153eSBarry Smith Consider trying `TSROSW` if the stiff part is linear or weakly nonlinear. 1553bcf0153eSBarry Smith 1554*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSARKIMEXSetType()`, `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()`, `TSARKIMEXGetFullyImplicit()`, 1555bcf0153eSBarry Smith `TSARKIMEX1BEE`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEX3`, `TSARKIMEXL2`, `TSARKIMEXA2`, `TSARKIMEXARS122`, 1556bcf0153eSBarry Smith `TSARKIMEX4`, `TSARKIMEX5`, `TSARKIMEXPRSSP2`, `TSARKIMEXARS443`, `TSARKIMEXBPR3`, `TSARKIMEXType`, `TSARKIMEXRegister()`, `TSType` 15578a381b04SJed Brown M*/ 1558d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 1559d71ae5a4SJacob Faibussowitsch { 156080ab5e31SHong Zhang TS_ARKIMEX *ark; 15618a381b04SJed Brown 15628a381b04SJed Brown PetscFunctionBegin; 15639566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 15648a381b04SJed Brown 15658a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 156680ab5e31SHong Zhang ts->ops->adjointreset = TSAdjointReset_ARKIMEX; 15678a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 15688a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1569f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 15708a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 157180ab5e31SHong Zhang ts->ops->adjointsetup = TSAdjointSetUp_ARKIMEX; 15728a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1573cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1574108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 157524655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 15768a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 15778a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 15788a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 157980ab5e31SHong Zhang ts->ops->getstages = TSGetStages_ARKIMEX; 158080ab5e31SHong Zhang ts->ops->adjointstep = TSAdjointStep_ARKIMEX; 15818a381b04SJed Brown 1582efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1583efd4aadfSBarry Smith 158480ab5e31SHong Zhang PetscCall(PetscNew(&ark)); 158580ab5e31SHong Zhang ts->data = (void *)ark; 158680ab5e31SHong Zhang ark->imex = PETSC_TRUE; 158780ab5e31SHong Zhang 158880ab5e31SHong Zhang ark->VecsDeltaLam = NULL; 158980ab5e31SHong Zhang ark->VecsSensiTemp = NULL; 159080ab5e31SHong Zhang ark->VecsSensiPTemp = NULL; 15918a381b04SJed Brown 15929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetType_C", TSARKIMEXGetType_ARKIMEX)); 15939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetType_C", TSARKIMEXSetType_ARKIMEX)); 15949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXSetFullyImplicit_C", TSARKIMEXSetFullyImplicit_ARKIMEX)); 15959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSARKIMEXGetFullyImplicit_C", TSARKIMEXGetFullyImplicit_ARKIMEX)); 159696400bd6SLisandro Dalcin 15979566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts, TSARKIMEXDefault)); 15983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15998a381b04SJed Brown } 1600