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 78a381b04SJed Brown F(t,X,Xdot) = G(t,X) 88a381b04SJed Brown 98a381b04SJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 108a381b04SJed Brown 118a381b04SJed Brown */ 12*b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 138a381b04SJed Brown 14108c343cSJed Brown static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3; 158a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled; 168a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized; 178a381b04SJed Brown 188a381b04SJed Brown typedef struct _ARKTableau *ARKTableau; 198a381b04SJed Brown struct _ARKTableau { 208a381b04SJed Brown char *name; 214f385281SJed Brown PetscInt order; /* Classical approximation order of the method */ 224f385281SJed Brown PetscInt s; /* Number of stages */ 234f385281SJed Brown PetscInt pinterp; /* Interpolation order */ 244f385281SJed Brown PetscReal *At,*bt,*ct; /* Stiff tableau */ 258a381b04SJed Brown PetscReal *A,*b,*c; /* Non-stiff tableau */ 26108c343cSJed Brown PetscReal *bembedt,*bembed; /* Embedded formula of order one less (order-1) */ 27cd652676SJed Brown PetscReal *binterpt,*binterp; /* Dense output formula */ 28108c343cSJed Brown PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 298a381b04SJed Brown }; 308a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 318a381b04SJed Brown struct _ARKTableauLink { 328a381b04SJed Brown struct _ARKTableau tab; 338a381b04SJed Brown ARKTableauLink next; 348a381b04SJed Brown }; 358a381b04SJed Brown static ARKTableauLink ARKTableauList; 368a381b04SJed Brown 378a381b04SJed Brown typedef struct { 388a381b04SJed Brown ARKTableau tableau; 398a381b04SJed Brown Vec *Y; /* States computed during the step */ 408a381b04SJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 418a381b04SJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 428a381b04SJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 438a381b04SJed Brown Vec Work; /* Generic work vector */ 448a381b04SJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 458a381b04SJed Brown PetscScalar *work; /* Scalar work */ 468a381b04SJed Brown PetscReal shift; 478a381b04SJed Brown PetscReal stage_time; 484cc180ffSJed Brown PetscBool imex; 49108c343cSJed Brown TSStepStatus status; 508a381b04SJed Brown } TS_ARKIMEX; 511f80e275SEmil Constantinescu /*MC 521f80e275SEmil Constantinescu TSARKIMEXARS122 - Second order ARK IMEX scheme. 538a381b04SJed Brown 541f80e275SEmil Constantinescu This method has one explicit stage and one implicit stage. 551f80e275SEmil Constantinescu 561f80e275SEmil Constantinescu References: 571f80e275SEmil Constantinescu U. Ascher, S. Ruuth, R. J. Spitheri, Implicit-explicit Runge-Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997), pp. 151–167. 581f80e275SEmil Constantinescu 591f80e275SEmil Constantinescu Level: advanced 601f80e275SEmil Constantinescu 611f80e275SEmil Constantinescu .seealso: TSARKIMEX 621f80e275SEmil Constantinescu M*/ 631f80e275SEmil Constantinescu /*MC 641f80e275SEmil Constantinescu TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part. 651f80e275SEmil Constantinescu 661f80e275SEmil 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. 671f80e275SEmil Constantinescu 681f80e275SEmil Constantinescu Level: advanced 691f80e275SEmil Constantinescu 701f80e275SEmil Constantinescu .seealso: TSARKIMEX 711f80e275SEmil Constantinescu M*/ 721f80e275SEmil Constantinescu /*MC 731f80e275SEmil Constantinescu TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part. 741f80e275SEmil Constantinescu 751f80e275SEmil Constantinescu This method has two implicit stages, and L-stable implicit scheme. 761f80e275SEmil Constantinescu 771f80e275SEmil Constantinescu References: 781f80e275SEmil Constantinescu 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, pp. 129-155 791f80e275SEmil Constantinescu 801f80e275SEmil Constantinescu Level: advanced 811f80e275SEmil Constantinescu 821f80e275SEmil Constantinescu .seealso: TSARKIMEX 831f80e275SEmil Constantinescu M*/ 841f80e275SEmil Constantinescu /*MC 851f80e275SEmil Constantinescu TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part. 861f80e275SEmil Constantinescu 871f80e275SEmil 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. 881f80e275SEmil Constantinescu 891f80e275SEmil Constantinescu Level: advanced 901f80e275SEmil Constantinescu 911f80e275SEmil Constantinescu .seealso: TSARKIMEX 921f80e275SEmil Constantinescu M*/ 9364f491ddSJed Brown /*MC 9464f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 9564f491ddSJed Brown 961f80e275SEmil Constantinescu This method has one explicit stage and two implicit stages. This method was provided by Emil Constantinescu. 9764f491ddSJed Brown 98b330ce4dSSatish Balay Level: advanced 99b330ce4dSSatish Balay 10064f491ddSJed Brown .seealso: TSARKIMEX 10164f491ddSJed Brown M*/ 10264f491ddSJed Brown /*MC 10364f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 10464f491ddSJed Brown 10564f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 10664f491ddSJed Brown 107b330ce4dSSatish Balay Level: advanced 108b330ce4dSSatish Balay 10964f491ddSJed Brown .seealso: TSARKIMEX 11064f491ddSJed Brown M*/ 11164f491ddSJed Brown /*MC 1126cf0794eSJed Brown TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme. 1136cf0794eSJed Brown 1146cf0794eSJed Brown This method has three implicit stages. 1156cf0794eSJed Brown 1166cf0794eSJed Brown References: 1176cf0794eSJed Brown 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, pp. 129-155 1186cf0794eSJed Brown 1196cf0794eSJed Brown This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375 1206cf0794eSJed Brown 1216cf0794eSJed Brown Level: advanced 1226cf0794eSJed Brown 1236cf0794eSJed Brown .seealso: TSARKIMEX 1246cf0794eSJed Brown M*/ 1256cf0794eSJed Brown /*MC 12664f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 12764f491ddSJed Brown 12864f491ddSJed Brown This method has one explicit stage and three implicit stages. 12964f491ddSJed Brown 13064f491ddSJed Brown References: 13164f491ddSJed Brown Kennedy and Carpenter 2003. 13264f491ddSJed Brown 133b330ce4dSSatish Balay Level: advanced 134b330ce4dSSatish Balay 13564f491ddSJed Brown .seealso: TSARKIMEX 13664f491ddSJed Brown M*/ 13764f491ddSJed Brown /*MC 1386cf0794eSJed Brown TSARKIMEXARS443 - Third order ARK IMEX scheme. 1396cf0794eSJed Brown 1406cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1416cf0794eSJed Brown 1426cf0794eSJed Brown References: 1436cf0794eSJed Brown U. Ascher, S. Ruuth, R. J. Spitheri, Implicit-explicit Runge-Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997), pp. 151–167. 1446cf0794eSJed Brown 1456cf0794eSJed Brown This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375 1466cf0794eSJed Brown 1476cf0794eSJed Brown Level: advanced 1486cf0794eSJed Brown 1496cf0794eSJed Brown .seealso: TSARKIMEX 1506cf0794eSJed Brown M*/ 1516cf0794eSJed Brown /*MC 1526cf0794eSJed Brown TSARKIMEXBPR3 - Third order ARK IMEX scheme. 1536cf0794eSJed Brown 1546cf0794eSJed Brown This method has one explicit stage and four implicit stages. 1556cf0794eSJed Brown 1566cf0794eSJed Brown References: 1576cf0794eSJed Brown This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375 1586cf0794eSJed Brown 1596cf0794eSJed Brown Level: advanced 1606cf0794eSJed Brown 1616cf0794eSJed Brown .seealso: TSARKIMEX 1626cf0794eSJed Brown M*/ 1636cf0794eSJed Brown /*MC 16464f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 16564f491ddSJed Brown 16664f491ddSJed Brown This method has one explicit stage and four implicit stages. 16764f491ddSJed Brown 16864f491ddSJed Brown References: 16964f491ddSJed Brown Kennedy and Carpenter 2003. 17064f491ddSJed Brown 171b330ce4dSSatish Balay Level: advanced 172b330ce4dSSatish Balay 17364f491ddSJed Brown .seealso: TSARKIMEX 17464f491ddSJed Brown M*/ 17564f491ddSJed Brown /*MC 17664f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 17764f491ddSJed Brown 17864f491ddSJed Brown This method has one explicit stage and five implicit stages. 17964f491ddSJed Brown 18064f491ddSJed Brown References: 18164f491ddSJed Brown Kennedy and Carpenter 2003. 18264f491ddSJed Brown 183b330ce4dSSatish Balay Level: advanced 184b330ce4dSSatish Balay 18564f491ddSJed Brown .seealso: TSARKIMEX 18664f491ddSJed Brown M*/ 18764f491ddSJed Brown 1888a381b04SJed Brown #undef __FUNCT__ 1898a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll" 1908a381b04SJed Brown /*@C 1918a381b04SJed Brown TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 1928a381b04SJed Brown 193fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 1948a381b04SJed Brown 1958a381b04SJed Brown Level: advanced 1968a381b04SJed Brown 1978a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all 1988a381b04SJed Brown 1998a381b04SJed Brown .seealso: TSARKIMEXRegisterDestroy() 2008a381b04SJed Brown @*/ 2018a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 2028a381b04SJed Brown { 2038a381b04SJed Brown PetscErrorCode ierr; 2048a381b04SJed Brown 2058a381b04SJed Brown PetscFunctionBegin; 2068a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 2078a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 2088a381b04SJed Brown { 2098a381b04SJed Brown const PetscReal 2101f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2111f80e275SEmil Constantinescu {0.5,0.0}}, 2121f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2131f80e275SEmil Constantinescu {0.0,0.5}}, 2141f80e275SEmil Constantinescu b[2] = {0.0,1.0}, 2151f80e275SEmil Constantinescu bembedt[2] = {0.5,0.5}; 2161f80e275SEmil 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*/ 2171f80e275SEmil Constantinescu ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr); 2181f80e275SEmil Constantinescu } 2191f80e275SEmil Constantinescu { 2201f80e275SEmil Constantinescu const PetscReal 2211f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2221f80e275SEmil Constantinescu {1.0,0.0}}, 2231f80e275SEmil Constantinescu At[2][2] = {{0.0,0.0}, 2241f80e275SEmil Constantinescu {0.5,0.5}}, 2251f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2261f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}; 2271f80e275SEmil 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*/ 2281f80e275SEmil Constantinescu ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr); 2291f80e275SEmil Constantinescu } 2301f80e275SEmil Constantinescu { 2311f80e275SEmil Constantinescu const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0); 2321f80e275SEmil Constantinescu const PetscReal 2331f80e275SEmil Constantinescu A[2][2] = {{0.0,0.0}, 2341f80e275SEmil Constantinescu {1.0,0.0}}, 2351f80e275SEmil Constantinescu At[2][2] = {{us2,0.0}, 2361f80e275SEmil Constantinescu {1.0-2.0*us2,us2}}, 2371f80e275SEmil Constantinescu b[2] = {0.5,0.5}, 2381f80e275SEmil Constantinescu bembedt[2] = {0.0,1.0}, 2391f80e275SEmil Constantinescu binterpt[2][2] = {{(us2-1.0)/(2.0*us2-1.0),-1/(2.0*(1.0-2.0*us2))},{1-(us2-1.0)/(2.0*us2-1.0),-1/(2.0*(1.0-2.0*us2))}}, 2401f80e275SEmil Constantinescu binterp[2][2] = {{1.0,-0.5},{0.0,0.5}}; 2411f80e275SEmil Constantinescu ierr = TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);CHKERRQ(ierr); 2421f80e275SEmil Constantinescu } 2431f80e275SEmil Constantinescu { 2441f80e275SEmil Constantinescu const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 2458a381b04SJed Brown A[3][3] = {{0,0,0}, 2461f80e275SEmil Constantinescu {2-s2,0,0}, 2471f80e275SEmil Constantinescu {0.55,0.45,0}}, 2481f80e275SEmil Constantinescu At[3][3] = {{0,0,0}, 2491f80e275SEmil Constantinescu {1-1/s2,1-1/s2,0}, 2501f80e275SEmil Constantinescu {1/(2*s2),1/(2*s2),1-1/s2}}, 251e93ac1c2SEmil Constantinescu bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)}, 252ce4a059fSEmil Constantinescu binterpt[3][2] = {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}}; 2531f80e275SEmil Constantinescu ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 2541f80e275SEmil Constantinescu } 2551f80e275SEmil Constantinescu { 2561f80e275SEmil Constantinescu const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 2571f80e275SEmil Constantinescu A[3][3] = {{0,0,0}, 2581f80e275SEmil Constantinescu {2-s2,0,0}, 2598a381b04SJed Brown {0.75,0.25,0}}, 2608a381b04SJed Brown At[3][3] = {{0,0,0}, 2611f80e275SEmil Constantinescu {1-1/s2,1-1/s2,0}, 2621f80e275SEmil Constantinescu {1/(2*s2),1/(2*s2),1-1/s2}}, 263e93ac1c2SEmil Constantinescu bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)}, 264ce4a059fSEmil Constantinescu binterpt[3][2] = {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}}; 265108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 2668a381b04SJed Brown } 26706db7b1cSJed Brown { /* Optimal for linear implicit part */ 26803403c7fSJed Brown const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 269a3a57f36SJed Brown A[3][3] = {{0,0,0}, 270a3a57f36SJed Brown {2-s2,0,0}, 271a3a57f36SJed Brown {(3-2*s2)/6,(3+2*s2)/6,0}}, 272a3a57f36SJed Brown At[3][3] = {{0,0,0}, 273a3a57f36SJed Brown {1-1/s2,1-1/s2,0}, 274cd652676SJed Brown {1/(2*s2),1/(2*s2),1-1/s2}}, 275e93ac1c2SEmil Constantinescu bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)}, 276ce4a059fSEmil Constantinescu binterpt[3][2] = {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}}; 2771f80e275SEmil Constantinescu ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 278a3a57f36SJed Brown } 2796cf0794eSJed Brown { /* Optimal for linear implicit part */ 2806cf0794eSJed Brown const PetscReal 2816cf0794eSJed Brown A[3][3] = {{0,0,0}, 2826cf0794eSJed Brown {0.5,0,0}, 2836cf0794eSJed Brown {0.5,0.5,0}}, 2846cf0794eSJed Brown At[3][3] = {{0.25,0,0}, 2856cf0794eSJed Brown {0,0.25,0}, 2866cf0794eSJed Brown {1./3,1./3,1./3}}; 287108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2886cf0794eSJed Brown } 289a3a57f36SJed Brown { 290a3a57f36SJed Brown const PetscReal 291a3a57f36SJed Brown A[4][4] = {{0,0,0,0}, 2924040e9f2SJed Brown {1767732205903./2027836641118.,0,0,0}, 2934040e9f2SJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 2944040e9f2SJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 295a3a57f36SJed Brown At[4][4] = {{0,0,0,0}, 2964040e9f2SJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 2974040e9f2SJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 2984040e9f2SJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 299cc46b9d1SJed Brown bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.}, 3004040e9f2SJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 3014040e9f2SJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 3024040e9f2SJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 3034040e9f2SJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 304108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 305a3a57f36SJed Brown } 306a3a57f36SJed Brown { 307a3a57f36SJed Brown const PetscReal 308e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3096cf0794eSJed Brown {1./2,0,0,0,0}, 3106cf0794eSJed Brown {11./18,1./18,0,0,0}, 3116cf0794eSJed Brown {5./6,-5./6,.5,0,0}, 3126cf0794eSJed Brown {1./4,7./4,3./4,-7./4,0}}, 3136cf0794eSJed Brown At[5][5] = {{0,0,0,0,0}, 3146cf0794eSJed Brown {0,1./2,0,0,0}, 3156cf0794eSJed Brown {0,1./6,1./2,0,0}, 3166cf0794eSJed Brown {0,-1./2,1./2,1./2,0}, 317108c343cSJed Brown {0,3./2,-3./2,1./2,1./2}}, 318108c343cSJed Brown *bembedt = PETSC_NULL; 319108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 3206cf0794eSJed Brown } 3216cf0794eSJed Brown { 3226cf0794eSJed Brown const PetscReal 323e74514c0SSatish Balay A[5][5] = {{0,0,0,0,0}, 3246cf0794eSJed Brown {1,0,0,0,0}, 3256cf0794eSJed Brown {4./9,2./9,0,0,0}, 3266cf0794eSJed Brown {1./4,0,3./4,0,0}, 3276cf0794eSJed Brown {1./4,0,3./5,0,0}}, 328e74514c0SSatish Balay At[5][5] = {{0,0,0,0,0}, 3296cf0794eSJed Brown {.5,.5,0,0,0}, 3306cf0794eSJed Brown {5./18,-1./9,.5,0,0}, 3316cf0794eSJed Brown {.5,0,0,.5,0}, 332108c343cSJed Brown {.25,0,.75,-.5,.5}}, 333108c343cSJed Brown *bembedt = PETSC_NULL; 334108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 3356cf0794eSJed Brown } 3366cf0794eSJed Brown { 3376cf0794eSJed Brown const PetscReal 338a3a57f36SJed Brown A[6][6] = {{0,0,0,0,0,0}, 339a3a57f36SJed Brown {1./2,0,0,0,0,0}, 3404040e9f2SJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 3414040e9f2SJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 3424040e9f2SJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 3434040e9f2SJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 344a3a57f36SJed Brown At[6][6] = {{0,0,0,0,0,0}, 345a3a57f36SJed Brown {1./4,1./4,0,0,0,0}, 3464040e9f2SJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 3474040e9f2SJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 3484040e9f2SJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 3494040e9f2SJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 350cc46b9d1SJed Brown bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.}, 3514040e9f2SJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 352cd652676SJed Brown {0,0,0}, 3534040e9f2SJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 3544040e9f2SJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 3554040e9f2SJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 3564040e9f2SJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 357108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 358a3a57f36SJed Brown } 359a3a57f36SJed Brown { 360a3a57f36SJed Brown const PetscReal 361a3a57f36SJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 362a3a57f36SJed Brown {41./100,0,0,0,0,0,0,0}, 3634040e9f2SJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 3644040e9f2SJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 3654040e9f2SJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 3664040e9f2SJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 3674040e9f2SJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 3684040e9f2SJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 369a3a57f36SJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 3704040e9f2SJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 3714040e9f2SJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 3724040e9f2SJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 3734040e9f2SJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 3744040e9f2SJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 3754040e9f2SJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 3764040e9f2SJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 377cc46b9d1SJed Brown bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.}, 3784040e9f2SJed Brown binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.}, 379cd652676SJed Brown {0 , 0 , 0 }, 380cd652676SJed Brown {0 , 0 , 0 }, 3814040e9f2SJed Brown {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.}, 3824040e9f2SJed Brown {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.}, 3834040e9f2SJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.}, 3844040e9f2SJed Brown {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.}, 3854040e9f2SJed Brown {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }}; 386108c343cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 387a3a57f36SJed Brown } 388a3a57f36SJed Brown 3898a381b04SJed Brown PetscFunctionReturn(0); 3908a381b04SJed Brown } 3918a381b04SJed Brown 3928a381b04SJed Brown #undef __FUNCT__ 3938a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy" 3948a381b04SJed Brown /*@C 3958a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 3968a381b04SJed Brown 3978a381b04SJed Brown Not Collective 3988a381b04SJed Brown 3998a381b04SJed Brown Level: advanced 4008a381b04SJed Brown 4018a381b04SJed Brown .keywords: TSARKIMEX, register, destroy 4028a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic() 4038a381b04SJed Brown @*/ 4048a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 4058a381b04SJed Brown { 4068a381b04SJed Brown PetscErrorCode ierr; 4078a381b04SJed Brown ARKTableauLink link; 4088a381b04SJed Brown 4098a381b04SJed Brown PetscFunctionBegin; 4108a381b04SJed Brown while ((link = ARKTableauList)) { 4118a381b04SJed Brown ARKTableau t = &link->tab; 4128a381b04SJed Brown ARKTableauList = link->next; 4138a381b04SJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 414108c343cSJed Brown ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr); 415cd652676SJed Brown ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 4168a381b04SJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 4178a381b04SJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 4188a381b04SJed Brown } 4198a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 4208a381b04SJed Brown PetscFunctionReturn(0); 4218a381b04SJed Brown } 4228a381b04SJed Brown 4238a381b04SJed Brown #undef __FUNCT__ 4248a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage" 4258a381b04SJed Brown /*@C 4268a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 4278a381b04SJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 4288a381b04SJed Brown when using static libraries. 4298a381b04SJed Brown 4308a381b04SJed Brown Input Parameter: 4318a381b04SJed Brown path - The dynamic library path, or PETSC_NULL 4328a381b04SJed Brown 4338a381b04SJed Brown Level: developer 4348a381b04SJed Brown 4358a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package 4368a381b04SJed Brown .seealso: PetscInitialize() 4378a381b04SJed Brown @*/ 4388a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[]) 4398a381b04SJed Brown { 4408a381b04SJed Brown PetscErrorCode ierr; 4418a381b04SJed Brown 4428a381b04SJed Brown PetscFunctionBegin; 4438a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 4448a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 4458a381b04SJed Brown ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 4468a381b04SJed Brown ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 4478a381b04SJed Brown PetscFunctionReturn(0); 4488a381b04SJed Brown } 4498a381b04SJed Brown 4508a381b04SJed Brown #undef __FUNCT__ 4518a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage" 4528a381b04SJed Brown /*@C 4538a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 4548a381b04SJed Brown called from PetscFinalize(). 4558a381b04SJed Brown 4568a381b04SJed Brown Level: developer 4578a381b04SJed Brown 4588a381b04SJed Brown .keywords: Petsc, destroy, package 4598a381b04SJed Brown .seealso: PetscFinalize() 4608a381b04SJed Brown @*/ 4618a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 4628a381b04SJed Brown { 4638a381b04SJed Brown PetscErrorCode ierr; 4648a381b04SJed Brown 4658a381b04SJed Brown PetscFunctionBegin; 4668a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 4678a381b04SJed Brown ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 4688a381b04SJed Brown PetscFunctionReturn(0); 4698a381b04SJed Brown } 4708a381b04SJed Brown 4718a381b04SJed Brown #undef __FUNCT__ 4728a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister" 473cd652676SJed Brown /*@C 474cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 475cd652676SJed Brown 476cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 477cd652676SJed Brown 478cd652676SJed Brown Input Parameters: 479cd652676SJed Brown + name - identifier for method 480cd652676SJed Brown . order - approximation order of method 481cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 482cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 483cd652676SJed Brown . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At) 484cd652676SJed Brown . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At) 485cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 486cd652676SJed Brown . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At) 487cd652676SJed Brown . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A) 488108c343cSJed Brown . bembedt - Stiff part of completion table for embedded method (dimension s; PETSC_NULL if not available) 489108c343cSJed Brown . bembed - Non-stiff part of completion table for embedded method (dimension s; PETSC_NULL to use bembedt if provided) 490cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 491cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 492cd652676SJed Brown - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt) 493cd652676SJed Brown 494cd652676SJed Brown Notes: 495cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 496cd652676SJed Brown 497cd652676SJed Brown Level: advanced 498cd652676SJed Brown 499cd652676SJed Brown .keywords: TS, register 500cd652676SJed Brown 501cd652676SJed Brown .seealso: TSARKIMEX 502cd652676SJed Brown @*/ 5038a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s, 5048a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 505cd652676SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 506108c343cSJed Brown const PetscReal bembedt[],const PetscReal bembed[], 507cd652676SJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 5088a381b04SJed Brown { 5098a381b04SJed Brown PetscErrorCode ierr; 5108a381b04SJed Brown ARKTableauLink link; 5118a381b04SJed Brown ARKTableau t; 5128a381b04SJed Brown PetscInt i,j; 5138a381b04SJed Brown 5148a381b04SJed Brown PetscFunctionBegin; 5158a381b04SJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 516cd652676SJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 5178a381b04SJed Brown t = &link->tab; 5188a381b04SJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 5198a381b04SJed Brown t->order = order; 5208a381b04SJed Brown t->s = s; 5218a381b04SJed Brown ierr = PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 5228a381b04SJed Brown ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 5238a381b04SJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 5248a381b04SJed Brown if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);} 5258a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 5268a381b04SJed Brown if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);} 5278a381b04SJed Brown else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i]; 5288a381b04SJed Brown if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);} 5298a381b04SJed 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]; 5308a381b04SJed Brown if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);} 5318a381b04SJed 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]; 532108c343cSJed Brown if (bembedt) { 533108c343cSJed Brown ierr = PetscMalloc2(s,PetscReal,&t->bembedt,s,PetscReal,&t->bembed);CHKERRQ(ierr); 534108c343cSJed Brown ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr); 535108c343cSJed Brown ierr = PetscMemcpy(t->bembed,bembed?bembed:bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr); 536108c343cSJed Brown } 537108c343cSJed Brown 5384f385281SJed Brown t->pinterp = pinterp; 539cd652676SJed Brown ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr); 540cd652676SJed Brown ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 541cd652676SJed Brown ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 5428a381b04SJed Brown link->next = ARKTableauList; 5438a381b04SJed Brown ARKTableauList = link; 5448a381b04SJed Brown PetscFunctionReturn(0); 5458a381b04SJed Brown } 5468a381b04SJed Brown 5478a381b04SJed Brown #undef __FUNCT__ 548108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX" 549108c343cSJed Brown /* 550108c343cSJed Brown The step completion formula is 551108c343cSJed Brown 552108c343cSJed Brown x1 = x0 - h bt^T YdotI + h b^T YdotRHS 553108c343cSJed Brown 554108c343cSJed Brown This function can be called before or after ts->vec_sol has been updated. 555108c343cSJed Brown Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order. 556108c343cSJed Brown We can write 557108c343cSJed Brown 558108c343cSJed Brown x1e = x0 - h bet^T YdotI + h be^T YdotRHS 559108c343cSJed Brown = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS 560108c343cSJed Brown = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS 561108c343cSJed Brown 562108c343cSJed Brown so we can evaluate the method with different order even after the step has been optimistically completed. 563108c343cSJed Brown */ 564108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 565108c343cSJed Brown { 566108c343cSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 567108c343cSJed Brown ARKTableau tab = ark->tableau; 568108c343cSJed Brown PetscScalar *w = ark->work; 569108c343cSJed Brown PetscReal h; 570108c343cSJed Brown PetscInt s = tab->s,j; 571108c343cSJed Brown PetscErrorCode ierr; 572108c343cSJed Brown 573108c343cSJed Brown PetscFunctionBegin; 574108c343cSJed Brown switch (ark->status) { 575108c343cSJed Brown case TS_STEP_INCOMPLETE: 576108c343cSJed Brown case TS_STEP_PENDING: 577108c343cSJed Brown h = ts->time_step; break; 578108c343cSJed Brown case TS_STEP_COMPLETE: 579108c343cSJed Brown h = ts->time_step_prev; break; 580b9ce6d65SJed Brown default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus"); 581108c343cSJed Brown } 582108c343cSJed Brown if (order == tab->order) { 583108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Use the standard completion formula (bt,b) */ 584108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 585108c343cSJed Brown for (j=0; j<s; j++) w[j] = -h*tab->bt[j]; 586108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 587108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->b[j]; 588108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 589108c343cSJed Brown } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 590108c343cSJed Brown if (done) *done = PETSC_TRUE; 591108c343cSJed Brown PetscFunctionReturn(0); 592108c343cSJed Brown } else if (order == tab->order-1) { 593108c343cSJed Brown if (!tab->bembedt) goto unavailable; 594108c343cSJed Brown if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */ 595108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 596108c343cSJed Brown for (j=0; j<s; j++) w[j] = -h*tab->bembedt[j]; 597108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr); 598108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 599108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 600108c343cSJed Brown } else { /* Rollback and re-complete using (bet-be,be-b) */ 601108c343cSJed Brown ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 602108c343cSJed Brown for (j=0; j<s; j++) w[j] = -h*(tab->bembedt[j] - tab->bt[j]); 603108c343cSJed Brown ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr); 604108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 605108c343cSJed Brown ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr); 606108c343cSJed Brown } 607108c343cSJed Brown if (done) *done = PETSC_TRUE; 608108c343cSJed Brown PetscFunctionReturn(0); 609108c343cSJed Brown } 610108c343cSJed Brown unavailable: 611108c343cSJed Brown if (done) *done = PETSC_FALSE; 612108c343cSJed Brown else SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 613108c343cSJed Brown PetscFunctionReturn(0); 614108c343cSJed Brown } 615108c343cSJed Brown 616108c343cSJed Brown #undef __FUNCT__ 6178a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX" 6188a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 6198a381b04SJed Brown { 6208a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 6218a381b04SJed Brown ARKTableau tab = ark->tableau; 6228a381b04SJed Brown const PetscInt s = tab->s; 6238a381b04SJed Brown const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c; 624406d0ec2SJed Brown PetscScalar *w = ark->work; 6258a381b04SJed Brown Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z; 626108c343cSJed Brown TSAdapt adapt; 6278a381b04SJed Brown SNES snes; 628108c343cSJed Brown PetscInt i,j,its,lits,reject,next_scheme; 629cdbf8f93SLisandro Dalcin PetscReal next_time_step; 630108c343cSJed Brown PetscReal t; 631108c343cSJed Brown PetscBool accept; 6328a381b04SJed Brown PetscErrorCode ierr; 6338a381b04SJed Brown 6348a381b04SJed Brown PetscFunctionBegin; 6358a381b04SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 636cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 6378a381b04SJed Brown t = ts->ptime; 638108c343cSJed Brown accept = PETSC_TRUE; 639108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 6408a381b04SJed Brown 64197335746SJed Brown for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 642108c343cSJed Brown PetscReal h = ts->time_step; 6438a381b04SJed Brown for (i=0; i<s; i++) { 6448a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 6458a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 6468a381b04SJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 6478a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 6488a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 6498a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 6508a381b04SJed Brown } else { 6518a381b04SJed Brown ark->stage_time = t + h*ct[i]; 6528a381b04SJed Brown ark->shift = 1./(h*At[i*s+i]); 6538a381b04SJed Brown /* Affine part */ 6548a381b04SJed Brown ierr = VecZeroEntries(W);CHKERRQ(ierr); 6558a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 6568a381b04SJed Brown ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 657f16577ceSEmil Constantinescu ierr = VecScale(W, ark->shift);CHKERRQ(ierr); 658f16577ceSEmil Constantinescu 6598a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 6608a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 6614f385281SJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 6624f385281SJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 663f16577ceSEmil Constantinescu 6648a381b04SJed Brown /* Initial guess taken from last stage */ 6658a381b04SJed Brown ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr); 6668a381b04SJed Brown ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 6678a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 6688a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 6698a381b04SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 67097335746SJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 67197335746SJed Brown ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr); 67297335746SJed Brown if (!accept) goto reject_step; 6738a381b04SJed Brown } 6748a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 6754cc180ffSJed Brown ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr); 6764cc180ffSJed Brown if (ark->imex) { 6778a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 6784cc180ffSJed Brown } else { 6794cc180ffSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 6804cc180ffSJed Brown } 6818a381b04SJed Brown } 682108c343cSJed Brown ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr); 683108c343cSJed Brown ark->status = TS_STEP_PENDING; 6848a381b04SJed Brown 685108c343cSJed Brown /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 686108c343cSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 687108c343cSJed Brown ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 688108c343cSJed Brown ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 689108c343cSJed Brown ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 690108c343cSJed Brown if (accept) { 691108c343cSJed Brown /* ignore next_scheme for now */ 6928a381b04SJed Brown ts->ptime += ts->time_step; 693cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 6948a381b04SJed Brown ts->steps++; 695108c343cSJed Brown ark->status = TS_STEP_COMPLETE; 696108c343cSJed Brown break; 697108c343cSJed Brown } else { /* Roll back the current step */ 698108c343cSJed Brown for (j=0; j<s; j++) w[j] = h*bt[j]; 699108c343cSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotI);CHKERRQ(ierr); 700108c343cSJed Brown for (j=0; j<s; j++) w[j] = -h*b[j]; 701108c343cSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotRHS);CHKERRQ(ierr); 702108c343cSJed Brown ts->time_step = next_time_step; 703108c343cSJed Brown ark->status = TS_STEP_INCOMPLETE; 704108c343cSJed Brown } 705476b6736SJed Brown reject_step: continue; 706108c343cSJed Brown } 707b2ce242eSJed Brown if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 7088a381b04SJed Brown PetscFunctionReturn(0); 7098a381b04SJed Brown } 7108a381b04SJed Brown 711cd652676SJed Brown #undef __FUNCT__ 712cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX" 713cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 714cd652676SJed Brown { 715cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7164f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 717108c343cSJed Brown PetscReal h; 718108c343cSJed Brown PetscReal tt,t; 719cd652676SJed Brown PetscScalar *bt,*b; 720cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 721cd652676SJed Brown PetscErrorCode ierr; 722cd652676SJed Brown 723cd652676SJed Brown PetscFunctionBegin; 724cd652676SJed Brown if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 725108c343cSJed Brown switch (ark->status) { 726108c343cSJed Brown case TS_STEP_INCOMPLETE: 727108c343cSJed Brown case TS_STEP_PENDING: 728108c343cSJed Brown h = ts->time_step; 729108c343cSJed Brown t = (itime - ts->ptime)/h; 730108c343cSJed Brown break; 731108c343cSJed Brown case TS_STEP_COMPLETE: 732108c343cSJed Brown h = ts->time_step_prev; 733108c343cSJed Brown t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 734108c343cSJed Brown break; 735b9ce6d65SJed Brown default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus"); 736108c343cSJed Brown } 737cd652676SJed Brown ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr); 738cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 7394f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 740cd652676SJed Brown for (i=0; i<s; i++) { 741108c343cSJed Brown bt[i] += h * Bt[i*pinterp+j] * tt * -1.0; 742108c343cSJed Brown b[i] += h * B[i*pinterp+j] * tt; 743cd652676SJed Brown } 744cd652676SJed Brown } 745cd652676SJed Brown if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved"); 746cd652676SJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 747cd652676SJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 748cd652676SJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 749cd652676SJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 750cd652676SJed Brown PetscFunctionReturn(0); 751cd652676SJed Brown } 752cd652676SJed Brown 7538a381b04SJed Brown /*------------------------------------------------------------*/ 7548a381b04SJed Brown #undef __FUNCT__ 7558a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX" 7568a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 7578a381b04SJed Brown { 7588a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7598a381b04SJed Brown PetscInt s; 7608a381b04SJed Brown PetscErrorCode ierr; 7618a381b04SJed Brown 7628a381b04SJed Brown PetscFunctionBegin; 7638a381b04SJed Brown if (!ark->tableau) PetscFunctionReturn(0); 7648a381b04SJed Brown s = ark->tableau->s; 7658a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 7668a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 7678a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 7688a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 7698a381b04SJed Brown ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 7708a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 7718a381b04SJed Brown ierr = PetscFree(ark->work);CHKERRQ(ierr); 7728a381b04SJed Brown PetscFunctionReturn(0); 7738a381b04SJed Brown } 7748a381b04SJed Brown 7758a381b04SJed Brown #undef __FUNCT__ 7768a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX" 7778a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 7788a381b04SJed Brown { 7798a381b04SJed Brown PetscErrorCode ierr; 7808a381b04SJed Brown 7818a381b04SJed Brown PetscFunctionBegin; 7828a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 7838a381b04SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 7848a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr); 7858a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr); 786995b3938SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",PETSC_NULL);CHKERRQ(ierr); 7878a381b04SJed Brown PetscFunctionReturn(0); 7888a381b04SJed Brown } 7898a381b04SJed Brown 7908a381b04SJed Brown /* 7918a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 7928a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 7938a381b04SJed Brown */ 7948a381b04SJed Brown #undef __FUNCT__ 7958a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 7968a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 7978a381b04SJed Brown { 7988a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7998a381b04SJed Brown PetscErrorCode ierr; 8008a381b04SJed Brown 8018a381b04SJed Brown PetscFunctionBegin; 8028a381b04SJed Brown ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 8034cc180ffSJed Brown ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr); 8048a381b04SJed Brown PetscFunctionReturn(0); 8058a381b04SJed Brown } 8068a381b04SJed Brown 8078a381b04SJed Brown #undef __FUNCT__ 8088a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 8098a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 8108a381b04SJed Brown { 8118a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 8128a381b04SJed Brown PetscErrorCode ierr; 8138a381b04SJed Brown 8148a381b04SJed Brown PetscFunctionBegin; 8158a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 8168a381b04SJed Brown ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 8178a381b04SJed Brown PetscFunctionReturn(0); 8188a381b04SJed Brown } 8198a381b04SJed Brown 8208a381b04SJed Brown #undef __FUNCT__ 8218a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX" 8228a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 8238a381b04SJed Brown { 8248a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 8258a381b04SJed Brown ARKTableau tab = ark->tableau; 8268a381b04SJed Brown PetscInt s = tab->s; 8278a381b04SJed Brown PetscErrorCode ierr; 8288a381b04SJed Brown 8298a381b04SJed Brown PetscFunctionBegin; 8308a381b04SJed Brown if (!ark->tableau) { 831e24355feSJed Brown ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 8328a381b04SJed Brown } 8338a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 8348a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 8358a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 8368a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 8378a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 8388a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 8398a381b04SJed Brown ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr); 8408a381b04SJed Brown PetscFunctionReturn(0); 8418a381b04SJed Brown } 8428a381b04SJed Brown /*------------------------------------------------------------*/ 8438a381b04SJed Brown 8448a381b04SJed Brown #undef __FUNCT__ 8458a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 8468a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts) 8478a381b04SJed Brown { 8484cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 8498a381b04SJed Brown PetscErrorCode ierr; 8508a381b04SJed Brown char arktype[256]; 8518a381b04SJed Brown 8528a381b04SJed Brown PetscFunctionBegin; 8538a381b04SJed Brown ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr); 8548a381b04SJed Brown { 8558a381b04SJed Brown ARKTableauLink link; 8568a381b04SJed Brown PetscInt count,choice; 8578a381b04SJed Brown PetscBool flg; 8588a381b04SJed Brown const char **namelist; 8598a381b04SJed Brown ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr); 8608a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 8618a381b04SJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 8628a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 8638a381b04SJed Brown ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 8648a381b04SJed Brown ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 8658a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 8664cc180ffSJed Brown flg = (PetscBool)!ark->imex; 8674cc180ffSJed Brown ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 8684cc180ffSJed Brown ark->imex = (PetscBool)!flg; 869d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 8708a381b04SJed Brown } 8718a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 8728a381b04SJed Brown PetscFunctionReturn(0); 8738a381b04SJed Brown } 8748a381b04SJed Brown 8758a381b04SJed Brown #undef __FUNCT__ 8768a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray" 8778a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 8788a381b04SJed Brown { 879257d2499SJed Brown PetscErrorCode ierr; 880f1d86077SJed Brown PetscInt i; 881f1d86077SJed Brown size_t left,count; 8828a381b04SJed Brown char *p; 8838a381b04SJed Brown 8848a381b04SJed Brown PetscFunctionBegin; 885f1d86077SJed Brown for (i=0,p=buf,left=len; i<n; i++) { 886f1d86077SJed Brown ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 8878a381b04SJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 8888a381b04SJed Brown left -= count; 8898a381b04SJed Brown p += count; 8908a381b04SJed Brown *p++ = ' '; 8918a381b04SJed Brown } 8928a381b04SJed Brown p[i ? 0 : -1] = 0; 8938a381b04SJed Brown PetscFunctionReturn(0); 8948a381b04SJed Brown } 8958a381b04SJed Brown 8968a381b04SJed Brown #undef __FUNCT__ 8978a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX" 8988a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 8998a381b04SJed Brown { 9008a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 9018a381b04SJed Brown ARKTableau tab = ark->tableau; 9028a381b04SJed Brown PetscBool iascii; 9038a381b04SJed Brown PetscErrorCode ierr; 9048a381b04SJed Brown 9058a381b04SJed Brown PetscFunctionBegin; 9068a381b04SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9078a381b04SJed Brown if (iascii) { 9088a381b04SJed Brown const TSARKIMEXType arktype; 9098a381b04SJed Brown char buf[512]; 9108a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 9118a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 9128a381b04SJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 91331f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 91431f6fcc0SJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 91531f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 9168a381b04SJed Brown } 917d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 9188a381b04SJed Brown PetscFunctionReturn(0); 9198a381b04SJed Brown } 9208a381b04SJed Brown 9218a381b04SJed Brown #undef __FUNCT__ 9228a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType" 9238a381b04SJed Brown /*@C 9248a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 9258a381b04SJed Brown 9268a381b04SJed Brown Logically collective 9278a381b04SJed Brown 9288a381b04SJed Brown Input Parameter: 9298a381b04SJed Brown + ts - timestepping context 9308a381b04SJed Brown - arktype - type of ARK-IMEX scheme 9318a381b04SJed Brown 9328a381b04SJed Brown Level: intermediate 9338a381b04SJed Brown 934020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5 9358a381b04SJed Brown @*/ 9368a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype) 9378a381b04SJed Brown { 9388a381b04SJed Brown PetscErrorCode ierr; 9398a381b04SJed Brown 9408a381b04SJed Brown PetscFunctionBegin; 9418a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9428a381b04SJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 9438a381b04SJed Brown PetscFunctionReturn(0); 9448a381b04SJed Brown } 9458a381b04SJed Brown 9468a381b04SJed Brown #undef __FUNCT__ 9478a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType" 9488a381b04SJed Brown /*@C 9498a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 9508a381b04SJed Brown 9518a381b04SJed Brown Logically collective 9528a381b04SJed Brown 9538a381b04SJed Brown Input Parameter: 9548a381b04SJed Brown . ts - timestepping context 9558a381b04SJed Brown 9568a381b04SJed Brown Output Parameter: 9578a381b04SJed Brown . arktype - type of ARK-IMEX scheme 9588a381b04SJed Brown 9598a381b04SJed Brown Level: intermediate 9608a381b04SJed Brown 9618a381b04SJed Brown .seealso: TSARKIMEXGetType() 9628a381b04SJed Brown @*/ 9638a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype) 9648a381b04SJed Brown { 9658a381b04SJed Brown PetscErrorCode ierr; 9668a381b04SJed Brown 9678a381b04SJed Brown PetscFunctionBegin; 9688a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9698a381b04SJed Brown ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 9708a381b04SJed Brown PetscFunctionReturn(0); 9718a381b04SJed Brown } 9728a381b04SJed Brown 9734cc180ffSJed Brown #undef __FUNCT__ 9744cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit" 9754cc180ffSJed Brown /*@C 9764cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 9774cc180ffSJed Brown 9784cc180ffSJed Brown Logically collective 9794cc180ffSJed Brown 9804cc180ffSJed Brown Input Parameter: 9814cc180ffSJed Brown + ts - timestepping context 9824cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 9834cc180ffSJed Brown 9844cc180ffSJed Brown Level: intermediate 9854cc180ffSJed Brown 9864cc180ffSJed Brown .seealso: TSARKIMEXGetType() 9874cc180ffSJed Brown @*/ 9884cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 9894cc180ffSJed Brown { 9904cc180ffSJed Brown PetscErrorCode ierr; 9914cc180ffSJed Brown 9924cc180ffSJed Brown PetscFunctionBegin; 9934cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 9944cc180ffSJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 9954cc180ffSJed Brown PetscFunctionReturn(0); 9964cc180ffSJed Brown } 9974cc180ffSJed Brown 9988a381b04SJed Brown EXTERN_C_BEGIN 9998a381b04SJed Brown #undef __FUNCT__ 10008a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 10018a381b04SJed Brown PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype) 10028a381b04SJed Brown { 10038a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 10048a381b04SJed Brown PetscErrorCode ierr; 10058a381b04SJed Brown 10068a381b04SJed Brown PetscFunctionBegin; 10078a381b04SJed Brown if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);} 10088a381b04SJed Brown *arktype = ark->tableau->name; 10098a381b04SJed Brown PetscFunctionReturn(0); 10108a381b04SJed Brown } 10118a381b04SJed Brown #undef __FUNCT__ 10128a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 10138a381b04SJed Brown PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype) 10148a381b04SJed Brown { 10158a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 10168a381b04SJed Brown PetscErrorCode ierr; 10178a381b04SJed Brown PetscBool match; 10188a381b04SJed Brown ARKTableauLink link; 10198a381b04SJed Brown 10208a381b04SJed Brown PetscFunctionBegin; 10218a381b04SJed Brown if (ark->tableau) { 10228a381b04SJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 10238a381b04SJed Brown if (match) PetscFunctionReturn(0); 10248a381b04SJed Brown } 10258a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 10268a381b04SJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 10278a381b04SJed Brown if (match) { 10288a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 10298a381b04SJed Brown ark->tableau = &link->tab; 10308a381b04SJed Brown PetscFunctionReturn(0); 10318a381b04SJed Brown } 10328a381b04SJed Brown } 10338a381b04SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 10348a381b04SJed Brown PetscFunctionReturn(0); 10358a381b04SJed Brown } 10364cc180ffSJed Brown #undef __FUNCT__ 10374cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX" 10384cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 10394cc180ffSJed Brown { 10404cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 10414cc180ffSJed Brown 10424cc180ffSJed Brown PetscFunctionBegin; 10434cc180ffSJed Brown ark->imex = (PetscBool)!flg; 10444cc180ffSJed Brown PetscFunctionReturn(0); 10454cc180ffSJed Brown } 10468a381b04SJed Brown EXTERN_C_END 10478a381b04SJed Brown 10488a381b04SJed Brown /* ------------------------------------------------------------ */ 10498a381b04SJed Brown /*MC 10508a381b04SJed Brown TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes 10518a381b04SJed Brown 1052fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 1053fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 1054fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 1055fca742c7SJed Brown 1056fca742c7SJed Brown Notes: 1057c8058688SBarry Smith The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 1058c8058688SBarry Smith 1059fca742c7SJed Brown This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 1060fca742c7SJed Brown 10618a381b04SJed Brown Level: beginner 10628a381b04SJed Brown 1063c8058688SBarry Smith .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3, 106426885f59SBarry Smith TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, TSARKIMEXRegister() 10658a381b04SJed Brown 10668a381b04SJed Brown M*/ 10678a381b04SJed Brown EXTERN_C_BEGIN 10688a381b04SJed Brown #undef __FUNCT__ 10698a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX" 10708a381b04SJed Brown PetscErrorCode TSCreate_ARKIMEX(TS ts) 10718a381b04SJed Brown { 10728a381b04SJed Brown TS_ARKIMEX *th; 10738a381b04SJed Brown PetscErrorCode ierr; 10748a381b04SJed Brown 10758a381b04SJed Brown PetscFunctionBegin; 10768a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 10778a381b04SJed Brown ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr); 10788a381b04SJed Brown #endif 10798a381b04SJed Brown 10808a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 10818a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 10828a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 10838a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 10848a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 1085cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 1086108c343cSJed Brown ts->ops->evaluatestep = TSEvaluateStep_ARKIMEX; 10878a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 10888a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 10898a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 10908a381b04SJed Brown 10918a381b04SJed Brown ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr); 10928a381b04SJed Brown ts->data = (void*)th; 10934cc180ffSJed Brown th->imex = PETSC_TRUE; 10948a381b04SJed Brown 10958a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 10968a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 10974cc180ffSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 10988a381b04SJed Brown PetscFunctionReturn(0); 10998a381b04SJed Brown } 11008a381b04SJed Brown EXTERN_C_END 1101