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: 669bd3e852SBarry Smith . -ts_arkimex_type ars122 679bd3e852SBarry Smith 681f80e275SEmil Constantinescu References: 6996a0c994SBarry Smith . 1. - 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 739bd3e852SBarry Smith .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: 819bd3e852SBarry Smith . -ts_arkimex_type a2 829bd3e852SBarry Smith 831f80e275SEmil Constantinescu Level: advanced 841f80e275SEmil Constantinescu 859bd3e852SBarry Smith .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: 939bd3e852SBarry Smith . -ts_arkimex_type l2 949bd3e852SBarry Smith 951f80e275SEmil Constantinescu References: 9696a0c994SBarry Smith . 1. - 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 1009bd3e852SBarry Smith .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: 1089bd3e852SBarry Smith . -ts_arkimex_type 1bee 1099bd3e852SBarry Smith 110e817cc15SEmil Constantinescu Level: advanced 111e817cc15SEmil Constantinescu 1129bd3e852SBarry Smith .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: 1209bd3e852SBarry Smith . -ts_arkimex_type 2c 1219bd3e852SBarry Smith 1221f80e275SEmil Constantinescu Level: advanced 1231f80e275SEmil Constantinescu 1249bd3e852SBarry Smith .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: 1329bd3e852SBarry Smith . -ts_arkimex_type 2d 1339bd3e852SBarry Smith 134b330ce4dSSatish Balay Level: advanced 135b330ce4dSSatish Balay 1369bd3e852SBarry Smith .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: 1449bd3e852SBarry Smith . -ts_arkimex_type 2e 1459bd3e852SBarry Smith 146b330ce4dSSatish Balay Level: advanced 147b330ce4dSSatish Balay 1489bd3e852SBarry Smith .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: 15696a0c994SBarry Smith . 1. - 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: 1619bd3e852SBarry Smith . -ts_arkimex_type prssp2 1629bd3e852SBarry Smith 1636cf0794eSJed Brown Level: advanced 1646cf0794eSJed Brown 1659bd3e852SBarry Smith .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: 1739bd3e852SBarry Smith . -ts_arkimex_type 3 1749bd3e852SBarry Smith 17564f491ddSJed Brown References: 17696a0c994SBarry Smith . 1. - Kennedy and Carpenter 2003. 17764f491ddSJed Brown 178b330ce4dSSatish Balay Level: advanced 179b330ce4dSSatish Balay 1809bd3e852SBarry Smith .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: 1889bd3e852SBarry Smith . -ts_arkimex_type ars443 1899bd3e852SBarry Smith 1906cf0794eSJed Brown References: 19196a0c994SBarry Smith + 1. - U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997). 192a8d69d7bSBarry Smith - 2. - 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 1969bd3e852SBarry Smith .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: 2049bd3e852SBarry Smith . -ts_arkimex_type bpr3 2059bd3e852SBarry Smith 2066cf0794eSJed Brown References: 207a8d69d7bSBarry Smith . This method is referred to as ARK3 in https://arxiv.org/abs/1110.4375 2086cf0794eSJed Brown 2096cf0794eSJed Brown Level: advanced 2106cf0794eSJed Brown 2119bd3e852SBarry Smith .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: 2199bd3e852SBarry Smith . -ts_arkimex_type 4 2209bd3e852SBarry Smith 22164f491ddSJed Brown References: 22296a0c994SBarry Smith . Kennedy and Carpenter 2003. 22364f491ddSJed Brown 224b330ce4dSSatish Balay Level: advanced 225b330ce4dSSatish Balay 2269bd3e852SBarry Smith .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: 2349bd3e852SBarry Smith . -ts_arkimex_type 5 2359bd3e852SBarry Smith 23664f491ddSJed Brown References: 23796a0c994SBarry Smith . Kennedy and Carpenter 2003. 23864f491ddSJed Brown 239b330ce4dSSatish Balay Level: advanced 240b330ce4dSSatish Balay 2419bd3e852SBarry Smith .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 2518a381b04SJed Brown .seealso: TSARKIMEXRegisterDestroy() 2528a381b04SJed Brown @*/ 2538a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 2548a381b04SJed Brown { 2558a381b04SJed Brown PetscErrorCode ierr; 2568a381b04SJed Brown 2578a381b04SJed Brown PetscFunctionBegin; 2588a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2598a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 260e817cc15SEmil Constantinescu 261e817cc15SEmil Constantinescu { 262e817cc15SEmil Constantinescu const PetscReal 263e817cc15SEmil Constantinescu A[3][3] = {{0.0,0.0,0.0}, 264e817cc15SEmil Constantinescu {0.0,0.0,0.0}, 265748ad121SEmil Constantinescu {0.0,0.5,0.0}}, 266e817cc15SEmil Constantinescu At[3][3] = {{1.0,0.0,0.0}, 267e817cc15SEmil Constantinescu {0.0,0.5,0.0}, 268e817cc15SEmil Constantinescu {0.0,0.5,0.5}}, 269e817cc15SEmil Constantinescu b[3] = {0.0,0.5,0.5}, 270e817cc15SEmil Constantinescu bembedt[3] = {1.0,0.0,0.0}; 2710298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 272e817cc15SEmil Constantinescu } 2738a381b04SJed Brown { 2748a381b04SJed Brown const PetscReal 2751f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2761f80e275SEmil Constantinescu {0.5,0.0}}, 2771f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2781f80e275SEmil Constantinescu {0.0,0.5}}, 2791f80e275SEmil Constantinescu b[2] = {0.0,1.0}, 2801f80e275SEmil Constantinescu 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*/ 2820298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 2831f80e275SEmil Constantinescu } 2841f80e275SEmil Constantinescu { 2851f80e275SEmil Constantinescu const PetscReal 2861f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2871f80e275SEmil Constantinescu {1.0,0.0}}, 2881f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2891f80e275SEmil Constantinescu {0.5,0.5}}, 2901f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2911f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}; 2921f80e275SEmil 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*/ 2930298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 2941f80e275SEmil Constantinescu } 2951f80e275SEmil Constantinescu { 296da80777bSKarl 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 */ 2971f80e275SEmil Constantinescu const PetscReal 2981f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2991f80e275SEmil Constantinescu {1.0,0.0}}, 300da80777bSKarl Rupp At[2][2] = {{0.2928932188134524755992,0.0}, 301da80777bSKarl Rupp {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}}, 3021f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 3031f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}, 304da80777bSKarl Rupp binterpt[2][2] = {{ (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}, 305da80777bSKarl Rupp {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}}, 3061f80e275SEmil Constantinescu binterp[2][2] = {{1.0,-0.5},{0.0,0.5}}; 3070298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);CHKERRQ(ierr); 3081f80e275SEmil Constantinescu } 3091f80e275SEmil Constantinescu { 310da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 311da80777bSKarl Rupp const PetscReal 3128a381b04SJed Brown A[3][3] = {{0,0,0}, 313da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 314617a39beSEmil Constantinescu {0.5,0.5,0}}, 315da80777bSKarl Rupp At[3][3] = {{0,0,0}, 316da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 317da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 318da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 319da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 320da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 321da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3220298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 3231f80e275SEmil Constantinescu } 3241f80e275SEmil Constantinescu { 325da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 326da80777bSKarl Rupp const PetscReal 3271f80e275SEmil Constantinescu A[3][3] = {{0,0,0}, 328da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 3298a381b04SJed Brown {0.75,0.25,0}}, 330da80777bSKarl Rupp At[3][3] = {{0,0,0}, 331da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 332da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 333da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 334da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 335da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 336da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3370298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 3388a381b04SJed Brown } 33906db7b1cSJed Brown { /* Optimal for linear implicit part */ 340da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 341da80777bSKarl Rupp const PetscReal 342da80777bSKarl Rupp A[3][3] = {{0,0,0}, 343da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 344da80777bSKarl Rupp {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}}, 345da80777bSKarl Rupp At[3][3] = {{0,0,0}, 346da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 347da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 348da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 349da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 350da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 351da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3520298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 353a3a57f36SJed Brown } 3546cf0794eSJed Brown { /* Optimal for linear implicit part */ 3556cf0794eSJed Brown const PetscReal 3566cf0794eSJed Brown A[3][3] = {{0,0,0}, 3576cf0794eSJed Brown {0.5,0,0}, 3586cf0794eSJed Brown {0.5,0.5,0}}, 3596cf0794eSJed Brown At[3][3] = {{0.25,0,0}, 3606cf0794eSJed Brown {0,0.25,0}, 3616cf0794eSJed Brown {1./3,1./3,1./3}}; 3620298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr); 3636cf0794eSJed Brown } 364a3a57f36SJed Brown { 365a3a57f36SJed Brown const PetscReal 366a3a57f36SJed Brown A[4][4] = {{0,0,0,0}, 3674040e9f2SJed Brown {1767732205903./2027836641118.,0,0,0}, 3684040e9f2SJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 3694040e9f2SJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 370a3a57f36SJed Brown At[4][4] = {{0,0,0,0}, 3714040e9f2SJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 3724040e9f2SJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 3734040e9f2SJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 374cc46b9d1SJed Brown bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.}, 3754040e9f2SJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 3764040e9f2SJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 3774040e9f2SJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 3784040e9f2SJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 3790298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 380a3a57f36SJed Brown } 381a3a57f36SJed Brown { 382a3a57f36SJed Brown const PetscReal 383e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3846cf0794eSJed Brown {1./2,0,0,0,0}, 3856cf0794eSJed Brown {11./18,1./18,0,0,0}, 3866cf0794eSJed Brown {5./6,-5./6,.5,0,0}, 3876cf0794eSJed Brown {1./4,7./4,3./4,-7./4,0}}, 3886cf0794eSJed Brown At[5][5] = {{0,0,0,0,0}, 3896cf0794eSJed Brown {0,1./2,0,0,0}, 3906cf0794eSJed Brown {0,1./6,1./2,0,0}, 3916cf0794eSJed Brown {0,-1./2,1./2,1./2,0}, 392108c343cSJed Brown {0,3./2,-3./2,1./2,1./2}}, 3930298fd71SBarry Smith *bembedt = NULL; 3940298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 3956cf0794eSJed Brown } 3966cf0794eSJed Brown { 3976cf0794eSJed Brown const PetscReal 398e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3996cf0794eSJed Brown {1,0,0,0,0}, 4006cf0794eSJed Brown {4./9,2./9,0,0,0}, 4016cf0794eSJed Brown {1./4,0,3./4,0,0}, 4026cf0794eSJed Brown {1./4,0,3./5,0,0}}, 403e74514c0SSatish Balay At[5][5] = {{0,0,0,0,0}, 4046cf0794eSJed Brown {.5,.5,0,0,0}, 4056cf0794eSJed Brown {5./18,-1./9,.5,0,0}, 4066cf0794eSJed Brown {.5,0,0,.5,0}, 407108c343cSJed Brown {.25,0,.75,-.5,.5}}, 4080298fd71SBarry Smith *bembedt = NULL; 4090298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 4106cf0794eSJed Brown } 4116cf0794eSJed Brown { 4126cf0794eSJed Brown const PetscReal 413a3a57f36SJed Brown A[6][6] = {{0,0,0,0,0,0}, 414a3a57f36SJed Brown {1./2,0,0,0,0,0}, 4154040e9f2SJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 4164040e9f2SJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 4174040e9f2SJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 4184040e9f2SJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 419a3a57f36SJed Brown At[6][6] = {{0,0,0,0,0,0}, 420a3a57f36SJed Brown {1./4,1./4,0,0,0,0}, 4214040e9f2SJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 4224040e9f2SJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 4234040e9f2SJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 4244040e9f2SJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 425cc46b9d1SJed Brown bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.}, 4264040e9f2SJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 427cd652676SJed Brown {0,0,0}, 4284040e9f2SJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 4294040e9f2SJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 4304040e9f2SJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 4314040e9f2SJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 4320298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 433a3a57f36SJed Brown } 434a3a57f36SJed Brown { 435a3a57f36SJed Brown const PetscReal 436a3a57f36SJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 437a3a57f36SJed Brown {41./100,0,0,0,0,0,0,0}, 4384040e9f2SJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 4394040e9f2SJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 4404040e9f2SJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 4414040e9f2SJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 4424040e9f2SJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 4434040e9f2SJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 444a3a57f36SJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 4454040e9f2SJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 4464040e9f2SJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 4474040e9f2SJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 4484040e9f2SJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 4494040e9f2SJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 4504040e9f2SJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 4514040e9f2SJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 452cc46b9d1SJed Brown bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.}, 4534040e9f2SJed Brown binterpt[8][3] = {{-17674230611817./10670229744614., 43486358583215./12773830924787., -9257016797708./5021505065439.}, 454cd652676SJed Brown {0, 0, 0 }, 455cd652676SJed Brown {0, 0, 0 }, 4564040e9f2SJed Brown {65168852399939./7868540260826., -91478233927265./11067650958493., 26096422576131./11239449250142.}, 4574040e9f2SJed Brown {15494834004392./5936557850923., -79368583304911./10890268929626., 92396832856987./20362823103730.}, 4584040e9f2SJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473., 30029262896817./10175596800299.}, 4594040e9f2SJed Brown {-19024464361622./5461577185407., 115839755401235./10719374521269., -26136350496073./3983972220547.}, 4604040e9f2SJed Brown {-6511271360970./6095937251113., 5843115559534./2180450260947., -5289405421727./3760307252460. }}; 4610298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 462a3a57f36SJed Brown } 4638a381b04SJed Brown PetscFunctionReturn(0); 4648a381b04SJed Brown } 4658a381b04SJed Brown 4668a381b04SJed Brown /*@C 4678a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 4688a381b04SJed Brown 4698a381b04SJed Brown Not Collective 4708a381b04SJed Brown 4718a381b04SJed Brown Level: advanced 4728a381b04SJed Brown 473607a6623SBarry Smith .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll() 4748a381b04SJed Brown @*/ 4758a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 4768a381b04SJed Brown { 4778a381b04SJed Brown PetscErrorCode ierr; 4788a381b04SJed Brown ARKTableauLink link; 4798a381b04SJed Brown 4808a381b04SJed Brown PetscFunctionBegin; 4818a381b04SJed Brown while ((link = ARKTableauList)) { 4828a381b04SJed Brown ARKTableau t = &link->tab; 4838a381b04SJed Brown ARKTableauList = link->next; 4848a381b04SJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 485108c343cSJed Brown ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr); 486cd652676SJed Brown ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 4878a381b04SJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 4888a381b04SJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 4898a381b04SJed Brown } 4908a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4918a381b04SJed Brown PetscFunctionReturn(0); 4928a381b04SJed Brown } 4938a381b04SJed Brown 4948a381b04SJed Brown /*@C 4958a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 4968a690491SBarry Smith from TSInitializePackage(). 4978a381b04SJed Brown 4988a381b04SJed Brown Level: developer 4998a381b04SJed Brown 5008a381b04SJed Brown .seealso: PetscInitialize() 5018a381b04SJed Brown @*/ 502607a6623SBarry Smith PetscErrorCode TSARKIMEXInitializePackage(void) 5038a381b04SJed Brown { 5048a381b04SJed Brown PetscErrorCode ierr; 5058a381b04SJed Brown 5068a381b04SJed Brown PetscFunctionBegin; 5078a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 5088a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 5098a381b04SJed Brown ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 5108a381b04SJed Brown ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 5118a381b04SJed Brown PetscFunctionReturn(0); 5128a381b04SJed Brown } 5138a381b04SJed Brown 5148a381b04SJed Brown /*@C 5158a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 5168a381b04SJed Brown called from PetscFinalize(). 5178a381b04SJed Brown 5188a381b04SJed Brown Level: developer 5198a381b04SJed Brown 5208a381b04SJed Brown .seealso: PetscFinalize() 5218a381b04SJed Brown @*/ 5228a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 5238a381b04SJed Brown { 5248a381b04SJed Brown PetscErrorCode ierr; 5258a381b04SJed Brown 5268a381b04SJed Brown PetscFunctionBegin; 5278a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 5288a381b04SJed Brown ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 5298a381b04SJed Brown PetscFunctionReturn(0); 5308a381b04SJed Brown } 5318a381b04SJed Brown 532cd652676SJed Brown /*@C 533cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 534cd652676SJed Brown 535cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 536cd652676SJed Brown 537cd652676SJed Brown Input Parameters: 538cd652676SJed Brown + name - identifier for method 539cd652676SJed Brown . order - approximation order of method 540cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 541cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 5420298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 5430298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 544cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 5450298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 5460298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 5470298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 5480298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 549cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 550cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 5510298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 552cd652676SJed Brown 553cd652676SJed Brown Notes: 554cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 555cd652676SJed Brown 556cd652676SJed Brown Level: advanced 557cd652676SJed Brown 558cd652676SJed Brown .seealso: TSARKIMEX 559cd652676SJed Brown @*/ 56019fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s, 5618a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 562cd652676SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 563108c343cSJed Brown const PetscReal bembedt[],const PetscReal bembed[], 564cd652676SJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 5658a381b04SJed Brown { 5668a381b04SJed Brown PetscErrorCode ierr; 5678a381b04SJed Brown ARKTableauLink link; 5688a381b04SJed Brown ARKTableau t; 5698a381b04SJed Brown PetscInt i,j; 5708a381b04SJed Brown 5718a381b04SJed Brown PetscFunctionBegin; 5721d36bdfdSBarry Smith ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 573071fcb05SBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 5748a381b04SJed Brown t = &link->tab; 5758a381b04SJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 5768a381b04SJed Brown t->order = order; 5778a381b04SJed Brown t->s = s; 578dcca6d9dSJed Brown ierr = PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 579580bdb30SBarry Smith ierr = PetscArraycpy(t->At,At,s*s);CHKERRQ(ierr); 580580bdb30SBarry Smith ierr = PetscArraycpy(t->A,A,s*s);CHKERRQ(ierr); 581580bdb30SBarry Smith if (bt) { ierr = PetscArraycpy(t->bt,bt,s);CHKERRQ(ierr); } 5828a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 583580bdb30SBarry Smith if (b) { ierr = PetscArraycpy(t->b,b,s);CHKERRQ(ierr); } 5845dceddf7SDebojyoti Ghosh else for (i=0; i<s; i++) t->b[i] = t->bt[i]; 585580bdb30SBarry Smith if (ct) { ierr = PetscArraycpy(t->ct,ct,s);CHKERRQ(ierr); } 5868a381b04SJed 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]; 587580bdb30SBarry Smith if (c) { ierr = PetscArraycpy(t->c,c,s);CHKERRQ(ierr); } 5888a381b04SJed 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]; 589e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 590e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 591e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 592e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 593e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5944e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 595108c343cSJed Brown if (bembedt) { 596dcca6d9dSJed Brown ierr = PetscMalloc2(s,&t->bembedt,s,&t->bembed);CHKERRQ(ierr); 597580bdb30SBarry Smith ierr = PetscArraycpy(t->bembedt,bembedt,s);CHKERRQ(ierr); 598580bdb30SBarry Smith ierr = PetscArraycpy(t->bembed,bembed ? bembed : bembedt,s);CHKERRQ(ierr); 599108c343cSJed Brown } 600108c343cSJed Brown 6014f385281SJed Brown t->pinterp = pinterp; 602dcca6d9dSJed Brown ierr = PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);CHKERRQ(ierr); 603580bdb30SBarry Smith ierr = PetscArraycpy(t->binterpt,binterpt,s*pinterp);CHKERRQ(ierr); 604580bdb30SBarry Smith ierr = PetscArraycpy(t->binterp,binterp ? binterp : binterpt,s*pinterp);CHKERRQ(ierr); 6058a381b04SJed Brown link->next = ARKTableauList; 6068a381b04SJed Brown ARKTableauList = link; 6078a381b04SJed Brown PetscFunctionReturn(0); 6088a381b04SJed Brown } 6098a381b04SJed Brown 610108c343cSJed Brown /* 611108c343cSJed Brown The step completion formula is 612108c343cSJed Brown 613108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 614108c343cSJed Brown 615108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 616108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 617108c343cSJed Brown We can write 618108c343cSJed Brown 619108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 620108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 621108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 622108c343cSJed Brown 623108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 624108c343cSJed Brown */ 625108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 626108c343cSJed Brown { 627108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 628108c343cSJed Brown ARKTableau tab = ark->tableau; 629108c343cSJed Brown PetscScalar *w = ark->work; 630108c343cSJed Brown PetscReal h; 631108c343cSJed Brown PetscInt s = tab->s,j; 632108c343cSJed Brown PetscErrorCode ierr; 633108c343cSJed Brown 634108c343cSJed Brown PetscFunctionBegin; 635108c343cSJed Brown switch (ark->status) { 636108c343cSJed Brown case TS_STEP_INCOMPLETE: 637108c343cSJed Brown case TS_STEP_PENDING: 638108c343cSJed Brown h = ts->time_step; break; 639108c343cSJed Brown case TS_STEP_COMPLETE: 640be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 641ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 642108c343cSJed Brown } 643108c343cSJed Brown if (order == tab->order) { 644e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 645740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 646e817cc15SEmil Constantinescu ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr); 647e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 648108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 649e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bt[j]; 650108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 651e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 652108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->b[j]; 653108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 654e817cc15SEmil Constantinescu } 655e817cc15SEmil Constantinescu } 656108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 657108c343cSJed Brown if (done) *done = PETSC_TRUE; 658108c343cSJed Brown PetscFunctionReturn(0); 659108c343cSJed Brown } else if (order == tab->order-1) { 660108c343cSJed Brown if (!tab->bembedt) goto unavailable; 661108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 662108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 663e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembedt[j]; 664108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 665108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 666108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 667108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 668108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 669e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]); 670108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr); 671108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 672108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 673108c343cSJed Brown } 674108c343cSJed Brown if (done) *done = PETSC_TRUE; 675108c343cSJed Brown PetscFunctionReturn(0); 676108c343cSJed Brown } 677108c343cSJed Brown unavailable: 678108c343cSJed Brown if (done) *done = PETSC_FALSE; 679a7fac7c2SEmil Constantinescu else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order); 680108c343cSJed Brown PetscFunctionReturn(0); 681108c343cSJed Brown } 682108c343cSJed Brown 683c79dcfbdSBarry Smith static PetscErrorCode TSARKIMEXTestMassIdentity(TS ts,PetscBool *id) 684c79dcfbdSBarry Smith { 685c79dcfbdSBarry Smith PetscErrorCode ierr; 686c79dcfbdSBarry Smith Vec Udot,Y1,Y2; 687c79dcfbdSBarry Smith TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 688c79dcfbdSBarry Smith PetscReal norm; 689c79dcfbdSBarry Smith 690c79dcfbdSBarry Smith PetscFunctionBegin; 691c79dcfbdSBarry Smith ierr = VecDuplicate(ts->vec_sol,&Udot);CHKERRQ(ierr); 692c79dcfbdSBarry Smith ierr = VecDuplicate(ts->vec_sol,&Y1);CHKERRQ(ierr); 693c79dcfbdSBarry Smith ierr = VecDuplicate(ts->vec_sol,&Y2);CHKERRQ(ierr); 694c79dcfbdSBarry Smith ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y1,ark->imex);CHKERRQ(ierr); 695c79dcfbdSBarry Smith ierr = VecSetRandom(Udot,NULL);CHKERRQ(ierr); 696c79dcfbdSBarry Smith ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,Udot,Y2,ark->imex);CHKERRQ(ierr); 697c79dcfbdSBarry Smith ierr = VecAXPY(Y2,-1.0,Y1);CHKERRQ(ierr); 698c79dcfbdSBarry Smith ierr = VecAXPY(Y2,-1.0,Udot);CHKERRQ(ierr); 699c79dcfbdSBarry Smith ierr = VecNorm(Y2,NORM_2,&norm);CHKERRQ(ierr); 700c79dcfbdSBarry Smith if (norm < 100.0*PETSC_MACHINE_EPSILON) { 701c79dcfbdSBarry Smith *id = PETSC_TRUE; 702c79dcfbdSBarry Smith } else { 703c79dcfbdSBarry Smith *id = PETSC_FALSE; 704c79dcfbdSBarry Smith ierr = PetscInfo1((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);CHKERRQ(ierr); 705c79dcfbdSBarry Smith } 706c79dcfbdSBarry Smith ierr = VecDestroy(&Udot);CHKERRQ(ierr); 707c79dcfbdSBarry Smith ierr = VecDestroy(&Y1);CHKERRQ(ierr); 708c79dcfbdSBarry Smith ierr = VecDestroy(&Y2);CHKERRQ(ierr); 709c79dcfbdSBarry Smith PetscFunctionReturn(0); 710c79dcfbdSBarry Smith } 711c79dcfbdSBarry Smith 71224655328SShri static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 71324655328SShri { 71424655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 71524655328SShri ARKTableau tab = ark->tableau; 71624655328SShri const PetscInt s = tab->s; 71724655328SShri const PetscReal *bt = tab->bt,*b = tab->b; 71824655328SShri PetscScalar *w = ark->work; 71924655328SShri Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 72024655328SShri PetscInt j; 721be5899b3SLisandro Dalcin PetscReal h; 72224655328SShri PetscErrorCode ierr; 72324655328SShri 72424655328SShri PetscFunctionBegin; 725be5899b3SLisandro Dalcin switch (ark->status) { 726be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 727be5899b3SLisandro Dalcin case TS_STEP_PENDING: 728be5899b3SLisandro Dalcin h = ts->time_step; break; 729be5899b3SLisandro Dalcin case TS_STEP_COMPLETE: 730be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 731be5899b3SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 732be5899b3SLisandro Dalcin } 73324655328SShri for (j=0; j<s; j++) w[j] = -h*bt[j]; 73424655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 73524655328SShri for (j=0; j<s; j++) w[j] = -h*b[j]; 73624655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 73724655328SShri PetscFunctionReturn(0); 73824655328SShri } 73924655328SShri 7408a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 7418a381b04SJed Brown { 7428a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7438a381b04SJed Brown ARKTableau tab = ark->tableau; 7448a381b04SJed Brown const PetscInt s = tab->s; 74524655328SShri const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 746406d0ec2SJed Brown PetscScalar *w = ark->work; 7471297b224SEmil Constantinescu Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z; 74896400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 749108c343cSJed Brown TSAdapt adapt; 7508a381b04SJed Brown SNES snes; 751fecfb714SLisandro Dalcin PetscInt i,j,its,lits; 752be5899b3SLisandro Dalcin PetscInt rejections = 0; 75396400bd6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 75496400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7558a381b04SJed Brown PetscErrorCode ierr; 7568a381b04SJed Brown 7578a381b04SJed Brown PetscFunctionBegin; 75896400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 75996400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 76096400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 76196400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 76296400bd6SLisandro Dalcin } 76396400bd6SLisandro Dalcin 76496400bd6SLisandro Dalcin if (!ts->steprollback) { 76596400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 76696400bd6SLisandro Dalcin ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr); 76796400bd6SLisandro Dalcin } 768fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 76996400bd6SLisandro Dalcin for (i = 0; i<s; i++) { 77096400bd6SLisandro Dalcin ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr); 77196400bd6SLisandro Dalcin ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr); 77296400bd6SLisandro Dalcin ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr); 77396400bd6SLisandro Dalcin } 77496400bd6SLisandro Dalcin } 77596400bd6SLisandro Dalcin } 77696400bd6SLisandro Dalcin 777fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 77896400bd6SLisandro Dalcin TS ts_start; 779c79dcfbdSBarry Smith if (PetscDefined(USE_DEBUG)) { 780c79dcfbdSBarry Smith PetscBool id = PETSC_FALSE; 781c79dcfbdSBarry Smith ierr = TSARKIMEXTestMassIdentity(ts,&id);CHKERRQ(ierr); 782c79dcfbdSBarry Smith if (!id) SETERRQ(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"); 783c79dcfbdSBarry Smith } 784baa10174SEmil Constantinescu ierr = TSClone(ts,&ts_start);CHKERRQ(ierr); 785e817cc15SEmil Constantinescu ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr); 786e817cc15SEmil Constantinescu ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr); 78719eac22cSLisandro Dalcin ierr = TSSetMaxSteps(ts_start,ts->steps+1);CHKERRQ(ierr); 78819eac22cSLisandro Dalcin ierr = TSSetMaxTime(ts_start,ts->ptime+ts->time_step);CHKERRQ(ierr); 789feed9e9dSBarry Smith ierr = TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 790740132f1SEmil Constantinescu ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr); 791e817cc15SEmil Constantinescu ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr); 792740132f1SEmil Constantinescu ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr); 793e817cc15SEmil Constantinescu ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr); 79434561852SEmil Constantinescu 795dcb233daSLisandro Dalcin ierr = TSRestartStep(ts_start);CHKERRQ(ierr); 796e817cc15SEmil Constantinescu ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr); 797e817cc15SEmil Constantinescu ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr); 79896400bd6SLisandro Dalcin ierr = TSGetTimeStep(ts_start,&ts->time_step);CHKERRQ(ierr); 799bbd56ea5SKarl Rupp 80085fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 80185fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data; 80285fc7851SLisandro Dalcin ierr = VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);CHKERRQ(ierr); 80385fc7851SLisandro Dalcin } 80496400bd6SLisandro Dalcin ts->steps++; 80534561852SEmil Constantinescu 806d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 807d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 80896400bd6SLisandro Dalcin { 80996400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 81096400bd6SLisandro Dalcin ierr = TSSetSNES(ts,snes);CHKERRQ(ierr); 81196400bd6SLisandro Dalcin } 812166a6834SEmil Constantinescu ierr = TSDestroy(&ts_start);CHKERRQ(ierr); 813e817cc15SEmil Constantinescu } 814e817cc15SEmil Constantinescu 815108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 81696400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 81796400bd6SLisandro Dalcin PetscReal t = ts->ptime; 818108c343cSJed Brown PetscReal h = ts->time_step; 8198a381b04SJed Brown for (i=0; i<s; i++) { 8209be3e283SDebojyoti Ghosh ark->stage_time = t + h*ct[i]; 82196400bd6SLisandro Dalcin ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr); 8228a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 8236c4ed002SBarry Smith if (i!=0 && ts->equation_type >= TS_EQ_IMPLICIT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Explicit stages other than the first one are not supported for implicit problems"); 8248a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 825e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 8268a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 8278a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 8288a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 8298a381b04SJed Brown } else { 830b296d7d5SJed Brown ark->scoeff = 1./At[i*s+i]; 8318a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 8328a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 833e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 8344f385281SJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 835c58d1302SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 836c58d1302SEmil Constantinescu ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr); 837fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 83856dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 83956dcabbaSDebojyoti Ghosh ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr); 84056dcabbaSDebojyoti Ghosh } else { 8418a381b04SJed Brown /* Initial guess taken from last stage */ 8428a381b04SJed Brown ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr); 84356dcabbaSDebojyoti Ghosh } 84496400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 845baa10174SEmil Constantinescu ierr = SNESSolve(snes,NULL,Y[i]);CHKERRQ(ierr); 8468a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 8478a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 8485ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 849552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 85096400bd6SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);CHKERRQ(ierr); 85196400bd6SLisandro Dalcin if (!stageok) { 8521be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8531be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 85496400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8551be93e3eSJed Brown goto reject_step; 8561be93e3eSJed Brown } 8578a381b04SJed Brown } 858e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 859e817cc15SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8606c4ed002SBarry Smith if (!tab->stiffly_accurate) SETERRQ1(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); 861df5e1e3dSEmil Constantinescu ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr); /* YdotI = YdotI(tn-1) */ 862e817cc15SEmil Constantinescu } else { 863df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 864e817cc15SEmil Constantinescu } 865e817cc15SEmil Constantinescu } else { 8665eca1a21SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8678a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 868df5e1e3dSEmil Constantinescu ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);/* YdotI = -G(t,Y,0) */ 869e817cc15SEmil Constantinescu ierr = VecScale(YdotI[i],-1.0);CHKERRQ(ierr); 8705eca1a21SEmil Constantinescu } else { 871df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 8725eca1a21SEmil Constantinescu } 8734cc180ffSJed Brown if (ark->imex) { 8748a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 8754cc180ffSJed Brown } else { 8764cc180ffSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 8774cc180ffSJed Brown } 8788a381b04SJed Brown } 87996400bd6SLisandro Dalcin ierr = TSPostStage(ts,ark->stage_time,i,Y);CHKERRQ(ierr); 880e817cc15SEmil Constantinescu } 88196400bd6SLisandro Dalcin 882be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 883fecfb714SLisandro Dalcin ierr = TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 884108c343cSJed Brown ark->status = TS_STEP_PENDING; 885552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 886108c343cSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 887fecfb714SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 888fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 88996400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 89096400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 89196400bd6SLisandro Dalcin ierr = TSRollBack_ARKIMEX(ts);CHKERRQ(ierr); 892be5899b3SLisandro Dalcin ts->time_step = next_time_step; 893be5899b3SLisandro Dalcin goto reject_step; 89496400bd6SLisandro Dalcin } 89596400bd6SLisandro Dalcin 8968a381b04SJed Brown ts->ptime += ts->time_step; 897cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 898108c343cSJed Brown break; 89996400bd6SLisandro Dalcin 90096400bd6SLisandro Dalcin reject_step: 901fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 902be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 90396400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 904be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 905108c343cSJed Brown } 906f85781f1SEmil Constantinescu } 9078a381b04SJed Brown PetscFunctionReturn(0); 9088a381b04SJed Brown } 9098a381b04SJed Brown 910cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 911cd652676SJed Brown { 912cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9134f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 914108c343cSJed Brown PetscReal h; 915108c343cSJed Brown PetscReal tt,t; 916cd652676SJed Brown PetscScalar *bt,*b; 917cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 918cd652676SJed Brown PetscErrorCode ierr; 919cd652676SJed Brown 920cd652676SJed Brown PetscFunctionBegin; 921ce94432eSBarry Smith if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 922108c343cSJed Brown switch (ark->status) { 923108c343cSJed Brown case TS_STEP_INCOMPLETE: 924108c343cSJed Brown case TS_STEP_PENDING: 925108c343cSJed Brown h = ts->time_step; 926108c343cSJed Brown t = (itime - ts->ptime)/h; 927108c343cSJed Brown break; 928108c343cSJed Brown case TS_STEP_COMPLETE: 929be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 930108c343cSJed Brown t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 931108c343cSJed Brown break; 932ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 933108c343cSJed Brown } 934dcca6d9dSJed Brown ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 935cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 9364f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 937cd652676SJed Brown for (i=0; i<s; i++) { 938c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 939108c343cSJed Brown b[i] += h * B[i*pinterp+j] * tt; 940cd652676SJed Brown } 941cd652676SJed Brown } 942cd652676SJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 943cd652676SJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 944cd652676SJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 945cd652676SJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 946cd652676SJed Brown PetscFunctionReturn(0); 947cd652676SJed Brown } 948cd652676SJed Brown 94956dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 95056dcabbaSDebojyoti Ghosh { 95156dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 95256dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 953be5899b3SLisandro Dalcin PetscReal h,h_prev,t,tt; 95456dcabbaSDebojyoti Ghosh PetscScalar *bt,*b; 95556dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 95656dcabbaSDebojyoti Ghosh PetscErrorCode ierr; 95756dcabbaSDebojyoti Ghosh 95856dcabbaSDebojyoti Ghosh PetscFunctionBegin; 95956dcabbaSDebojyoti Ghosh if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 960be5899b3SLisandro Dalcin ierr = PetscCalloc2(s,&bt,s,&b);CHKERRQ(ierr); 96181d12688SDebojyoti Ghosh h = ts->time_step; 962be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 963be5899b3SLisandro Dalcin t = 1 + h/h_prev*c; 96456dcabbaSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 96556dcabbaSDebojyoti Ghosh for (i=0; i<s; i++) { 96681d12688SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 96756dcabbaSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 96856dcabbaSDebojyoti Ghosh } 96956dcabbaSDebojyoti Ghosh } 97096400bd6SLisandro Dalcin if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 97156dcabbaSDebojyoti Ghosh ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr); 97256dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr); 97356dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr); 97456dcabbaSDebojyoti Ghosh ierr = PetscFree2(bt,b);CHKERRQ(ierr); 97556dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 97656dcabbaSDebojyoti Ghosh } 97756dcabbaSDebojyoti Ghosh 9788a381b04SJed Brown /*------------------------------------------------------------*/ 97996400bd6SLisandro Dalcin 98096400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauReset(TS ts) 98196400bd6SLisandro Dalcin { 98296400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 98396400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 98496400bd6SLisandro Dalcin PetscErrorCode ierr; 98596400bd6SLisandro Dalcin 98696400bd6SLisandro Dalcin PetscFunctionBegin; 98796400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 98896400bd6SLisandro Dalcin ierr = PetscFree(ark->work);CHKERRQ(ierr); 98996400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y);CHKERRQ(ierr); 99096400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI);CHKERRQ(ierr); 99196400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS);CHKERRQ(ierr); 99296400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y_prev);CHKERRQ(ierr); 99396400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 99496400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 99596400bd6SLisandro Dalcin PetscFunctionReturn(0); 99696400bd6SLisandro Dalcin } 99796400bd6SLisandro Dalcin 9988a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 9998a381b04SJed Brown { 10008a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 10018a381b04SJed Brown PetscErrorCode ierr; 10028a381b04SJed Brown 10038a381b04SJed Brown PetscFunctionBegin; 100496400bd6SLisandro Dalcin ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr); 10058a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 1006e817cc15SEmil Constantinescu ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr); 10078a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 10088a381b04SJed Brown PetscFunctionReturn(0); 10098a381b04SJed Brown } 10108a381b04SJed Brown 1011d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1012d5e6173cSPeter Brune { 1013d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 1014d5e6173cSPeter Brune PetscErrorCode ierr; 1015d5e6173cSPeter Brune 1016d5e6173cSPeter Brune PetscFunctionBegin; 1017d5e6173cSPeter Brune if (Z) { 1018d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1019d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1020d5e6173cSPeter Brune } else *Z = ax->Z; 1021d5e6173cSPeter Brune } 1022d5e6173cSPeter Brune if (Ydot) { 1023d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1024d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1025d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 1026d5e6173cSPeter Brune } 1027d5e6173cSPeter Brune PetscFunctionReturn(0); 1028d5e6173cSPeter Brune } 1029d5e6173cSPeter Brune 1030d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1031d5e6173cSPeter Brune { 1032d5e6173cSPeter Brune PetscErrorCode ierr; 1033d5e6173cSPeter Brune 1034d5e6173cSPeter Brune PetscFunctionBegin; 1035d5e6173cSPeter Brune if (Z) { 1036d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1037d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1038d5e6173cSPeter Brune } 1039d5e6173cSPeter Brune } 1040d5e6173cSPeter Brune if (Ydot) { 1041d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1042d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1043d5e6173cSPeter Brune } 1044d5e6173cSPeter Brune } 1045d5e6173cSPeter Brune PetscFunctionReturn(0); 1046d5e6173cSPeter Brune } 1047d5e6173cSPeter Brune 10488a381b04SJed Brown /* 10498a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 10508a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 10518a381b04SJed Brown */ 10528a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 10538a381b04SJed Brown { 10548a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1055d5e6173cSPeter Brune DM dm,dmsave; 1056d5e6173cSPeter Brune Vec Z,Ydot; 1057b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10588a381b04SJed Brown PetscErrorCode ierr; 10598a381b04SJed Brown 10608a381b04SJed Brown PetscFunctionBegin; 1061d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1062d5e6173cSPeter Brune ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1063b296d7d5SJed Brown ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 1064d5e6173cSPeter Brune dmsave = ts->dm; 1065d5e6173cSPeter Brune ts->dm = dm; 1066740132f1SEmil Constantinescu 1067d5e6173cSPeter Brune ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr); 1068e817cc15SEmil Constantinescu 1069d5e6173cSPeter Brune ts->dm = dmsave; 1070d5e6173cSPeter Brune ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 10718a381b04SJed Brown PetscFunctionReturn(0); 10728a381b04SJed Brown } 10738a381b04SJed Brown 1074d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 10758a381b04SJed Brown { 10768a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1077d5e6173cSPeter Brune DM dm,dmsave; 1078d5e6173cSPeter Brune Vec Ydot; 1079b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10808a381b04SJed Brown PetscErrorCode ierr; 10818a381b04SJed Brown 10828a381b04SJed Brown PetscFunctionBegin; 1083d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 10840298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 10858a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1086d5e6173cSPeter Brune dmsave = ts->dm; 1087d5e6173cSPeter Brune ts->dm = dm; 1088740132f1SEmil Constantinescu 1089d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr); 1090740132f1SEmil Constantinescu 1091d5e6173cSPeter Brune ts->dm = dmsave; 10920298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1093d5e6173cSPeter Brune PetscFunctionReturn(0); 1094d5e6173cSPeter Brune } 1095d5e6173cSPeter Brune 1096d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1097d5e6173cSPeter Brune { 1098d5e6173cSPeter Brune PetscFunctionBegin; 1099d5e6173cSPeter Brune PetscFunctionReturn(0); 1100d5e6173cSPeter Brune } 1101d5e6173cSPeter Brune 1102d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1103d5e6173cSPeter Brune { 1104d5e6173cSPeter Brune TS ts = (TS)ctx; 1105d5e6173cSPeter Brune PetscErrorCode ierr; 1106d5e6173cSPeter Brune Vec Z,Z_c; 1107d5e6173cSPeter Brune 1108d5e6173cSPeter Brune PetscFunctionBegin; 11090298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 11100298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1111d5e6173cSPeter Brune ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 1112d5e6173cSPeter Brune ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 11130298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 11140298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 11158a381b04SJed Brown PetscFunctionReturn(0); 11168a381b04SJed Brown } 11178a381b04SJed Brown 1118cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1119cdb298fcSPeter Brune { 1120cdb298fcSPeter Brune PetscFunctionBegin; 1121cdb298fcSPeter Brune PetscFunctionReturn(0); 1122cdb298fcSPeter Brune } 1123cdb298fcSPeter Brune 1124cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1125cdb298fcSPeter Brune { 1126cdb298fcSPeter Brune TS ts = (TS)ctx; 1127cdb298fcSPeter Brune PetscErrorCode ierr; 1128cdb298fcSPeter Brune Vec Z,Z_c; 1129cdb298fcSPeter Brune 1130cdb298fcSPeter Brune PetscFunctionBegin; 11310298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11320298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1133cdb298fcSPeter Brune 1134cdb298fcSPeter Brune ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1135cdb298fcSPeter Brune ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1136cdb298fcSPeter Brune 11370298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11380298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1139cdb298fcSPeter Brune PetscFunctionReturn(0); 1140cdb298fcSPeter Brune } 1141cdb298fcSPeter Brune 114296400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 114396400bd6SLisandro Dalcin { 114496400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 114596400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 114696400bd6SLisandro Dalcin PetscErrorCode ierr; 114796400bd6SLisandro Dalcin 114896400bd6SLisandro Dalcin PetscFunctionBegin; 114996400bd6SLisandro Dalcin ierr = PetscMalloc1(tab->s,&ark->work);CHKERRQ(ierr); 115096400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);CHKERRQ(ierr); 115196400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);CHKERRQ(ierr); 115296400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);CHKERRQ(ierr); 115396400bd6SLisandro Dalcin if (ark->extrapolate) { 115496400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 115596400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 115696400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 115796400bd6SLisandro Dalcin } 115896400bd6SLisandro Dalcin PetscFunctionReturn(0); 115996400bd6SLisandro Dalcin } 116096400bd6SLisandro Dalcin 11618a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 11628a381b04SJed Brown { 11638a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11648a381b04SJed Brown PetscErrorCode ierr; 1165d5e6173cSPeter Brune DM dm; 116696400bd6SLisandro Dalcin SNES snes; 1167f9c1d6abSBarry Smith 11688a381b04SJed Brown PetscFunctionBegin; 116996400bd6SLisandro Dalcin ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr); 11708a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 1171e817cc15SEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr); 11728a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 1173d5e6173cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1174d5e6173cSPeter Brune ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1175cdb298fcSPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 117696400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 11778a381b04SJed Brown PetscFunctionReturn(0); 11788a381b04SJed Brown } 11798a381b04SJed Brown /*------------------------------------------------------------*/ 11808a381b04SJed Brown 11814416b707SBarry Smith static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 11828a381b04SJed Brown { 11834cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11848a381b04SJed Brown PetscErrorCode ierr; 11858a381b04SJed Brown 11868a381b04SJed Brown PetscFunctionBegin; 1187e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");CHKERRQ(ierr); 11888a381b04SJed Brown { 11898a381b04SJed Brown ARKTableauLink link; 11908a381b04SJed Brown PetscInt count,choice; 11918a381b04SJed Brown PetscBool flg; 11928a381b04SJed Brown const char **namelist; 11938a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 1194f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 11958a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 119696400bd6SLisandro Dalcin ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);CHKERRQ(ierr); 119796400bd6SLisandro Dalcin if (flg) {ierr = TSARKIMEXSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11988a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 119996400bd6SLisandro Dalcin 12004cc180ffSJed Brown flg = (PetscBool) !ark->imex; 12010298fd71SBarry Smith ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr); 12024cc180ffSJed Brown ark->imex = (PetscBool) !flg; 120303842d09SLisandro Dalcin ierr = 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);CHKERRQ(ierr); 12048a381b04SJed Brown } 12058a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 12068a381b04SJed Brown PetscFunctionReturn(0); 12078a381b04SJed Brown } 12088a381b04SJed Brown 12098a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 12108a381b04SJed Brown { 12118a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 12128a381b04SJed Brown PetscBool iascii; 12138a381b04SJed Brown PetscErrorCode ierr; 12148a381b04SJed Brown 12158a381b04SJed Brown PetscFunctionBegin; 1216251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12178a381b04SJed Brown if (iascii) { 12189c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 121919fd82e9SBarry Smith TSARKIMEXType arktype; 12208a381b04SJed Brown char buf[512]; 12213a28c0e5SStefano Zampini PetscBool flg; 12223a28c0e5SStefano Zampini 12238a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 12243a28c0e5SStefano Zampini ierr = TSARKIMEXGetFullyImplicit(ts,&flg);CHKERRQ(ierr); 12258a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 12268caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 122731f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 12288caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 12293a28c0e5SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Fully implicit: %s\n",flg ? "yes" : "no");CHKERRQ(ierr); 1230e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 1231e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr); 1232e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr); 123331f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 12348a381b04SJed Brown } 12358a381b04SJed Brown PetscFunctionReturn(0); 12368a381b04SJed Brown } 12378a381b04SJed Brown 1238f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1239f2c2a1b9SBarry Smith { 1240f2c2a1b9SBarry Smith PetscErrorCode ierr; 1241f2c2a1b9SBarry Smith SNES snes; 12429c334d8fSLisandro Dalcin TSAdapt adapt; 1243f2c2a1b9SBarry Smith 1244f2c2a1b9SBarry Smith PetscFunctionBegin; 12459c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12469c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1247f2c2a1b9SBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1248f2c2a1b9SBarry Smith ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 1249ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 12500298fd71SBarry Smith ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 12510298fd71SBarry Smith ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 1252f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1253f2c2a1b9SBarry Smith } 1254f2c2a1b9SBarry Smith 12558a381b04SJed Brown /*@C 12568a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 12578a381b04SJed Brown 12588a381b04SJed Brown Logically collective 12598a381b04SJed Brown 1260*d8d19677SJose E. Roman Input Parameters: 12618a381b04SJed Brown + ts - timestepping context 12628a381b04SJed Brown - arktype - type of ARK-IMEX scheme 12638a381b04SJed Brown 12649bd3e852SBarry Smith Options Database: 12659bd3e852SBarry Smith . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> 12669bd3e852SBarry Smith 12678a381b04SJed Brown Level: intermediate 12688a381b04SJed Brown 12699bd3e852SBarry Smith .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, 12709bd3e852SBarry Smith TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 12718a381b04SJed Brown @*/ 127219fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 12738a381b04SJed Brown { 12748a381b04SJed Brown PetscErrorCode ierr; 12758a381b04SJed Brown 12768a381b04SJed Brown PetscFunctionBegin; 12778a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1278b92453a8SLisandro Dalcin PetscValidCharPointer(arktype,2); 127919fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 12808a381b04SJed Brown PetscFunctionReturn(0); 12818a381b04SJed Brown } 12828a381b04SJed Brown 12838a381b04SJed Brown /*@C 12848a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 12858a381b04SJed Brown 12868a381b04SJed Brown Logically collective 12878a381b04SJed Brown 12888a381b04SJed Brown Input Parameter: 12898a381b04SJed Brown . ts - timestepping context 12908a381b04SJed Brown 12918a381b04SJed Brown Output Parameter: 12928a381b04SJed Brown . arktype - type of ARK-IMEX scheme 12938a381b04SJed Brown 12948a381b04SJed Brown Level: intermediate 12958a381b04SJed Brown 12968a381b04SJed Brown .seealso: TSARKIMEXGetType() 12978a381b04SJed Brown @*/ 129819fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 12998a381b04SJed Brown { 13008a381b04SJed Brown PetscErrorCode ierr; 13018a381b04SJed Brown 13028a381b04SJed Brown PetscFunctionBegin; 13038a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 130419fd82e9SBarry Smith ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 13058a381b04SJed Brown PetscFunctionReturn(0); 13068a381b04SJed Brown } 13078a381b04SJed Brown 130816353aafSBarry Smith /*@ 13094cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 13104cc180ffSJed Brown 13114cc180ffSJed Brown Logically collective 13124cc180ffSJed Brown 1313*d8d19677SJose E. Roman Input Parameters: 13144cc180ffSJed Brown + ts - timestepping context 13154cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 13164cc180ffSJed Brown 13174cc180ffSJed Brown Level: intermediate 13184cc180ffSJed Brown 13193a28c0e5SStefano Zampini .seealso: TSARKIMEXGetType(), TSARKIMEXGetFullyImplicit() 13204cc180ffSJed Brown @*/ 13214cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 13224cc180ffSJed Brown { 13234cc180ffSJed Brown PetscErrorCode ierr; 13244cc180ffSJed Brown 13254cc180ffSJed Brown PetscFunctionBegin; 13264cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13273a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts,flg,2); 13284cc180ffSJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 13294cc180ffSJed Brown PetscFunctionReturn(0); 13304cc180ffSJed Brown } 13314cc180ffSJed Brown 13323a28c0e5SStefano Zampini /*@ 13333a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 13343a28c0e5SStefano Zampini 13353a28c0e5SStefano Zampini Logically collective 13363a28c0e5SStefano Zampini 13373a28c0e5SStefano Zampini Input Parameter: 13383a28c0e5SStefano Zampini . ts - timestepping context 13393a28c0e5SStefano Zampini 13407a7aea1fSJed Brown Output Parameter: 13413a28c0e5SStefano Zampini . flg - PETSC_TRUE for fully implicit 13423a28c0e5SStefano Zampini 13433a28c0e5SStefano Zampini Level: intermediate 13443a28c0e5SStefano Zampini 13453a28c0e5SStefano Zampini .seealso: TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit() 13463a28c0e5SStefano Zampini @*/ 13473a28c0e5SStefano Zampini PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts,PetscBool *flg) 13483a28c0e5SStefano Zampini { 13493a28c0e5SStefano Zampini PetscErrorCode ierr; 13503a28c0e5SStefano Zampini 13513a28c0e5SStefano Zampini PetscFunctionBegin; 13523a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13533a28c0e5SStefano Zampini PetscValidPointer(flg,2); 13543a28c0e5SStefano Zampini ierr = PetscUseMethod(ts,"TSARKIMEXGetFullyImplicit_C",(TS,PetscBool*),(ts,flg));CHKERRQ(ierr); 13553a28c0e5SStefano Zampini PetscFunctionReturn(0); 13563a28c0e5SStefano Zampini } 13573a28c0e5SStefano Zampini 1358e0877f53SBarry Smith static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 13598a381b04SJed Brown { 13608a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13618a381b04SJed Brown 13628a381b04SJed Brown PetscFunctionBegin; 13638a381b04SJed Brown *arktype = ark->tableau->name; 13648a381b04SJed Brown PetscFunctionReturn(0); 13658a381b04SJed Brown } 1366e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype) 13678a381b04SJed Brown { 13688a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13698a381b04SJed Brown PetscErrorCode ierr; 13708a381b04SJed Brown PetscBool match; 13718a381b04SJed Brown ARKTableauLink link; 13728a381b04SJed Brown 13738a381b04SJed Brown PetscFunctionBegin; 13748a381b04SJed Brown if (ark->tableau) { 13758a381b04SJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 13768a381b04SJed Brown if (match) PetscFunctionReturn(0); 13778a381b04SJed Brown } 13788a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 13798a381b04SJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 13808a381b04SJed Brown if (match) { 138196400bd6SLisandro Dalcin if (ts->setupcalled) {ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr);} 13828a381b04SJed Brown ark->tableau = &link->tab; 138396400bd6SLisandro Dalcin if (ts->setupcalled) {ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr);} 13842ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13858a381b04SJed Brown PetscFunctionReturn(0); 13868a381b04SJed Brown } 13878a381b04SJed Brown } 1388ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 13898a381b04SJed Brown } 1390e0877f53SBarry Smith 1391e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 13924cc180ffSJed Brown { 13934cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13944cc180ffSJed Brown 13954cc180ffSJed Brown PetscFunctionBegin; 13964cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13974cc180ffSJed Brown PetscFunctionReturn(0); 13984cc180ffSJed Brown } 13998a381b04SJed Brown 14003a28c0e5SStefano Zampini static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts,PetscBool *flg) 14013a28c0e5SStefano Zampini { 14023a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 14033a28c0e5SStefano Zampini 14043a28c0e5SStefano Zampini PetscFunctionBegin; 14053a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 14063a28c0e5SStefano Zampini PetscFunctionReturn(0); 14073a28c0e5SStefano Zampini } 14083a28c0e5SStefano Zampini 1409b3a6b972SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1410b3a6b972SJed Brown { 1411b3a6b972SJed Brown PetscErrorCode ierr; 1412b3a6b972SJed Brown 1413b3a6b972SJed Brown PetscFunctionBegin; 1414b3a6b972SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 1415b3a6b972SJed Brown if (ts->dm) { 1416b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1417b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1418b3a6b972SJed Brown } 1419b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1420b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr); 1421b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr); 1422b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 14233a28c0e5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 1424b3a6b972SJed Brown PetscFunctionReturn(0); 1425b3a6b972SJed Brown } 1426b3a6b972SJed Brown 14278a381b04SJed Brown /* ------------------------------------------------------------ */ 14288a381b04SJed Brown /*MC 1429c79dcfbdSBarry Smith TSARKIMEX - ODE and DAE solver using additive Runge-Kutta IMEX schemes 14308a381b04SJed Brown 1431fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1432fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1433fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1434fca742c7SJed Brown 1435fca742c7SJed Brown Notes: 1436a4386c9eSJed Brown The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1437c8058688SBarry Smith 14385eca1a21SEmil Constantinescu If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 14395eca1a21SEmil Constantinescu 1440a4386c9eSJed 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). 1441fca742c7SJed Brown 1442d0685a90SJed Brown Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1443d0685a90SJed Brown 14448a381b04SJed Brown Level: beginner 14458a381b04SJed Brown 14463a28c0e5SStefano Zampini .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEXGetFullyImplicit(), 14473a28c0e5SStefano Zampini TSARKIMEX1BEE, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122, 1448d0685a90SJed Brown TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister() 14498a381b04SJed Brown 14508a381b04SJed Brown M*/ 14518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 14528a381b04SJed Brown { 14538a381b04SJed Brown TS_ARKIMEX *th; 14548a381b04SJed Brown PetscErrorCode ierr; 14558a381b04SJed Brown 14568a381b04SJed Brown PetscFunctionBegin; 1457607a6623SBarry Smith ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 14588a381b04SJed Brown 14598a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 14608a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 14618a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1462f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 14638a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 14648a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1465cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1466108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 146724655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 14688a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 14698a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 14708a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 14718a381b04SJed Brown 1472efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1473efd4aadfSBarry Smith 1474b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 14758a381b04SJed Brown ts->data = (void*)th; 14764cc180ffSJed Brown th->imex = PETSC_TRUE; 14778a381b04SJed Brown 1478bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 1479bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 1480bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 14813a28c0e5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetFullyImplicit_C",TSARKIMEXGetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 148296400bd6SLisandro Dalcin 148396400bd6SLisandro Dalcin ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 14848a381b04SJed Brown PetscFunctionReturn(0); 14858a381b04SJed Brown } 1486