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 651f80e275SEmil Constantinescu References: 6696a0c994SBarry 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). 671f80e275SEmil Constantinescu 681f80e275SEmil Constantinescu Level: advanced 691f80e275SEmil Constantinescu 701f80e275SEmil Constantinescu .seealso: TSARKIMEX 711f80e275SEmil Constantinescu M*/ 721f80e275SEmil Constantinescu /*MC 731f80e275SEmil Constantinescu TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 741f80e275SEmil Constantinescu 751f80e275SEmil 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. 761f80e275SEmil Constantinescu 771f80e275SEmil Constantinescu Level: advanced 781f80e275SEmil Constantinescu 791f80e275SEmil Constantinescu .seealso: TSARKIMEX 801f80e275SEmil Constantinescu M*/ 811f80e275SEmil Constantinescu /*MC 821f80e275SEmil Constantinescu TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 831f80e275SEmil Constantinescu 841f80e275SEmil Constantinescu This method has two implicit stages, and L-stable implicit scheme. 851f80e275SEmil Constantinescu 861f80e275SEmil Constantinescu References: 8796a0c994SBarry 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. 881f80e275SEmil Constantinescu 891f80e275SEmil Constantinescu Level: advanced 901f80e275SEmil Constantinescu 911f80e275SEmil Constantinescu .seealso: TSARKIMEX 921f80e275SEmil Constantinescu M*/ 931f80e275SEmil Constantinescu /*MC 94e817cc15SEmil Constantinescu TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method. 95e817cc15SEmil Constantinescu 96e817cc15SEmil Constantinescu This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used. 97e817cc15SEmil Constantinescu 98e817cc15SEmil Constantinescu Level: advanced 99e817cc15SEmil Constantinescu 100e817cc15SEmil Constantinescu .seealso: TSARKIMEX 101e817cc15SEmil Constantinescu M*/ 102e817cc15SEmil Constantinescu /*MC 1031f80e275SEmil Constantinescu TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 1041f80e275SEmil Constantinescu 1051f80e275SEmil 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. 1061f80e275SEmil Constantinescu 1071f80e275SEmil Constantinescu Level: advanced 1081f80e275SEmil Constantinescu 1091f80e275SEmil Constantinescu .seealso: TSARKIMEX 1101f80e275SEmil Constantinescu M*/ 11164f491ddSJed Brown /*MC 11264f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 11364f491ddSJed Brown 114617a39beSEmil 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. 11564f491ddSJed Brown 116b330ce4dSSatish Balay Level: advanced 117b330ce4dSSatish Balay 11864f491ddSJed Brown .seealso: TSARKIMEX 11964f491ddSJed Brown M*/ 12064f491ddSJed Brown /*MC 12164f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 12264f491ddSJed Brown 12364f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 12464f491ddSJed Brown 125b330ce4dSSatish Balay Level: advanced 126b330ce4dSSatish Balay 12764f491ddSJed Brown .seealso: TSARKIMEX 12864f491ddSJed Brown M*/ 12964f491ddSJed Brown /*MC 1306cf0794eSJed Brown TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 1316cf0794eSJed Brown 1326cf0794eSJed Brown This method has three implicit stages. 1336cf0794eSJed Brown 1346cf0794eSJed Brown References: 13596a0c994SBarry 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. 1366cf0794eSJed Brown 1376cf0794eSJed Brown This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375 1386cf0794eSJed Brown 1396cf0794eSJed Brown Level: advanced 1406cf0794eSJed Brown 1416cf0794eSJed Brown .seealso: TSARKIMEX 1426cf0794eSJed Brown M*/ 1436cf0794eSJed Brown /*MC 14464f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 14564f491ddSJed Brown 14664f491ddSJed Brown This method has one explicit stage and three implicit stages. 14764f491ddSJed Brown 14864f491ddSJed Brown References: 14996a0c994SBarry Smith . 1. - Kennedy and Carpenter 2003. 15064f491ddSJed Brown 151b330ce4dSSatish Balay Level: advanced 152b330ce4dSSatish Balay 15364f491ddSJed Brown .seealso: TSARKIMEX 15464f491ddSJed Brown M*/ 15564f491ddSJed Brown /*MC 1566cf0794eSJed Brown TSARKIMEXARS443 - Third order ARK IMEX scheme. 1576cf0794eSJed Brown 1586cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1596cf0794eSJed Brown 1606cf0794eSJed Brown References: 16196a0c994SBarry 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). 16296a0c994SBarry Smith - 2. - This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375 1636cf0794eSJed Brown 1646cf0794eSJed Brown Level: advanced 1656cf0794eSJed Brown 1666cf0794eSJed Brown .seealso: TSARKIMEX 1676cf0794eSJed Brown M*/ 1686cf0794eSJed Brown /*MC 1696cf0794eSJed Brown TSARKIMEXBPR3 - Third order ARK IMEX scheme. 1706cf0794eSJed Brown 1716cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1726cf0794eSJed Brown 1736cf0794eSJed Brown References: 17496a0c994SBarry Smith . This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375 1756cf0794eSJed Brown 1766cf0794eSJed Brown Level: advanced 1776cf0794eSJed Brown 1786cf0794eSJed Brown .seealso: TSARKIMEX 1796cf0794eSJed Brown M*/ 1806cf0794eSJed Brown /*MC 18164f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 18264f491ddSJed Brown 18364f491ddSJed Brown This method has one explicit stage and four implicit stages. 18464f491ddSJed Brown 18564f491ddSJed Brown References: 18696a0c994SBarry Smith . Kennedy and Carpenter 2003. 18764f491ddSJed Brown 188b330ce4dSSatish Balay Level: advanced 189b330ce4dSSatish Balay 19064f491ddSJed Brown .seealso: TSARKIMEX 19164f491ddSJed Brown M*/ 19264f491ddSJed Brown /*MC 19364f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 19464f491ddSJed Brown 19564f491ddSJed Brown This method has one explicit stage and five implicit stages. 19664f491ddSJed Brown 19764f491ddSJed Brown References: 19896a0c994SBarry Smith . Kennedy and Carpenter 2003. 19964f491ddSJed Brown 200b330ce4dSSatish Balay Level: advanced 201b330ce4dSSatish Balay 20264f491ddSJed Brown .seealso: TSARKIMEX 20364f491ddSJed Brown M*/ 20464f491ddSJed Brown 2058a381b04SJed Brown #undef __FUNCT__ 2068a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll" 2078a381b04SJed Brown /*@C 2088a381b04SJed Brown TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 2098a381b04SJed Brown 210fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 2118a381b04SJed Brown 2128a381b04SJed Brown Level: advanced 2138a381b04SJed Brown 2148a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all 2158a381b04SJed Brown 2168a381b04SJed Brown .seealso: TSARKIMEXRegisterDestroy() 2178a381b04SJed Brown @*/ 2188a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 2198a381b04SJed Brown { 2208a381b04SJed Brown PetscErrorCode ierr; 2218a381b04SJed Brown 2228a381b04SJed Brown PetscFunctionBegin; 2238a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2248a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 225e817cc15SEmil Constantinescu 226e817cc15SEmil Constantinescu { 227e817cc15SEmil Constantinescu const PetscReal 228e817cc15SEmil Constantinescu A[3][3] = {{0.0,0.0,0.0}, 229e817cc15SEmil Constantinescu {0.0,0.0,0.0}, 230748ad121SEmil Constantinescu {0.0,0.5,0.0}}, 231e817cc15SEmil Constantinescu At[3][3] = {{1.0,0.0,0.0}, 232e817cc15SEmil Constantinescu {0.0,0.5,0.0}, 233e817cc15SEmil Constantinescu {0.0,0.5,0.5}}, 234e817cc15SEmil Constantinescu b[3] = {0.0,0.5,0.5}, 235e817cc15SEmil Constantinescu bembedt[3] = {1.0,0.0,0.0}; 2360298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 237e817cc15SEmil Constantinescu } 2388a381b04SJed Brown { 2398a381b04SJed Brown const PetscReal 2401f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2411f80e275SEmil Constantinescu {0.5,0.0}}, 2421f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2431f80e275SEmil Constantinescu {0.0,0.5}}, 2441f80e275SEmil Constantinescu b[2] = {0.0,1.0}, 2451f80e275SEmil Constantinescu bembedt[2] = {0.5,0.5}; 2461f80e275SEmil 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*/ 2470298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 2481f80e275SEmil Constantinescu } 2491f80e275SEmil Constantinescu { 2501f80e275SEmil Constantinescu const PetscReal 2511f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2521f80e275SEmil Constantinescu {1.0,0.0}}, 2531f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2541f80e275SEmil Constantinescu {0.5,0.5}}, 2551f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2561f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}; 2571f80e275SEmil 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*/ 2580298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr); 2591f80e275SEmil Constantinescu } 2601f80e275SEmil Constantinescu { 261da80777bSKarl 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 */ 2621f80e275SEmil Constantinescu const PetscReal 2631f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2641f80e275SEmil Constantinescu {1.0,0.0}}, 265da80777bSKarl Rupp At[2][2] = {{0.2928932188134524755992,0.0}, 266da80777bSKarl Rupp {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}}, 2671f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2681f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}, 269da80777bSKarl Rupp binterpt[2][2] = {{ (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}, 270da80777bSKarl Rupp {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}}, 2711f80e275SEmil Constantinescu binterp[2][2] = {{1.0,-0.5},{0.0,0.5}}; 2720298fd71SBarry 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); 2731f80e275SEmil Constantinescu } 2741f80e275SEmil Constantinescu { 275da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 276da80777bSKarl Rupp const PetscReal 2778a381b04SJed Brown A[3][3] = {{0,0,0}, 278da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 279617a39beSEmil Constantinescu {0.5,0.5,0}}, 280da80777bSKarl Rupp At[3][3] = {{0,0,0}, 281da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 282da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 283da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 284da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 285da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 286da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 2870298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 2881f80e275SEmil Constantinescu } 2891f80e275SEmil Constantinescu { 290da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 291da80777bSKarl Rupp const PetscReal 2921f80e275SEmil Constantinescu A[3][3] = {{0,0,0}, 293da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 2948a381b04SJed Brown {0.75,0.25,0}}, 295da80777bSKarl Rupp At[3][3] = {{0,0,0}, 296da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 297da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 298da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 299da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 300da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 301da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3020298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 3038a381b04SJed Brown } 30406db7b1cSJed Brown { /* Optimal for linear implicit part */ 305da80777bSKarl Rupp /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time */ 306da80777bSKarl Rupp const PetscReal 307da80777bSKarl Rupp A[3][3] = {{0,0,0}, 308da80777bSKarl Rupp {2-1.414213562373095048802,0,0}, 309da80777bSKarl Rupp {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}}, 310da80777bSKarl Rupp At[3][3] = {{0,0,0}, 311da80777bSKarl Rupp {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0}, 312da80777bSKarl Rupp {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}}, 313da80777bSKarl Rupp bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)}, 314da80777bSKarl Rupp binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 315da80777bSKarl Rupp {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)}, 316da80777bSKarl Rupp {1.0-1.414213562373095048802,1.0/1.414213562373095048802}}; 3170298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 318a3a57f36SJed Brown } 3196cf0794eSJed Brown { /* Optimal for linear implicit part */ 3206cf0794eSJed Brown const PetscReal 3216cf0794eSJed Brown A[3][3] = {{0,0,0}, 3226cf0794eSJed Brown {0.5,0,0}, 3236cf0794eSJed Brown {0.5,0.5,0}}, 3246cf0794eSJed Brown At[3][3] = {{0.25,0,0}, 3256cf0794eSJed Brown {0,0.25,0}, 3266cf0794eSJed Brown {1./3,1./3,1./3}}; 3270298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr); 3286cf0794eSJed Brown } 329a3a57f36SJed Brown { 330a3a57f36SJed Brown const PetscReal 331a3a57f36SJed Brown A[4][4] = {{0,0,0,0}, 3324040e9f2SJed Brown {1767732205903./2027836641118.,0,0,0}, 3334040e9f2SJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 3344040e9f2SJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 335a3a57f36SJed Brown At[4][4] = {{0,0,0,0}, 3364040e9f2SJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 3374040e9f2SJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 3384040e9f2SJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 339cc46b9d1SJed Brown bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.}, 3404040e9f2SJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 3414040e9f2SJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 3424040e9f2SJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 3434040e9f2SJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 3440298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr); 345a3a57f36SJed Brown } 346a3a57f36SJed Brown { 347a3a57f36SJed Brown const PetscReal 348e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3496cf0794eSJed Brown {1./2,0,0,0,0}, 3506cf0794eSJed Brown {11./18,1./18,0,0,0}, 3516cf0794eSJed Brown {5./6,-5./6,.5,0,0}, 3526cf0794eSJed Brown {1./4,7./4,3./4,-7./4,0}}, 3536cf0794eSJed Brown At[5][5] = {{0,0,0,0,0}, 3546cf0794eSJed Brown {0,1./2,0,0,0}, 3556cf0794eSJed Brown {0,1./6,1./2,0,0}, 3566cf0794eSJed Brown {0,-1./2,1./2,1./2,0}, 357108c343cSJed Brown {0,3./2,-3./2,1./2,1./2}}, 3580298fd71SBarry Smith *bembedt = NULL; 3590298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 3606cf0794eSJed Brown } 3616cf0794eSJed Brown { 3626cf0794eSJed Brown const PetscReal 363e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3646cf0794eSJed Brown {1,0,0,0,0}, 3656cf0794eSJed Brown {4./9,2./9,0,0,0}, 3666cf0794eSJed Brown {1./4,0,3./4,0,0}, 3676cf0794eSJed Brown {1./4,0,3./5,0,0}}, 368e74514c0SSatish Balay At[5][5] = {{0,0,0,0,0}, 3696cf0794eSJed Brown {.5,.5,0,0,0}, 3706cf0794eSJed Brown {5./18,-1./9,.5,0,0}, 3716cf0794eSJed Brown {.5,0,0,.5,0}, 372108c343cSJed Brown {.25,0,.75,-.5,.5}}, 3730298fd71SBarry Smith *bembedt = NULL; 3740298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr); 3756cf0794eSJed Brown } 3766cf0794eSJed Brown { 3776cf0794eSJed Brown const PetscReal 378a3a57f36SJed Brown A[6][6] = {{0,0,0,0,0,0}, 379a3a57f36SJed Brown {1./2,0,0,0,0,0}, 3804040e9f2SJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 3814040e9f2SJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 3824040e9f2SJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 3834040e9f2SJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 384a3a57f36SJed Brown At[6][6] = {{0,0,0,0,0,0}, 385a3a57f36SJed Brown {1./4,1./4,0,0,0,0}, 3864040e9f2SJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 3874040e9f2SJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 3884040e9f2SJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 3894040e9f2SJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 390cc46b9d1SJed Brown bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.}, 3914040e9f2SJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 392cd652676SJed Brown {0,0,0}, 3934040e9f2SJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 3944040e9f2SJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 3954040e9f2SJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 3964040e9f2SJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 3970298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 398a3a57f36SJed Brown } 399a3a57f36SJed Brown { 400a3a57f36SJed Brown const PetscReal 401a3a57f36SJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 402a3a57f36SJed Brown {41./100,0,0,0,0,0,0,0}, 4034040e9f2SJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 4044040e9f2SJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 4054040e9f2SJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 4064040e9f2SJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 4074040e9f2SJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 4084040e9f2SJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 409a3a57f36SJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 4104040e9f2SJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 4114040e9f2SJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 4124040e9f2SJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 4134040e9f2SJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 4144040e9f2SJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 4154040e9f2SJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 4164040e9f2SJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 417cc46b9d1SJed Brown bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.}, 4184040e9f2SJed Brown binterpt[8][3] = {{-17674230611817./10670229744614., 43486358583215./12773830924787., -9257016797708./5021505065439.}, 419cd652676SJed Brown {0, 0, 0 }, 420cd652676SJed Brown {0, 0, 0 }, 4214040e9f2SJed Brown {65168852399939./7868540260826., -91478233927265./11067650958493., 26096422576131./11239449250142.}, 4224040e9f2SJed Brown {15494834004392./5936557850923., -79368583304911./10890268929626., 92396832856987./20362823103730.}, 4234040e9f2SJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473., 30029262896817./10175596800299.}, 4244040e9f2SJed Brown {-19024464361622./5461577185407., 115839755401235./10719374521269., -26136350496073./3983972220547.}, 4254040e9f2SJed Brown {-6511271360970./6095937251113., 5843115559534./2180450260947., -5289405421727./3760307252460. }}; 4260298fd71SBarry Smith ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr); 427a3a57f36SJed Brown } 4288a381b04SJed Brown PetscFunctionReturn(0); 4298a381b04SJed Brown } 4308a381b04SJed Brown 4318a381b04SJed Brown #undef __FUNCT__ 4328a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy" 4338a381b04SJed Brown /*@C 4348a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 4358a381b04SJed Brown 4368a381b04SJed Brown Not Collective 4378a381b04SJed Brown 4388a381b04SJed Brown Level: advanced 4398a381b04SJed Brown 4408a381b04SJed Brown .keywords: TSARKIMEX, register, destroy 441607a6623SBarry Smith .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll() 4428a381b04SJed Brown @*/ 4438a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 4448a381b04SJed Brown { 4458a381b04SJed Brown PetscErrorCode ierr; 4468a381b04SJed Brown ARKTableauLink link; 4478a381b04SJed Brown 4488a381b04SJed Brown PetscFunctionBegin; 4498a381b04SJed Brown while ((link = ARKTableauList)) { 4508a381b04SJed Brown ARKTableau t = &link->tab; 4518a381b04SJed Brown ARKTableauList = link->next; 4528a381b04SJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 453108c343cSJed Brown ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr); 454cd652676SJed Brown ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 4558a381b04SJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 4568a381b04SJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 4578a381b04SJed Brown } 4588a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4598a381b04SJed Brown PetscFunctionReturn(0); 4608a381b04SJed Brown } 4618a381b04SJed Brown 4628a381b04SJed Brown #undef __FUNCT__ 4638a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage" 4648a381b04SJed Brown /*@C 4658a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 4668a381b04SJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 4678a381b04SJed Brown when using static libraries. 4688a381b04SJed Brown 4698a381b04SJed Brown Level: developer 4708a381b04SJed Brown 4718a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package 4728a381b04SJed Brown .seealso: PetscInitialize() 4738a381b04SJed Brown @*/ 474607a6623SBarry Smith PetscErrorCode TSARKIMEXInitializePackage(void) 4758a381b04SJed Brown { 4768a381b04SJed Brown PetscErrorCode ierr; 4778a381b04SJed Brown 4788a381b04SJed Brown PetscFunctionBegin; 4798a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 4808a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 4818a381b04SJed Brown ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 4828a381b04SJed Brown ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 4838a381b04SJed Brown PetscFunctionReturn(0); 4848a381b04SJed Brown } 4858a381b04SJed Brown 4868a381b04SJed Brown #undef __FUNCT__ 4878a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage" 4888a381b04SJed Brown /*@C 4898a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 4908a381b04SJed Brown called from PetscFinalize(). 4918a381b04SJed Brown 4928a381b04SJed Brown Level: developer 4938a381b04SJed Brown 4948a381b04SJed Brown .keywords: Petsc, destroy, package 4958a381b04SJed Brown .seealso: PetscFinalize() 4968a381b04SJed Brown @*/ 4978a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 4988a381b04SJed Brown { 4998a381b04SJed Brown PetscErrorCode ierr; 5008a381b04SJed Brown 5018a381b04SJed Brown PetscFunctionBegin; 5028a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 5038a381b04SJed Brown ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 5048a381b04SJed Brown PetscFunctionReturn(0); 5058a381b04SJed Brown } 5068a381b04SJed Brown 5078a381b04SJed Brown #undef __FUNCT__ 5088a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister" 509cd652676SJed Brown /*@C 510cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 511cd652676SJed Brown 512cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 513cd652676SJed Brown 514cd652676SJed Brown Input Parameters: 515cd652676SJed Brown + name - identifier for method 516cd652676SJed Brown . order - approximation order of method 517cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 518cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 5190298fd71SBarry Smith . bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At) 5200298fd71SBarry Smith . ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At) 521cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 5220298fd71SBarry Smith . b - Non-stiff step completion table (dimension s; NULL to use last row of At) 5230298fd71SBarry Smith . c - Non-stiff abscissa (dimension s; NULL to use row sums of A) 5240298fd71SBarry Smith . bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available) 5250298fd71SBarry Smith . bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided) 526cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 527cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 5280298fd71SBarry Smith - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt) 529cd652676SJed Brown 530cd652676SJed Brown Notes: 531cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 532cd652676SJed Brown 533cd652676SJed Brown Level: advanced 534cd652676SJed Brown 535cd652676SJed Brown .keywords: TS, register 536cd652676SJed Brown 537cd652676SJed Brown .seealso: TSARKIMEX 538cd652676SJed Brown @*/ 53919fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s, 5408a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 541cd652676SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 542108c343cSJed Brown const PetscReal bembedt[],const PetscReal bembed[], 543cd652676SJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 5448a381b04SJed Brown { 5458a381b04SJed Brown PetscErrorCode ierr; 5468a381b04SJed Brown ARKTableauLink link; 5478a381b04SJed Brown ARKTableau t; 5488a381b04SJed Brown PetscInt i,j; 5498a381b04SJed Brown 5508a381b04SJed Brown PetscFunctionBegin; 5511795a4d1SJed Brown ierr = PetscCalloc1(1,&link);CHKERRQ(ierr); 5528a381b04SJed Brown t = &link->tab; 5538a381b04SJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 5548a381b04SJed Brown t->order = order; 5558a381b04SJed Brown t->s = s; 556dcca6d9dSJed 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); 5578a381b04SJed Brown ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 5588a381b04SJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 5598a381b04SJed Brown if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); } 5608a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 5618a381b04SJed Brown if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 5625dceddf7SDebojyoti Ghosh else for (i=0; i<s; i++) t->b[i] = t->bt[i]; 5638a381b04SJed Brown if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); } 5648a381b04SJed 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]; 5658a381b04SJed Brown if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 5668a381b04SJed 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]; 567e817cc15SEmil Constantinescu t->stiffly_accurate = PETSC_TRUE; 568e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE; 569e817cc15SEmil Constantinescu t->explicit_first_stage = PETSC_TRUE; 570e817cc15SEmil Constantinescu for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE; 571e817cc15SEmil Constantinescu /*def of FSAL can be made more precise*/ 5724e9d4bf5SJed Brown t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate); 573108c343cSJed Brown if (bembedt) { 574dcca6d9dSJed Brown ierr = PetscMalloc2(s,&t->bembedt,s,&t->bembed);CHKERRQ(ierr); 575108c343cSJed Brown ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr); 576108c343cSJed Brown ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr); 577108c343cSJed Brown } 578108c343cSJed Brown 5794f385281SJed Brown t->pinterp = pinterp; 580dcca6d9dSJed Brown ierr = PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);CHKERRQ(ierr); 581cd652676SJed Brown ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 582cd652676SJed Brown ierr = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 5838a381b04SJed Brown link->next = ARKTableauList; 5848a381b04SJed Brown ARKTableauList = link; 5858a381b04SJed Brown PetscFunctionReturn(0); 5868a381b04SJed Brown } 5878a381b04SJed Brown 5888a381b04SJed Brown #undef __FUNCT__ 589108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX" 590108c343cSJed Brown /* 591108c343cSJed Brown The step completion formula is 592108c343cSJed Brown 593108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 594108c343cSJed Brown 595108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 596108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 597108c343cSJed Brown We can write 598108c343cSJed Brown 599108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 600108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 601108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 602108c343cSJed Brown 603108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 604108c343cSJed Brown */ 605108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 606108c343cSJed Brown { 607108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 608108c343cSJed Brown ARKTableau tab = ark->tableau; 609108c343cSJed Brown PetscScalar *w = ark->work; 610108c343cSJed Brown PetscReal h; 611108c343cSJed Brown PetscInt s = tab->s,j; 612108c343cSJed Brown PetscErrorCode ierr; 613108c343cSJed Brown 614108c343cSJed Brown PetscFunctionBegin; 615108c343cSJed Brown switch (ark->status) { 616108c343cSJed Brown case TS_STEP_INCOMPLETE: 617108c343cSJed Brown case TS_STEP_PENDING: 618108c343cSJed Brown h = ts->time_step; break; 619108c343cSJed Brown case TS_STEP_COMPLETE: 620be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 621ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 622108c343cSJed Brown } 623108c343cSJed Brown if (order == tab->order) { 624e817cc15SEmil Constantinescu if (ark->status == TS_STEP_INCOMPLETE) { 625740132f1SEmil Constantinescu if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */ 626e817cc15SEmil Constantinescu ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr); 627e817cc15SEmil Constantinescu } else { /* Use the standard completion formula (bt,b) */ 628108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 629e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bt[j]; 630108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 631e817cc15SEmil Constantinescu if (ark->imex) { /* Method is IMEX, complete the explicit formula */ 632108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->b[j]; 633108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 634e817cc15SEmil Constantinescu } 635e817cc15SEmil Constantinescu } 636108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 637108c343cSJed Brown if (done) *done = PETSC_TRUE; 638108c343cSJed Brown PetscFunctionReturn(0); 639108c343cSJed Brown } else if (order == tab->order-1) { 640108c343cSJed Brown if (!tab->bembedt) goto unavailable; 641108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 642108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 643e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembedt[j]; 644108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 645108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 646108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 647108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 648108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 649e817cc15SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]); 650108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr); 651108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 652108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 653108c343cSJed Brown } 654108c343cSJed Brown if (done) *done = PETSC_TRUE; 655108c343cSJed Brown PetscFunctionReturn(0); 656108c343cSJed Brown } 657108c343cSJed Brown unavailable: 658108c343cSJed Brown if (done) *done = PETSC_FALSE; 659a7fac7c2SEmil 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); 660108c343cSJed Brown PetscFunctionReturn(0); 661108c343cSJed Brown } 662108c343cSJed Brown 663108c343cSJed Brown #undef __FUNCT__ 66424655328SShri #define __FUNCT__ "TSRollBack_ARKIMEX" 66524655328SShri static PetscErrorCode TSRollBack_ARKIMEX(TS ts) 66624655328SShri { 66724655328SShri TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 66824655328SShri ARKTableau tab = ark->tableau; 66924655328SShri const PetscInt s = tab->s; 67024655328SShri const PetscReal *bt = tab->bt,*b = tab->b; 67124655328SShri PetscScalar *w = ark->work; 67224655328SShri Vec *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS; 67324655328SShri PetscInt j; 674be5899b3SLisandro Dalcin PetscReal h; 67524655328SShri PetscErrorCode ierr; 67624655328SShri 67724655328SShri PetscFunctionBegin; 678be5899b3SLisandro Dalcin switch (ark->status) { 679be5899b3SLisandro Dalcin case TS_STEP_INCOMPLETE: 680be5899b3SLisandro Dalcin case TS_STEP_PENDING: 681be5899b3SLisandro Dalcin h = ts->time_step; break; 682be5899b3SLisandro Dalcin case TS_STEP_COMPLETE: 683be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 684be5899b3SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 685be5899b3SLisandro Dalcin } 68624655328SShri for (j=0; j<s; j++) w[j] = -h*bt[j]; 68724655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 68824655328SShri for (j=0; j<s; j++) w[j] = -h*b[j]; 68924655328SShri ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 69024655328SShri PetscFunctionReturn(0); 69124655328SShri } 69224655328SShri 69324655328SShri #undef __FUNCT__ 6948a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX" 6958a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 6968a381b04SJed Brown { 6978a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 6988a381b04SJed Brown ARKTableau tab = ark->tableau; 6998a381b04SJed Brown const PetscInt s = tab->s; 70024655328SShri const PetscReal *At = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c; 701406d0ec2SJed Brown PetscScalar *w = ark->work; 7021297b224SEmil Constantinescu Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z; 70396400bd6SLisandro Dalcin PetscBool extrapolate = ark->extrapolate; 704108c343cSJed Brown TSAdapt adapt; 7058a381b04SJed Brown SNES snes; 706fecfb714SLisandro Dalcin PetscInt i,j,its,lits; 707be5899b3SLisandro Dalcin PetscInt rejections = 0; 70896400bd6SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 70996400bd6SLisandro Dalcin PetscReal next_time_step = ts->time_step; 7108a381b04SJed Brown PetscErrorCode ierr; 7118a381b04SJed Brown 7128a381b04SJed Brown PetscFunctionBegin; 71396400bd6SLisandro Dalcin if (ark->extrapolate && !ark->Y_prev) { 71496400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 71596400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 71696400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 71796400bd6SLisandro Dalcin } 71896400bd6SLisandro Dalcin 71996400bd6SLisandro Dalcin if (!ts->steprollback) { 72096400bd6SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */ 72196400bd6SLisandro Dalcin ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr); 72296400bd6SLisandro Dalcin } 723fecfb714SLisandro Dalcin if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */ 72496400bd6SLisandro Dalcin for (i = 0; i<s; i++) { 72596400bd6SLisandro Dalcin ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr); 72696400bd6SLisandro Dalcin ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr); 72796400bd6SLisandro Dalcin ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr); 72896400bd6SLisandro Dalcin } 72996400bd6SLisandro Dalcin } 73096400bd6SLisandro Dalcin } 73196400bd6SLisandro Dalcin 732fecfb714SLisandro Dalcin if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) { 73396400bd6SLisandro Dalcin TS ts_start; 734baa10174SEmil Constantinescu ierr = TSClone(ts,&ts_start);CHKERRQ(ierr); 735e817cc15SEmil Constantinescu ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr); 736e817cc15SEmil Constantinescu ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr); 737eb082435SEmil Constantinescu ierr = TSSetDuration(ts_start,1,ts->ptime+ts->time_step);CHKERRQ(ierr); 738feed9e9dSBarry Smith ierr = TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 739740132f1SEmil Constantinescu ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr); 740e817cc15SEmil Constantinescu ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr); 741740132f1SEmil Constantinescu ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr); 742e817cc15SEmil Constantinescu ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr); 74334561852SEmil Constantinescu 744e817cc15SEmil Constantinescu ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr); 745e817cc15SEmil Constantinescu ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr); 74696400bd6SLisandro Dalcin ierr = TSGetTimeStep(ts_start,&ts->time_step);CHKERRQ(ierr); 747bbd56ea5SKarl Rupp 74885fc7851SLisandro Dalcin { /* Save the initial slope for the next step */ 74985fc7851SLisandro Dalcin TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data; 75085fc7851SLisandro Dalcin ierr = VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);CHKERRQ(ierr); 75185fc7851SLisandro Dalcin } 75296400bd6SLisandro Dalcin ts->steps++; 753be5899b3SLisandro Dalcin ts->total_steps++; 75434561852SEmil Constantinescu 755d15a3a53SEmil Constantinescu /* Set the correct TS in SNES */ 756d15a3a53SEmil Constantinescu /* We'll try to bypass this by changing the method on the fly */ 75796400bd6SLisandro Dalcin { 75896400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 75996400bd6SLisandro Dalcin ierr = TSSetSNES(ts,snes);CHKERRQ(ierr); 76096400bd6SLisandro Dalcin } 761166a6834SEmil Constantinescu ierr = TSDestroy(&ts_start);CHKERRQ(ierr); 762e817cc15SEmil Constantinescu } 763e817cc15SEmil Constantinescu 764108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 76596400bd6SLisandro Dalcin while (!ts->reason && ark->status != TS_STEP_COMPLETE) { 76696400bd6SLisandro Dalcin PetscReal t = ts->ptime; 767108c343cSJed Brown PetscReal h = ts->time_step; 7688a381b04SJed Brown for (i=0; i<s; i++) { 7699be3e283SDebojyoti Ghosh ark->stage_time = t + h*ct[i]; 77096400bd6SLisandro Dalcin ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr); 7718a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 7726c4ed002SBarry 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"); 7738a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 774e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 7758a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 7768a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 7778a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 7788a381b04SJed Brown } else { 779b296d7d5SJed Brown ark->scoeff = 1./At[i*s+i]; 7808a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 7818a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 782e817cc15SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 7834f385281SJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 784c58d1302SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 785c58d1302SEmil Constantinescu ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr); 786fecfb714SLisandro Dalcin if (extrapolate && !ts->steprestart) { 78756dcabbaSDebojyoti Ghosh /* Initial guess extrapolated from previous time step stage values */ 78856dcabbaSDebojyoti Ghosh ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr); 78956dcabbaSDebojyoti Ghosh } else { 7908a381b04SJed Brown /* Initial guess taken from last stage */ 7918a381b04SJed Brown ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr); 79256dcabbaSDebojyoti Ghosh } 79396400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 794baa10174SEmil Constantinescu ierr = SNESSolve(snes,NULL,Y[i]);CHKERRQ(ierr); 7958a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 7968a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 7975ef26d82SJed Brown ts->snes_its += its; ts->ksp_its += lits; 798552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 79996400bd6SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);CHKERRQ(ierr); 80096400bd6SLisandro Dalcin if (!stageok) { 8011be93e3eSJed Brown /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to 8021be93e3eSJed Brown * use extrapolation to initialize the solves on the next attempt. */ 80396400bd6SLisandro Dalcin extrapolate = PETSC_FALSE; 8041be93e3eSJed Brown goto reject_step; 8051be93e3eSJed Brown } 8068a381b04SJed Brown } 807e817cc15SEmil Constantinescu if (ts->equation_type >= TS_EQ_IMPLICIT) { 808e817cc15SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8096c4ed002SBarry 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); 810df5e1e3dSEmil Constantinescu ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr); /* YdotI = YdotI(tn-1) */ 811e817cc15SEmil Constantinescu } else { 812df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 813e817cc15SEmil Constantinescu } 814e817cc15SEmil Constantinescu } else { 8155eca1a21SEmil Constantinescu if (i==0 && tab->explicit_first_stage) { 8168a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 817df5e1e3dSEmil Constantinescu ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);/* YdotI = -G(t,Y,0) */ 818e817cc15SEmil Constantinescu ierr = VecScale(YdotI[i],-1.0);CHKERRQ(ierr); 8195eca1a21SEmil Constantinescu } else { 820df5e1e3dSEmil Constantinescu ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* YdotI = shift*(X-Z) */ 8215eca1a21SEmil Constantinescu } 8224cc180ffSJed Brown if (ark->imex) { 8238a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 8244cc180ffSJed Brown } else { 8254cc180ffSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 8264cc180ffSJed Brown } 8278a381b04SJed Brown } 82896400bd6SLisandro Dalcin ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr); 829e817cc15SEmil Constantinescu } 83096400bd6SLisandro Dalcin 831be5899b3SLisandro Dalcin ark->status = TS_STEP_INCOMPLETE; 832fecfb714SLisandro Dalcin ierr = TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 833108c343cSJed Brown ark->status = TS_STEP_PENDING; 834552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 835108c343cSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 836fecfb714SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 837fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 83896400bd6SLisandro Dalcin ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 83996400bd6SLisandro Dalcin if (!accept) { /* Roll back the current step */ 84096400bd6SLisandro Dalcin ierr = TSRollBack_ARKIMEX(ts);CHKERRQ(ierr); 841be5899b3SLisandro Dalcin ts->time_step = next_time_step; 842be5899b3SLisandro Dalcin goto reject_step; 84396400bd6SLisandro Dalcin } 84496400bd6SLisandro Dalcin 8458a381b04SJed Brown ts->ptime += ts->time_step; 846cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 847108c343cSJed Brown break; 84896400bd6SLisandro Dalcin 84996400bd6SLisandro Dalcin reject_step: 850fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 851be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 85296400bd6SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 853be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 854108c343cSJed Brown } 855f85781f1SEmil Constantinescu } 8568a381b04SJed Brown PetscFunctionReturn(0); 8578a381b04SJed Brown } 8588a381b04SJed Brown 859cd652676SJed Brown #undef __FUNCT__ 860cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX" 861cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 862cd652676SJed Brown { 863cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 8644f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 865108c343cSJed Brown PetscReal h; 866108c343cSJed Brown PetscReal tt,t; 867cd652676SJed Brown PetscScalar *bt,*b; 868cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 869cd652676SJed Brown PetscErrorCode ierr; 870cd652676SJed Brown 871cd652676SJed Brown PetscFunctionBegin; 872ce94432eSBarry Smith if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 873108c343cSJed Brown switch (ark->status) { 874108c343cSJed Brown case TS_STEP_INCOMPLETE: 875108c343cSJed Brown case TS_STEP_PENDING: 876108c343cSJed Brown h = ts->time_step; 877108c343cSJed Brown t = (itime - ts->ptime)/h; 878108c343cSJed Brown break; 879108c343cSJed Brown case TS_STEP_COMPLETE: 880be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 881108c343cSJed Brown t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 882108c343cSJed Brown break; 883ce94432eSBarry Smith default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 884108c343cSJed Brown } 885dcca6d9dSJed Brown ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr); 886cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 8874f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 888cd652676SJed Brown for (i=0; i<s; i++) { 889c1758d98SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 890108c343cSJed Brown b[i] += h * B[i*pinterp+j] * tt; 891cd652676SJed Brown } 892cd652676SJed Brown } 893cd652676SJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 894cd652676SJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 895cd652676SJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 896cd652676SJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 897cd652676SJed Brown PetscFunctionReturn(0); 898cd652676SJed Brown } 899cd652676SJed Brown 90056dcabbaSDebojyoti Ghosh #undef __FUNCT__ 90156dcabbaSDebojyoti Ghosh #define __FUNCT__ "TSExtrapolate_ARKIMEX" 90256dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X) 90356dcabbaSDebojyoti Ghosh { 90456dcabbaSDebojyoti Ghosh TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 90556dcabbaSDebojyoti Ghosh PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 906be5899b3SLisandro Dalcin PetscReal h,h_prev,t,tt; 90756dcabbaSDebojyoti Ghosh PetscScalar *bt,*b; 90856dcabbaSDebojyoti Ghosh const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 90956dcabbaSDebojyoti Ghosh PetscErrorCode ierr; 91056dcabbaSDebojyoti Ghosh 91156dcabbaSDebojyoti Ghosh PetscFunctionBegin; 91256dcabbaSDebojyoti Ghosh if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 913be5899b3SLisandro Dalcin ierr = PetscCalloc2(s,&bt,s,&b);CHKERRQ(ierr); 91481d12688SDebojyoti Ghosh h = ts->time_step; 915be5899b3SLisandro Dalcin h_prev = ts->ptime - ts->ptime_prev; 916be5899b3SLisandro Dalcin t = 1 + h/h_prev*c; 91756dcabbaSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 91856dcabbaSDebojyoti Ghosh for (i=0; i<s; i++) { 91981d12688SDebojyoti Ghosh bt[i] += h * Bt[i*pinterp+j] * tt; 92056dcabbaSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 92156dcabbaSDebojyoti Ghosh } 92256dcabbaSDebojyoti Ghosh } 92396400bd6SLisandro Dalcin if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored"); 92456dcabbaSDebojyoti Ghosh ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr); 92556dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr); 92656dcabbaSDebojyoti Ghosh ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr); 92756dcabbaSDebojyoti Ghosh ierr = PetscFree2(bt,b);CHKERRQ(ierr); 92856dcabbaSDebojyoti Ghosh PetscFunctionReturn(0); 92956dcabbaSDebojyoti Ghosh } 93056dcabbaSDebojyoti Ghosh 9318a381b04SJed Brown /*------------------------------------------------------------*/ 93296400bd6SLisandro Dalcin 93396400bd6SLisandro Dalcin #undef __FUNCT__ 93496400bd6SLisandro Dalcin #define __FUNCT__ "TSARKIMEXTableauReset" 93596400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauReset(TS ts) 93696400bd6SLisandro Dalcin { 93796400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 93896400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 93996400bd6SLisandro Dalcin PetscErrorCode ierr; 94096400bd6SLisandro Dalcin 94196400bd6SLisandro Dalcin PetscFunctionBegin; 94296400bd6SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 94396400bd6SLisandro Dalcin ierr = PetscFree(ark->work);CHKERRQ(ierr); 94496400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y);CHKERRQ(ierr); 94596400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI);CHKERRQ(ierr); 94696400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS);CHKERRQ(ierr); 94796400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->Y_prev);CHKERRQ(ierr); 94896400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 94996400bd6SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 95096400bd6SLisandro Dalcin PetscFunctionReturn(0); 95196400bd6SLisandro Dalcin } 95296400bd6SLisandro Dalcin 9538a381b04SJed Brown #undef __FUNCT__ 9548a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX" 9558a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 9568a381b04SJed Brown { 9578a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9588a381b04SJed Brown PetscErrorCode ierr; 9598a381b04SJed Brown 9608a381b04SJed Brown PetscFunctionBegin; 96196400bd6SLisandro Dalcin ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr); 9628a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 963e817cc15SEmil Constantinescu ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr); 9648a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 9658a381b04SJed Brown PetscFunctionReturn(0); 9668a381b04SJed Brown } 9678a381b04SJed Brown 9688a381b04SJed Brown #undef __FUNCT__ 9698a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX" 9708a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 9718a381b04SJed Brown { 9728a381b04SJed Brown PetscErrorCode ierr; 9738a381b04SJed Brown 9748a381b04SJed Brown PetscFunctionBegin; 9758a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 9768a381b04SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 977bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr); 978bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr); 979bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr); 9808a381b04SJed Brown PetscFunctionReturn(0); 9818a381b04SJed Brown } 9828a381b04SJed Brown 983d5e6173cSPeter Brune 984d5e6173cSPeter Brune #undef __FUNCT__ 985d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXGetVecs" 986d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 987d5e6173cSPeter Brune { 988d5e6173cSPeter Brune TS_ARKIMEX *ax = (TS_ARKIMEX*)ts->data; 989d5e6173cSPeter Brune PetscErrorCode ierr; 990d5e6173cSPeter Brune 991d5e6173cSPeter Brune PetscFunctionBegin; 992d5e6173cSPeter Brune if (Z) { 993d5e6173cSPeter Brune if (dm && dm != ts->dm) { 994d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 995d5e6173cSPeter Brune } else *Z = ax->Z; 996d5e6173cSPeter Brune } 997d5e6173cSPeter Brune if (Ydot) { 998d5e6173cSPeter Brune if (dm && dm != ts->dm) { 999d5e6173cSPeter Brune ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1000d5e6173cSPeter Brune } else *Ydot = ax->Ydot; 1001d5e6173cSPeter Brune } 1002d5e6173cSPeter Brune PetscFunctionReturn(0); 1003d5e6173cSPeter Brune } 1004d5e6173cSPeter Brune 1005d5e6173cSPeter Brune 1006d5e6173cSPeter Brune #undef __FUNCT__ 1007d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXRestoreVecs" 1008d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot) 1009d5e6173cSPeter Brune { 1010d5e6173cSPeter Brune PetscErrorCode ierr; 1011d5e6173cSPeter Brune 1012d5e6173cSPeter Brune PetscFunctionBegin; 1013d5e6173cSPeter Brune if (Z) { 1014d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1015d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr); 1016d5e6173cSPeter Brune } 1017d5e6173cSPeter Brune } 1018d5e6173cSPeter Brune if (Ydot) { 1019d5e6173cSPeter Brune if (dm && dm != ts->dm) { 1020d5e6173cSPeter Brune ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr); 1021d5e6173cSPeter Brune } 1022d5e6173cSPeter Brune } 1023d5e6173cSPeter Brune PetscFunctionReturn(0); 1024d5e6173cSPeter Brune } 1025d5e6173cSPeter Brune 10268a381b04SJed Brown /* 10278a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 10288a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 10298a381b04SJed Brown */ 10308a381b04SJed Brown #undef __FUNCT__ 10318a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 10328a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 10338a381b04SJed Brown { 10348a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1035d5e6173cSPeter Brune DM dm,dmsave; 1036d5e6173cSPeter Brune Vec Z,Ydot; 1037b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10388a381b04SJed Brown PetscErrorCode ierr; 10398a381b04SJed Brown 10408a381b04SJed Brown PetscFunctionBegin; 1041d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1042d5e6173cSPeter Brune ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1043b296d7d5SJed Brown ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 1044d5e6173cSPeter Brune dmsave = ts->dm; 1045d5e6173cSPeter Brune ts->dm = dm; 1046740132f1SEmil Constantinescu 1047d5e6173cSPeter Brune ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr); 1048e817cc15SEmil Constantinescu 1049d5e6173cSPeter Brune ts->dm = dmsave; 1050d5e6173cSPeter Brune ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 10518a381b04SJed Brown PetscFunctionReturn(0); 10528a381b04SJed Brown } 10538a381b04SJed Brown 10548a381b04SJed Brown #undef __FUNCT__ 10558a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 1056d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 10578a381b04SJed Brown { 10588a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 1059d5e6173cSPeter Brune DM dm,dmsave; 1060d5e6173cSPeter Brune Vec Ydot; 1061b296d7d5SJed Brown PetscReal shift = ark->scoeff / ts->time_step; 10628a381b04SJed Brown PetscErrorCode ierr; 10638a381b04SJed Brown 10648a381b04SJed Brown PetscFunctionBegin; 1065d5e6173cSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 10660298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 10678a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 1068d5e6173cSPeter Brune dmsave = ts->dm; 1069d5e6173cSPeter Brune ts->dm = dm; 1070740132f1SEmil Constantinescu 1071d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr); 1072740132f1SEmil Constantinescu 1073d5e6173cSPeter Brune ts->dm = dmsave; 10740298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 1075d5e6173cSPeter Brune PetscFunctionReturn(0); 1076d5e6173cSPeter Brune } 1077d5e6173cSPeter Brune 1078d5e6173cSPeter Brune #undef __FUNCT__ 1079d5e6173cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSARKIMEX" 1080d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx) 1081d5e6173cSPeter Brune { 1082d5e6173cSPeter Brune PetscFunctionBegin; 1083d5e6173cSPeter Brune PetscFunctionReturn(0); 1084d5e6173cSPeter Brune } 1085d5e6173cSPeter Brune 1086d5e6173cSPeter Brune #undef __FUNCT__ 1087d5e6173cSPeter Brune #define __FUNCT__ "DMRestrictHook_TSARKIMEX" 1088d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 1089d5e6173cSPeter Brune { 1090d5e6173cSPeter Brune TS ts = (TS)ctx; 1091d5e6173cSPeter Brune PetscErrorCode ierr; 1092d5e6173cSPeter Brune Vec Z,Z_c; 1093d5e6173cSPeter Brune 1094d5e6173cSPeter Brune PetscFunctionBegin; 10950298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 10960298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 1097d5e6173cSPeter Brune ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 1098d5e6173cSPeter Brune ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 10990298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr); 11000298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr); 11018a381b04SJed Brown PetscFunctionReturn(0); 11028a381b04SJed Brown } 11038a381b04SJed Brown 1104cdb298fcSPeter Brune 1105cdb298fcSPeter Brune #undef __FUNCT__ 1106cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainHook_TSARKIMEX" 1107cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx) 1108cdb298fcSPeter Brune { 1109cdb298fcSPeter Brune PetscFunctionBegin; 1110cdb298fcSPeter Brune PetscFunctionReturn(0); 1111cdb298fcSPeter Brune } 1112cdb298fcSPeter Brune 1113cdb298fcSPeter Brune #undef __FUNCT__ 1114cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX" 1115cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 1116cdb298fcSPeter Brune { 1117cdb298fcSPeter Brune TS ts = (TS)ctx; 1118cdb298fcSPeter Brune PetscErrorCode ierr; 1119cdb298fcSPeter Brune Vec Z,Z_c; 1120cdb298fcSPeter Brune 1121cdb298fcSPeter Brune PetscFunctionBegin; 11220298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11230298fd71SBarry Smith ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1124cdb298fcSPeter Brune 1125cdb298fcSPeter Brune ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1126cdb298fcSPeter Brune ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1127cdb298fcSPeter Brune 11280298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr); 11290298fd71SBarry Smith ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr); 1130cdb298fcSPeter Brune PetscFunctionReturn(0); 1131cdb298fcSPeter Brune } 1132cdb298fcSPeter Brune 11338a381b04SJed Brown #undef __FUNCT__ 113496400bd6SLisandro Dalcin #define __FUNCT__ "TSARKIMEXTableauSetUp" 113596400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauSetUp(TS ts) 113696400bd6SLisandro Dalcin { 113796400bd6SLisandro Dalcin TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 113896400bd6SLisandro Dalcin ARKTableau tab = ark->tableau; 113996400bd6SLisandro Dalcin PetscErrorCode ierr; 114096400bd6SLisandro Dalcin 114196400bd6SLisandro Dalcin PetscFunctionBegin; 114296400bd6SLisandro Dalcin ierr = PetscMalloc1(tab->s,&ark->work);CHKERRQ(ierr); 114396400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);CHKERRQ(ierr); 114496400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);CHKERRQ(ierr); 114596400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);CHKERRQ(ierr); 114696400bd6SLisandro Dalcin if (ark->extrapolate) { 114796400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr); 114896400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr); 114996400bd6SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr); 115096400bd6SLisandro Dalcin } 115196400bd6SLisandro Dalcin PetscFunctionReturn(0); 115296400bd6SLisandro Dalcin } 115396400bd6SLisandro Dalcin 115496400bd6SLisandro Dalcin #undef __FUNCT__ 11558a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX" 11568a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 11578a381b04SJed Brown { 11588a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11598a381b04SJed Brown PetscErrorCode ierr; 1160d5e6173cSPeter Brune DM dm; 116196400bd6SLisandro Dalcin SNES snes; 1162f9c1d6abSBarry Smith 11638a381b04SJed Brown PetscFunctionBegin; 116496400bd6SLisandro Dalcin ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr); 11658a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 1166e817cc15SEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr); 11678a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 1168d5e6173cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1169d5e6173cSPeter Brune if (dm) { 1170d5e6173cSPeter Brune ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1171cdb298fcSPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr); 1172d5e6173cSPeter Brune } 117396400bd6SLisandro Dalcin ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 11748a381b04SJed Brown PetscFunctionReturn(0); 11758a381b04SJed Brown } 11768a381b04SJed Brown /*------------------------------------------------------------*/ 11778a381b04SJed Brown 11788a381b04SJed Brown #undef __FUNCT__ 11798a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 11804416b707SBarry Smith static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 11818a381b04SJed Brown { 11824cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 11838a381b04SJed Brown PetscErrorCode ierr; 11848a381b04SJed Brown 11858a381b04SJed Brown PetscFunctionBegin; 1186e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");CHKERRQ(ierr); 11878a381b04SJed Brown { 11888a381b04SJed Brown ARKTableauLink link; 11898a381b04SJed Brown PetscInt count,choice; 11908a381b04SJed Brown PetscBool flg; 11918a381b04SJed Brown const char **namelist; 11928a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 1193785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 11948a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 119596400bd6SLisandro Dalcin ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);CHKERRQ(ierr); 119696400bd6SLisandro Dalcin if (flg) {ierr = TSARKIMEXSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11978a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 119896400bd6SLisandro Dalcin 11994cc180ffSJed Brown flg = (PetscBool) !ark->imex; 12000298fd71SBarry Smith ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr); 12014cc180ffSJed Brown ark->imex = (PetscBool) !flg; 120203842d09SLisandro 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); 12038a381b04SJed Brown } 12048a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 12058a381b04SJed Brown PetscFunctionReturn(0); 12068a381b04SJed Brown } 12078a381b04SJed Brown 12088a381b04SJed Brown #undef __FUNCT__ 12098a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray" 12108a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 12118a381b04SJed Brown { 1212257d2499SJed Brown PetscErrorCode ierr; 1213f1d86077SJed Brown PetscInt i; 1214f1d86077SJed Brown size_t left,count; 12158a381b04SJed Brown char *p; 12168a381b04SJed Brown 12178a381b04SJed Brown PetscFunctionBegin; 1218f1d86077SJed Brown for (i=0,p=buf,left=len; i<n; i++) { 1219*da649d3eSBarry Smith ierr = PetscSNPrintfCount(p,left,fmt,&count,(double)x[i]);CHKERRQ(ierr); 12208a381b04SJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 12218a381b04SJed Brown left -= count; 12228a381b04SJed Brown p += count; 12238a381b04SJed Brown *p++ = ' '; 12248a381b04SJed Brown } 12258a381b04SJed Brown p[i ? 0 : -1] = 0; 12268a381b04SJed Brown PetscFunctionReturn(0); 12278a381b04SJed Brown } 12288a381b04SJed Brown 12298a381b04SJed Brown #undef __FUNCT__ 12308a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX" 12318a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 12328a381b04SJed Brown { 12338a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 12348a381b04SJed Brown PetscBool iascii; 12358a381b04SJed Brown PetscErrorCode ierr; 12368a381b04SJed Brown 12378a381b04SJed Brown PetscFunctionBegin; 1238251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 12398a381b04SJed Brown if (iascii) { 12409c334d8fSLisandro Dalcin ARKTableau tab = ark->tableau; 124119fd82e9SBarry Smith TSARKIMEXType arktype; 12428a381b04SJed Brown char buf[512]; 12438a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 12448a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 12458caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 124631f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 12478caf3d72SBarry Smith ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 1248e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 1249e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr); 1250e817cc15SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr); 125131f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 12528a381b04SJed Brown } 1253be5899b3SLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 125496400bd6SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 12558a381b04SJed Brown PetscFunctionReturn(0); 12568a381b04SJed Brown } 12578a381b04SJed Brown 12588a381b04SJed Brown #undef __FUNCT__ 1259f2c2a1b9SBarry Smith #define __FUNCT__ "TSLoad_ARKIMEX" 1260f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer) 1261f2c2a1b9SBarry Smith { 1262f2c2a1b9SBarry Smith PetscErrorCode ierr; 1263f2c2a1b9SBarry Smith SNES snes; 12649c334d8fSLisandro Dalcin TSAdapt adapt; 1265f2c2a1b9SBarry Smith 1266f2c2a1b9SBarry Smith PetscFunctionBegin; 12679c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12689c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 1269f2c2a1b9SBarry Smith ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1270f2c2a1b9SBarry Smith ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 1271ad6bc421SBarry Smith /* function and Jacobian context for SNES when used with TS is always ts object */ 12720298fd71SBarry Smith ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 12730298fd71SBarry Smith ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 1274f2c2a1b9SBarry Smith PetscFunctionReturn(0); 1275f2c2a1b9SBarry Smith } 1276f2c2a1b9SBarry Smith 1277f2c2a1b9SBarry Smith #undef __FUNCT__ 12788a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType" 12798a381b04SJed Brown /*@C 12808a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 12818a381b04SJed Brown 12828a381b04SJed Brown Logically collective 12838a381b04SJed Brown 12848a381b04SJed Brown Input Parameter: 12858a381b04SJed Brown + ts - timestepping context 12868a381b04SJed Brown - arktype - type of ARK-IMEX scheme 12878a381b04SJed Brown 12888a381b04SJed Brown Level: intermediate 12898a381b04SJed Brown 1290020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 12918a381b04SJed Brown @*/ 129219fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype) 12938a381b04SJed Brown { 12948a381b04SJed Brown PetscErrorCode ierr; 12958a381b04SJed Brown 12968a381b04SJed Brown PetscFunctionBegin; 12978a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 129819fd82e9SBarry Smith ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 12998a381b04SJed Brown PetscFunctionReturn(0); 13008a381b04SJed Brown } 13018a381b04SJed Brown 13028a381b04SJed Brown #undef __FUNCT__ 13038a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType" 13048a381b04SJed Brown /*@C 13058a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 13068a381b04SJed Brown 13078a381b04SJed Brown Logically collective 13088a381b04SJed Brown 13098a381b04SJed Brown Input Parameter: 13108a381b04SJed Brown . ts - timestepping context 13118a381b04SJed Brown 13128a381b04SJed Brown Output Parameter: 13138a381b04SJed Brown . arktype - type of ARK-IMEX scheme 13148a381b04SJed Brown 13158a381b04SJed Brown Level: intermediate 13168a381b04SJed Brown 13178a381b04SJed Brown .seealso: TSARKIMEXGetType() 13188a381b04SJed Brown @*/ 131919fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype) 13208a381b04SJed Brown { 13218a381b04SJed Brown PetscErrorCode ierr; 13228a381b04SJed Brown 13238a381b04SJed Brown PetscFunctionBegin; 13248a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 132519fd82e9SBarry Smith ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 13268a381b04SJed Brown PetscFunctionReturn(0); 13278a381b04SJed Brown } 13288a381b04SJed Brown 13294cc180ffSJed Brown #undef __FUNCT__ 13304cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit" 133116353aafSBarry Smith /*@ 13324cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 13334cc180ffSJed Brown 13344cc180ffSJed Brown Logically collective 13354cc180ffSJed Brown 13364cc180ffSJed Brown Input Parameter: 13374cc180ffSJed Brown + ts - timestepping context 13384cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 13394cc180ffSJed Brown 13404cc180ffSJed Brown Level: intermediate 13414cc180ffSJed Brown 13424cc180ffSJed Brown .seealso: TSARKIMEXGetType() 13434cc180ffSJed Brown @*/ 13444cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 13454cc180ffSJed Brown { 13464cc180ffSJed Brown PetscErrorCode ierr; 13474cc180ffSJed Brown 13484cc180ffSJed Brown PetscFunctionBegin; 13494cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 13504cc180ffSJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 13514cc180ffSJed Brown PetscFunctionReturn(0); 13524cc180ffSJed Brown } 13534cc180ffSJed Brown 13548a381b04SJed Brown #undef __FUNCT__ 13558a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 1356e0877f53SBarry Smith static PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype) 13578a381b04SJed Brown { 13588a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13598a381b04SJed Brown 13608a381b04SJed Brown PetscFunctionBegin; 13618a381b04SJed Brown *arktype = ark->tableau->name; 13628a381b04SJed Brown PetscFunctionReturn(0); 13638a381b04SJed Brown } 13648a381b04SJed Brown #undef __FUNCT__ 13658a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 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);} 13848a381b04SJed Brown PetscFunctionReturn(0); 13858a381b04SJed Brown } 13868a381b04SJed Brown } 1387ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 13888a381b04SJed Brown PetscFunctionReturn(0); 13898a381b04SJed Brown } 1390e0877f53SBarry Smith 13914cc180ffSJed Brown #undef __FUNCT__ 13924cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX" 1393e0877f53SBarry Smith static PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 13944cc180ffSJed Brown { 13954cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 13964cc180ffSJed Brown 13974cc180ffSJed Brown PetscFunctionBegin; 13984cc180ffSJed Brown ark->imex = (PetscBool)!flg; 13994cc180ffSJed Brown PetscFunctionReturn(0); 14004cc180ffSJed Brown } 14018a381b04SJed Brown 14028a381b04SJed Brown /* ------------------------------------------------------------ */ 14038a381b04SJed Brown /*MC 1404a4386c9eSJed Brown TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes 14058a381b04SJed Brown 1406fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1407fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1408fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1409fca742c7SJed Brown 1410fca742c7SJed Brown Notes: 1411a4386c9eSJed Brown The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1412c8058688SBarry Smith 14135eca1a21SEmil Constantinescu If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information. 14145eca1a21SEmil Constantinescu 1415a4386c9eSJed 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). 1416fca742c7SJed Brown 1417d0685a90SJed Brown Consider trying TSROSW if the stiff part is linear or weakly nonlinear. 1418d0685a90SJed Brown 14198a381b04SJed Brown Level: beginner 14208a381b04SJed Brown 1421d0685a90SJed Brown .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX1BEE, 1422d0685a90SJed Brown TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122, 1423d0685a90SJed Brown TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister() 14248a381b04SJed Brown 14258a381b04SJed Brown M*/ 14268a381b04SJed Brown #undef __FUNCT__ 14278a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX" 14288cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts) 14298a381b04SJed Brown { 14308a381b04SJed Brown TS_ARKIMEX *th; 14318a381b04SJed Brown PetscErrorCode ierr; 14328a381b04SJed Brown 14338a381b04SJed Brown PetscFunctionBegin; 1434607a6623SBarry Smith ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr); 14358a381b04SJed Brown 14368a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 14378a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 14388a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 1439f2c2a1b9SBarry Smith ts->ops->load = TSLoad_ARKIMEX; 14408a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 14418a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1442cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1443108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 144424655328SShri ts->ops->rollback = TSRollBack_ARKIMEX; 14458a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 14468a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 14478a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 14488a381b04SJed Brown 1449b00a9115SJed Brown ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 14508a381b04SJed Brown ts->data = (void*)th; 14514cc180ffSJed Brown th->imex = PETSC_TRUE; 14528a381b04SJed Brown 1453bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 1454bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 1455bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 145696400bd6SLisandro Dalcin 145796400bd6SLisandro Dalcin ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 14588a381b04SJed Brown PetscFunctionReturn(0); 14598a381b04SJed Brown } 1460