xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision da80777bb9afa786145217fc25b54e33c3ebf541)
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 */
12b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
138a381b04SJed Brown 
1419fd82e9SBarry Smith static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3;
158a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled;
168a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized;
17e817cc15SEmil Constantinescu static PetscInt explicit_stage_time_id;
188a381b04SJed Brown 
198a381b04SJed Brown typedef struct _ARKTableau *ARKTableau;
208a381b04SJed Brown struct _ARKTableau {
218a381b04SJed Brown   char      *name;
224f385281SJed Brown   PetscInt  order;               /* Classical approximation order of the method */
234f385281SJed Brown   PetscInt  s;                   /* Number of stages */
24e817cc15SEmil Constantinescu   PetscBool stiffly_accurate;    /* The implicit part is stiffly accurate*/
25e817cc15SEmil Constantinescu   PetscBool FSAL_implicit;       /* The implicit part is FSAL*/
26e817cc15SEmil Constantinescu   PetscBool explicit_first_stage;/* The implicit part has an explicit first stage*/
274f385281SJed Brown   PetscInt  pinterp;             /* Interpolation order */
284f385281SJed Brown   PetscReal *At,*bt,*ct;        /* Stiff tableau */
298a381b04SJed Brown   PetscReal *A,*b,*c;           /* Non-stiff tableau */
30108c343cSJed Brown   PetscReal *bembedt,*bembed;   /* Embedded formula of order one less (order-1) */
31cd652676SJed Brown   PetscReal *binterpt,*binterp; /* Dense output formula */
32108c343cSJed Brown   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
338a381b04SJed Brown };
348a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
358a381b04SJed Brown struct _ARKTableauLink {
368a381b04SJed Brown   struct _ARKTableau tab;
378a381b04SJed Brown   ARKTableauLink next;
388a381b04SJed Brown };
398a381b04SJed Brown static ARKTableauLink ARKTableauList;
408a381b04SJed Brown 
418a381b04SJed Brown typedef struct {
428a381b04SJed Brown   ARKTableau   tableau;
438a381b04SJed Brown   Vec          *Y;               /* States computed during the step */
448a381b04SJed Brown   Vec          *YdotI;           /* Time derivatives for the stiff part */
458a381b04SJed Brown   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
46e817cc15SEmil Constantinescu   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
478a381b04SJed Brown   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
488a381b04SJed Brown   Vec          Work;             /* Generic work vector */
498a381b04SJed Brown   Vec          Z;                /* Ydot = shift(Y-Z) */
508a381b04SJed Brown   PetscScalar  *work;            /* Scalar work */
51b296d7d5SJed Brown   PetscReal    scoeff;           /* shift = scoeff/dt */
528a381b04SJed Brown   PetscReal    stage_time;
534cc180ffSJed Brown   PetscBool    imex;
54108c343cSJed Brown   TSStepStatus status;
558a381b04SJed Brown } TS_ARKIMEX;
561f80e275SEmil Constantinescu /*MC
571f80e275SEmil Constantinescu      TSARKIMEXARS122 - Second order ARK IMEX scheme.
588a381b04SJed Brown 
591f80e275SEmil Constantinescu      This method has one explicit stage and one implicit stage.
601f80e275SEmil Constantinescu 
611f80e275SEmil Constantinescu      References:
621f80e275SEmil 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.
631f80e275SEmil Constantinescu 
641f80e275SEmil Constantinescu      Level: advanced
651f80e275SEmil Constantinescu 
661f80e275SEmil Constantinescu .seealso: TSARKIMEX
671f80e275SEmil Constantinescu M*/
681f80e275SEmil Constantinescu /*MC
691f80e275SEmil Constantinescu      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.
701f80e275SEmil Constantinescu 
711f80e275SEmil 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.
721f80e275SEmil Constantinescu 
731f80e275SEmil Constantinescu      Level: advanced
741f80e275SEmil Constantinescu 
751f80e275SEmil Constantinescu .seealso: TSARKIMEX
761f80e275SEmil Constantinescu M*/
771f80e275SEmil Constantinescu /*MC
781f80e275SEmil Constantinescu      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.
791f80e275SEmil Constantinescu 
801f80e275SEmil Constantinescu      This method has two implicit stages, and L-stable implicit scheme.
811f80e275SEmil Constantinescu 
821f80e275SEmil Constantinescu     References:
831f80e275SEmil 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
841f80e275SEmil Constantinescu 
851f80e275SEmil Constantinescu      Level: advanced
861f80e275SEmil Constantinescu 
871f80e275SEmil Constantinescu .seealso: TSARKIMEX
881f80e275SEmil Constantinescu M*/
891f80e275SEmil Constantinescu /*MC
90e817cc15SEmil Constantinescu      TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.
91e817cc15SEmil Constantinescu 
92e817cc15SEmil Constantinescu      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.
93e817cc15SEmil Constantinescu 
94e817cc15SEmil Constantinescu      Level: advanced
95e817cc15SEmil Constantinescu 
96e817cc15SEmil Constantinescu .seealso: TSARKIMEX
97e817cc15SEmil Constantinescu M*/
98e817cc15SEmil Constantinescu /*MC
991f80e275SEmil Constantinescu      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.
1001f80e275SEmil Constantinescu 
1011f80e275SEmil 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.
1021f80e275SEmil Constantinescu 
1031f80e275SEmil Constantinescu      Level: advanced
1041f80e275SEmil Constantinescu 
1051f80e275SEmil Constantinescu .seealso: TSARKIMEX
1061f80e275SEmil Constantinescu M*/
10764f491ddSJed Brown /*MC
10864f491ddSJed Brown      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
10964f491ddSJed Brown 
110617a39beSEmil 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.
11164f491ddSJed Brown 
112b330ce4dSSatish Balay      Level: advanced
113b330ce4dSSatish Balay 
11464f491ddSJed Brown .seealso: TSARKIMEX
11564f491ddSJed Brown M*/
11664f491ddSJed Brown /*MC
11764f491ddSJed Brown      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
11864f491ddSJed Brown 
11964f491ddSJed Brown      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
12064f491ddSJed Brown 
121b330ce4dSSatish Balay      Level: advanced
122b330ce4dSSatish Balay 
12364f491ddSJed Brown .seealso: TSARKIMEX
12464f491ddSJed Brown M*/
12564f491ddSJed Brown /*MC
1266cf0794eSJed Brown      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.
1276cf0794eSJed Brown 
1286cf0794eSJed Brown      This method has three implicit stages.
1296cf0794eSJed Brown 
1306cf0794eSJed Brown      References:
1316cf0794eSJed 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
1326cf0794eSJed Brown 
1336cf0794eSJed Brown      This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375
1346cf0794eSJed Brown 
1356cf0794eSJed Brown      Level: advanced
1366cf0794eSJed Brown 
1376cf0794eSJed Brown .seealso: TSARKIMEX
1386cf0794eSJed Brown M*/
1396cf0794eSJed Brown /*MC
14064f491ddSJed Brown      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
14164f491ddSJed Brown 
14264f491ddSJed Brown      This method has one explicit stage and three implicit stages.
14364f491ddSJed Brown 
14464f491ddSJed Brown      References:
14564f491ddSJed Brown      Kennedy and Carpenter 2003.
14664f491ddSJed Brown 
147b330ce4dSSatish Balay      Level: advanced
148b330ce4dSSatish Balay 
14964f491ddSJed Brown .seealso: TSARKIMEX
15064f491ddSJed Brown M*/
15164f491ddSJed Brown /*MC
1526cf0794eSJed Brown      TSARKIMEXARS443 - 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      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.
1586cf0794eSJed Brown 
1596cf0794eSJed Brown      This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375
1606cf0794eSJed Brown 
1616cf0794eSJed Brown      Level: advanced
1626cf0794eSJed Brown 
1636cf0794eSJed Brown .seealso: TSARKIMEX
1646cf0794eSJed Brown M*/
1656cf0794eSJed Brown /*MC
1666cf0794eSJed Brown      TSARKIMEXBPR3 - Third order ARK IMEX scheme.
1676cf0794eSJed Brown 
1686cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1696cf0794eSJed Brown 
1706cf0794eSJed Brown      References:
1716cf0794eSJed Brown      This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375
1726cf0794eSJed Brown 
1736cf0794eSJed Brown      Level: advanced
1746cf0794eSJed Brown 
1756cf0794eSJed Brown .seealso: TSARKIMEX
1766cf0794eSJed Brown M*/
1776cf0794eSJed Brown /*MC
17864f491ddSJed Brown      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
17964f491ddSJed Brown 
18064f491ddSJed Brown      This method has one explicit stage and four implicit stages.
18164f491ddSJed Brown 
18264f491ddSJed Brown      References:
18364f491ddSJed Brown      Kennedy and Carpenter 2003.
18464f491ddSJed Brown 
185b330ce4dSSatish Balay      Level: advanced
186b330ce4dSSatish Balay 
18764f491ddSJed Brown .seealso: TSARKIMEX
18864f491ddSJed Brown M*/
18964f491ddSJed Brown /*MC
19064f491ddSJed Brown      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
19164f491ddSJed Brown 
19264f491ddSJed Brown      This method has one explicit stage and five implicit stages.
19364f491ddSJed Brown 
19464f491ddSJed Brown      References:
19564f491ddSJed Brown      Kennedy and Carpenter 2003.
19664f491ddSJed Brown 
197b330ce4dSSatish Balay      Level: advanced
198b330ce4dSSatish Balay 
19964f491ddSJed Brown .seealso: TSARKIMEX
20064f491ddSJed Brown M*/
20164f491ddSJed Brown 
2028a381b04SJed Brown #undef __FUNCT__
2038a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
2048a381b04SJed Brown /*@C
2058a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
2068a381b04SJed Brown 
207fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
2088a381b04SJed Brown 
2098a381b04SJed Brown   Level: advanced
2108a381b04SJed Brown 
2118a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
2128a381b04SJed Brown 
2138a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
2148a381b04SJed Brown @*/
2158a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
2168a381b04SJed Brown {
2178a381b04SJed Brown   PetscErrorCode ierr;
2188a381b04SJed Brown 
2198a381b04SJed Brown   PetscFunctionBegin;
2208a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
2218a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
222e817cc15SEmil Constantinescu 
223e817cc15SEmil Constantinescu   {
224e817cc15SEmil Constantinescu     const PetscReal
225e817cc15SEmil Constantinescu       A[3][3] = {{0.0,0.0,0.0},
226e817cc15SEmil Constantinescu                  {0.0,0.0,0.0},
227748ad121SEmil Constantinescu                  {0.0,0.5,0.0}},
228e817cc15SEmil Constantinescu       At[3][3] = {{1.0,0.0,0.0},
229e817cc15SEmil Constantinescu                   {0.0,0.5,0.0},
230e817cc15SEmil Constantinescu                   {0.0,0.5,0.5}},
231e817cc15SEmil Constantinescu         b[3] = {0.0,0.5,0.5},
232e817cc15SEmil Constantinescu           bembedt[3] = {1.0,0.0,0.0};
2332c0c504eSEmil Constantinescu           ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr);
234e817cc15SEmil Constantinescu   }
2358a381b04SJed Brown   {
2368a381b04SJed Brown     const PetscReal
2371f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2381f80e275SEmil Constantinescu                  {0.5,0.0}},
2391f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2401f80e275SEmil Constantinescu                   {0.0,0.5}},
2411f80e275SEmil Constantinescu         b[2] = {0.0,1.0},
2421f80e275SEmil Constantinescu           bembedt[2] = {0.5,0.5};
2431f80e275SEmil 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*/
2441f80e275SEmil 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);
2451f80e275SEmil Constantinescu   }
2461f80e275SEmil Constantinescu   {
2471f80e275SEmil Constantinescu     const PetscReal
2481f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2491f80e275SEmil Constantinescu                  {1.0,0.0}},
2501f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2511f80e275SEmil Constantinescu                   {0.5,0.5}},
2521f80e275SEmil Constantinescu         b[2] = {0.5,0.5},
2531f80e275SEmil Constantinescu           bembedt[2] = {0.0,1.0};
2541f80e275SEmil 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*/
2551f80e275SEmil 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);
2561f80e275SEmil Constantinescu   }
2571f80e275SEmil Constantinescu   {
258*da80777bSKarl 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   */
2591f80e275SEmil Constantinescu     const PetscReal
2601f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2611f80e275SEmil Constantinescu                  {1.0,0.0}},
262*da80777bSKarl Rupp       /*At[2][2] = {{us2,0.0},
263*da80777bSKarl Rupp                   {1.0-2.0*us2,us2}},*/
264*da80777bSKarl Rupp       At[2][2] = {{0.2928932188134524755992,0.0},
265*da80777bSKarl Rupp                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
2661f80e275SEmil Constantinescu       b[2] = {0.5,0.5},
2671f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0},
268*da80777bSKarl Rupp       /*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))}},*/
269*da80777bSKarl Rupp       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
270*da80777bSKarl 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}};
2721f80e275SEmil 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);
2731f80e275SEmil Constantinescu   }
2741f80e275SEmil Constantinescu   {
275*da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
276*da80777bSKarl Rupp     const PetscReal
2778a381b04SJed Brown       A[3][3] = {{0,0,0},
278*da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
279617a39beSEmil Constantinescu                  {0.5,0.5,0}},
280*da80777bSKarl Rupp       /*At[3][3] = {{0,0,0},
2811f80e275SEmil Constantinescu                   {1-1/s2,1-1/s2,0},
282*da80777bSKarl Rupp                   {1/(2*s2),1/(2*s2),1-1/s2}},*/
283*da80777bSKarl Rupp       At[3][3] = {{0,0,0},
284*da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
285*da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
286*da80777bSKarl Rupp       /*bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)}, */
287*da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
288*da80777bSKarl Rupp       /*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}};*/
289*da80777bSKarl Rupp       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
290*da80777bSKarl Rupp                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
291*da80777bSKarl Rupp                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
2921f80e275SEmil 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);
2931f80e275SEmil Constantinescu   }
2941f80e275SEmil Constantinescu   {
295*da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
296*da80777bSKarl Rupp     const PetscReal
2971f80e275SEmil Constantinescu       A[3][3] = {{0,0,0},
298*da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
2998a381b04SJed Brown                  {0.75,0.25,0}},
300*da80777bSKarl Rupp       /*At[3][3] = {{0,0,0},
3011f80e275SEmil Constantinescu                   {1-1/s2,1-1/s2,0},
302*da80777bSKarl Rupp                   {1/(2*s2),1/(2*s2),1-1/s2}},*/
303*da80777bSKarl Rupp       At[3][3] = {{0,0,0},
304*da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
305*da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
306*da80777bSKarl Rupp       /*bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)},*/
307*da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
308*da80777bSKarl Rupp       /*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}};*/
309*da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
310*da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
311*da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
312108c343cSJed 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);
3138a381b04SJed Brown   }
31406db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
315*da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
316*da80777bSKarl Rupp     const PetscReal
317*da80777bSKarl Rupp       /*A[3][3] = {{0,0,0},
318a3a57f36SJed Brown                  {2-s2,0,0},
319*da80777bSKarl Rupp                  {(3-2*s2)/6,(3+2*s2)/6,0}},*/
320*da80777bSKarl Rupp       A[3][3] = {{0,0,0},
321*da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
322*da80777bSKarl Rupp                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
323*da80777bSKarl Rupp       /*At[3][3] = {{0,0,0},
324a3a57f36SJed Brown                   {1-1/s2,1-1/s2,0},
325*da80777bSKarl Rupp                   {1/(2*s2),1/(2*s2),1-1/s2}},*/
326*da80777bSKarl Rupp       At[3][3] = {{0,0,0},
327*da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
328*da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
329*da80777bSKarl Rupp       /*bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)},*/
330*da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
331*da80777bSKarl Rupp       /*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}};*/
332*da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
333*da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
334*da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3351f80e275SEmil 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);
336a3a57f36SJed Brown   }
3376cf0794eSJed Brown   {                             /* Optimal for linear implicit part */
3386cf0794eSJed Brown     const PetscReal
3396cf0794eSJed Brown       A[3][3] = {{0,0,0},
3406cf0794eSJed Brown                  {0.5,0,0},
3416cf0794eSJed Brown                  {0.5,0.5,0}},
3426cf0794eSJed Brown       At[3][3] = {{0.25,0,0},
3436cf0794eSJed Brown                   {0,0.25,0},
3446cf0794eSJed Brown                   {1./3,1./3,1./3}};
345108c343cSJed 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);
3466cf0794eSJed Brown   }
347a3a57f36SJed Brown   {
348a3a57f36SJed Brown     const PetscReal
349a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
3504040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
3514040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
3524040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
353a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
3544040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
3554040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
3564040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
357cc46b9d1SJed Brown       bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
3584040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
3594040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
3604040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
3614040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
362108c343cSJed 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);
363a3a57f36SJed Brown   }
364a3a57f36SJed Brown   {
365a3a57f36SJed Brown     const PetscReal
366e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3676cf0794eSJed Brown                  {1./2,0,0,0,0},
3686cf0794eSJed Brown                  {11./18,1./18,0,0,0},
3696cf0794eSJed Brown                  {5./6,-5./6,.5,0,0},
3706cf0794eSJed Brown                  {1./4,7./4,3./4,-7./4,0}},
3716cf0794eSJed Brown       At[5][5] = {{0,0,0,0,0},
3726cf0794eSJed Brown                   {0,1./2,0,0,0},
3736cf0794eSJed Brown                   {0,1./6,1./2,0,0},
3746cf0794eSJed Brown                   {0,-1./2,1./2,1./2,0},
375108c343cSJed Brown                   {0,3./2,-3./2,1./2,1./2}},
376108c343cSJed Brown       *bembedt = PETSC_NULL;
377108c343cSJed 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);
3786cf0794eSJed Brown   }
3796cf0794eSJed Brown   {
3806cf0794eSJed Brown     const PetscReal
381e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3826cf0794eSJed Brown                  {1,0,0,0,0},
3836cf0794eSJed Brown                  {4./9,2./9,0,0,0},
3846cf0794eSJed Brown                  {1./4,0,3./4,0,0},
3856cf0794eSJed Brown                  {1./4,0,3./5,0,0}},
386e74514c0SSatish Balay       At[5][5] = {{0,0,0,0,0},
3876cf0794eSJed Brown                   {.5,.5,0,0,0},
3886cf0794eSJed Brown                   {5./18,-1./9,.5,0,0},
3896cf0794eSJed Brown                   {.5,0,0,.5,0},
390108c343cSJed Brown                   {.25,0,.75,-.5,.5}},
391108c343cSJed Brown       *bembedt = PETSC_NULL;
392108c343cSJed 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);
3936cf0794eSJed Brown   }
3946cf0794eSJed Brown   {
3956cf0794eSJed Brown     const PetscReal
396a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
397a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
3984040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
3994040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
4004040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
4014040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
402a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
403a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
4044040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
4054040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
4064040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
4074040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
408cc46b9d1SJed Brown       bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
4094040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
410cd652676SJed Brown                         {0,0,0},
4114040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
4124040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
4134040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
4144040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
415108c343cSJed 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);
416a3a57f36SJed Brown   }
417a3a57f36SJed Brown   {
418a3a57f36SJed Brown     const PetscReal
419a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
420a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
4214040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
4224040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
4234040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
4244040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
4254040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
4264040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
427a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
4284040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
4294040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
4304040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
4314040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
4324040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
4334040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
4344040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
435cc46b9d1SJed Brown       bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
4364040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614. ,  43486358583215./12773830924787. , -9257016797708./5021505065439.},
437cd652676SJed Brown                         {0                               ,  0                              , 0                            },
438cd652676SJed Brown                         {0                               ,  0                              , 0                            },
4394040e9f2SJed Brown                         {65168852399939./7868540260826.  ,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
4404040e9f2SJed Brown                         {15494834004392./5936557850923.  ,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
4414040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473. , 30029262896817./10175596800299.},
4424040e9f2SJed Brown                         {-19024464361622./5461577185407. ,  115839755401235./10719374521269., -26136350496073./3983972220547.},
4434040e9f2SJed Brown                         {-6511271360970./6095937251113.  ,  5843115559534./2180450260947.   , -5289405421727./3760307252460. }};
444108c343cSJed 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);
445a3a57f36SJed Brown   }
446a3a57f36SJed Brown 
4478a381b04SJed Brown   PetscFunctionReturn(0);
4488a381b04SJed Brown }
4498a381b04SJed Brown 
4508a381b04SJed Brown #undef __FUNCT__
4518a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
4528a381b04SJed Brown /*@C
4538a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
4548a381b04SJed Brown 
4558a381b04SJed Brown    Not Collective
4568a381b04SJed Brown 
4578a381b04SJed Brown    Level: advanced
4588a381b04SJed Brown 
4598a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
4608a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
4618a381b04SJed Brown @*/
4628a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
4638a381b04SJed Brown {
4648a381b04SJed Brown   PetscErrorCode ierr;
4658a381b04SJed Brown   ARKTableauLink link;
4668a381b04SJed Brown 
4678a381b04SJed Brown   PetscFunctionBegin;
4688a381b04SJed Brown   while ((link = ARKTableauList)) {
4698a381b04SJed Brown     ARKTableau t = &link->tab;
4708a381b04SJed Brown     ARKTableauList = link->next;
4718a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
472108c343cSJed Brown     ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr);
473cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
4748a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
4758a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
4768a381b04SJed Brown   }
4778a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
4788a381b04SJed Brown   PetscFunctionReturn(0);
4798a381b04SJed Brown }
4808a381b04SJed Brown 
4818a381b04SJed Brown #undef __FUNCT__
4828a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
4838a381b04SJed Brown /*@C
4848a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
4858a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
4868a381b04SJed Brown   when using static libraries.
4878a381b04SJed Brown 
4888a381b04SJed Brown   Input Parameter:
4898a381b04SJed Brown   path - The dynamic library path, or PETSC_NULL
4908a381b04SJed Brown 
4918a381b04SJed Brown   Level: developer
4928a381b04SJed Brown 
4938a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
4948a381b04SJed Brown .seealso: PetscInitialize()
4958a381b04SJed Brown @*/
4968a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
4978a381b04SJed Brown {
4988a381b04SJed Brown   PetscErrorCode ierr;
4998a381b04SJed Brown 
5008a381b04SJed Brown   PetscFunctionBegin;
5018a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
5028a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
5038a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
504e817cc15SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
5058a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
5068a381b04SJed Brown   PetscFunctionReturn(0);
5078a381b04SJed Brown }
5088a381b04SJed Brown 
5098a381b04SJed Brown #undef __FUNCT__
5108a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
5118a381b04SJed Brown /*@C
5128a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
5138a381b04SJed Brown   called from PetscFinalize().
5148a381b04SJed Brown 
5158a381b04SJed Brown   Level: developer
5168a381b04SJed Brown 
5178a381b04SJed Brown .keywords: Petsc, destroy, package
5188a381b04SJed Brown .seealso: PetscFinalize()
5198a381b04SJed Brown @*/
5208a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
5218a381b04SJed Brown {
5228a381b04SJed Brown   PetscErrorCode ierr;
5238a381b04SJed Brown 
5248a381b04SJed Brown   PetscFunctionBegin;
5258a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
5268a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
5278a381b04SJed Brown   PetscFunctionReturn(0);
5288a381b04SJed Brown }
5298a381b04SJed Brown 
5308a381b04SJed Brown #undef __FUNCT__
5318a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
532cd652676SJed Brown /*@C
533cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
534cd652676SJed Brown 
535cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
536cd652676SJed Brown 
537cd652676SJed Brown    Input Parameters:
538cd652676SJed Brown +  name - identifier for method
539cd652676SJed Brown .  order - approximation order of method
540cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
541cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
542cd652676SJed Brown .  bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
543cd652676SJed Brown .  ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
544cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
545cd652676SJed Brown .  b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
546cd652676SJed Brown .  c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
547108c343cSJed Brown .  bembedt - Stiff part of completion table for embedded method (dimension s; PETSC_NULL if not available)
548108c343cSJed Brown .  bembed - Non-stiff part of completion table for embedded method (dimension s; PETSC_NULL to use bembedt if provided)
549cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
550cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
551cd652676SJed Brown -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
552cd652676SJed Brown 
553cd652676SJed Brown    Notes:
554cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
555cd652676SJed Brown 
556cd652676SJed Brown    Level: advanced
557cd652676SJed Brown 
558cd652676SJed Brown .keywords: TS, register
559cd652676SJed Brown 
560cd652676SJed Brown .seealso: TSARKIMEX
561cd652676SJed Brown @*/
56219fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
5638a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
564cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
565108c343cSJed Brown                                  const PetscReal bembedt[],const PetscReal bembed[],
566cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
5678a381b04SJed Brown {
5688a381b04SJed Brown   PetscErrorCode ierr;
5698a381b04SJed Brown   ARKTableauLink link;
5708a381b04SJed Brown   ARKTableau     t;
5718a381b04SJed Brown   PetscInt       i,j;
5728a381b04SJed Brown 
5738a381b04SJed Brown   PetscFunctionBegin;
5748a381b04SJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
575cd652676SJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
5768a381b04SJed Brown   t = &link->tab;
5778a381b04SJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5788a381b04SJed Brown   t->order = order;
5798a381b04SJed Brown   t->s = s;
5808a381b04SJed 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);
5818a381b04SJed Brown   ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
5828a381b04SJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5838a381b04SJed Brown   if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);}
5848a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
5858a381b04SJed Brown   if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);}
5868a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
5878a381b04SJed Brown   if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);}
5888a381b04SJed 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];
5898a381b04SJed Brown   if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);}
5908a381b04SJed 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];
591e817cc15SEmil Constantinescu   t->stiffly_accurate = PETSC_TRUE;
592e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
593e817cc15SEmil Constantinescu   t->explicit_first_stage = PETSC_TRUE;
594e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
595e817cc15SEmil Constantinescu   /*def of FSAL can be made more precise*/
5964e9d4bf5SJed Brown   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
597108c343cSJed Brown   if (bembedt) {
598108c343cSJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembedt,s,PetscReal,&t->bembed);CHKERRQ(ierr);
599108c343cSJed Brown     ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr);
600108c343cSJed Brown     ierr = PetscMemcpy(t->bembed,bembed?bembed:bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr);
601108c343cSJed Brown   }
602108c343cSJed Brown 
6034f385281SJed Brown   t->pinterp = pinterp;
604cd652676SJed Brown   ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
605cd652676SJed Brown   ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
606cd652676SJed Brown   ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
6078a381b04SJed Brown   link->next = ARKTableauList;
6088a381b04SJed Brown   ARKTableauList = link;
6098a381b04SJed Brown   PetscFunctionReturn(0);
6108a381b04SJed Brown }
6118a381b04SJed Brown 
6128a381b04SJed Brown #undef __FUNCT__
613108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX"
614108c343cSJed Brown /*
615108c343cSJed Brown  The step completion formula is
616108c343cSJed Brown 
617108c343cSJed Brown  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
618108c343cSJed Brown 
619108c343cSJed Brown  This function can be called before or after ts->vec_sol has been updated.
620108c343cSJed Brown  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
621108c343cSJed Brown  We can write
622108c343cSJed Brown 
623108c343cSJed Brown  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
624108c343cSJed Brown      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
625108c343cSJed Brown      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
626108c343cSJed Brown 
627108c343cSJed Brown  so we can evaluate the method with different order even after the step has been optimistically completed.
628108c343cSJed Brown */
629108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
630108c343cSJed Brown {
631108c343cSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
632108c343cSJed Brown   ARKTableau     tab  = ark->tableau;
633108c343cSJed Brown   PetscScalar    *w = ark->work;
634108c343cSJed Brown   PetscReal      h;
635108c343cSJed Brown   PetscInt       s = tab->s,j;
636108c343cSJed Brown   PetscErrorCode ierr;
637108c343cSJed Brown 
638108c343cSJed Brown   PetscFunctionBegin;
639108c343cSJed Brown   switch (ark->status) {
640108c343cSJed Brown   case TS_STEP_INCOMPLETE:
641108c343cSJed Brown   case TS_STEP_PENDING:
642108c343cSJed Brown     h = ts->time_step; break;
643108c343cSJed Brown   case TS_STEP_COMPLETE:
644108c343cSJed Brown     h = ts->time_step_prev; break;
645b9ce6d65SJed Brown   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
646108c343cSJed Brown   }
647108c343cSJed Brown   if (order == tab->order) {
648e817cc15SEmil Constantinescu     if (ark->status == TS_STEP_INCOMPLETE) {
649740132f1SEmil Constantinescu       if (!ark->imex && tab->stiffly_accurate) {/* Only the stiffly accurate implicit formula is used */
650e817cc15SEmil Constantinescu         ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr);
651e817cc15SEmil Constantinescu       } else { /* Use the standard completion formula (bt,b) */
652108c343cSJed Brown         ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
653e817cc15SEmil Constantinescu         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
654108c343cSJed Brown         ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
655e817cc15SEmil Constantinescu         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
656108c343cSJed Brown           for (j=0; j<s; j++) w[j] = h*tab->b[j];
657108c343cSJed Brown           ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
658e817cc15SEmil Constantinescu         }
659e817cc15SEmil Constantinescu       }
660108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
661108c343cSJed Brown     if (done) *done = PETSC_TRUE;
662108c343cSJed Brown     PetscFunctionReturn(0);
663108c343cSJed Brown   } else if (order == tab->order-1) {
664108c343cSJed Brown     if (!tab->bembedt) goto unavailable;
665108c343cSJed Brown     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
666108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
667e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
668108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
669108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
670108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
671108c343cSJed Brown     } else {                    /* Rollback and re-complete using (bet-be,be-b) */
672108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
673e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
674108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr);
675108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
676108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
677108c343cSJed Brown     }
678108c343cSJed Brown     if (done) *done = PETSC_TRUE;
679108c343cSJed Brown     PetscFunctionReturn(0);
680108c343cSJed Brown   }
681108c343cSJed Brown   unavailable:
682108c343cSJed Brown   if (done) *done = PETSC_FALSE;
683108c343cSJed 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);
684108c343cSJed Brown   PetscFunctionReturn(0);
685108c343cSJed Brown }
686108c343cSJed Brown 
687108c343cSJed Brown #undef __FUNCT__
6888a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
6898a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
6908a381b04SJed Brown {
6918a381b04SJed Brown   TS_ARKIMEX          *ark = (TS_ARKIMEX*)ts->data;
6928a381b04SJed Brown   ARKTableau          tab  = ark->tableau;
6938a381b04SJed Brown   const PetscInt      s    = tab->s;
6948a381b04SJed Brown   const PetscReal     *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
695406d0ec2SJed Brown   PetscScalar         *w   = ark->work;
696e817cc15SEmil Constantinescu   Vec                 *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,W = ark->Work,Z = ark->Z;
697108c343cSJed Brown   TSAdapt             adapt;
6988a381b04SJed Brown   SNES                snes;
699108c343cSJed Brown   PetscInt            i,j,its,lits,reject,next_scheme;
700cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
701108c343cSJed Brown   PetscReal           t;
702108c343cSJed Brown   PetscBool           accept;
7038a381b04SJed Brown   PetscErrorCode      ierr;
7048a381b04SJed Brown 
7058a381b04SJed Brown   PetscFunctionBegin;
706e817cc15SEmil Constantinescu   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage) {
707e817cc15SEmil Constantinescu     PetscReal valid_time;
708e817cc15SEmil Constantinescu     PetscBool isvalid;
709e817cc15SEmil Constantinescu     ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,
710e817cc15SEmil Constantinescu                                           explicit_stage_time_id,
711e817cc15SEmil Constantinescu                                           valid_time,
712e817cc15SEmil Constantinescu                                           isvalid);
713e817cc15SEmil Constantinescu     CHKERRQ(ierr);
714e817cc15SEmil Constantinescu     if (!isvalid || valid_time != ts->ptime) {
715e817cc15SEmil Constantinescu       TS             ts_start;
716e817cc15SEmil Constantinescu       SNES           snes_start;
717740132f1SEmil Constantinescu       DM dm;
718740132f1SEmil Constantinescu       PetscReal atol;
719740132f1SEmil Constantinescu       Vec vatol;
720740132f1SEmil Constantinescu       PetscReal rtol;
721740132f1SEmil Constantinescu       Vec vrtol;
72219436ca2SJed Brown 
72319436ca2SJed Brown       ierr = TSCreate(PETSC_COMM_WORLD,&ts_start);CHKERRQ(ierr);
72419436ca2SJed Brown       ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr);
72519436ca2SJed Brown       ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr);
726e817cc15SEmil Constantinescu       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
727740132f1SEmil Constantinescu       ierr = TSSetDM(ts_start,dm);CHKERRQ(ierr);
728e817cc15SEmil Constantinescu       ts_start->adapt=ts->adapt;
729740132f1SEmil Constantinescu       PetscObjectReference((PetscObject)ts_start->adapt);
730e817cc15SEmil Constantinescu       ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
731e817cc15SEmil Constantinescu       ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr);
732e817cc15SEmil Constantinescu       ierr = TSSetDuration(ts_start,1,ts->time_step);CHKERRQ(ierr);
733740132f1SEmil Constantinescu       ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr);
734e817cc15SEmil Constantinescu       ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
735740132f1SEmil Constantinescu       ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
736e817cc15SEmil Constantinescu       ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
737e817cc15SEmil Constantinescu       ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr);
738740132f1SEmil Constantinescu       ierr = TSGetTolerances(ts,&atol,&vatol,&rtol,&vrtol);CHKERRQ(ierr);
739740132f1SEmil Constantinescu       ierr = TSSetTolerances(ts_start,atol,vatol,rtol,vrtol);CHKERRQ(ierr);
740e817cc15SEmil Constantinescu       ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
741e817cc15SEmil Constantinescu       ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr);
742740132f1SEmil Constantinescu       ts->time_step = ts_start->time_step;
743740132f1SEmil Constantinescu       ts->steps++;
744e817cc15SEmil Constantinescu       ierr = VecCopy(((TS_ARKIMEX *)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr);
745740132f1SEmil Constantinescu       ierr = TSDestroy(&ts_start);CHKERRQ(ierr);
746740132f1SEmil Constantinescu       ierr = TSSetSNES(ts,snes_start);CHKERRQ(ierr);
747e817cc15SEmil Constantinescu     }
748e817cc15SEmil Constantinescu   }
749e817cc15SEmil Constantinescu 
7508a381b04SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
751cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
7528a381b04SJed Brown   t = ts->ptime;
753108c343cSJed Brown   accept = PETSC_TRUE;
754108c343cSJed Brown   ark->status = TS_STEP_INCOMPLETE;
7558a381b04SJed Brown 
756e817cc15SEmil Constantinescu 
75797335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
758108c343cSJed Brown     PetscReal h = ts->time_step;
759b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
7608a381b04SJed Brown     for (i=0; i<s; i++) {
7618a381b04SJed Brown       if (At[i*s+i] == 0) {           /* This stage is explicit */
7628a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
763e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7648a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
7658a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7668a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7678a381b04SJed Brown       } else {
7688a381b04SJed Brown         ark->stage_time = t + h*ct[i];
769b296d7d5SJed Brown         ark->scoeff = 1./At[i*s+i];
770b8123daeSJed Brown         ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
7718a381b04SJed Brown         /* Affine part */
7728a381b04SJed Brown         ierr = VecZeroEntries(W);CHKERRQ(ierr);
7738a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7748a381b04SJed Brown         ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
775b296d7d5SJed Brown         ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr);
776f16577ceSEmil Constantinescu 
7778a381b04SJed Brown         /* Ydot = shift*(Y-Z) */
7788a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
779e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7804f385281SJed Brown         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
781f16577ceSEmil Constantinescu 
7828a381b04SJed Brown         /* Initial guess taken from last stage */
7838a381b04SJed Brown         ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr);
7848a381b04SJed Brown         ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
785e817cc15SEmil Constantinescu         ierr = (ts->ops->snesfunction)(snes,Y[i],W,ts);CHKERRQ(ierr);
7868a381b04SJed Brown         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
7878a381b04SJed Brown         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
7885ef26d82SJed Brown         ts->snes_its += its; ts->ksp_its += lits;
789ad6bc421SBarry Smith         ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
79097335746SJed Brown         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
79197335746SJed Brown         if (!accept) goto reject_step;
7928a381b04SJed Brown       }
793e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) {
794e817cc15SEmil Constantinescu         if (i==0 && tab->explicit_first_stage) {
795e817cc15SEmil Constantinescu           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);
796e817cc15SEmil Constantinescu         } else {
797e817cc15SEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
798e817cc15SEmil Constantinescu         }
799e817cc15SEmil Constantinescu       } else {
8008a381b04SJed Brown         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
8014cc180ffSJed Brown         ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
802e817cc15SEmil Constantinescu         ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr);
8034cc180ffSJed Brown         if (ark->imex) {
8048a381b04SJed Brown           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
8054cc180ffSJed Brown         } else {
8064cc180ffSJed Brown           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
8074cc180ffSJed Brown         }
8088a381b04SJed Brown       }
809e817cc15SEmil Constantinescu     }
810108c343cSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
811108c343cSJed Brown     ark->status = TS_STEP_PENDING;
8128a381b04SJed Brown 
813108c343cSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
814ad6bc421SBarry Smith     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
815108c343cSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
816108c343cSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
817108c343cSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
818108c343cSJed Brown     if (accept) {
819108c343cSJed Brown       /* ignore next_scheme for now */
8208a381b04SJed Brown       ts->ptime += ts->time_step;
821cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
8228a381b04SJed Brown       ts->steps++;
823e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) {/* save the initial slope for the next step*/
824e817cc15SEmil Constantinescu         ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
825e817cc15SEmil Constantinescu       }
826108c343cSJed Brown       ark->status = TS_STEP_COMPLETE;
827e817cc15SEmil Constantinescu       if (tab->explicit_first_stage) {
828e817cc15SEmil Constantinescu         ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
829e817cc15SEmil Constantinescu       }
830e817cc15SEmil Constantinescu 
831108c343cSJed Brown       break;
832108c343cSJed Brown     } else {                    /* Roll back the current step */
8332c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*bt[j];
834108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotI);CHKERRQ(ierr);
8352c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
836108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotRHS);CHKERRQ(ierr);
837108c343cSJed Brown       ts->time_step = next_time_step;
838108c343cSJed Brown       ark->status = TS_STEP_INCOMPLETE;
839108c343cSJed Brown     }
840476b6736SJed Brown     reject_step: continue;
841108c343cSJed Brown   }
842b2ce242eSJed Brown   if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
8438a381b04SJed Brown   PetscFunctionReturn(0);
8448a381b04SJed Brown }
8458a381b04SJed Brown 
846cd652676SJed Brown #undef __FUNCT__
847cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
848cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
849cd652676SJed Brown {
850cd652676SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8514f385281SJed Brown   PetscInt        s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
852108c343cSJed Brown   PetscReal       h;
853108c343cSJed Brown   PetscReal       tt,t;
854cd652676SJed Brown   PetscScalar     *bt,*b;
855cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
856cd652676SJed Brown   PetscErrorCode  ierr;
857cd652676SJed Brown 
858cd652676SJed Brown   PetscFunctionBegin;
859cd652676SJed Brown   if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
860108c343cSJed Brown   switch (ark->status) {
861108c343cSJed Brown   case TS_STEP_INCOMPLETE:
862108c343cSJed Brown   case TS_STEP_PENDING:
863108c343cSJed Brown     h = ts->time_step;
864108c343cSJed Brown     t = (itime - ts->ptime)/h;
865108c343cSJed Brown     break;
866108c343cSJed Brown   case TS_STEP_COMPLETE:
867108c343cSJed Brown     h = ts->time_step_prev;
868108c343cSJed Brown     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
869108c343cSJed Brown     break;
870b9ce6d65SJed Brown   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
871108c343cSJed Brown   }
872cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
873cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
8744f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
875cd652676SJed Brown     for (i=0; i<s; i++) {
876108c343cSJed Brown       bt[i] += h * Bt[i*pinterp+j] * tt * -1.0;
877108c343cSJed Brown       b[i]  += h * B[i*pinterp+j] * tt;
878cd652676SJed Brown     }
879cd652676SJed Brown   }
880cd652676SJed 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");
881cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
882cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
883cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
884cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
885cd652676SJed Brown   PetscFunctionReturn(0);
886cd652676SJed Brown }
887cd652676SJed Brown 
8888a381b04SJed Brown /*------------------------------------------------------------*/
8898a381b04SJed Brown #undef __FUNCT__
8908a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
8918a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
8928a381b04SJed Brown {
8938a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8948a381b04SJed Brown   PetscInt        s;
8958a381b04SJed Brown   PetscErrorCode  ierr;
8968a381b04SJed Brown 
8978a381b04SJed Brown   PetscFunctionBegin;
8988a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
8998a381b04SJed Brown   s = ark->tableau->s;
9008a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
9018a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
9028a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
9038a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
9048a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
905e817cc15SEmil Constantinescu   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
9068a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
9078a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
9088a381b04SJed Brown   PetscFunctionReturn(0);
9098a381b04SJed Brown }
9108a381b04SJed Brown 
9118a381b04SJed Brown #undef __FUNCT__
9128a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
9138a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
9148a381b04SJed Brown {
9158a381b04SJed Brown   PetscErrorCode  ierr;
9168a381b04SJed Brown 
9178a381b04SJed Brown   PetscFunctionBegin;
9188a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
9198a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
9208a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr);
9218a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr);
922995b3938SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",PETSC_NULL);CHKERRQ(ierr);
9238a381b04SJed Brown   PetscFunctionReturn(0);
9248a381b04SJed Brown }
9258a381b04SJed Brown 
926d5e6173cSPeter Brune 
927d5e6173cSPeter Brune #undef __FUNCT__
928d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXGetVecs"
929d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
930d5e6173cSPeter Brune {
931d5e6173cSPeter Brune   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
932d5e6173cSPeter Brune   PetscErrorCode ierr;
933d5e6173cSPeter Brune 
934d5e6173cSPeter Brune   PetscFunctionBegin;
935d5e6173cSPeter Brune   if (Z) {
936d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
937d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
938d5e6173cSPeter Brune     } else *Z = ax->Z;
939d5e6173cSPeter Brune   }
940d5e6173cSPeter Brune   if (Ydot) {
941d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
942d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
943d5e6173cSPeter Brune     } else *Ydot = ax->Ydot;
944d5e6173cSPeter Brune   }
945d5e6173cSPeter Brune   PetscFunctionReturn(0);
946d5e6173cSPeter Brune }
947d5e6173cSPeter Brune 
948d5e6173cSPeter Brune 
949d5e6173cSPeter Brune #undef __FUNCT__
950d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXRestoreVecs"
951d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
952d5e6173cSPeter Brune {
953d5e6173cSPeter Brune   PetscErrorCode ierr;
954d5e6173cSPeter Brune 
955d5e6173cSPeter Brune   PetscFunctionBegin;
956d5e6173cSPeter Brune   if (Z) {
957d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
958d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
959d5e6173cSPeter Brune     }
960d5e6173cSPeter Brune   }
961d5e6173cSPeter Brune   if (Ydot) {
962d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
963d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
964d5e6173cSPeter Brune     }
965d5e6173cSPeter Brune   }
966d5e6173cSPeter Brune   PetscFunctionReturn(0);
967d5e6173cSPeter Brune }
968d5e6173cSPeter Brune 
9698a381b04SJed Brown /*
9708a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9718a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9728a381b04SJed Brown */
9738a381b04SJed Brown #undef __FUNCT__
9748a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
9758a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
9768a381b04SJed Brown {
9778a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
978d5e6173cSPeter Brune   DM             dm,dmsave;
979d5e6173cSPeter Brune   Vec            Z,Ydot;
980b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
9818a381b04SJed Brown   PetscErrorCode ierr;
9828a381b04SJed Brown 
9838a381b04SJed Brown   PetscFunctionBegin;
984d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
985d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
986b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
987d5e6173cSPeter Brune   dmsave = ts->dm;
988d5e6173cSPeter Brune   ts->dm = dm;
989740132f1SEmil Constantinescu 
990d5e6173cSPeter Brune   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
991e817cc15SEmil Constantinescu 
992d5e6173cSPeter Brune   ts->dm = dmsave;
993d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
9948a381b04SJed Brown   PetscFunctionReturn(0);
9958a381b04SJed Brown }
9968a381b04SJed Brown 
9978a381b04SJed Brown #undef __FUNCT__
9988a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
9998a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
10008a381b04SJed Brown {
10018a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1002d5e6173cSPeter Brune   DM             dm,dmsave;
1003d5e6173cSPeter Brune   Vec            Ydot;
1004b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10058a381b04SJed Brown   PetscErrorCode ierr;
10068a381b04SJed Brown 
10078a381b04SJed Brown   PetscFunctionBegin;
1008d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1009d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,PETSC_NULL,&Ydot);CHKERRQ(ierr);
10108a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1011d5e6173cSPeter Brune   dmsave = ts->dm;
1012d5e6173cSPeter Brune   ts->dm = dm;
1013740132f1SEmil Constantinescu 
1014b296d7d5SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,str,ark->imex);CHKERRQ(ierr);
1015740132f1SEmil Constantinescu 
1016d5e6173cSPeter Brune   ts->dm = dmsave;
1017d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,PETSC_NULL,&Ydot);CHKERRQ(ierr);
1018d5e6173cSPeter Brune   PetscFunctionReturn(0);
1019d5e6173cSPeter Brune }
1020d5e6173cSPeter Brune 
1021d5e6173cSPeter Brune #undef __FUNCT__
1022d5e6173cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSARKIMEX"
1023d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1024d5e6173cSPeter Brune {
1025d5e6173cSPeter Brune 
1026d5e6173cSPeter Brune   PetscFunctionBegin;
1027d5e6173cSPeter Brune   PetscFunctionReturn(0);
1028d5e6173cSPeter Brune }
1029d5e6173cSPeter Brune 
1030d5e6173cSPeter Brune #undef __FUNCT__
1031d5e6173cSPeter Brune #define __FUNCT__ "DMRestrictHook_TSARKIMEX"
1032d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1033d5e6173cSPeter Brune {
1034d5e6173cSPeter Brune   TS             ts = (TS)ctx;
1035d5e6173cSPeter Brune   PetscErrorCode ierr;
1036d5e6173cSPeter Brune   Vec            Z,Z_c;
1037d5e6173cSPeter Brune 
1038d5e6173cSPeter Brune   PetscFunctionBegin;
1039d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,fine,&Z,PETSC_NULL);CHKERRQ(ierr);
1040d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1041d5e6173cSPeter Brune   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1042d5e6173cSPeter Brune   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
1043d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,PETSC_NULL);CHKERRQ(ierr);
1044d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,PETSC_NULL);CHKERRQ(ierr);
10458a381b04SJed Brown   PetscFunctionReturn(0);
10468a381b04SJed Brown }
10478a381b04SJed Brown 
1048cdb298fcSPeter Brune 
1049cdb298fcSPeter Brune #undef __FUNCT__
1050cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainHook_TSARKIMEX"
1051cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1052cdb298fcSPeter Brune {
1053cdb298fcSPeter Brune 
1054cdb298fcSPeter Brune   PetscFunctionBegin;
1055cdb298fcSPeter Brune   PetscFunctionReturn(0);
1056cdb298fcSPeter Brune }
1057cdb298fcSPeter Brune 
1058cdb298fcSPeter Brune #undef __FUNCT__
1059cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX"
1060cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1061cdb298fcSPeter Brune {
1062cdb298fcSPeter Brune   TS             ts = (TS)ctx;
1063cdb298fcSPeter Brune   PetscErrorCode ierr;
1064cdb298fcSPeter Brune   Vec            Z,Z_c;
1065cdb298fcSPeter Brune 
1066cdb298fcSPeter Brune   PetscFunctionBegin;
1067cdb298fcSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,&Z,PETSC_NULL);CHKERRQ(ierr);
1068cdb298fcSPeter Brune   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1069cdb298fcSPeter Brune 
1070cdb298fcSPeter Brune   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1071cdb298fcSPeter Brune   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1072cdb298fcSPeter Brune 
1073cdb298fcSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,PETSC_NULL);CHKERRQ(ierr);
1074cdb298fcSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1075cdb298fcSPeter Brune   PetscFunctionReturn(0);
1076cdb298fcSPeter Brune }
1077cdb298fcSPeter Brune 
10788a381b04SJed Brown #undef __FUNCT__
10798a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
10808a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
10818a381b04SJed Brown {
10828a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1083f2c2a1b9SBarry Smith   ARKTableau     tab;
1084f2c2a1b9SBarry Smith   PetscInt       s;
10858a381b04SJed Brown   PetscErrorCode ierr;
1086d5e6173cSPeter Brune   DM             dm;
1087f9c1d6abSBarry Smith 
10888a381b04SJed Brown   PetscFunctionBegin;
10898a381b04SJed Brown   if (!ark->tableau) {
1090e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
10918a381b04SJed Brown   }
1092f2c2a1b9SBarry Smith   tab  = ark->tableau;
1093f2c2a1b9SBarry Smith   s    = tab->s;
10948a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
10958a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
10968a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
10978a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
10988a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
1099e817cc15SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
11008a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
11018a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
1102d5e6173cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1103d5e6173cSPeter Brune   if (dm) {
1104d5e6173cSPeter Brune     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1105cdb298fcSPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1106d5e6173cSPeter Brune   }
11078a381b04SJed Brown   PetscFunctionReturn(0);
11088a381b04SJed Brown }
11098a381b04SJed Brown /*------------------------------------------------------------*/
11108a381b04SJed Brown 
11118a381b04SJed Brown #undef __FUNCT__
11128a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
11138a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
11148a381b04SJed Brown {
11154cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11168a381b04SJed Brown   PetscErrorCode ierr;
11178a381b04SJed Brown   char           arktype[256];
11188a381b04SJed Brown 
11198a381b04SJed Brown   PetscFunctionBegin;
11208a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
11218a381b04SJed Brown   {
11228a381b04SJed Brown     ARKTableauLink link;
11238a381b04SJed Brown     PetscInt       count,choice;
11248a381b04SJed Brown     PetscBool      flg;
11258a381b04SJed Brown     const char     **namelist;
11268caf3d72SBarry Smith     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr);
11278a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
11288a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
11298a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11308a381b04SJed Brown     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
11318a381b04SJed Brown     ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
11328a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
11334cc180ffSJed Brown     flg = (PetscBool)!ark->imex;
11344cc180ffSJed Brown     ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
11354cc180ffSJed Brown     ark->imex = (PetscBool)!flg;
1136d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
11378a381b04SJed Brown   }
11388a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
11398a381b04SJed Brown   PetscFunctionReturn(0);
11408a381b04SJed Brown }
11418a381b04SJed Brown 
11428a381b04SJed Brown #undef __FUNCT__
11438a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
11448a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
11458a381b04SJed Brown {
1146257d2499SJed Brown   PetscErrorCode ierr;
1147f1d86077SJed Brown   PetscInt       i;
1148f1d86077SJed Brown   size_t         left,count;
11498a381b04SJed Brown   char           *p;
11508a381b04SJed Brown 
11518a381b04SJed Brown   PetscFunctionBegin;
1152f1d86077SJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
1153f1d86077SJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
11548a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
11558a381b04SJed Brown     left -= count;
11568a381b04SJed Brown     p += count;
11578a381b04SJed Brown     *p++ = ' ';
11588a381b04SJed Brown   }
11598a381b04SJed Brown   p[i ? 0 : -1] = 0;
11608a381b04SJed Brown   PetscFunctionReturn(0);
11618a381b04SJed Brown }
11628a381b04SJed Brown 
11638a381b04SJed Brown #undef __FUNCT__
11648a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
11658a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
11668a381b04SJed Brown {
11678a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11688a381b04SJed Brown   ARKTableau     tab = ark->tableau;
11698a381b04SJed Brown   PetscBool      iascii;
11708a381b04SJed Brown   PetscErrorCode ierr;
1171559eea31SJed Brown   TSAdapt        adapt;
11728a381b04SJed Brown 
11738a381b04SJed Brown   PetscFunctionBegin;
1174251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11758a381b04SJed Brown   if (iascii) {
117619fd82e9SBarry Smith     TSARKIMEXType arktype;
11778a381b04SJed Brown     char buf[512];
11788a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
11798a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
11808caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
118131f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
11828caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1183e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate?"yes":"no");CHKERRQ(ierr);
1184e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage?"yes":"no");CHKERRQ(ierr);
1185e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit?"yes":"no");CHKERRQ(ierr);
118631f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
11878a381b04SJed Brown   }
1188ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
1189559eea31SJed Brown   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
1190d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
11918a381b04SJed Brown   PetscFunctionReturn(0);
11928a381b04SJed Brown }
11938a381b04SJed Brown 
11948a381b04SJed Brown #undef __FUNCT__
1195f2c2a1b9SBarry Smith #define __FUNCT__ "TSLoad_ARKIMEX"
1196f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1197f2c2a1b9SBarry Smith {
1198f2c2a1b9SBarry Smith   PetscErrorCode ierr;
1199f2c2a1b9SBarry Smith   SNES           snes;
1200ad6bc421SBarry Smith   TSAdapt        tsadapt;
1201f2c2a1b9SBarry Smith 
1202f2c2a1b9SBarry Smith   PetscFunctionBegin;
1203ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&tsadapt);CHKERRQ(ierr);
1204ad6bc421SBarry Smith   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
1205f2c2a1b9SBarry Smith   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1206f2c2a1b9SBarry Smith   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1207ad6bc421SBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
1208ad6bc421SBarry Smith   ierr = SNESSetFunction(snes,PETSC_NULL,PETSC_NULL,ts);CHKERRQ(ierr);
1209ad6bc421SBarry Smith   ierr = SNESSetJacobian(snes,PETSC_NULL,PETSC_NULL,PETSC_NULL,ts);CHKERRQ(ierr);
1210f2c2a1b9SBarry Smith   PetscFunctionReturn(0);
1211f2c2a1b9SBarry Smith }
1212f2c2a1b9SBarry Smith 
1213f2c2a1b9SBarry Smith #undef __FUNCT__
12148a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
12158a381b04SJed Brown /*@C
12168a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
12178a381b04SJed Brown 
12188a381b04SJed Brown   Logically collective
12198a381b04SJed Brown 
12208a381b04SJed Brown   Input Parameter:
12218a381b04SJed Brown +  ts - timestepping context
12228a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
12238a381b04SJed Brown 
12248a381b04SJed Brown   Level: intermediate
12258a381b04SJed Brown 
1226020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
12278a381b04SJed Brown @*/
122819fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
12298a381b04SJed Brown {
12308a381b04SJed Brown   PetscErrorCode ierr;
12318a381b04SJed Brown 
12328a381b04SJed Brown   PetscFunctionBegin;
12338a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
123419fd82e9SBarry Smith   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
12358a381b04SJed Brown   PetscFunctionReturn(0);
12368a381b04SJed Brown }
12378a381b04SJed Brown 
12388a381b04SJed Brown #undef __FUNCT__
12398a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
12408a381b04SJed Brown /*@C
12418a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
12428a381b04SJed Brown 
12438a381b04SJed Brown   Logically collective
12448a381b04SJed Brown 
12458a381b04SJed Brown   Input Parameter:
12468a381b04SJed Brown .  ts - timestepping context
12478a381b04SJed Brown 
12488a381b04SJed Brown   Output Parameter:
12498a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
12508a381b04SJed Brown 
12518a381b04SJed Brown   Level: intermediate
12528a381b04SJed Brown 
12538a381b04SJed Brown .seealso: TSARKIMEXGetType()
12548a381b04SJed Brown @*/
125519fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
12568a381b04SJed Brown {
12578a381b04SJed Brown   PetscErrorCode ierr;
12588a381b04SJed Brown 
12598a381b04SJed Brown   PetscFunctionBegin;
12608a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
126119fd82e9SBarry Smith   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
12628a381b04SJed Brown   PetscFunctionReturn(0);
12638a381b04SJed Brown }
12648a381b04SJed Brown 
12654cc180ffSJed Brown #undef __FUNCT__
12664cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
12674cc180ffSJed Brown /*@C
12684cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
12694cc180ffSJed Brown 
12704cc180ffSJed Brown   Logically collective
12714cc180ffSJed Brown 
12724cc180ffSJed Brown   Input Parameter:
12734cc180ffSJed Brown +  ts - timestepping context
12744cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
12754cc180ffSJed Brown 
12764cc180ffSJed Brown   Level: intermediate
12774cc180ffSJed Brown 
12784cc180ffSJed Brown .seealso: TSARKIMEXGetType()
12794cc180ffSJed Brown @*/
12804cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
12814cc180ffSJed Brown {
12824cc180ffSJed Brown   PetscErrorCode ierr;
12834cc180ffSJed Brown 
12844cc180ffSJed Brown   PetscFunctionBegin;
12854cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12864cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
12874cc180ffSJed Brown   PetscFunctionReturn(0);
12884cc180ffSJed Brown }
12894cc180ffSJed Brown 
12908a381b04SJed Brown EXTERN_C_BEGIN
12918a381b04SJed Brown #undef __FUNCT__
12928a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
129319fd82e9SBarry Smith PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
12948a381b04SJed Brown {
12958a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
12968a381b04SJed Brown   PetscErrorCode ierr;
12978a381b04SJed Brown 
12988a381b04SJed Brown   PetscFunctionBegin;
1299f2c2a1b9SBarry Smith   if (!ark->tableau) {
1300f2c2a1b9SBarry Smith     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1301f2c2a1b9SBarry Smith   }
13028a381b04SJed Brown   *arktype = ark->tableau->name;
13038a381b04SJed Brown   PetscFunctionReturn(0);
13048a381b04SJed Brown }
13058a381b04SJed Brown #undef __FUNCT__
13068a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
130719fd82e9SBarry Smith PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
13088a381b04SJed Brown {
13098a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13108a381b04SJed Brown   PetscErrorCode ierr;
13118a381b04SJed Brown   PetscBool match;
13128a381b04SJed Brown   ARKTableauLink link;
13138a381b04SJed Brown 
13148a381b04SJed Brown   PetscFunctionBegin;
13158a381b04SJed Brown   if (ark->tableau) {
13168a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
13178a381b04SJed Brown     if (match) PetscFunctionReturn(0);
13188a381b04SJed Brown   }
13198a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
13208a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
13218a381b04SJed Brown     if (match) {
13228a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
13238a381b04SJed Brown       ark->tableau = &link->tab;
13248a381b04SJed Brown       PetscFunctionReturn(0);
13258a381b04SJed Brown     }
13268a381b04SJed Brown   }
13278a381b04SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
13288a381b04SJed Brown   PetscFunctionReturn(0);
13298a381b04SJed Brown }
13304cc180ffSJed Brown #undef __FUNCT__
13314cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
13324cc180ffSJed Brown PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
13334cc180ffSJed Brown {
13344cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13354cc180ffSJed Brown 
13364cc180ffSJed Brown   PetscFunctionBegin;
13374cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
13384cc180ffSJed Brown   PetscFunctionReturn(0);
13394cc180ffSJed Brown }
13408a381b04SJed Brown EXTERN_C_END
13418a381b04SJed Brown 
13428a381b04SJed Brown /* ------------------------------------------------------------ */
13438a381b04SJed Brown /*MC
1344a4386c9eSJed Brown       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
13458a381b04SJed Brown 
1346fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1347fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1348fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1349fca742c7SJed Brown 
1350fca742c7SJed Brown   Notes:
1351a4386c9eSJed Brown   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1352c8058688SBarry Smith 
1353a4386c9eSJed 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).
1354fca742c7SJed Brown 
13558a381b04SJed Brown   Level: beginner
13568a381b04SJed Brown 
1357c8058688SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
1358a4386c9eSJed Brown            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
13598a381b04SJed Brown 
13608a381b04SJed Brown M*/
13618a381b04SJed Brown EXTERN_C_BEGIN
13628a381b04SJed Brown #undef __FUNCT__
13638a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
13648a381b04SJed Brown PetscErrorCode  TSCreate_ARKIMEX(TS ts)
13658a381b04SJed Brown {
13668a381b04SJed Brown   TS_ARKIMEX     *th;
13678a381b04SJed Brown   PetscErrorCode ierr;
13688a381b04SJed Brown 
13698a381b04SJed Brown   PetscFunctionBegin;
13708a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
13718a381b04SJed Brown   ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr);
13728a381b04SJed Brown #endif
13738a381b04SJed Brown 
13748a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
13758a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
13768a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
1377f2c2a1b9SBarry Smith   ts->ops->load           = TSLoad_ARKIMEX;
13788a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
13798a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
1380cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1381108c343cSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
13828a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
13838a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
13848a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
13858a381b04SJed Brown 
13868a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
13878a381b04SJed Brown   ts->data = (void*)th;
13884cc180ffSJed Brown   th->imex = PETSC_TRUE;
13898a381b04SJed Brown 
13908a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
13918a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
13924cc180ffSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
13938a381b04SJed Brown   PetscFunctionReturn(0);
13948a381b04SJed Brown }
13958a381b04SJed Brown EXTERN_C_END
1396