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 103e817cc15SEmil Constantinescu 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 68324655328SShri static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 68424655328SShri { 68524655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 68624655328SShri ARKTableau tab = ark->tableau; 68724655328SShri const PetscInt s = tab->s; 68824655328SShri const PetscReal *bt = tab->bt,*b = tab->b; 68924655328SShri PetscScalar *w = ark->work; 69024655328SShri Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 69124655328SShri PetscInt j; 692be5899b3SLisandro Dalcin PetscReal h; 69324655328SShri PetscErrorCode ierr; 69424655328SShri 69524655328SShri PetscFunctionBegin; 696be5899b3SLisandro Dalcin switch (ark->status) { 697be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 698be5899b3SLisandro Dalcin case TS_STEP_PENDING: 699be5899b3SLisandro Dalcin h = ts->time_step; break; 700be5899b3SLisandro Dalcin case TS_STEP_COMPLETE: 701be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 702be5899b3SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 703be5899b3SLisandro Dalcin } 70424655328SShri for (j=0; j<s; j++) w[j] = -h*bt[j]; 70524655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 70624655328SShri for (j=0; j<s; j++) w[j] = -h*b[j]; 70724655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 70824655328SShri PetscFunctionReturn(0); 70924655328SShri } 71024655328SShri 7118a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 7128a381b04SJed Brown { 7138a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7148a381b04SJed Brown ARKTableau tab = ark->tableau; 7158a381b04SJed Brown const PetscInt s = tab->s; 71624655328SShri const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 717406d0ec2SJed Brown PetscScalar *w = ark->work; 7181297b224SEmil Constantinescu Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z; 71996400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 720108c343cSJed Brown TSAdapt adapt; 7218a381b04SJed Brown SNES snes; 722fecfb714SLisandro Dalcin PetscInt i,j,its,lits; 723be5899b3SLisandro Dalcin PetscInt rejections = 0; 72496400bd6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 72596400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7268a381b04SJed Brown PetscErrorCode ierr; 7278a381b04SJed Brown 7288a381b04SJed Brown PetscFunctionBegin; 72996400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 73096400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 73196400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 73296400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 73396400bd6SLisandro Dalcin } 73496400bd6SLisandro Dalcin 73596400bd6SLisandro Dalcin if (!ts->steprollback) { 73696400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 73796400bd6SLisandro Dalcin ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr); 73896400bd6SLisandro Dalcin } 739fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 74096400bd6SLisandro Dalcin for (i = 0; i<s; i++) { 74196400bd6SLisandro Dalcin ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr); 74296400bd6SLisandro Dalcin ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr); 74396400bd6SLisandro Dalcin ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr); 74496400bd6SLisandro Dalcin } 74596400bd6SLisandro Dalcin } 74696400bd6SLisandro Dalcin } 74796400bd6SLisandro Dalcin 748fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 74996400bd6SLisandro Dalcin TS ts_start; 750baa10174SEmil Constantinescu ierr = TSClone(ts,&ts_start);CHKERRQ(ierr); 751e817cc15SEmil Constantinescu ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr); 752e817cc15SEmil Constantinescu ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr); 75319eac22cSLisandro Dalcin ierr = TSSetMaxSteps(ts_start,ts->steps+1);CHKERRQ(ierr); 75419eac22cSLisandro Dalcin ierr = TSSetMaxTime(ts_start,ts->ptime+ts->time_step);CHKERRQ(ierr); 755feed9e9dSBarry Smith ierr = TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 756740132f1SEmil Constantinescu ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr); 757e817cc15SEmil Constantinescu ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr); 758740132f1SEmil Constantinescu ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr); 759e817cc15SEmil Constantinescu ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr); 76034561852SEmil Constantinescu 761dcb233daSLisandro Dalcin ierr = TSRestartStep(ts_start);CHKERRQ(ierr); 762e817cc15SEmil Constantinescu ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr); 763e817cc15SEmil Constantinescu ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr); 76496400bd6SLisandro Dalcin ierr = TSGetTimeStep(ts_start,&ts->time_step);CHKERRQ(ierr); 765bbd56ea5SKarl Rupp 76685fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 76785fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data; 76885fc7851SLisandro Dalcin ierr = VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);CHKERRQ(ierr); 76985fc7851SLisandro Dalcin } 77096400bd6SLisandro Dalcin ts->steps++; 77134561852SEmil Constantinescu 772d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 773d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 77496400bd6SLisandro Dalcin { 77596400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 77696400bd6SLisandro Dalcin ierr = TSSetSNES(ts,snes);CHKERRQ(ierr); 77796400bd6SLisandro Dalcin } 778166a6834SEmil Constantinescu ierr = TSDestroy(&ts_start);CHKERRQ(ierr); 779e817cc15SEmil Constantinescu } 780e817cc15SEmil Constantinescu 781108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 78296400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 78396400bd6SLisandro Dalcin PetscReal t = ts->ptime; 784108c343cSJed Brown PetscReal h = ts->time_step; 7858a381b04SJed Brown for (i=0; i<s; i++) { 7869be3e283SDebojyoti Ghosh ark->stage_time = t + h*ct[i]; 78796400bd6SLisandro Dalcin ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr); 7888a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 7896c4ed002SBarry 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"); 7908a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 791e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 7928a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 7938a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 7948a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 7958a381b04SJed Brown } else { 796b296d7d5SJed Brown ark->scoeff = 1./At[i*s+i]; 7978a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 7988a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 799e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 8004f385281SJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 801c58d1302SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 802c58d1302SEmil Constantinescu ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr); 803fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 80456dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 80556dcabbaSDebojyoti Ghosh ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr); 80656dcabbaSDebojyoti Ghosh } else { 8078a381b04SJed Brown /* Initial guess taken from last stage */ 8088a381b04SJed Brown ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr); 80956dcabbaSDebojyoti Ghosh } 81096400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 811baa10174SEmil Constantinescu ierr = SNESSolve(snes,NULL,Y[i]);CHKERRQ(ierr); 8128a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 8138a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 8145ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 815552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 81696400bd6SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);CHKERRQ(ierr); 81796400bd6SLisandro Dalcin if (!stageok) { 8181be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8191be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 82096400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8211be93e3eSJed Brown goto reject_step; 8221be93e3eSJed Brown } 8238a381b04SJed Brown } 824e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 825e817cc15SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8266c4ed002SBarry 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); 827df5e1e3dSEmil Constantinescu ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr); /* YdotI = YdotI(tn-1) */ 828e817cc15SEmil Constantinescu } else { 829df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 830e817cc15SEmil Constantinescu } 831e817cc15SEmil Constantinescu } else { 8325eca1a21SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8338a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 834df5e1e3dSEmil Constantinescu ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);/* YdotI = -G(t,Y,0) */ 835e817cc15SEmil Constantinescu ierr = VecScale(YdotI[i],-1.0);CHKERRQ(ierr); 8365eca1a21SEmil Constantinescu } else { 837df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 8385eca1a21SEmil Constantinescu } 8394cc180ffSJed Brown if (ark->imex) { 8408a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 8414cc180ffSJed Brown } else { 8424cc180ffSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 8434cc180ffSJed Brown } 8448a381b04SJed Brown } 84596400bd6SLisandro Dalcin ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr); 846e817cc15SEmil Constantinescu } 84796400bd6SLisandro Dalcin 848be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 849fecfb714SLisandro Dalcin ierr = TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 850108c343cSJed Brown ark->status = TS_STEP_PENDING; 851552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 852108c343cSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 853fecfb714SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 854fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 85596400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 85696400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 85796400bd6SLisandro Dalcin ierr = TSRollBack_ARKIMEX(ts);CHKERRQ(ierr); 858be5899b3SLisandro Dalcin ts->time_step = next_time_step; 859be5899b3SLisandro Dalcin goto reject_step; 86096400bd6SLisandro Dalcin } 86196400bd6SLisandro Dalcin 8628a381b04SJed Brown ts->ptime += ts->time_step; 863cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 864108c343cSJed Brown break; 86596400bd6SLisandro Dalcin 86696400bd6SLisandro Dalcin reject_step: 867fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 868be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 86996400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 870be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 871108c343cSJed Brown } 872f85781f1SEmil Constantinescu } 8738a381b04SJed Brown PetscFunctionReturn(0); 8748a381b04SJed Brown } 8758a381b04SJed Brown 876cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 877cd652676SJed Brown { 878cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 8794f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 880108c343cSJed Brown PetscReal h; 881108c343cSJed Brown PetscReal tt,t; 882cd652676SJed Brown PetscScalar *bt,*b; 883cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 884cd652676SJed Brown PetscErrorCode ierr; 885cd652676SJed Brown 886cd652676SJed Brown PetscFunctionBegin; 887ce94432eSBarry Smith if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 888108c343cSJed Brown switch (ark->status) { 889108c343cSJed Brown case TS_STEP_INCOMPLETE: 890108c343cSJed Brown case TS_STEP_PENDING: 891108c343cSJed Brown h = ts->time_step; 892108c343cSJed Brown t = (itime - ts->ptime)/h; 893108c343cSJed Brown break; 894108c343cSJed Brown case TS_STEP_COMPLETE: 895be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 896108c343cSJed Brown t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 897108c343cSJed Brown break; 898ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 899108c343cSJed Brown } 900dcca6d9dSJed Brown ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 901cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 9024f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 903cd652676SJed Brown for (i=0; i<s; i++) { 904c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 905108c343cSJed Brown b[i] += h * B[i*pinterp+j] * tt; 906cd652676SJed Brown } 907cd652676SJed Brown } 908cd652676SJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 909cd652676SJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 910cd652676SJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 911cd652676SJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 912cd652676SJed Brown PetscFunctionReturn(0); 913cd652676SJed Brown } 914cd652676SJed Brown 91556dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 91656dcabbaSDebojyoti Ghosh { 91756dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 91856dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 919be5899b3SLisandro Dalcin PetscReal h,h_prev,t,tt; 92056dcabbaSDebojyoti Ghosh PetscScalar *bt,*b; 92156dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 92256dcabbaSDebojyoti Ghosh PetscErrorCode ierr; 92356dcabbaSDebojyoti Ghosh 92456dcabbaSDebojyoti Ghosh PetscFunctionBegin; 92556dcabbaSDebojyoti Ghosh if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 926be5899b3SLisandro Dalcin ierr = PetscCalloc2(s,&bt,s,&b);CHKERRQ(ierr); 92781d12688SDebojyoti Ghosh h = ts->time_step; 928be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 929be5899b3SLisandro Dalcin t = 1 + h/h_prev*c; 93056dcabbaSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 93156dcabbaSDebojyoti Ghosh for (i=0; i<s; i++) { 93281d12688SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 93356dcabbaSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 93456dcabbaSDebojyoti Ghosh } 93556dcabbaSDebojyoti Ghosh } 93696400bd6SLisandro Dalcin if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 93756dcabbaSDebojyoti Ghosh ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr); 93856dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr); 93956dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr); 94056dcabbaSDebojyoti Ghosh ierr = PetscFree2(bt,b);CHKERRQ(ierr); 94156dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 94256dcabbaSDebojyoti Ghosh } 94356dcabbaSDebojyoti Ghosh 9448a381b04SJed Brown /*------------------------------------------------------------*/ 94596400bd6SLisandro Dalcin 94696400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauReset(TS ts) 94796400bd6SLisandro Dalcin { 94896400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 94996400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 95096400bd6SLisandro Dalcin PetscErrorCode ierr; 95196400bd6SLisandro Dalcin 95296400bd6SLisandro Dalcin PetscFunctionBegin; 95396400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 95496400bd6SLisandro Dalcin ierr = PetscFree(ark->work);CHKERRQ(ierr); 95596400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y);CHKERRQ(ierr); 95696400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI);CHKERRQ(ierr); 95796400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS);CHKERRQ(ierr); 95896400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y_prev);CHKERRQ(ierr); 95996400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 96096400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 96196400bd6SLisandro Dalcin PetscFunctionReturn(0); 96296400bd6SLisandro Dalcin } 96396400bd6SLisandro Dalcin 9648a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 9658a381b04SJed Brown { 9668a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9678a381b04SJed Brown PetscErrorCode ierr; 9688a381b04SJed Brown 9698a381b04SJed Brown PetscFunctionBegin; 97096400bd6SLisandro Dalcin ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr); 9718a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 972e817cc15SEmil Constantinescu ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr); 9738a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 9748a381b04SJed Brown PetscFunctionReturn(0); 9758a381b04SJed Brown } 9768a381b04SJed Brown 977d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 978d5e6173cSPeter Brune { 979d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 980d5e6173cSPeter Brune PetscErrorCode ierr; 981d5e6173cSPeter Brune 982d5e6173cSPeter Brune PetscFunctionBegin; 983d5e6173cSPeter Brune if (Z) { 984d5e6173cSPeter Brune if (dm && dm != ts->dm) { 985d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 986d5e6173cSPeter Brune } else *Z = ax->Z; 987d5e6173cSPeter Brune } 988d5e6173cSPeter Brune if (Ydot) { 989d5e6173cSPeter Brune if (dm && dm != ts->dm) { 990d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 991d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 992d5e6173cSPeter Brune } 993d5e6173cSPeter Brune PetscFunctionReturn(0); 994d5e6173cSPeter Brune } 995d5e6173cSPeter Brune 996d5e6173cSPeter Brune 997d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 998d5e6173cSPeter Brune { 999d5e6173cSPeter Brune PetscErrorCode ierr; 1000d5e6173cSPeter Brune 1001d5e6173cSPeter Brune PetscFunctionBegin; 1002d5e6173cSPeter Brune if (Z) { 1003d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1004d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1005d5e6173cSPeter Brune } 1006d5e6173cSPeter Brune } 1007d5e6173cSPeter Brune if (Ydot) { 1008d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1009d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1010d5e6173cSPeter Brune } 1011d5e6173cSPeter Brune } 1012d5e6173cSPeter Brune PetscFunctionReturn(0); 1013d5e6173cSPeter Brune } 1014d5e6173cSPeter Brune 10158a381b04SJed Brown /* 10168a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 10178a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 10188a381b04SJed Brown */ 10198a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 10208a381b04SJed Brown { 10218a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1022d5e6173cSPeter Brune DM dm,dmsave; 1023d5e6173cSPeter Brune Vec Z,Ydot; 1024b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10258a381b04SJed Brown PetscErrorCode ierr; 10268a381b04SJed Brown 10278a381b04SJed Brown PetscFunctionBegin; 1028d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1029d5e6173cSPeter Brune ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1030b296d7d5SJed Brown ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 1031d5e6173cSPeter Brune dmsave = ts->dm; 1032d5e6173cSPeter Brune ts->dm = dm; 1033740132f1SEmil Constantinescu 1034d5e6173cSPeter Brune ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr); 1035e817cc15SEmil Constantinescu 1036d5e6173cSPeter Brune ts->dm = dmsave; 1037d5e6173cSPeter Brune ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 10388a381b04SJed Brown PetscFunctionReturn(0); 10398a381b04SJed Brown } 10408a381b04SJed Brown 1041d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 10428a381b04SJed Brown { 10438a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1044d5e6173cSPeter Brune DM dm,dmsave; 1045d5e6173cSPeter Brune Vec Ydot; 1046b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10478a381b04SJed Brown PetscErrorCode ierr; 10488a381b04SJed Brown 10498a381b04SJed Brown PetscFunctionBegin; 1050d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 10510298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 10528a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1053d5e6173cSPeter Brune dmsave = ts->dm; 1054d5e6173cSPeter Brune ts->dm = dm; 1055740132f1SEmil Constantinescu 1056d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr); 1057740132f1SEmil Constantinescu 1058d5e6173cSPeter Brune ts->dm = dmsave; 10590298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1060d5e6173cSPeter Brune PetscFunctionReturn(0); 1061d5e6173cSPeter Brune } 1062d5e6173cSPeter Brune 1063d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1064d5e6173cSPeter Brune { 1065d5e6173cSPeter Brune PetscFunctionBegin; 1066d5e6173cSPeter Brune PetscFunctionReturn(0); 1067d5e6173cSPeter Brune } 1068d5e6173cSPeter Brune 1069d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1070d5e6173cSPeter Brune { 1071d5e6173cSPeter Brune TS ts = (TS)ctx; 1072d5e6173cSPeter Brune PetscErrorCode ierr; 1073d5e6173cSPeter Brune Vec Z,Z_c; 1074d5e6173cSPeter Brune 1075d5e6173cSPeter Brune PetscFunctionBegin; 10760298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 10770298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1078d5e6173cSPeter Brune ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 1079d5e6173cSPeter Brune ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 10800298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 10810298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 10828a381b04SJed Brown PetscFunctionReturn(0); 10838a381b04SJed Brown } 10848a381b04SJed Brown 1085cdb298fcSPeter Brune 1086cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1087cdb298fcSPeter Brune { 1088cdb298fcSPeter Brune PetscFunctionBegin; 1089cdb298fcSPeter Brune PetscFunctionReturn(0); 1090cdb298fcSPeter Brune } 1091cdb298fcSPeter Brune 1092cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1093cdb298fcSPeter Brune { 1094cdb298fcSPeter Brune TS ts = (TS)ctx; 1095cdb298fcSPeter Brune PetscErrorCode ierr; 1096cdb298fcSPeter Brune Vec Z,Z_c; 1097cdb298fcSPeter Brune 1098cdb298fcSPeter Brune PetscFunctionBegin; 10990298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11000298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1101cdb298fcSPeter Brune 1102cdb298fcSPeter Brune ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1103cdb298fcSPeter Brune ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1104cdb298fcSPeter Brune 11050298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11060298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1107cdb298fcSPeter Brune PetscFunctionReturn(0); 1108cdb298fcSPeter Brune } 1109cdb298fcSPeter Brune 111096400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 111196400bd6SLisandro Dalcin { 111296400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 111396400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 111496400bd6SLisandro Dalcin PetscErrorCode ierr; 111596400bd6SLisandro Dalcin 111696400bd6SLisandro Dalcin PetscFunctionBegin; 111796400bd6SLisandro Dalcin ierr = PetscMalloc1(tab->s,&ark->work);CHKERRQ(ierr); 111896400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);CHKERRQ(ierr); 111996400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);CHKERRQ(ierr); 112096400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);CHKERRQ(ierr); 112196400bd6SLisandro Dalcin if (ark->extrapolate) { 112296400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 112396400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 112496400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 112596400bd6SLisandro Dalcin } 112696400bd6SLisandro Dalcin PetscFunctionReturn(0); 112796400bd6SLisandro Dalcin } 112896400bd6SLisandro Dalcin 11298a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 11308a381b04SJed Brown { 11318a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11328a381b04SJed Brown PetscErrorCode ierr; 1133d5e6173cSPeter Brune DM dm; 113496400bd6SLisandro Dalcin SNES snes; 1135f9c1d6abSBarry Smith 11368a381b04SJed Brown PetscFunctionBegin; 113796400bd6SLisandro Dalcin ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr); 11388a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 1139e817cc15SEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr); 11408a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 1141d5e6173cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1142d5e6173cSPeter Brune ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1143cdb298fcSPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 114496400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 11458a381b04SJed Brown PetscFunctionReturn(0); 11468a381b04SJed Brown } 11478a381b04SJed Brown /*------------------------------------------------------------*/ 11488a381b04SJed Brown 11494416b707SBarry Smith static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 11508a381b04SJed Brown { 11514cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11528a381b04SJed Brown PetscErrorCode ierr; 11538a381b04SJed Brown 11548a381b04SJed Brown PetscFunctionBegin; 1155e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");CHKERRQ(ierr); 11568a381b04SJed Brown { 11578a381b04SJed Brown ARKTableauLink link; 11588a381b04SJed Brown PetscInt count,choice; 11598a381b04SJed Brown PetscBool flg; 11608a381b04SJed Brown const char **namelist; 11618a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 1162f489ac74SBarry Smith ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 11638a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 116496400bd6SLisandro Dalcin ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);CHKERRQ(ierr); 116596400bd6SLisandro Dalcin if (flg) {ierr = TSARKIMEXSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11668a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 116796400bd6SLisandro Dalcin 11684cc180ffSJed Brown flg = (PetscBool) !ark->imex; 11690298fd71SBarry Smith ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr); 11704cc180ffSJed Brown ark->imex = (PetscBool) !flg; 117103842d09SLisandro 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); 11728a381b04SJed Brown } 11738a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 11748a381b04SJed Brown PetscFunctionReturn(0); 11758a381b04SJed Brown } 11768a381b04SJed Brown 11778a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 11788a381b04SJed Brown { 11798a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11808a381b04SJed Brown PetscBool iascii; 11818a381b04SJed Brown PetscErrorCode ierr; 11828a381b04SJed Brown 11838a381b04SJed Brown PetscFunctionBegin; 1184251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11858a381b04SJed Brown if (iascii) { 11869c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 118719fd82e9SBarry Smith TSARKIMEXType arktype; 11888a381b04SJed Brown char buf[512]; 1189*3a28c0e5SStefano Zampini PetscBool flg; 1190*3a28c0e5SStefano Zampini 11918a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 1192*3a28c0e5SStefano Zampini ierr = TSARKIMEXGetFullyImplicit(ts,&flg);CHKERRQ(ierr); 11938a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 11948caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 119531f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 11968caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1197*3a28c0e5SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Fully implicit: %s\n",flg ? "yes" : "no");CHKERRQ(ierr); 1198e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 1199e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr); 1200e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr); 120131f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 12028a381b04SJed Brown } 12038a381b04SJed Brown PetscFunctionReturn(0); 12048a381b04SJed Brown } 12058a381b04SJed Brown 1206f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1207f2c2a1b9SBarry Smith { 1208f2c2a1b9SBarry Smith PetscErrorCode ierr; 1209f2c2a1b9SBarry Smith SNES snes; 12109c334d8fSLisandro Dalcin TSAdapt adapt; 1211f2c2a1b9SBarry Smith 1212f2c2a1b9SBarry Smith PetscFunctionBegin; 12139c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12149c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1215f2c2a1b9SBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1216f2c2a1b9SBarry Smith ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 1217ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 12180298fd71SBarry Smith ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 12190298fd71SBarry Smith ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 1220f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1221f2c2a1b9SBarry Smith } 1222f2c2a1b9SBarry Smith 12238a381b04SJed Brown /*@C 12248a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 12258a381b04SJed Brown 12268a381b04SJed Brown Logically collective 12278a381b04SJed Brown 12288a381b04SJed Brown Input Parameter: 12298a381b04SJed Brown + ts - timestepping context 12308a381b04SJed Brown - arktype - type of ARK-IMEX scheme 12318a381b04SJed Brown 12329bd3e852SBarry Smith Options Database: 12339bd3e852SBarry Smith . -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5> 12349bd3e852SBarry Smith 12358a381b04SJed Brown Level: intermediate 12368a381b04SJed Brown 12379bd3e852SBarry Smith .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, 12389bd3e852SBarry Smith TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 12398a381b04SJed Brown @*/ 124019fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 12418a381b04SJed Brown { 12428a381b04SJed Brown PetscErrorCode ierr; 12438a381b04SJed Brown 12448a381b04SJed Brown PetscFunctionBegin; 12458a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1246b92453a8SLisandro Dalcin PetscValidCharPointer(arktype,2); 124719fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 12488a381b04SJed Brown PetscFunctionReturn(0); 12498a381b04SJed Brown } 12508a381b04SJed Brown 12518a381b04SJed Brown /*@C 12528a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 12538a381b04SJed Brown 12548a381b04SJed Brown Logically collective 12558a381b04SJed Brown 12568a381b04SJed Brown Input Parameter: 12578a381b04SJed Brown . ts - timestepping context 12588a381b04SJed Brown 12598a381b04SJed Brown Output Parameter: 12608a381b04SJed Brown . arktype - type of ARK-IMEX scheme 12618a381b04SJed Brown 12628a381b04SJed Brown Level: intermediate 12638a381b04SJed Brown 12648a381b04SJed Brown .seealso: TSARKIMEXGetType() 12658a381b04SJed Brown @*/ 126619fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 12678a381b04SJed Brown { 12688a381b04SJed Brown PetscErrorCode ierr; 12698a381b04SJed Brown 12708a381b04SJed Brown PetscFunctionBegin; 12718a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 127219fd82e9SBarry Smith ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 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 12814cc180ffSJed Brown Input Parameter: 12824cc180ffSJed Brown + ts - timestepping context 12834cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 12844cc180ffSJed Brown 12854cc180ffSJed Brown Level: intermediate 12864cc180ffSJed Brown 1287*3a28c0e5SStefano Zampini .seealso: TSARKIMEXGetType(), TSARKIMEXGetFullyImplicit() 12884cc180ffSJed Brown @*/ 12894cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 12904cc180ffSJed Brown { 12914cc180ffSJed Brown PetscErrorCode ierr; 12924cc180ffSJed Brown 12934cc180ffSJed Brown PetscFunctionBegin; 12944cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1295*3a28c0e5SStefano Zampini PetscValidLogicalCollectiveBool(ts,flg,2); 12964cc180ffSJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 12974cc180ffSJed Brown PetscFunctionReturn(0); 12984cc180ffSJed Brown } 12994cc180ffSJed Brown 1300*3a28c0e5SStefano Zampini /*@ 1301*3a28c0e5SStefano Zampini TSARKIMEXGetFullyImplicit - Inquires if both parts of the equation are solved implicitly 1302*3a28c0e5SStefano Zampini 1303*3a28c0e5SStefano Zampini Logically collective 1304*3a28c0e5SStefano Zampini 1305*3a28c0e5SStefano Zampini Input Parameter: 1306*3a28c0e5SStefano Zampini . ts - timestepping context 1307*3a28c0e5SStefano Zampini 1308*3a28c0e5SStefano Zampini Output Parameter 1309*3a28c0e5SStefano Zampini . flg - PETSC_TRUE for fully implicit 1310*3a28c0e5SStefano Zampini 1311*3a28c0e5SStefano Zampini Level: intermediate 1312*3a28c0e5SStefano Zampini 1313*3a28c0e5SStefano Zampini .seealso: TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit() 1314*3a28c0e5SStefano Zampini @*/ 1315*3a28c0e5SStefano Zampini PetscErrorCode TSARKIMEXGetFullyImplicit(TS ts,PetscBool *flg) 1316*3a28c0e5SStefano Zampini { 1317*3a28c0e5SStefano Zampini PetscErrorCode ierr; 1318*3a28c0e5SStefano Zampini 1319*3a28c0e5SStefano Zampini PetscFunctionBegin; 1320*3a28c0e5SStefano Zampini PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1321*3a28c0e5SStefano Zampini PetscValidPointer(flg,2); 1322*3a28c0e5SStefano Zampini ierr = PetscUseMethod(ts,"TSARKIMEXGetFullyImplicit_C",(TS,PetscBool*),(ts,flg));CHKERRQ(ierr); 1323*3a28c0e5SStefano Zampini PetscFunctionReturn(0); 1324*3a28c0e5SStefano Zampini } 1325*3a28c0e5SStefano Zampini 1326e0877f53SBarry Smith static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 13278a381b04SJed Brown { 13288a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13298a381b04SJed Brown 13308a381b04SJed Brown PetscFunctionBegin; 13318a381b04SJed Brown *arktype = ark->tableau->name; 13328a381b04SJed Brown PetscFunctionReturn(0); 13338a381b04SJed Brown } 1334e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype) 13358a381b04SJed Brown { 13368a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13378a381b04SJed Brown PetscErrorCode ierr; 13388a381b04SJed Brown PetscBool match; 13398a381b04SJed Brown ARKTableauLink link; 13408a381b04SJed Brown 13418a381b04SJed Brown PetscFunctionBegin; 13428a381b04SJed Brown if (ark->tableau) { 13438a381b04SJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 13448a381b04SJed Brown if (match) PetscFunctionReturn(0); 13458a381b04SJed Brown } 13468a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 13478a381b04SJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 13488a381b04SJed Brown if (match) { 134996400bd6SLisandro Dalcin if (ts->setupcalled) {ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr);} 13508a381b04SJed Brown ark->tableau = &link->tab; 135196400bd6SLisandro Dalcin if (ts->setupcalled) {ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr);} 13522ffb9264SLisandro Dalcin ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 13538a381b04SJed Brown PetscFunctionReturn(0); 13548a381b04SJed Brown } 13558a381b04SJed Brown } 1356ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 13578a381b04SJed Brown PetscFunctionReturn(0); 13588a381b04SJed Brown } 1359e0877f53SBarry Smith 1360e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 13614cc180ffSJed Brown { 13624cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13634cc180ffSJed Brown 13644cc180ffSJed Brown PetscFunctionBegin; 13654cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13664cc180ffSJed Brown PetscFunctionReturn(0); 13674cc180ffSJed Brown } 13688a381b04SJed Brown 1369*3a28c0e5SStefano Zampini static PetscErrorCode TSARKIMEXGetFullyImplicit_ARKIMEX(TS ts,PetscBool *flg) 1370*3a28c0e5SStefano Zampini { 1371*3a28c0e5SStefano Zampini TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1372*3a28c0e5SStefano Zampini 1373*3a28c0e5SStefano Zampini PetscFunctionBegin; 1374*3a28c0e5SStefano Zampini *flg = (PetscBool)!ark->imex; 1375*3a28c0e5SStefano Zampini PetscFunctionReturn(0); 1376*3a28c0e5SStefano Zampini } 1377*3a28c0e5SStefano Zampini 1378b3a6b972SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 1379b3a6b972SJed Brown { 1380b3a6b972SJed Brown PetscErrorCode ierr; 1381b3a6b972SJed Brown 1382b3a6b972SJed Brown PetscFunctionBegin; 1383b3a6b972SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 1384b3a6b972SJed Brown if (ts->dm) { 1385b3a6b972SJed Brown ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1386b3a6b972SJed Brown ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1387b3a6b972SJed Brown } 1388b3a6b972SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 1389b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr); 1390b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr); 1391b3a6b972SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 1392*3a28c0e5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 1393b3a6b972SJed Brown PetscFunctionReturn(0); 1394b3a6b972SJed Brown } 1395b3a6b972SJed Brown 13968a381b04SJed Brown /* ------------------------------------------------------------ */ 13978a381b04SJed Brown /*MC 1398a4386c9eSJed Brown TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes 13998a381b04SJed Brown 1400fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1401fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1402fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1403fca742c7SJed Brown 1404fca742c7SJed Brown Notes: 1405a4386c9eSJed Brown The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1406c8058688SBarry Smith 14075eca1a21SEmil Constantinescu If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 14085eca1a21SEmil Constantinescu 1409a4386c9eSJed 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). 1410fca742c7SJed Brown 1411d0685a90SJed Brown Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1412d0685a90SJed Brown 14138a381b04SJed Brown Level: beginner 14148a381b04SJed Brown 1415*3a28c0e5SStefano Zampini .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEXGetFullyImplicit(), 1416*3a28c0e5SStefano Zampini TSARKIMEX1BEE, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122, 1417d0685a90SJed Brown TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister() 14188a381b04SJed Brown 14198a381b04SJed Brown M*/ 14208cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 14218a381b04SJed Brown { 14228a381b04SJed Brown TS_ARKIMEX *th; 14238a381b04SJed Brown PetscErrorCode ierr; 14248a381b04SJed Brown 14258a381b04SJed Brown PetscFunctionBegin; 1426607a6623SBarry Smith ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 14278a381b04SJed Brown 14288a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 14298a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 14308a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1431f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 14328a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 14338a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1434cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1435108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 143624655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 14378a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 14388a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 14398a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 14408a381b04SJed Brown 1441efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1442efd4aadfSBarry Smith 1443b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 14448a381b04SJed Brown ts->data = (void*)th; 14454cc180ffSJed Brown th->imex = PETSC_TRUE; 14468a381b04SJed Brown 1447bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 1448bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 1449bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 1450*3a28c0e5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetFullyImplicit_C",TSARKIMEXGetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 145196400bd6SLisandro Dalcin 145296400bd6SLisandro Dalcin ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 14538a381b04SJed Brown PetscFunctionReturn(0); 14548a381b04SJed Brown } 1455