18a381b04SJed Brown /* 28a381b04SJed Brown Code for timestepping with additive Runge-Kutta IMEX method 38a381b04SJed Brown 48a381b04SJed Brown Notes: 58a381b04SJed Brown The general system is written as 68a381b04SJed Brown 7f9c1d6abSBarry Smith F(t,U,Udot) = G(t,U) 88a381b04SJed Brown 98a381b04SJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 108a381b04SJed Brown 118a381b04SJed Brown */ 12af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 131e25c274SJed Brown #include <petscdm.h> 148a381b04SJed Brown 1519fd82e9SBarry Smith static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3; 168a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled; 178a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized; 1856dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec); 198a381b04SJed Brown 208a381b04SJed Brown typedef struct _ARKTableau *ARKTableau; 218a381b04SJed Brown struct _ARKTableau { 228a381b04SJed Brown char *name; 234f385281SJed Brown PetscInt order; /* Classical approximation order of the method */ 244f385281SJed Brown PetscInt s; /* Number of stages */ 25e817cc15SEmil Constantinescu PetscBool stiffly_accurate; /* The implicit part is stiffly accurate*/ 26e817cc15SEmil Constantinescu PetscBool FSAL_implicit; /* The implicit part is FSAL*/ 27e817cc15SEmil Constantinescu PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/ 284f385281SJed Brown PetscInt pinterp; /* Interpolation order */ 294f385281SJed Brown PetscReal *At,*bt,*ct; /* Stiff tableau */ 308a381b04SJed Brown PetscReal *A,*b,*c; /* Non-stiff tableau */ 31108c343cSJed Brown PetscReal *bembedt,*bembed; /* Embedded formula of order one less (order-1) */ 32cd652676SJed Brown PetscReal *binterpt,*binterp; /* Dense output formula */ 33108c343cSJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 348a381b04SJed Brown }; 358a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 368a381b04SJed Brown struct _ARKTableauLink { 378a381b04SJed Brown struct _ARKTableau tab; 388a381b04SJed Brown ARKTableauLink next; 398a381b04SJed Brown }; 408a381b04SJed Brown static ARKTableauLink ARKTableauList; 418a381b04SJed Brown 428a381b04SJed Brown typedef struct { 438a381b04SJed Brown ARKTableau tableau; 448a381b04SJed Brown Vec *Y; /* States computed during the step */ 458a381b04SJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 468a381b04SJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 4756dcabbaSDebojyoti Ghosh Vec *Y_prev; /* States computed during the previous time step */ 4856dcabbaSDebojyoti Ghosh Vec *YdotI_prev; /* Time derivatives for the stiff part for the previous time step*/ 4956dcabbaSDebojyoti Ghosh Vec *YdotRHS_prev; /* Function evaluations for the non-stiff part for the previous time step*/ 50e817cc15SEmil Constantinescu Vec Ydot0; /* Holds the slope from the previous step in FSAL case */ 518a381b04SJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 528a381b04SJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 538a381b04SJed Brown PetscScalar *work; /* Scalar work */ 54b296d7d5SJed Brown PetscReal scoeff; /* shift = scoeff/dt */ 558a381b04SJed Brown PetscReal stage_time; 564cc180ffSJed Brown PetscBool imex; 5796400bd6SLisandro Dalcin PetscBool extrapolate; /* Extrapolate initial guess from previous time-step stage values */ 58108c343cSJed Brown TSStepStatus status; 598a381b04SJed Brown } TS_ARKIMEX; 601f80e275SEmil Constantinescu /*MC 611f80e275SEmil Constantinescu TSARKIMEXARS122 - Second order ARK IMEX scheme. 628a381b04SJed Brown 631f80e275SEmil Constantinescu This method has one explicit stage and one implicit stage. 641f80e275SEmil Constantinescu 659bd3e852SBarry Smith Options Database: 6667b8a455SSatish Balay . -ts_arkimex_type ars122 - set arkimex type to ars122 679bd3e852SBarry Smith 681f80e275SEmil Constantinescu References: 69606c0280SSatish Balay . * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 701f80e275SEmil Constantinescu 711f80e275SEmil Constantinescu Level: advanced 721f80e275SEmil Constantinescu 73db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 741f80e275SEmil Constantinescu M*/ 751f80e275SEmil Constantinescu /*MC 761f80e275SEmil Constantinescu TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 771f80e275SEmil Constantinescu 781f80e275SEmil Constantinescu This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu. 791f80e275SEmil Constantinescu 809bd3e852SBarry Smith Options Database: 8167b8a455SSatish Balay . -ts_arkimex_type a2 - set arkimex type to a2 829bd3e852SBarry Smith 831f80e275SEmil Constantinescu Level: advanced 841f80e275SEmil Constantinescu 85db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 861f80e275SEmil Constantinescu M*/ 871f80e275SEmil Constantinescu /*MC 881f80e275SEmil Constantinescu TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 891f80e275SEmil Constantinescu 901f80e275SEmil Constantinescu This method has two implicit stages, and L-stable implicit scheme. 911f80e275SEmil Constantinescu 929bd3e852SBarry Smith Options Database: 9367b8a455SSatish Balay . -ts_arkimex_type l2 - set arkimex type to l2 949bd3e852SBarry Smith 951f80e275SEmil Constantinescu References: 96606c0280SSatish Balay . * - L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005. 971f80e275SEmil Constantinescu 981f80e275SEmil Constantinescu Level: advanced 991f80e275SEmil Constantinescu 100db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1011f80e275SEmil Constantinescu M*/ 1021f80e275SEmil Constantinescu /*MC 103c79dcfbdSBarry Smith TSARKIMEX1BEE - First order backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method. 104e817cc15SEmil Constantinescu 105e817cc15SEmil Constantinescu This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used. 106e817cc15SEmil Constantinescu 1079bd3e852SBarry Smith Options Database: 10867b8a455SSatish Balay . -ts_arkimex_type 1bee - set arkimex type to 1bee 1099bd3e852SBarry Smith 110e817cc15SEmil Constantinescu Level: advanced 111e817cc15SEmil Constantinescu 112db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 113e817cc15SEmil Constantinescu M*/ 114e817cc15SEmil Constantinescu /*MC 1151f80e275SEmil Constantinescu TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 1161f80e275SEmil Constantinescu 1171f80e275SEmil Constantinescu This method has one explicit stage and two implicit stages. The implicit part is the same as in TSARKIMEX2D and TSARKIMEX2E, but the explicit part has a larger stability region on the negative real axis. This method was provided by Emil Constantinescu. 1181f80e275SEmil Constantinescu 1199bd3e852SBarry Smith Options Database: 12067b8a455SSatish Balay . -ts_arkimex_type 2c - set arkimex type to 2c 1219bd3e852SBarry Smith 1221f80e275SEmil Constantinescu Level: advanced 1231f80e275SEmil Constantinescu 124db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1251f80e275SEmil Constantinescu M*/ 12664f491ddSJed Brown /*MC 12764f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 12864f491ddSJed Brown 129617a39beSEmil Constantinescu This method has one explicit stage and two implicit stages. The stability function is independent of the explicit part in the infinity limit of the implict component. This method was provided by Emil Constantinescu. 13064f491ddSJed Brown 1319bd3e852SBarry Smith Options Database: 13267b8a455SSatish Balay . -ts_arkimex_type 2d - set arkimex type to 2d 1339bd3e852SBarry Smith 134b330ce4dSSatish Balay Level: advanced 135b330ce4dSSatish Balay 136db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 13764f491ddSJed Brown M*/ 13864f491ddSJed Brown /*MC 13964f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 14064f491ddSJed Brown 14164f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 14264f491ddSJed Brown 1439bd3e852SBarry Smith Options Database: 14467b8a455SSatish Balay . -ts_arkimex_type 2e - set arkimex type to 2e 1459bd3e852SBarry Smith 146b330ce4dSSatish Balay Level: advanced 147b330ce4dSSatish Balay 148db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 14964f491ddSJed Brown M*/ 15064f491ddSJed Brown /*MC 1516cf0794eSJed Brown TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 1526cf0794eSJed Brown 1536cf0794eSJed Brown This method has three implicit stages. 1546cf0794eSJed Brown 1556cf0794eSJed Brown References: 156606c0280SSatish Balay . * - L. Pareschi, G. Russo, Implicit Explicit Runge Kutta schemes and applications to hyperbolic systems with relaxations. Journal of Scientific Computing Volume: 25, Issue: 1, October, 2005. 1576cf0794eSJed Brown 158a8d69d7bSBarry Smith This method is referred to as SSP2-(3,3,2) in https://arxiv.org/abs/1110.4375 1596cf0794eSJed Brown 1609bd3e852SBarry Smith Options Database: 16167b8a455SSatish Balay . -ts_arkimex_type prssp2 - set arkimex type to prssp2 1629bd3e852SBarry Smith 1636cf0794eSJed Brown Level: advanced 1646cf0794eSJed Brown 165db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1666cf0794eSJed Brown M*/ 1676cf0794eSJed Brown /*MC 16864f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 16964f491ddSJed Brown 17064f491ddSJed Brown This method has one explicit stage and three implicit stages. 17164f491ddSJed Brown 1729bd3e852SBarry Smith Options Database: 17367b8a455SSatish Balay . -ts_arkimex_type 3 - set arkimex type to 3 1749bd3e852SBarry Smith 17564f491ddSJed Brown References: 176606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 17764f491ddSJed Brown 178b330ce4dSSatish Balay Level: advanced 179b330ce4dSSatish Balay 180db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 18164f491ddSJed Brown M*/ 18264f491ddSJed Brown /*MC 1836cf0794eSJed Brown TSARKIMEXARS443 - Third order ARK IMEX scheme. 1846cf0794eSJed Brown 1856cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1866cf0794eSJed Brown 1879bd3e852SBarry Smith Options Database: 18867b8a455SSatish Balay . -ts_arkimex_type ars443 - set arkimex type to ars443 1899bd3e852SBarry Smith 1906cf0794eSJed Brown References: 191606c0280SSatish Balay + * - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 192606c0280SSatish Balay - * - This method is referred to as ARS(4,4,3) in https://arxiv.org/abs/1110.4375 1936cf0794eSJed Brown 1946cf0794eSJed Brown Level: advanced 1956cf0794eSJed Brown 196db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 1976cf0794eSJed Brown M*/ 1986cf0794eSJed Brown /*MC 1996cf0794eSJed Brown TSARKIMEXBPR3 - Third order ARK IMEX scheme. 2006cf0794eSJed Brown 2016cf0794eSJed Brown This method has one explicit stage and four implicit stages. 2026cf0794eSJed Brown 2039bd3e852SBarry Smith Options Database: 20467b8a455SSatish Balay . -ts_arkimex_type bpr3 - set arkimex type to bpr3 2059bd3e852SBarry Smith 2066cf0794eSJed Brown References: 207606c0280SSatish Balay . * - This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375 2086cf0794eSJed Brown 2096cf0794eSJed Brown Level: advanced 2106cf0794eSJed Brown 211db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 2126cf0794eSJed Brown M*/ 2136cf0794eSJed Brown /*MC 21464f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 21564f491ddSJed Brown 21664f491ddSJed Brown This method has one explicit stage and four implicit stages. 21764f491ddSJed Brown 2189bd3e852SBarry Smith Options Database: 21967b8a455SSatish Balay . -ts_arkimex_type 4 - set arkimex type to4 2209bd3e852SBarry Smith 22164f491ddSJed Brown References: 222606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 22364f491ddSJed Brown 224b330ce4dSSatish Balay Level: advanced 225b330ce4dSSatish Balay 226db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 22764f491ddSJed Brown M*/ 22864f491ddSJed Brown /*MC 22964f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 23064f491ddSJed Brown 23164f491ddSJed Brown This method has one explicit stage and five implicit stages. 23264f491ddSJed Brown 2339bd3e852SBarry Smith Options Database: 23467b8a455SSatish Balay . -ts_arkimex_type 5 - set arkimex type to 5 2359bd3e852SBarry Smith 23664f491ddSJed Brown References: 237606c0280SSatish Balay . * - Kennedy and Carpenter 2003. 23864f491ddSJed Brown 239b330ce4dSSatish Balay Level: advanced 240b330ce4dSSatish Balay 241db781477SPatrick Sanan .seealso: `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEXSetType()` 24264f491ddSJed Brown M*/ 24364f491ddSJed Brown 2448a381b04SJed Brown /*@C 2458a381b04SJed Brown TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 2468a381b04SJed Brown 247fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 2488a381b04SJed Brown 2498a381b04SJed Brown Level: advanced 2508a381b04SJed Brown 251db781477SPatrick Sanan .seealso: `TSARKIMEXRegisterDestroy()` 2528a381b04SJed Brown @*/ 2538a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 2548a381b04SJed Brown { 2558a381b04SJed Brown PetscFunctionBegin; 2568a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2578a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 258e817cc15SEmil Constantinescu 259e817cc15SEmil Constantinescu { 260e817cc15SEmil Constantinescu const PetscReal 261e817cc15SEmil Constantinescu A[3][3] = {{0.0,0.0,0.0}, 262e817cc15SEmil Constantinescu {0.0,0.0,0.0}, 263748ad121SEmil Constantinescu {0.0,0.5,0.0}}, 264e817cc15SEmil Constantinescu At[3][3] = {{1.0,0.0,0.0}, 265e817cc15SEmil Constantinescu {0.0,0.5,0.0}, 266e817cc15SEmil Constantinescu {0.0,0.5,0.5}}, 267e817cc15SEmil Constantinescu b[3] = {0.0,0.5,0.5}, 268e817cc15SEmil Constantinescu bembedt[3] = {1.0,0.0,0.0}; 2699566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL)); 270e817cc15SEmil Constantinescu } 2718a381b04SJed Brown { 2728a381b04SJed Brown const PetscReal 2731f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2741f80e275SEmil Constantinescu {0.5,0.0}}, 2751f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2761f80e275SEmil Constantinescu {0.0,0.5}}, 2771f80e275SEmil Constantinescu b[2] = {0.0,1.0}, 2781f80e275SEmil Constantinescu bembedt[2] = {0.5,0.5}; 2791f80e275SEmil 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*/ 2809566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL)); 2811f80e275SEmil Constantinescu } 2821f80e275SEmil Constantinescu { 2831f80e275SEmil Constantinescu const PetscReal 2841f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2851f80e275SEmil Constantinescu {1.0,0.0}}, 2861f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2871f80e275SEmil Constantinescu {0.5,0.5}}, 2881f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2891f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}; 2901f80e275SEmil 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*/ 2919566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL)); 2921f80e275SEmil Constantinescu } 2931f80e275SEmil Constantinescu { 294da80777bSKarl 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 */ 2951f80e275SEmil Constantinescu const PetscReal 2961f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2971f80e275SEmil Constantinescu {1.0,0.0}}, 298da80777bSKarl Rupp At[2][2] = {{0.2928932188134524755992,0.0}, 299da80777bSKarl Rupp {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}}, 3001f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 3011f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}, 302da80777bSKarl Rupp binterpt[2][2] = {{ (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}, 303da80777bSKarl Rupp {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}}, 3041f80e275SEmil Constantinescu binterp[2][2] = {{1.0,-0.5},{0.0,0.5}}; 3059566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0])); 3061f80e275SEmil Constantinescu } 3071f80e275SEmil Constantinescu { 308da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 309da80777bSKarl Rupp const PetscReal 3108a381b04SJed Brown A[3][3] = {{0,0,0}, 311da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 312617a39beSEmil Constantinescu {0.5,0.5,0}}, 313da80777bSKarl Rupp At[3][3] = {{0,0,0}, 314da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 315da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 316da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 317da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 318da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 319da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3209566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL)); 3211f80e275SEmil Constantinescu } 3221f80e275SEmil Constantinescu { 323da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 324da80777bSKarl Rupp const PetscReal 3251f80e275SEmil Constantinescu A[3][3] = {{0,0,0}, 326da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 3278a381b04SJed Brown {0.75,0.25,0}}, 328da80777bSKarl Rupp At[3][3] = {{0,0,0}, 329da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 330da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 331da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 332da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 333da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 334da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3359566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL)); 3368a381b04SJed Brown } 33706db7b1cSJed Brown { /* Optimal for linear implicit part */ 338da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 339da80777bSKarl Rupp const PetscReal 340da80777bSKarl Rupp A[3][3] = {{0,0,0}, 341da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 342da80777bSKarl Rupp {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}}, 343da80777bSKarl Rupp At[3][3] = {{0,0,0}, 344da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 345da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 346da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 347da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 348da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 349da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL)); 351a3a57f36SJed Brown } 3526cf0794eSJed Brown { /* Optimal for linear implicit part */ 3536cf0794eSJed Brown const PetscReal 3546cf0794eSJed Brown A[3][3] = {{0,0,0}, 3556cf0794eSJed Brown {0.5,0,0}, 3566cf0794eSJed Brown {0.5,0.5,0}}, 3576cf0794eSJed Brown At[3][3] = {{0.25,0,0}, 3586cf0794eSJed Brown {0,0.25,0}, 3596cf0794eSJed Brown {1./3,1./3,1./3}}; 3609566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL)); 3616cf0794eSJed Brown } 362a3a57f36SJed Brown { 363a3a57f36SJed Brown const PetscReal 364a3a57f36SJed Brown A[4][4] = {{0,0,0,0}, 3654040e9f2SJed Brown {1767732205903./2027836641118.,0,0,0}, 3664040e9f2SJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 3674040e9f2SJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 368a3a57f36SJed Brown At[4][4] = {{0,0,0,0}, 3694040e9f2SJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 3704040e9f2SJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 3714040e9f2SJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 372cc46b9d1SJed Brown bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.}, 3734040e9f2SJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 3744040e9f2SJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 3754040e9f2SJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 3764040e9f2SJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 3779566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL)); 378a3a57f36SJed Brown } 379a3a57f36SJed Brown { 380a3a57f36SJed Brown const PetscReal 381e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3826cf0794eSJed Brown {1./2,0,0,0,0}, 3836cf0794eSJed Brown {11./18,1./18,0,0,0}, 3846cf0794eSJed Brown {5./6,-5./6,.5,0,0}, 3856cf0794eSJed Brown {1./4,7./4,3./4,-7./4,0}}, 3866cf0794eSJed Brown At[5][5] = {{0,0,0,0,0}, 3876cf0794eSJed Brown {0,1./2,0,0,0}, 3886cf0794eSJed Brown {0,1./6,1./2,0,0}, 3896cf0794eSJed Brown {0,-1./2,1./2,1./2,0}, 390108c343cSJed Brown {0,3./2,-3./2,1./2,1./2}}, 3910298fd71SBarry Smith *bembedt = NULL; 3929566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL)); 3936cf0794eSJed Brown } 3946cf0794eSJed Brown { 3956cf0794eSJed Brown const PetscReal 396e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3976cf0794eSJed Brown {1,0,0,0,0}, 3986cf0794eSJed Brown {4./9,2./9,0,0,0}, 3996cf0794eSJed Brown {1./4,0,3./4,0,0}, 4006cf0794eSJed Brown {1./4,0,3./5,0,0}}, 401e74514c0SSatish Balay At[5][5] = {{0,0,0,0,0}, 4026cf0794eSJed Brown {.5,.5,0,0,0}, 4036cf0794eSJed Brown {5./18,-1./9,.5,0,0}, 4046cf0794eSJed Brown {.5,0,0,.5,0}, 405108c343cSJed Brown {.25,0,.75,-.5,.5}}, 4060298fd71SBarry Smith *bembedt = NULL; 4079566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL)); 4086cf0794eSJed Brown } 4096cf0794eSJed Brown { 4106cf0794eSJed Brown const PetscReal 411a3a57f36SJed Brown A[6][6] = {{0,0,0,0,0,0}, 412a3a57f36SJed Brown {1./2,0,0,0,0,0}, 4134040e9f2SJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 4144040e9f2SJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 4154040e9f2SJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 4164040e9f2SJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 417a3a57f36SJed Brown At[6][6] = {{0,0,0,0,0,0}, 418a3a57f36SJed Brown {1./4,1./4,0,0,0,0}, 4194040e9f2SJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 4204040e9f2SJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 4214040e9f2SJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 4224040e9f2SJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 423cc46b9d1SJed Brown bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.}, 4244040e9f2SJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 425cd652676SJed Brown {0,0,0}, 4264040e9f2SJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 4274040e9f2SJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 4284040e9f2SJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 4294040e9f2SJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 4309566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL)); 431a3a57f36SJed Brown } 432a3a57f36SJed Brown { 433a3a57f36SJed Brown const PetscReal 434a3a57f36SJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 435a3a57f36SJed Brown {41./100,0,0,0,0,0,0,0}, 4364040e9f2SJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 4374040e9f2SJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 4384040e9f2SJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 4394040e9f2SJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 4404040e9f2SJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 4414040e9f2SJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 442a3a57f36SJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 4434040e9f2SJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 4444040e9f2SJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 4454040e9f2SJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 4464040e9f2SJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 4474040e9f2SJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 4484040e9f2SJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 4494040e9f2SJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 450cc46b9d1SJed Brown bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.}, 4514040e9f2SJed Brown binterpt[8][3] = {{-17674230611817./10670229744614., 43486358583215./12773830924787., -9257016797708./5021505065439.}, 452cd652676SJed Brown {0, 0, 0 }, 453cd652676SJed Brown {0, 0, 0 }, 4544040e9f2SJed Brown {65168852399939./7868540260826., -91478233927265./11067650958493., 26096422576131./11239449250142.}, 4554040e9f2SJed Brown {15494834004392./5936557850923., -79368583304911./10890268929626., 92396832856987./20362823103730.}, 4564040e9f2SJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473., 30029262896817./10175596800299.}, 4574040e9f2SJed Brown {-19024464361622./5461577185407., 115839755401235./10719374521269., -26136350496073./3983972220547.}, 4584040e9f2SJed Brown {-6511271360970./6095937251113., 5843115559534./2180450260947., -5289405421727./3760307252460. }}; 4599566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL)); 460a3a57f36SJed Brown } 4618a381b04SJed Brown PetscFunctionReturn(0); 4628a381b04SJed Brown } 4638a381b04SJed Brown 4648a381b04SJed Brown /*@C 4658a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 4668a381b04SJed Brown 4678a381b04SJed Brown Not Collective 4688a381b04SJed Brown 4698a381b04SJed Brown Level: advanced 4708a381b04SJed Brown 471db781477SPatrick Sanan .seealso: `TSARKIMEXRegister()`, `TSARKIMEXRegisterAll()` 4728a381b04SJed Brown @*/ 4738a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 4748a381b04SJed Brown { 4758a381b04SJed Brown ARKTableauLink link; 4768a381b04SJed Brown 4778a381b04SJed Brown PetscFunctionBegin; 4788a381b04SJed Brown while ((link = ARKTableauList)) { 4798a381b04SJed Brown ARKTableau t = &link->tab; 4808a381b04SJed Brown ARKTableauList = link->next; 4819566063dSJacob Faibussowitsch PetscCall(PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c)); 4829566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->bembedt,t->bembed)); 4839566063dSJacob Faibussowitsch PetscCall(PetscFree2(t->binterpt,t->binterp)); 4849566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 4859566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 4868a381b04SJed Brown } 4878a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4888a381b04SJed Brown PetscFunctionReturn(0); 4898a381b04SJed Brown } 4908a381b04SJed Brown 4918a381b04SJed Brown /*@C 4928a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 4938a690491SBarry Smith from TSInitializePackage(). 4948a381b04SJed Brown 4958a381b04SJed Brown Level: developer 4968a381b04SJed Brown 497db781477SPatrick Sanan .seealso: `PetscInitialize()` 4988a381b04SJed Brown @*/ 499607a6623SBarry Smith PetscErrorCode TSARKIMEXInitializePackage(void) 5008a381b04SJed Brown { 5018a381b04SJed Brown PetscFunctionBegin; 5028a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 5038a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 5049566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterAll()); 5059566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSARKIMEXFinalizePackage)); 5068a381b04SJed Brown PetscFunctionReturn(0); 5078a381b04SJed Brown } 5088a381b04SJed Brown 5098a381b04SJed Brown /*@C 5108a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 5118a381b04SJed Brown called from PetscFinalize(). 5128a381b04SJed Brown 5138a381b04SJed Brown Level: developer 5148a381b04SJed Brown 515db781477SPatrick Sanan .seealso: `PetscFinalize()` 5168a381b04SJed Brown @*/ 5178a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 5188a381b04SJed Brown { 5198a381b04SJed Brown PetscFunctionBegin; 5208a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 5219566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegisterDestroy()); 5228a381b04SJed Brown PetscFunctionReturn(0); 5238a381b04SJed Brown } 5248a381b04SJed Brown 525cd652676SJed Brown /*@C 526cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 527cd652676SJed Brown 528cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 529cd652676SJed Brown 530cd652676SJed Brown Input Parameters: 531cd652676SJed Brown + name - identifier for method 532cd652676SJed Brown . order - approximation order of method 533cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 534cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 5350298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 5360298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 537cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 5380298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 5390298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 5400298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 5410298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 542cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 543cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 5440298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 545cd652676SJed Brown 546cd652676SJed Brown Notes: 547cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 548cd652676SJed Brown 549cd652676SJed Brown Level: advanced 550cd652676SJed Brown 551db781477SPatrick Sanan .seealso: `TSARKIMEX` 552cd652676SJed Brown @*/ 55319fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s, 5548a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 555cd652676SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 556108c343cSJed Brown const PetscReal bembedt[],const PetscReal bembed[], 557cd652676SJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 5588a381b04SJed Brown { 5598a381b04SJed Brown ARKTableauLink link; 5608a381b04SJed Brown ARKTableau t; 5618a381b04SJed Brown PetscInt i,j; 5628a381b04SJed Brown 5638a381b04SJed Brown PetscFunctionBegin; 5649566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 5659566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 5668a381b04SJed Brown t = &link->tab; 5679566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&t->name)); 5688a381b04SJed Brown t->order = order; 5698a381b04SJed Brown t->s = s; 5709566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c)); 5719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->At,At,s*s)); 5729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A,A,s*s)); 5739566063dSJacob Faibussowitsch if (bt) PetscCall(PetscArraycpy(t->bt,bt,s)); 5748a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 5759566063dSJacob Faibussowitsch if (b) PetscCall(PetscArraycpy(t->b,b,s)); 5765dceddf7SDebojyoti Ghosh else for (i=0; i<s; i++) t->b[i] = t->bt[i]; 5779566063dSJacob Faibussowitsch if (ct) PetscCall(PetscArraycpy(t->ct,ct,s)); 5788a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 5799566063dSJacob Faibussowitsch if (c) PetscCall(PetscArraycpy(t->c,c,s)); 5808a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 581e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 582e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 583e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 584e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 585e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5864e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 587108c343cSJed Brown if (bembedt) { 5889566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s,&t->bembedt,s,&t->bembed)); 5899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembedt,bembedt,s)); 5909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bembed,bembed ? bembed : bembedt,s)); 591108c343cSJed Brown } 592108c343cSJed Brown 5934f385281SJed Brown t->pinterp = pinterp; 5949566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp)); 5959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterpt,binterpt,s*pinterp)); 5969566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp,binterp ? binterp : binterpt,s*pinterp)); 5978a381b04SJed Brown link->next = ARKTableauList; 5988a381b04SJed Brown ARKTableauList = link; 5998a381b04SJed Brown PetscFunctionReturn(0); 6008a381b04SJed Brown } 6018a381b04SJed Brown 602108c343cSJed Brown /* 603108c343cSJed Brown The step completion formula is 604108c343cSJed Brown 605108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 606108c343cSJed Brown 607108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 608108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 609108c343cSJed Brown We can write 610108c343cSJed Brown 611108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 612108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 613108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 614108c343cSJed Brown 615108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 616108c343cSJed Brown */ 617108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 618108c343cSJed Brown { 619108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 620108c343cSJed Brown ARKTableau tab = ark->tableau; 621108c343cSJed Brown PetscScalar *w = ark->work; 622108c343cSJed Brown PetscReal h; 623108c343cSJed Brown PetscInt s = tab->s,j; 624108c343cSJed Brown 625108c343cSJed Brown PetscFunctionBegin; 626108c343cSJed Brown switch (ark->status) { 627108c343cSJed Brown case TS_STEP_INCOMPLETE: 628108c343cSJed Brown case TS_STEP_PENDING: 629108c343cSJed Brown h = ts->time_step; break; 630108c343cSJed Brown case TS_STEP_COMPLETE: 631be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 632ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 633108c343cSJed Brown } 634108c343cSJed Brown if (order == tab->order) { 635e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 636740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 6379566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[s-1],X)); 638e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 6399566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 640e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bt[j]; 6419566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,w,ark->YdotI)); 642e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 643108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->b[j]; 6449566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,w,ark->YdotRHS)); 645e817cc15SEmil Constantinescu } 646e817cc15SEmil Constantinescu } 6479566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol,X)); 648108c343cSJed Brown if (done) *done = PETSC_TRUE; 649108c343cSJed Brown PetscFunctionReturn(0); 650108c343cSJed Brown } else if (order == tab->order-1) { 651108c343cSJed Brown if (!tab->bembedt) goto unavailable; 652108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 6539566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 654e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembedt[j]; 6559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,w,ark->YdotI)); 656108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 6579566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,w,ark->YdotRHS)); 658108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 6599566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 660e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]); 6619566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,tab->s,w,ark->YdotI)); 662108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 6639566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,w,ark->YdotRHS)); 664108c343cSJed Brown } 665108c343cSJed Brown if (done) *done = PETSC_TRUE; 666108c343cSJed Brown PetscFunctionReturn(0); 667108c343cSJed Brown } 668108c343cSJed Brown unavailable: 669108c343cSJed Brown if (done) *done = PETSC_FALSE; 67063a3b9bcSJacob Faibussowitsch else 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,tab->order,order); 671108c343cSJed Brown PetscFunctionReturn(0); 672108c343cSJed Brown } 673108c343cSJed Brown 674c79dcfbdSBarry Smith static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts,PetscBool *id) 675c79dcfbdSBarry Smith { 676c79dcfbdSBarry Smith Vec Udot,Y1,Y2; 677c79dcfbdSBarry Smith TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 678c79dcfbdSBarry Smith PetscReal norm; 679c79dcfbdSBarry Smith 680c79dcfbdSBarry Smith PetscFunctionBegin; 6819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&Udot)); 6829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&Y1)); 6839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&Y2)); 6849566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y1,ark->imex)); 6859566063dSJacob Faibussowitsch PetscCall(VecSetRandom(Udot,NULL)); 6869566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y2,ark->imex)); 6879566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2,-1.0,Y1)); 6889566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y2,-1.0,Udot)); 6899566063dSJacob Faibussowitsch PetscCall(VecNorm(Y2,NORM_2,&norm)); 690c79dcfbdSBarry Smith if (norm < 100.0*PETSC_MACHINE_EPSILON) { 691c79dcfbdSBarry Smith *id = PETSC_TRUE; 692c79dcfbdSBarry Smith } else { 693c79dcfbdSBarry Smith *id = PETSC_FALSE; 6949566063dSJacob 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)); 695c79dcfbdSBarry Smith } 6969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Udot)); 6979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y1)); 6989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y2)); 699c79dcfbdSBarry Smith PetscFunctionReturn(0); 700c79dcfbdSBarry Smith } 701c79dcfbdSBarry Smith 70224655328SShri static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 70324655328SShri { 70424655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 70524655328SShri ARKTableau tab = ark->tableau; 70624655328SShri const PetscInt s = tab->s; 70724655328SShri const PetscReal *bt = tab->bt,*b = tab->b; 70824655328SShri PetscScalar *w = ark->work; 70924655328SShri Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 71024655328SShri PetscInt j; 711be5899b3SLisandro Dalcin PetscReal h; 71224655328SShri 71324655328SShri PetscFunctionBegin; 714be5899b3SLisandro Dalcin switch (ark->status) { 715be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 716be5899b3SLisandro Dalcin case TS_STEP_PENDING: 717be5899b3SLisandro Dalcin h = ts->time_step; break; 718be5899b3SLisandro Dalcin case TS_STEP_COMPLETE: 719be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 720be5899b3SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 721be5899b3SLisandro Dalcin } 72224655328SShri for (j=0; j<s; j++) w[j] = -h*bt[j]; 7239566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotI)); 72424655328SShri for (j=0; j<s; j++) w[j] = -h*b[j]; 7259566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol,s,w,YdotRHS)); 72624655328SShri PetscFunctionReturn(0); 72724655328SShri } 72824655328SShri 7298a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 7308a381b04SJed Brown { 7318a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7328a381b04SJed Brown ARKTableau tab = ark->tableau; 7338a381b04SJed Brown const PetscInt s = tab->s; 73424655328SShri const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 735406d0ec2SJed Brown PetscScalar *w = ark->work; 7361297b224SEmil Constantinescu Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z; 73796400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 738108c343cSJed Brown TSAdapt adapt; 7398a381b04SJed Brown SNES snes; 740fecfb714SLisandro Dalcin PetscInt i,j,its,lits; 741be5899b3SLisandro Dalcin PetscInt rejections = 0; 74296400bd6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 74396400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7448a381b04SJed Brown 7458a381b04SJed Brown PetscFunctionBegin; 74696400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 7479566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev)); 7489566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev)); 7499566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev)); 75096400bd6SLisandro Dalcin } 75196400bd6SLisandro Dalcin 75296400bd6SLisandro Dalcin if (!ts->steprollback) { 75396400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 7549566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[s-1],Ydot0)); 75596400bd6SLisandro Dalcin } 756fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 75796400bd6SLisandro Dalcin for (i = 0; i<s; i++) { 7589566063dSJacob Faibussowitsch PetscCall(VecCopy(Y[i],ark->Y_prev[i])); 7599566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotRHS[i],ark->YdotRHS_prev[i])); 7609566063dSJacob Faibussowitsch PetscCall(VecCopy(YdotI[i],ark->YdotI_prev[i])); 76196400bd6SLisandro Dalcin } 76296400bd6SLisandro Dalcin } 76396400bd6SLisandro Dalcin } 76496400bd6SLisandro Dalcin 765fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 76696400bd6SLisandro Dalcin TS ts_start; 767c79dcfbdSBarry Smith if (PetscDefined(USE_DEBUG)) { 768c79dcfbdSBarry Smith PetscBool id = PETSC_FALSE; 7699566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTestMassIdentity(ts,&id)); 7703c633725SBarry 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"); 771c79dcfbdSBarry Smith } 7729566063dSJacob Faibussowitsch PetscCall(TSClone(ts,&ts_start)); 7739566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts_start,ts->vec_sol)); 7749566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts_start,ts->ptime)); 7759566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts_start,ts->steps+1)); 7769566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts_start,ts->ptime+ts->time_step)); 7779566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER)); 7789566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts_start,ts->time_step)); 7799566063dSJacob Faibussowitsch PetscCall(TSSetType(ts_start,TSARKIMEX)); 7809566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE)); 7819566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts_start,TSARKIMEX1BEE)); 78234561852SEmil Constantinescu 7839566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts_start)); 7849566063dSJacob Faibussowitsch PetscCall(TSSolve(ts_start,ts->vec_sol)); 7859566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts_start,&ts->ptime)); 7869566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts_start,&ts->time_step)); 787bbd56ea5SKarl Rupp 78885fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 78985fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data; 7909566063dSJacob Faibussowitsch PetscCall(VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0)); 79185fc7851SLisandro Dalcin } 79296400bd6SLisandro Dalcin ts->steps++; 79334561852SEmil Constantinescu 794d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 795d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 79696400bd6SLisandro Dalcin { 7979566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 7989566063dSJacob Faibussowitsch PetscCall(TSSetSNES(ts,snes)); 79996400bd6SLisandro Dalcin } 8009566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts_start)); 801e817cc15SEmil Constantinescu } 802e817cc15SEmil Constantinescu 803108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 80496400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 80596400bd6SLisandro Dalcin PetscReal t = ts->ptime; 806108c343cSJed Brown PetscReal h = ts->time_step; 8078a381b04SJed Brown for (i=0; i<s; i++) { 8089be3e283SDebojyoti Ghosh ark->stage_time = t + h*ct[i]; 8099566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,ark->stage_time)); 8108a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 8113c633725SBarry 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"); 8129566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,Y[i])); 813e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 8149566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],i,w,YdotI)); 8158a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 8169566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],i,w,YdotRHS)); 8178a381b04SJed Brown } else { 818b296d7d5SJed Brown ark->scoeff = 1./At[i*s+i]; 8198a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 8209566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,Z)); 821e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 8229566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z,i,w,YdotI)); 823c58d1302SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 8249566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z,i,w,YdotRHS)); 825fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 82656dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 8279566063dSJacob Faibussowitsch PetscCall(TSExtrapolate_ARKIMEX(ts,c[i],Y[i])); 82856dcabbaSDebojyoti Ghosh } else { 8298a381b04SJed Brown /* Initial guess taken from last stage */ 8309566063dSJacob Faibussowitsch PetscCall(VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i])); 83156dcabbaSDebojyoti Ghosh } 8329566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 8339566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,Y[i])); 8349566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 8359566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes,&lits)); 8365ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 8379566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 8389566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok)); 83996400bd6SLisandro Dalcin if (!stageok) { 8401be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8411be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 84296400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8431be93e3eSJed Brown goto reject_step; 8441be93e3eSJed Brown } 8458a381b04SJed Brown } 846e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 847e817cc15SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8483c633725SBarry Smith 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",ark->tableau->name); 8499566063dSJacob Faibussowitsch PetscCall(VecCopy(Ydot0,YdotI[0])); /* YdotI = YdotI(tn-1) */ 850e817cc15SEmil Constantinescu } else { 8519566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i])); /* YdotI = shift*(X-Z) */ 852e817cc15SEmil Constantinescu } 853e817cc15SEmil Constantinescu } else { 8545eca1a21SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8559566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Ydot)); 8569566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex));/* YdotI = -G(t,Y,0) */ 8579566063dSJacob Faibussowitsch PetscCall(VecScale(YdotI[i],-1.0)); 8585eca1a21SEmil Constantinescu } else { 8599566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i])); /* YdotI = shift*(X-Z) */ 8605eca1a21SEmil Constantinescu } 8614cc180ffSJed Brown if (ark->imex) { 8629566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i])); 8634cc180ffSJed Brown } else { 8649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(YdotRHS[i])); 8654cc180ffSJed Brown } 8668a381b04SJed Brown } 8679566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,ark->stage_time,i,Y)); 868e817cc15SEmil Constantinescu } 86996400bd6SLisandro Dalcin 870be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 8719566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL)); 872108c343cSJed Brown ark->status = TS_STEP_PENDING; 8739566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 8749566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 8759566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE)); 8769566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept)); 87796400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 87896400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 8799566063dSJacob Faibussowitsch PetscCall(TSRollBack_ARKIMEX(ts)); 880be5899b3SLisandro Dalcin ts->time_step = next_time_step; 881be5899b3SLisandro Dalcin goto reject_step; 88296400bd6SLisandro Dalcin } 88396400bd6SLisandro Dalcin 8848a381b04SJed Brown ts->ptime += ts->time_step; 885cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 886108c343cSJed Brown break; 88796400bd6SLisandro Dalcin 88896400bd6SLisandro Dalcin reject_step: 889fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 890be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 89196400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 89263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n",ts->steps,rejections)); 893108c343cSJed Brown } 894f85781f1SEmil Constantinescu } 8958a381b04SJed Brown PetscFunctionReturn(0); 8968a381b04SJed Brown } 8978a381b04SJed Brown 898cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 899cd652676SJed Brown { 900cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9014f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 902108c343cSJed Brown PetscReal h; 903108c343cSJed Brown PetscReal tt,t; 904cd652676SJed Brown PetscScalar *bt,*b; 905cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 906cd652676SJed Brown 907cd652676SJed Brown PetscFunctionBegin; 9083c633725SBarry Smith PetscCheck(Bt && B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 909108c343cSJed Brown switch (ark->status) { 910108c343cSJed Brown case TS_STEP_INCOMPLETE: 911108c343cSJed Brown case TS_STEP_PENDING: 912108c343cSJed Brown h = ts->time_step; 913108c343cSJed Brown t = (itime - ts->ptime)/h; 914108c343cSJed Brown break; 915108c343cSJed Brown case TS_STEP_COMPLETE: 916be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 917108c343cSJed Brown t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 918108c343cSJed Brown break; 919ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 920108c343cSJed Brown } 9219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s,&bt,s,&b)); 922cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 9234f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 924cd652676SJed Brown for (i=0; i<s; i++) { 925c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 926108c343cSJed Brown b[i] += h * B[i*pinterp+j] * tt; 927cd652676SJed Brown } 928cd652676SJed Brown } 9299566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y[0],X)); 9309566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,bt,ark->YdotI)); 9319566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,b,ark->YdotRHS)); 9329566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt,b)); 933cd652676SJed Brown PetscFunctionReturn(0); 934cd652676SJed Brown } 935cd652676SJed Brown 93656dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 93756dcabbaSDebojyoti Ghosh { 93856dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 93956dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 940be5899b3SLisandro Dalcin PetscReal h,h_prev,t,tt; 94156dcabbaSDebojyoti Ghosh PetscScalar *bt,*b; 94256dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 94356dcabbaSDebojyoti Ghosh 94456dcabbaSDebojyoti Ghosh PetscFunctionBegin; 9453c633725SBarry Smith PetscCheck(Bt && B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 9469566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(s,&bt,s,&b)); 94781d12688SDebojyoti Ghosh h = ts->time_step; 948be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 949be5899b3SLisandro Dalcin t = 1 + h/h_prev*c; 95056dcabbaSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 95156dcabbaSDebojyoti Ghosh for (i=0; i<s; i++) { 95281d12688SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 95356dcabbaSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 95456dcabbaSDebojyoti Ghosh } 95556dcabbaSDebojyoti Ghosh } 9563c633725SBarry Smith PetscCheck(ark->Y_prev,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 9579566063dSJacob Faibussowitsch PetscCall(VecCopy(ark->Y_prev[0],X)); 9589566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,bt,ark->YdotI_prev)); 9599566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,b,ark->YdotRHS_prev)); 9609566063dSJacob Faibussowitsch PetscCall(PetscFree2(bt,b)); 96156dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 96256dcabbaSDebojyoti Ghosh } 96356dcabbaSDebojyoti Ghosh 9648a381b04SJed Brown /*------------------------------------------------------------*/ 96596400bd6SLisandro Dalcin 96696400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauReset(TS ts) 96796400bd6SLisandro Dalcin { 96896400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 96996400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 97096400bd6SLisandro Dalcin 97196400bd6SLisandro Dalcin PetscFunctionBegin; 97296400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 9739566063dSJacob Faibussowitsch PetscCall(PetscFree(ark->work)); 9749566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->Y)); 9759566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->YdotI)); 9769566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->YdotRHS)); 9779566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->Y_prev)); 9789566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->YdotI_prev)); 9799566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&ark->YdotRHS_prev)); 98096400bd6SLisandro Dalcin PetscFunctionReturn(0); 98196400bd6SLisandro Dalcin } 98296400bd6SLisandro Dalcin 9838a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 9848a381b04SJed Brown { 9858a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9868a381b04SJed Brown 9878a381b04SJed Brown PetscFunctionBegin; 9889566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauReset(ts)); 9899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot)); 9909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Ydot0)); 9919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ark->Z)); 9928a381b04SJed Brown PetscFunctionReturn(0); 9938a381b04SJed Brown } 9948a381b04SJed Brown 995d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 996d5e6173cSPeter Brune { 997d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 998d5e6173cSPeter Brune 999d5e6173cSPeter Brune PetscFunctionBegin; 1000d5e6173cSPeter Brune if (Z) { 1001d5e6173cSPeter Brune if (dm && dm != ts->dm) { 10029566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z)); 1003d5e6173cSPeter Brune } else *Z = ax->Z; 1004d5e6173cSPeter Brune } 1005d5e6173cSPeter Brune if (Ydot) { 1006d5e6173cSPeter Brune if (dm && dm != ts->dm) { 10079566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot)); 1008d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 1009d5e6173cSPeter Brune } 1010d5e6173cSPeter Brune PetscFunctionReturn(0); 1011d5e6173cSPeter Brune } 1012d5e6173cSPeter Brune 1013d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1014d5e6173cSPeter Brune { 1015d5e6173cSPeter Brune PetscFunctionBegin; 1016d5e6173cSPeter Brune if (Z) { 1017d5e6173cSPeter Brune if (dm && dm != ts->dm) { 10189566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z)); 1019d5e6173cSPeter Brune } 1020d5e6173cSPeter Brune } 1021d5e6173cSPeter Brune if (Ydot) { 1022d5e6173cSPeter Brune if (dm && dm != ts->dm) { 10239566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot)); 1024d5e6173cSPeter Brune } 1025d5e6173cSPeter Brune } 1026d5e6173cSPeter Brune PetscFunctionReturn(0); 1027d5e6173cSPeter Brune } 1028d5e6173cSPeter Brune 10298a381b04SJed Brown /* 10308a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 10318a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 10328a381b04SJed Brown */ 10338a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 10348a381b04SJed Brown { 10358a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1036d5e6173cSPeter Brune DM dm,dmsave; 1037d5e6173cSPeter Brune Vec Z,Ydot; 1038b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10398a381b04SJed Brown 10408a381b04SJed Brown PetscFunctionBegin; 10419566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 10429566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,dm,&Z,&Ydot)); 10439566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X)); /* Ydot = shift*(X-Z) */ 1044d5e6173cSPeter Brune dmsave = ts->dm; 1045d5e6173cSPeter Brune ts->dm = dm; 1046740132f1SEmil Constantinescu 10479566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex)); 1048e817cc15SEmil Constantinescu 1049d5e6173cSPeter Brune ts->dm = dmsave; 10509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot)); 10518a381b04SJed Brown PetscFunctionReturn(0); 10528a381b04SJed Brown } 10538a381b04SJed Brown 1054d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 10558a381b04SJed Brown { 10568a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1057d5e6173cSPeter Brune DM dm,dmsave; 1058d5e6173cSPeter Brune Vec Ydot; 1059b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10608a381b04SJed Brown 10618a381b04SJed Brown PetscFunctionBegin; 10629566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 10639566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,dm,NULL,&Ydot)); 10648a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1065d5e6173cSPeter Brune dmsave = ts->dm; 1066d5e6173cSPeter Brune ts->dm = dm; 1067740132f1SEmil Constantinescu 10689566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex)); 1069740132f1SEmil Constantinescu 1070d5e6173cSPeter Brune ts->dm = dmsave; 10719566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot)); 1072d5e6173cSPeter Brune PetscFunctionReturn(0); 1073d5e6173cSPeter Brune } 1074d5e6173cSPeter Brune 1075d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1076d5e6173cSPeter Brune { 1077d5e6173cSPeter Brune PetscFunctionBegin; 1078d5e6173cSPeter Brune PetscFunctionReturn(0); 1079d5e6173cSPeter Brune } 1080d5e6173cSPeter Brune 1081d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1082d5e6173cSPeter Brune { 1083d5e6173cSPeter Brune TS ts = (TS)ctx; 1084d5e6173cSPeter Brune Vec Z,Z_c; 1085d5e6173cSPeter Brune 1086d5e6173cSPeter Brune PetscFunctionBegin; 10879566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,fine,&Z,NULL)); 10889566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL)); 10899566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct,Z,Z_c)); 10909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Z_c,rscale,Z_c)); 10919566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,fine,&Z,NULL)); 10929566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL)); 10938a381b04SJed Brown PetscFunctionReturn(0); 10948a381b04SJed Brown } 10958a381b04SJed Brown 1096cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1097cdb298fcSPeter Brune { 1098cdb298fcSPeter Brune PetscFunctionBegin; 1099cdb298fcSPeter Brune PetscFunctionReturn(0); 1100cdb298fcSPeter Brune } 1101cdb298fcSPeter Brune 1102cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1103cdb298fcSPeter Brune { 1104cdb298fcSPeter Brune TS ts = (TS)ctx; 1105cdb298fcSPeter Brune Vec Z,Z_c; 1106cdb298fcSPeter Brune 1107cdb298fcSPeter Brune PetscFunctionBegin; 11089566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,dm,&Z,NULL)); 11099566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL)); 1110cdb298fcSPeter Brune 11119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD)); 11129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD)); 1113cdb298fcSPeter Brune 11149566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,dm,&Z,NULL)); 11159566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL)); 1116cdb298fcSPeter Brune PetscFunctionReturn(0); 1117cdb298fcSPeter Brune } 1118cdb298fcSPeter Brune 111996400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 112096400bd6SLisandro Dalcin { 112196400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 112296400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 112396400bd6SLisandro Dalcin 112496400bd6SLisandro Dalcin PetscFunctionBegin; 11259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&ark->work)); 11269566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y)); 11279566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI)); 11289566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS)); 112996400bd6SLisandro Dalcin if (ark->extrapolate) { 11309566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev)); 11319566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev)); 11329566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev)); 113396400bd6SLisandro Dalcin } 113496400bd6SLisandro Dalcin PetscFunctionReturn(0); 113596400bd6SLisandro Dalcin } 113696400bd6SLisandro Dalcin 11378a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 11388a381b04SJed Brown { 11398a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1140d5e6173cSPeter Brune DM dm; 114196400bd6SLisandro Dalcin SNES snes; 1142f9c1d6abSBarry Smith 11438a381b04SJed Brown PetscFunctionBegin; 11449566063dSJacob Faibussowitsch PetscCall(TSARKIMEXTableauSetUp(ts)); 11459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&ark->Ydot)); 11469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&ark->Ydot0)); 11479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&ark->Z)); 11489566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 11499566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts)); 11509566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts)); 11519566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 11528a381b04SJed Brown PetscFunctionReturn(0); 11538a381b04SJed Brown } 11548a381b04SJed Brown /*------------------------------------------------------------*/ 11558a381b04SJed Brown 11564416b707SBarry Smith static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 11578a381b04SJed Brown { 11584cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11598a381b04SJed Brown 11608a381b04SJed Brown PetscFunctionBegin; 1161d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"ARKIMEX ODE solver options"); 11628a381b04SJed Brown { 11638a381b04SJed Brown ARKTableauLink link; 11648a381b04SJed Brown PetscInt count,choice; 11658a381b04SJed Brown PetscBool flg; 11668a381b04SJed Brown const char **namelist; 11678a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 11689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count,(char***)&namelist)); 11698a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11709566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg)); 11719566063dSJacob Faibussowitsch if (flg) PetscCall(TSARKIMEXSetType(ts,namelist[choice])); 11729566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 117396400bd6SLisandro Dalcin 11744cc180ffSJed Brown flg = (PetscBool) !ark->imex; 11759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL)); 11764cc180ffSJed Brown ark->imex = (PetscBool) !flg; 11779566063dSJacob 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)); 11788a381b04SJed Brown } 1179d0609cedSBarry Smith PetscOptionsHeadEnd(); 11808a381b04SJed Brown PetscFunctionReturn(0); 11818a381b04SJed Brown } 11828a381b04SJed Brown 11838a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 11848a381b04SJed Brown { 11858a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11868a381b04SJed Brown PetscBool iascii; 11878a381b04SJed Brown 11888a381b04SJed Brown PetscFunctionBegin; 11899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 11908a381b04SJed Brown if (iascii) { 11919c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 119219fd82e9SBarry Smith TSARKIMEXType arktype; 11938a381b04SJed Brown char buf[512]; 11943a28c0e5SStefano Zampini PetscBool flg; 11953a28c0e5SStefano Zampini 11969566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetType(ts,&arktype)); 11979566063dSJacob Faibussowitsch PetscCall(TSARKIMEXGetFullyImplicit(ts,&flg)); 11989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype)); 11999566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct)); 12009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf)); 12019566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c)); 12029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Fully implicit: %s\n",flg ? "yes" : "no")); 12039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no")); 12049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no")); 12059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no")); 12069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf)); 12078a381b04SJed Brown } 12088a381b04SJed Brown PetscFunctionReturn(0); 12098a381b04SJed Brown } 12108a381b04SJed Brown 1211f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1212f2c2a1b9SBarry Smith { 1213f2c2a1b9SBarry Smith SNES snes; 12149c334d8fSLisandro Dalcin TSAdapt adapt; 1215f2c2a1b9SBarry Smith 1216f2c2a1b9SBarry Smith PetscFunctionBegin; 12179566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 12189566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt,viewer)); 12199566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 12209566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes,viewer)); 1221ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 12229566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,NULL,NULL,ts)); 12239566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,NULL,NULL,NULL,ts)); 1224f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1225f2c2a1b9SBarry Smith } 1226f2c2a1b9SBarry Smith 12278a381b04SJed Brown /*@C 12288a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 12298a381b04SJed Brown 12308a381b04SJed Brown Logically collective 12318a381b04SJed Brown 1232d8d19677SJose E. Roman Input Parameters: 12338a381b04SJed Brown + ts - timestepping context 12348a381b04SJed Brown - arktype - type of ARK-IMEX scheme 12358a381b04SJed Brown 12369bd3e852SBarry Smith Options Database: 123767b8a455SSatish Balay . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> - set ARK IMEX scheme type 12389bd3e852SBarry Smith 12398a381b04SJed Brown Level: intermediate 12408a381b04SJed Brown 1241db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEX`, `TSARKIMEXType`, `TSARKIMEX1BEE`, `TSARKIMEXA2`, `TSARKIMEXL2`, `TSARKIMEXARS122`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEXPRSSP2`, 1242db781477SPatrick Sanan `TSARKIMEX3`, `TSARKIMEXBPR3`, `TSARKIMEXARS443`, `TSARKIMEX4`, `TSARKIMEX5` 12438a381b04SJed Brown @*/ 124419fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 12458a381b04SJed Brown { 12468a381b04SJed Brown PetscFunctionBegin; 12478a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1248b92453a8SLisandro Dalcin PetscValidCharPointer(arktype,2); 1249cac4c232SBarry Smith PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype)); 12508a381b04SJed Brown PetscFunctionReturn(0); 12518a381b04SJed Brown } 12528a381b04SJed Brown 12538a381b04SJed Brown /*@C 12548a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 12558a381b04SJed Brown 12568a381b04SJed Brown Logically collective 12578a381b04SJed Brown 12588a381b04SJed Brown Input Parameter: 12598a381b04SJed Brown . ts - timestepping context 12608a381b04SJed Brown 12618a381b04SJed Brown Output Parameter: 12628a381b04SJed Brown . arktype - type of ARK-IMEX scheme 12638a381b04SJed Brown 12648a381b04SJed Brown Level: intermediate 12658a381b04SJed Brown 1266db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()` 12678a381b04SJed Brown @*/ 126819fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 12698a381b04SJed Brown { 12708a381b04SJed Brown PetscFunctionBegin; 12718a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1272cac4c232SBarry Smith PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype)); 12738a381b04SJed Brown PetscFunctionReturn(0); 12748a381b04SJed Brown } 12758a381b04SJed Brown 127616353aafSBarry Smith /*@ 12774cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 12784cc180ffSJed Brown 12794cc180ffSJed Brown Logically collective 12804cc180ffSJed Brown 1281d8d19677SJose E. Roman Input Parameters: 12824cc180ffSJed Brown + ts - timestepping context 12834cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 12844cc180ffSJed Brown 12854cc180ffSJed Brown Level: intermediate 12864cc180ffSJed Brown 1287db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEXGetFullyImplicit()` 12884cc180ffSJed Brown @*/ 12894cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 12904cc180ffSJed Brown { 12914cc180ffSJed Brown PetscFunctionBegin; 12924cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12933a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts,flg,2); 1294cac4c232SBarry Smith PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg)); 12954cc180ffSJed Brown PetscFunctionReturn(0); 12964cc180ffSJed Brown } 12974cc180ffSJed Brown 12983a28c0e5SStefano Zampini /*@ 12993a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 13003a28c0e5SStefano Zampini 13013a28c0e5SStefano Zampini Logically collective 13023a28c0e5SStefano Zampini 13033a28c0e5SStefano Zampini Input Parameter: 13043a28c0e5SStefano Zampini . ts - timestepping context 13053a28c0e5SStefano Zampini 13067a7aea1fSJed Brown Output Parameter: 13073a28c0e5SStefano Zampini . flg - PETSC_TRUE for fully implicit 13083a28c0e5SStefano Zampini 13093a28c0e5SStefano Zampini Level: intermediate 13103a28c0e5SStefano Zampini 1311db781477SPatrick Sanan .seealso: `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()` 13123a28c0e5SStefano Zampini @*/ 13133a28c0e5SStefano Zampini PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts,PetscBool *flg) 13143a28c0e5SStefano Zampini { 13153a28c0e5SStefano Zampini PetscFunctionBegin; 13163a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1317dadcf809SJacob Faibussowitsch PetscValidBoolPointer(flg,2); 1318cac4c232SBarry Smith PetscUseMethod(ts,"TSARKIMEXGetFullyImplicit_C",(TS,PetscBool*),(ts,flg)); 13193a28c0e5SStefano Zampini PetscFunctionReturn(0); 13203a28c0e5SStefano Zampini } 13213a28c0e5SStefano Zampini 1322e0877f53SBarry Smith static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 13238a381b04SJed Brown { 13248a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13258a381b04SJed Brown 13268a381b04SJed Brown PetscFunctionBegin; 13278a381b04SJed Brown *arktype = ark->tableau->name; 13288a381b04SJed Brown PetscFunctionReturn(0); 13298a381b04SJed Brown } 1330e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype) 13318a381b04SJed Brown { 13328a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13338a381b04SJed Brown PetscBool match; 13348a381b04SJed Brown ARKTableauLink link; 13358a381b04SJed Brown 13368a381b04SJed Brown PetscFunctionBegin; 13378a381b04SJed Brown if (ark->tableau) { 13389566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(ark->tableau->name,arktype,&match)); 13398a381b04SJed Brown if (match) PetscFunctionReturn(0); 13408a381b04SJed Brown } 13418a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 13429566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name,arktype,&match)); 13438a381b04SJed Brown if (match) { 13449566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauReset(ts)); 13458a381b04SJed Brown ark->tableau = &link->tab; 13469566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSARKIMEXTableauSetUp(ts)); 13472ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13488a381b04SJed Brown PetscFunctionReturn(0); 13498a381b04SJed Brown } 13508a381b04SJed Brown } 135198921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 13528a381b04SJed Brown } 1353e0877f53SBarry Smith 1354e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 13554cc180ffSJed Brown { 13564cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13574cc180ffSJed Brown 13584cc180ffSJed Brown PetscFunctionBegin; 13594cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13604cc180ffSJed Brown PetscFunctionReturn(0); 13614cc180ffSJed Brown } 13628a381b04SJed Brown 13633a28c0e5SStefano Zampini static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts,PetscBool *flg) 13643a28c0e5SStefano Zampini { 13653a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13663a28c0e5SStefano Zampini 13673a28c0e5SStefano Zampini PetscFunctionBegin; 13683a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 13693a28c0e5SStefano Zampini PetscFunctionReturn(0); 13703a28c0e5SStefano Zampini } 13713a28c0e5SStefano Zampini 1372b3a6b972SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1373b3a6b972SJed Brown { 1374b3a6b972SJed Brown PetscFunctionBegin; 13759566063dSJacob Faibussowitsch PetscCall(TSReset_ARKIMEX(ts)); 1376b3a6b972SJed Brown if (ts->dm) { 13779566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts)); 13789566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts)); 1379b3a6b972SJed Brown } 13809566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 13819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL)); 13829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL)); 13839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL)); 1384*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetFullyImplicit_C",NULL)); 1385b3a6b972SJed Brown PetscFunctionReturn(0); 1386b3a6b972SJed Brown } 1387b3a6b972SJed Brown 13888a381b04SJed Brown /* ------------------------------------------------------------ */ 13898a381b04SJed Brown /*MC 1390c79dcfbdSBarry Smith TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes 13918a381b04SJed Brown 1392fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1393fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1394fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1395fca742c7SJed Brown 1396fca742c7SJed Brown Notes: 1397a4386c9eSJed Brown The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1398c8058688SBarry Smith 13995eca1a21SEmil Constantinescu If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 14005eca1a21SEmil Constantinescu 1401a4386c9eSJed Brown Methods with an explicit stage can only be used with ODE in which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 1402fca742c7SJed Brown 1403d0685a90SJed Brown Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1404d0685a90SJed Brown 14058a381b04SJed Brown Level: beginner 14068a381b04SJed Brown 1407db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSARKIMEXSetType()`, `TSARKIMEXGetType()`, `TSARKIMEXSetFullyImplicit()`, `TSARKIMEXGetFullyImplicit()`, 1408db781477SPatrick Sanan `TSARKIMEX1BEE`, `TSARKIMEX2C`, `TSARKIMEX2D`, `TSARKIMEX2E`, `TSARKIMEX3`, `TSARKIMEXL2`, `TSARKIMEXA2`, `TSARKIMEXARS122`, 1409db781477SPatrick Sanan `TSARKIMEX4`, `TSARKIMEX5`, `TSARKIMEXPRSSP2`, `TSARKIMEXARS443`, `TSARKIMEXBPR3`, `TSARKIMEXType`, `TSARKIMEXRegister()` 14108a381b04SJed Brown 14118a381b04SJed Brown M*/ 14128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 14138a381b04SJed Brown { 14148a381b04SJed Brown TS_ARKIMEX *th; 14158a381b04SJed Brown 14168a381b04SJed Brown PetscFunctionBegin; 14179566063dSJacob Faibussowitsch PetscCall(TSARKIMEXInitializePackage()); 14188a381b04SJed Brown 14198a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 14208a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 14218a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1422f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 14238a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 14248a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1425cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1426108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 142724655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 14288a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 14298a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 14308a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 14318a381b04SJed Brown 1432efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1433efd4aadfSBarry Smith 14349566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&th)); 14358a381b04SJed Brown ts->data = (void*)th; 14364cc180ffSJed Brown th->imex = PETSC_TRUE; 14378a381b04SJed Brown 14389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX)); 14399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX)); 14409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX)); 14419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetFullyImplicit_C",TSARKIMEXGetFullyImplicit_ARKIMEX)); 144296400bd6SLisandro Dalcin 14439566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetType(ts,TSARKIMEXDefault)); 14448a381b04SJed Brown PetscFunctionReturn(0); 14458a381b04SJed Brown } 1446