xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 9eef816d3d04fdaa28de0101cfc8616bb284fa90)
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*/
131e25c274SJed Brown #include <petscdm.h>
148a381b04SJed Brown 
1519fd82e9SBarry Smith static TSARKIMEXType  TSARKIMEXDefault = TSARKIMEX3;
168a381b04SJed Brown static PetscBool      TSARKIMEXRegisterAllCalled;
178a381b04SJed Brown static PetscBool      TSARKIMEXPackageInitialized;
18e817cc15SEmil Constantinescu static PetscInt       explicit_stage_time_id;
1956dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec);
208a381b04SJed Brown 
218a381b04SJed Brown typedef struct _ARKTableau *ARKTableau;
228a381b04SJed Brown struct _ARKTableau {
238a381b04SJed Brown   char      *name;
244f385281SJed Brown   PetscInt  order;                /* Classical approximation order of the method */
254f385281SJed Brown   PetscInt  s;                    /* Number of stages */
26e817cc15SEmil Constantinescu   PetscBool stiffly_accurate;     /* The implicit part is stiffly accurate*/
27e817cc15SEmil Constantinescu   PetscBool FSAL_implicit;        /* The implicit part is FSAL*/
28e817cc15SEmil Constantinescu   PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/
294f385281SJed Brown   PetscInt  pinterp;              /* Interpolation order */
304f385281SJed Brown   PetscReal *At,*bt,*ct;          /* Stiff tableau */
318a381b04SJed Brown   PetscReal *A,*b,*c;             /* Non-stiff tableau */
32108c343cSJed Brown   PetscReal *bembedt,*bembed;     /* Embedded formula of order one less (order-1) */
33cd652676SJed Brown   PetscReal *binterpt,*binterp;   /* Dense output formula */
34108c343cSJed Brown   PetscReal ccfl;                 /* Placeholder for CFL coefficient relative to forward Euler */
358a381b04SJed Brown };
368a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
378a381b04SJed Brown struct _ARKTableauLink {
388a381b04SJed Brown   struct _ARKTableau tab;
398a381b04SJed Brown   ARKTableauLink     next;
408a381b04SJed Brown };
418a381b04SJed Brown static ARKTableauLink ARKTableauList;
428a381b04SJed Brown 
438a381b04SJed Brown typedef struct {
448a381b04SJed Brown   ARKTableau   tableau;
458a381b04SJed Brown   Vec          *Y;               /* States computed during the step */
468a381b04SJed Brown   Vec          *YdotI;           /* Time derivatives for the stiff part */
478a381b04SJed Brown   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
48*9eef816dSJed Brown   PetscBool    prev_step_valid;  /* Stored previous step (Y_prev, YdotI_prev, YdotRHS_prev) is valid */
4956dcabbaSDebojyoti Ghosh   Vec          *Y_prev;          /* States computed during the previous time step */
5056dcabbaSDebojyoti Ghosh   Vec          *YdotI_prev;      /* Time derivatives for the stiff part for the previous time step*/
5156dcabbaSDebojyoti Ghosh   Vec          *YdotRHS_prev;    /* Function evaluations for the non-stiff part for the previous time step*/
52e817cc15SEmil Constantinescu   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
538a381b04SJed Brown   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
548a381b04SJed Brown   Vec          Work;             /* Generic work vector */
558a381b04SJed Brown   Vec          Z;                /* Ydot = shift(Y-Z) */
568a381b04SJed Brown   PetscScalar  *work;            /* Scalar work */
57b296d7d5SJed Brown   PetscReal    scoeff;           /* shift = scoeff/dt */
588a381b04SJed Brown   PetscReal    stage_time;
594cc180ffSJed Brown   PetscBool    imex;
6056dcabbaSDebojyoti Ghosh   PetscBool    init_guess_extrp; /* Extrapolate initial guess from previous time-step stage values */
61108c343cSJed Brown   TSStepStatus status;
628a381b04SJed Brown } TS_ARKIMEX;
631f80e275SEmil Constantinescu /*MC
641f80e275SEmil Constantinescu      TSARKIMEXARS122 - Second order ARK IMEX scheme.
658a381b04SJed Brown 
661f80e275SEmil Constantinescu      This method has one explicit stage and one implicit stage.
671f80e275SEmil Constantinescu 
681f80e275SEmil Constantinescu      References:
691f80e275SEmil 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.
701f80e275SEmil Constantinescu 
711f80e275SEmil Constantinescu      Level: advanced
721f80e275SEmil Constantinescu 
731f80e275SEmil Constantinescu .seealso: TSARKIMEX
741f80e275SEmil Constantinescu M*/
751f80e275SEmil Constantinescu /*MC
761f80e275SEmil Constantinescu      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.
771f80e275SEmil Constantinescu 
781f80e275SEmil Constantinescu      This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu.
791f80e275SEmil Constantinescu 
801f80e275SEmil Constantinescu      Level: advanced
811f80e275SEmil Constantinescu 
821f80e275SEmil Constantinescu .seealso: TSARKIMEX
831f80e275SEmil Constantinescu M*/
841f80e275SEmil Constantinescu /*MC
851f80e275SEmil Constantinescu      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.
861f80e275SEmil Constantinescu 
871f80e275SEmil Constantinescu      This method has two implicit stages, and L-stable implicit scheme.
881f80e275SEmil Constantinescu 
891f80e275SEmil Constantinescu     References:
901f80e275SEmil 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
911f80e275SEmil Constantinescu 
921f80e275SEmil Constantinescu      Level: advanced
931f80e275SEmil Constantinescu 
941f80e275SEmil Constantinescu .seealso: TSARKIMEX
951f80e275SEmil Constantinescu M*/
961f80e275SEmil Constantinescu /*MC
97e817cc15SEmil Constantinescu      TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.
98e817cc15SEmil Constantinescu 
99e817cc15SEmil Constantinescu      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.
100e817cc15SEmil Constantinescu 
101e817cc15SEmil Constantinescu      Level: advanced
102e817cc15SEmil Constantinescu 
103e817cc15SEmil Constantinescu .seealso: TSARKIMEX
104e817cc15SEmil Constantinescu M*/
105e817cc15SEmil Constantinescu /*MC
1061f80e275SEmil Constantinescu      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.
1071f80e275SEmil Constantinescu 
1081f80e275SEmil 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.
1091f80e275SEmil Constantinescu 
1101f80e275SEmil Constantinescu      Level: advanced
1111f80e275SEmil Constantinescu 
1121f80e275SEmil Constantinescu .seealso: TSARKIMEX
1131f80e275SEmil Constantinescu M*/
11464f491ddSJed Brown /*MC
11564f491ddSJed Brown      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
11664f491ddSJed Brown 
117617a39beSEmil 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.
11864f491ddSJed Brown 
119b330ce4dSSatish Balay      Level: advanced
120b330ce4dSSatish Balay 
12164f491ddSJed Brown .seealso: TSARKIMEX
12264f491ddSJed Brown M*/
12364f491ddSJed Brown /*MC
12464f491ddSJed Brown      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
12564f491ddSJed Brown 
12664f491ddSJed Brown      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
12764f491ddSJed Brown 
128b330ce4dSSatish Balay      Level: advanced
129b330ce4dSSatish Balay 
13064f491ddSJed Brown .seealso: TSARKIMEX
13164f491ddSJed Brown M*/
13264f491ddSJed Brown /*MC
1336cf0794eSJed Brown      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.
1346cf0794eSJed Brown 
1356cf0794eSJed Brown      This method has three implicit stages.
1366cf0794eSJed Brown 
1376cf0794eSJed Brown      References:
1386cf0794eSJed 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
1396cf0794eSJed Brown 
1406cf0794eSJed Brown      This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375
1416cf0794eSJed Brown 
1426cf0794eSJed Brown      Level: advanced
1436cf0794eSJed Brown 
1446cf0794eSJed Brown .seealso: TSARKIMEX
1456cf0794eSJed Brown M*/
1466cf0794eSJed Brown /*MC
14764f491ddSJed Brown      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
14864f491ddSJed Brown 
14964f491ddSJed Brown      This method has one explicit stage and three implicit stages.
15064f491ddSJed Brown 
15164f491ddSJed Brown      References:
15264f491ddSJed Brown      Kennedy and Carpenter 2003.
15364f491ddSJed Brown 
154b330ce4dSSatish Balay      Level: advanced
155b330ce4dSSatish Balay 
15664f491ddSJed Brown .seealso: TSARKIMEX
15764f491ddSJed Brown M*/
15864f491ddSJed Brown /*MC
1596cf0794eSJed Brown      TSARKIMEXARS443 - Third order ARK IMEX scheme.
1606cf0794eSJed Brown 
1616cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1626cf0794eSJed Brown 
1636cf0794eSJed Brown      References:
1646cf0794eSJed 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.
1656cf0794eSJed Brown 
1666cf0794eSJed Brown      This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375
1676cf0794eSJed Brown 
1686cf0794eSJed Brown      Level: advanced
1696cf0794eSJed Brown 
1706cf0794eSJed Brown .seealso: TSARKIMEX
1716cf0794eSJed Brown M*/
1726cf0794eSJed Brown /*MC
1736cf0794eSJed Brown      TSARKIMEXBPR3 - Third order ARK IMEX scheme.
1746cf0794eSJed Brown 
1756cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1766cf0794eSJed Brown 
1776cf0794eSJed Brown      References:
1786cf0794eSJed Brown      This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375
1796cf0794eSJed Brown 
1806cf0794eSJed Brown      Level: advanced
1816cf0794eSJed Brown 
1826cf0794eSJed Brown .seealso: TSARKIMEX
1836cf0794eSJed Brown M*/
1846cf0794eSJed Brown /*MC
18564f491ddSJed Brown      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
18664f491ddSJed Brown 
18764f491ddSJed Brown      This method has one explicit stage and four implicit stages.
18864f491ddSJed Brown 
18964f491ddSJed Brown      References:
19064f491ddSJed Brown      Kennedy and Carpenter 2003.
19164f491ddSJed Brown 
192b330ce4dSSatish Balay      Level: advanced
193b330ce4dSSatish Balay 
19464f491ddSJed Brown .seealso: TSARKIMEX
19564f491ddSJed Brown M*/
19664f491ddSJed Brown /*MC
19764f491ddSJed Brown      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
19864f491ddSJed Brown 
19964f491ddSJed Brown      This method has one explicit stage and five implicit stages.
20064f491ddSJed Brown 
20164f491ddSJed Brown      References:
20264f491ddSJed Brown      Kennedy and Carpenter 2003.
20364f491ddSJed Brown 
204b330ce4dSSatish Balay      Level: advanced
205b330ce4dSSatish Balay 
20664f491ddSJed Brown .seealso: TSARKIMEX
20764f491ddSJed Brown M*/
20864f491ddSJed Brown 
2098a381b04SJed Brown #undef __FUNCT__
2108a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
2118a381b04SJed Brown /*@C
2128a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
2138a381b04SJed Brown 
214fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
2158a381b04SJed Brown 
2168a381b04SJed Brown   Level: advanced
2178a381b04SJed Brown 
2188a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
2198a381b04SJed Brown 
2208a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
2218a381b04SJed Brown @*/
2228a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
2238a381b04SJed Brown {
2248a381b04SJed Brown   PetscErrorCode ierr;
2258a381b04SJed Brown 
2268a381b04SJed Brown   PetscFunctionBegin;
2278a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
2288a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
229e817cc15SEmil Constantinescu 
230e817cc15SEmil Constantinescu   {
231e817cc15SEmil Constantinescu     const PetscReal
232e817cc15SEmil Constantinescu       A[3][3] = {{0.0,0.0,0.0},
233e817cc15SEmil Constantinescu                  {0.0,0.0,0.0},
234748ad121SEmil Constantinescu                  {0.0,0.5,0.0}},
235e817cc15SEmil Constantinescu       At[3][3] = {{1.0,0.0,0.0},
236e817cc15SEmil Constantinescu                   {0.0,0.5,0.0},
237e817cc15SEmil Constantinescu                   {0.0,0.5,0.5}},
238e817cc15SEmil Constantinescu       b[3]       = {0.0,0.5,0.5},
239e817cc15SEmil Constantinescu       bembedt[3] = {1.0,0.0,0.0};
2400298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
241e817cc15SEmil Constantinescu   }
2428a381b04SJed Brown   {
2438a381b04SJed Brown     const PetscReal
2441f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2451f80e275SEmil Constantinescu                  {0.5,0.0}},
2461f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2471f80e275SEmil Constantinescu                   {0.0,0.5}},
2481f80e275SEmil Constantinescu       b[2]       = {0.0,1.0},
2491f80e275SEmil Constantinescu       bembedt[2] = {0.5,0.5};
2501f80e275SEmil 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*/
2510298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2521f80e275SEmil Constantinescu   }
2531f80e275SEmil Constantinescu   {
2541f80e275SEmil Constantinescu     const PetscReal
2551f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2561f80e275SEmil Constantinescu                  {1.0,0.0}},
2571f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2581f80e275SEmil Constantinescu                   {0.5,0.5}},
2591f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
2601f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0};
2611f80e275SEmil 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*/
2620298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2631f80e275SEmil Constantinescu   }
2641f80e275SEmil Constantinescu   {
265da80777bSKarl 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   */
2661f80e275SEmil Constantinescu     const PetscReal
2671f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2681f80e275SEmil Constantinescu                  {1.0,0.0}},
269da80777bSKarl Rupp       At[2][2] = {{0.2928932188134524755992,0.0},
270da80777bSKarl Rupp                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
2711f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
2721f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0},
273da80777bSKarl Rupp       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
274da80777bSKarl Rupp                         {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}},
2751f80e275SEmil Constantinescu       binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
2760298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);CHKERRQ(ierr);
2771f80e275SEmil Constantinescu   }
2781f80e275SEmil Constantinescu   {
279da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
280da80777bSKarl Rupp     const PetscReal
2818a381b04SJed Brown       A[3][3] = {{0,0,0},
282da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
283617a39beSEmil Constantinescu                  {0.5,0.5,0}},
284da80777bSKarl Rupp       At[3][3] = {{0,0,0},
285da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
286da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
287da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
288da80777bSKarl Rupp       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
289da80777bSKarl Rupp                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
290da80777bSKarl Rupp                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
2910298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
2921f80e275SEmil Constantinescu   }
2931f80e275SEmil Constantinescu   {
294da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
295da80777bSKarl Rupp     const PetscReal
2961f80e275SEmil Constantinescu       A[3][3] = {{0,0,0},
297da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
2988a381b04SJed Brown                  {0.75,0.25,0}},
299da80777bSKarl Rupp       At[3][3] = {{0,0,0},
300da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
301da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
302da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
303da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
304da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
305da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3060298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
3078a381b04SJed Brown   }
30806db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
309da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
310da80777bSKarl Rupp     const PetscReal
311da80777bSKarl Rupp       A[3][3] = {{0,0,0},
312da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
313da80777bSKarl Rupp                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
314da80777bSKarl Rupp       At[3][3] = {{0,0,0},
315da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
316da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
317da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
318da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
319da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
320da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3210298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
322a3a57f36SJed Brown   }
3236cf0794eSJed Brown   {                             /* Optimal for linear implicit part */
3246cf0794eSJed Brown     const PetscReal
3256cf0794eSJed Brown       A[3][3] = {{0,0,0},
3266cf0794eSJed Brown                  {0.5,0,0},
3276cf0794eSJed Brown                  {0.5,0.5,0}},
3286cf0794eSJed Brown       At[3][3] = {{0.25,0,0},
3296cf0794eSJed Brown                   {0,0.25,0},
3306cf0794eSJed Brown                   {1./3,1./3,1./3}};
3310298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr);
3326cf0794eSJed Brown   }
333a3a57f36SJed Brown   {
334a3a57f36SJed Brown     const PetscReal
335a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
3364040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
3374040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
3384040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
339a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
3404040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
3414040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
3424040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
343cc46b9d1SJed Brown       bembedt[4]     = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
3444040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
3454040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
3464040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
3474040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
3480298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
349a3a57f36SJed Brown   }
350a3a57f36SJed Brown   {
351a3a57f36SJed Brown     const PetscReal
352e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3536cf0794eSJed Brown                  {1./2,0,0,0,0},
3546cf0794eSJed Brown                  {11./18,1./18,0,0,0},
3556cf0794eSJed Brown                  {5./6,-5./6,.5,0,0},
3566cf0794eSJed Brown                  {1./4,7./4,3./4,-7./4,0}},
3576cf0794eSJed Brown       At[5][5] = {{0,0,0,0,0},
3586cf0794eSJed Brown                   {0,1./2,0,0,0},
3596cf0794eSJed Brown                   {0,1./6,1./2,0,0},
3606cf0794eSJed Brown                   {0,-1./2,1./2,1./2,0},
361108c343cSJed Brown                   {0,3./2,-3./2,1./2,1./2}},
3620298fd71SBarry Smith     *bembedt = NULL;
3630298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
3646cf0794eSJed Brown   }
3656cf0794eSJed Brown   {
3666cf0794eSJed Brown     const PetscReal
367e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3686cf0794eSJed Brown                  {1,0,0,0,0},
3696cf0794eSJed Brown                  {4./9,2./9,0,0,0},
3706cf0794eSJed Brown                  {1./4,0,3./4,0,0},
3716cf0794eSJed Brown                  {1./4,0,3./5,0,0}},
372e74514c0SSatish Balay       At[5][5] = {{0,0,0,0,0},
3736cf0794eSJed Brown                   {.5,.5,0,0,0},
3746cf0794eSJed Brown                   {5./18,-1./9,.5,0,0},
3756cf0794eSJed Brown                   {.5,0,0,.5,0},
376108c343cSJed Brown                   {.25,0,.75,-.5,.5}},
3770298fd71SBarry Smith     *bembedt = NULL;
3780298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
3796cf0794eSJed Brown   }
3806cf0794eSJed Brown   {
3816cf0794eSJed Brown     const PetscReal
382a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
383a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
3844040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
3854040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
3864040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
3874040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
388a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
389a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
3904040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
3914040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
3924040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
3934040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
394cc46b9d1SJed Brown       bembedt[6]     = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
3954040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
396cd652676SJed Brown                         {0,0,0},
3974040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
3984040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
3994040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
4004040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
4010298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
402a3a57f36SJed Brown   }
403a3a57f36SJed Brown   {
404a3a57f36SJed Brown     const PetscReal
405a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
406a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
4074040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
4084040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
4094040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
4104040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
4114040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
4124040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
413a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
4144040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
4154040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
4164040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
4174040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
4184040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
4194040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
4204040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
421cc46b9d1SJed Brown       bembedt[8]     = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
4224040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614.,  43486358583215./12773830924787., -9257016797708./5021505065439.},
423cd652676SJed Brown                         {0,  0, 0                            },
424cd652676SJed Brown                         {0,  0, 0                            },
4254040e9f2SJed Brown                         {65168852399939./7868540260826.,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
4264040e9f2SJed Brown                         {15494834004392./5936557850923.,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
4274040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473., 30029262896817./10175596800299.},
4284040e9f2SJed Brown                         {-19024464361622./5461577185407.,  115839755401235./10719374521269., -26136350496073./3983972220547.},
4294040e9f2SJed Brown                         {-6511271360970./6095937251113.,  5843115559534./2180450260947., -5289405421727./3760307252460. }};
4300298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
431a3a57f36SJed Brown   }
4328a381b04SJed Brown   PetscFunctionReturn(0);
4338a381b04SJed Brown }
4348a381b04SJed Brown 
4358a381b04SJed Brown #undef __FUNCT__
4368a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
4378a381b04SJed Brown /*@C
4388a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
4398a381b04SJed Brown 
4408a381b04SJed Brown    Not Collective
4418a381b04SJed Brown 
4428a381b04SJed Brown    Level: advanced
4438a381b04SJed Brown 
4448a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
445607a6623SBarry Smith .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll()
4468a381b04SJed Brown @*/
4478a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
4488a381b04SJed Brown {
4498a381b04SJed Brown   PetscErrorCode ierr;
4508a381b04SJed Brown   ARKTableauLink link;
4518a381b04SJed Brown 
4528a381b04SJed Brown   PetscFunctionBegin;
4538a381b04SJed Brown   while ((link = ARKTableauList)) {
4548a381b04SJed Brown     ARKTableau t = &link->tab;
4558a381b04SJed Brown     ARKTableauList = link->next;
4568a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
457108c343cSJed Brown     ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr);
458cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
4598a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
4608a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
4618a381b04SJed Brown   }
4628a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
4638a381b04SJed Brown   PetscFunctionReturn(0);
4648a381b04SJed Brown }
4658a381b04SJed Brown 
4668a381b04SJed Brown #undef __FUNCT__
4678a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
4688a381b04SJed Brown /*@C
4698a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
4708a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
4718a381b04SJed Brown   when using static libraries.
4728a381b04SJed Brown 
4738a381b04SJed Brown   Level: developer
4748a381b04SJed Brown 
4758a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
4768a381b04SJed Brown .seealso: PetscInitialize()
4778a381b04SJed Brown @*/
478607a6623SBarry Smith PetscErrorCode TSARKIMEXInitializePackage(void)
4798a381b04SJed Brown {
4808a381b04SJed Brown   PetscErrorCode ierr;
4818a381b04SJed Brown 
4828a381b04SJed Brown   PetscFunctionBegin;
4838a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
4848a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
4858a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
486e817cc15SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
4878a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
4888a381b04SJed Brown   PetscFunctionReturn(0);
4898a381b04SJed Brown }
4908a381b04SJed Brown 
4918a381b04SJed Brown #undef __FUNCT__
4928a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
4938a381b04SJed Brown /*@C
4948a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
4958a381b04SJed Brown   called from PetscFinalize().
4968a381b04SJed Brown 
4978a381b04SJed Brown   Level: developer
4988a381b04SJed Brown 
4998a381b04SJed Brown .keywords: Petsc, destroy, package
5008a381b04SJed Brown .seealso: PetscFinalize()
5018a381b04SJed Brown @*/
5028a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
5038a381b04SJed Brown {
5048a381b04SJed Brown   PetscErrorCode ierr;
5058a381b04SJed Brown 
5068a381b04SJed Brown   PetscFunctionBegin;
5078a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
5088a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
5098a381b04SJed Brown   PetscFunctionReturn(0);
5108a381b04SJed Brown }
5118a381b04SJed Brown 
5128a381b04SJed Brown #undef __FUNCT__
5138a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
514cd652676SJed Brown /*@C
515cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
516cd652676SJed Brown 
517cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
518cd652676SJed Brown 
519cd652676SJed Brown    Input Parameters:
520cd652676SJed Brown +  name - identifier for method
521cd652676SJed Brown .  order - approximation order of method
522cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
523cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
5240298fd71SBarry Smith .  bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At)
5250298fd71SBarry Smith .  ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At)
526cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
5270298fd71SBarry Smith .  b - Non-stiff step completion table (dimension s; NULL to use last row of At)
5280298fd71SBarry Smith .  c - Non-stiff abscissa (dimension s; NULL to use row sums of A)
5290298fd71SBarry Smith .  bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available)
5300298fd71SBarry Smith .  bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided)
531cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
532cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
5330298fd71SBarry Smith -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt)
534cd652676SJed Brown 
535cd652676SJed Brown    Notes:
536cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
537cd652676SJed Brown 
538cd652676SJed Brown    Level: advanced
539cd652676SJed Brown 
540cd652676SJed Brown .keywords: TS, register
541cd652676SJed Brown 
542cd652676SJed Brown .seealso: TSARKIMEX
543cd652676SJed Brown @*/
54419fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
5458a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
546cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
547108c343cSJed Brown                                  const PetscReal bembedt[],const PetscReal bembed[],
548cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
5498a381b04SJed Brown {
5508a381b04SJed Brown   PetscErrorCode ierr;
5518a381b04SJed Brown   ARKTableauLink link;
5528a381b04SJed Brown   ARKTableau     t;
5538a381b04SJed Brown   PetscInt       i,j;
5548a381b04SJed Brown 
5558a381b04SJed Brown   PetscFunctionBegin;
5568a381b04SJed Brown   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
557cd652676SJed Brown   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
5588a381b04SJed Brown   t        = &link->tab;
5598a381b04SJed Brown   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5608a381b04SJed Brown   t->order = order;
5618a381b04SJed Brown   t->s     = s;
5628a381b04SJed 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);
5638a381b04SJed Brown   ierr     = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
5648a381b04SJed Brown   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5658a381b04SJed Brown   if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); }
5668a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
5678a381b04SJed Brown   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
5685dceddf7SDebojyoti Ghosh   else for (i=0; i<s; i++) t->b[i] = t->bt[i];
5698a381b04SJed Brown   if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); }
5708a381b04SJed 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];
5718a381b04SJed Brown   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
5728a381b04SJed 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];
573e817cc15SEmil Constantinescu   t->stiffly_accurate = PETSC_TRUE;
574e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
575e817cc15SEmil Constantinescu   t->explicit_first_stage = PETSC_TRUE;
576e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
577e817cc15SEmil Constantinescu   /*def of FSAL can be made more precise*/
5784e9d4bf5SJed Brown   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
579108c343cSJed Brown   if (bembedt) {
580108c343cSJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembedt,s,PetscReal,&t->bembed);CHKERRQ(ierr);
581108c343cSJed Brown     ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr);
582108c343cSJed Brown     ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr);
583108c343cSJed Brown   }
584108c343cSJed Brown 
5854f385281SJed Brown   t->pinterp     = pinterp;
586cd652676SJed Brown   ierr           = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
587cd652676SJed Brown   ierr           = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
588cd652676SJed Brown   ierr           = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
5898a381b04SJed Brown   link->next     = ARKTableauList;
5908a381b04SJed Brown   ARKTableauList = link;
5918a381b04SJed Brown   PetscFunctionReturn(0);
5928a381b04SJed Brown }
5938a381b04SJed Brown 
5948a381b04SJed Brown #undef __FUNCT__
595108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX"
596108c343cSJed Brown /*
597108c343cSJed Brown  The step completion formula is
598108c343cSJed Brown 
599108c343cSJed Brown  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
600108c343cSJed Brown 
601108c343cSJed Brown  This function can be called before or after ts->vec_sol has been updated.
602108c343cSJed Brown  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
603108c343cSJed Brown  We can write
604108c343cSJed Brown 
605108c343cSJed Brown  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
606108c343cSJed Brown      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
607108c343cSJed Brown      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
608108c343cSJed Brown 
609108c343cSJed Brown  so we can evaluate the method with different order even after the step has been optimistically completed.
610108c343cSJed Brown */
611108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
612108c343cSJed Brown {
613108c343cSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
614108c343cSJed Brown   ARKTableau     tab  = ark->tableau;
615108c343cSJed Brown   PetscScalar    *w   = ark->work;
616108c343cSJed Brown   PetscReal      h;
617108c343cSJed Brown   PetscInt       s = tab->s,j;
618108c343cSJed Brown   PetscErrorCode ierr;
619108c343cSJed Brown 
620108c343cSJed Brown   PetscFunctionBegin;
621108c343cSJed Brown   switch (ark->status) {
622108c343cSJed Brown   case TS_STEP_INCOMPLETE:
623108c343cSJed Brown   case TS_STEP_PENDING:
624108c343cSJed Brown     h = ts->time_step; break;
625108c343cSJed Brown   case TS_STEP_COMPLETE:
626108c343cSJed Brown     h = ts->time_step_prev; break;
627ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
628108c343cSJed Brown   }
629108c343cSJed Brown   if (order == tab->order) {
630e817cc15SEmil Constantinescu     if (ark->status == TS_STEP_INCOMPLETE) {
631740132f1SEmil Constantinescu       if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */
632e817cc15SEmil Constantinescu         ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr);
633e817cc15SEmil Constantinescu       } else { /* Use the standard completion formula (bt,b) */
634108c343cSJed Brown         ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
635e817cc15SEmil Constantinescu         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
636108c343cSJed Brown         ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
637e817cc15SEmil Constantinescu         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
638108c343cSJed Brown           for (j=0; j<s; j++) w[j] = h*tab->b[j];
639108c343cSJed Brown           ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
640e817cc15SEmil Constantinescu         }
641e817cc15SEmil Constantinescu       }
642108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
643108c343cSJed Brown     if (done) *done = PETSC_TRUE;
644108c343cSJed Brown     PetscFunctionReturn(0);
645108c343cSJed Brown   } else if (order == tab->order-1) {
646108c343cSJed Brown     if (!tab->bembedt) goto unavailable;
647108c343cSJed Brown     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
648108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
649e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
650108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
651108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
652108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
653108c343cSJed Brown     } else {                    /* Rollback and re-complete using (bet-be,be-b) */
654108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
655e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
656108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr);
657108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
658108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
659108c343cSJed Brown     }
660108c343cSJed Brown     if (done) *done = PETSC_TRUE;
661108c343cSJed Brown     PetscFunctionReturn(0);
662108c343cSJed Brown   }
663108c343cSJed Brown unavailable:
664108c343cSJed Brown   if (done) *done = PETSC_FALSE;
665ce94432eSBarry Smith   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
666108c343cSJed Brown   PetscFunctionReturn(0);
667108c343cSJed Brown }
668108c343cSJed Brown 
669108c343cSJed Brown #undef __FUNCT__
6708a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
6718a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
6728a381b04SJed Brown {
6738a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
6748a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
6758a381b04SJed Brown   const PetscInt  s    = tab->s;
6768a381b04SJed Brown   const PetscReal *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
677406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
678e817cc15SEmil Constantinescu   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,W = ark->Work,Z = ark->Z;
67956dcabbaSDebojyoti Ghosh   PetscBool       init_guess_extrp = ark->init_guess_extrp;
680108c343cSJed Brown   TSAdapt         adapt;
6818a381b04SJed Brown   SNES            snes;
682108c343cSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
683cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
684108c343cSJed Brown   PetscReal       t;
685108c343cSJed Brown   PetscBool       accept;
6868a381b04SJed Brown   PetscErrorCode  ierr;
6878a381b04SJed Brown 
6888a381b04SJed Brown   PetscFunctionBegin;
689e817cc15SEmil Constantinescu   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage) {
690e817cc15SEmil Constantinescu     PetscReal valid_time;
691e817cc15SEmil Constantinescu     PetscBool isvalid;
692e817cc15SEmil Constantinescu     ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,
693e817cc15SEmil Constantinescu                                           explicit_stage_time_id,
694e817cc15SEmil Constantinescu                                           valid_time,
695e817cc15SEmil Constantinescu                                           isvalid);
696e817cc15SEmil Constantinescu     CHKERRQ(ierr);
697e817cc15SEmil Constantinescu     if (!isvalid || valid_time != ts->ptime) {
698e817cc15SEmil Constantinescu       TS        ts_start;
699e817cc15SEmil Constantinescu       SNES      snes_start;
700740132f1SEmil Constantinescu       DM        dm;
701740132f1SEmil Constantinescu       PetscReal atol;
702740132f1SEmil Constantinescu       Vec       vatol;
703740132f1SEmil Constantinescu       PetscReal rtol;
704740132f1SEmil Constantinescu       Vec       vrtol;
70519436ca2SJed Brown 
70619436ca2SJed Brown       ierr = TSCreate(PETSC_COMM_WORLD,&ts_start);CHKERRQ(ierr);
70719436ca2SJed Brown       ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr);
70819436ca2SJed Brown       ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr);
709e817cc15SEmil Constantinescu       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
710740132f1SEmil Constantinescu       ierr = TSSetDM(ts_start,dm);CHKERRQ(ierr);
711bbd56ea5SKarl Rupp 
712e817cc15SEmil Constantinescu       ts_start->adapt=ts->adapt;
713740132f1SEmil Constantinescu       PetscObjectReference((PetscObject)ts_start->adapt);
714bbd56ea5SKarl Rupp 
715e817cc15SEmil Constantinescu       ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
716e817cc15SEmil Constantinescu       ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr);
717eb082435SEmil Constantinescu       ierr = TSSetDuration(ts_start,1,ts->ptime+ts->time_step);CHKERRQ(ierr);
718740132f1SEmil Constantinescu       ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr);
719e817cc15SEmil Constantinescu       ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
720740132f1SEmil Constantinescu       ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
721e817cc15SEmil Constantinescu       ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
722e817cc15SEmil Constantinescu       ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr);
723740132f1SEmil Constantinescu       ierr = TSGetTolerances(ts,&atol,&vatol,&rtol,&vrtol);CHKERRQ(ierr);
724740132f1SEmil Constantinescu       ierr = TSSetTolerances(ts_start,atol,vatol,rtol,vrtol);CHKERRQ(ierr);
725e817cc15SEmil Constantinescu       ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
726e817cc15SEmil Constantinescu       ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr);
727bbd56ea5SKarl Rupp 
728740132f1SEmil Constantinescu       ts->time_step = ts_start->time_step;
729740132f1SEmil Constantinescu       ts->steps++;
730e817cc15SEmil Constantinescu       ierr = VecCopy(((TS_ARKIMEX*)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr);
731166a6834SEmil Constantinescu       ts_start->snes=NULL;
732740132f1SEmil Constantinescu       ierr = TSSetSNES(ts,snes_start);CHKERRQ(ierr);
733166a6834SEmil Constantinescu       ierr = SNESDestroy(&snes_start);CHKERRQ(ierr);
734166a6834SEmil Constantinescu       ierr = TSDestroy(&ts_start);CHKERRQ(ierr);
735e817cc15SEmil Constantinescu     }
736e817cc15SEmil Constantinescu   }
737e817cc15SEmil Constantinescu 
7388a381b04SJed Brown   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
739cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
7408a381b04SJed Brown   t              = ts->ptime;
741108c343cSJed Brown   accept         = PETSC_TRUE;
742108c343cSJed Brown   ark->status    = TS_STEP_INCOMPLETE;
7438a381b04SJed Brown 
744e817cc15SEmil Constantinescu 
74597335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
746108c343cSJed Brown     PetscReal h = ts->time_step;
747b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
7488a381b04SJed Brown     for (i=0; i<s; i++) {
7498a381b04SJed Brown       if (At[i*s+i] == 0) {           /* This stage is explicit */
7508a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
751e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7528a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
7538a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7548a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7558a381b04SJed Brown       } else {
7568a381b04SJed Brown         ark->stage_time = t + h*ct[i];
757b296d7d5SJed Brown         ark->scoeff     = 1./At[i*s+i];
758b8123daeSJed Brown         ierr            = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
7598a381b04SJed Brown         /* Affine part */
7608a381b04SJed Brown         ierr = VecZeroEntries(W);CHKERRQ(ierr);
7618a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7628a381b04SJed Brown         ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
763b296d7d5SJed Brown         ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr);
764f16577ceSEmil Constantinescu 
7658a381b04SJed Brown         /* Ydot = shift*(Y-Z) */
7668a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
767e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7684f385281SJed Brown         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
769f16577ceSEmil Constantinescu 
770*9eef816dSJed Brown         if (init_guess_extrp && ark->prev_step_valid) {
77156dcabbaSDebojyoti Ghosh           /* Initial guess extrapolated from previous time step stage values */
77256dcabbaSDebojyoti Ghosh           ierr        = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr);
77356dcabbaSDebojyoti Ghosh         } else {
7748a381b04SJed Brown           /* Initial guess taken from last stage */
7758a381b04SJed Brown           ierr        = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr);
77656dcabbaSDebojyoti Ghosh         }
7778a381b04SJed Brown         ierr          = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
778e817cc15SEmil Constantinescu         ierr          = (ts->ops->snesfunction)(snes,Y[i],W,ts);CHKERRQ(ierr);
7798a381b04SJed Brown         ierr          = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
7808a381b04SJed Brown         ierr          = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
7815ef26d82SJed Brown         ts->snes_its += its; ts->ksp_its += lits;
782552698daSJed Brown         ierr          = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
78397335746SJed Brown         ierr          = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
78497335746SJed Brown         if (!accept) goto reject_step;
7858a381b04SJed Brown       }
786e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) {
787e817cc15SEmil Constantinescu         if (i==0 && tab->explicit_first_stage) {
788e817cc15SEmil Constantinescu           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);
789e817cc15SEmil Constantinescu         } else {
790e817cc15SEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
791e817cc15SEmil Constantinescu         }
792e817cc15SEmil Constantinescu       } else {
7938a381b04SJed Brown         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
7944cc180ffSJed Brown         ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
795e817cc15SEmil Constantinescu         ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr);
7964cc180ffSJed Brown         if (ark->imex) {
7978a381b04SJed Brown           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
7984cc180ffSJed Brown         } else {
7994cc180ffSJed Brown           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
8004cc180ffSJed Brown         }
8018a381b04SJed Brown       }
802e817cc15SEmil Constantinescu     }
8030298fd71SBarry Smith     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
804108c343cSJed Brown     ark->status = TS_STEP_PENDING;
8058a381b04SJed Brown 
806108c343cSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
807552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
808108c343cSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
809108c343cSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
810108c343cSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
811108c343cSJed Brown     if (accept) {
812108c343cSJed Brown       /* ignore next_scheme for now */
8138a381b04SJed Brown       ts->ptime    += ts->time_step;
814cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
8158a381b04SJed Brown       ts->steps++;
816e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) { /* save the initial slope for the next step*/
817e817cc15SEmil Constantinescu         ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
818e817cc15SEmil Constantinescu       }
819108c343cSJed Brown       ark->status = TS_STEP_COMPLETE;
820e817cc15SEmil Constantinescu       if (tab->explicit_first_stage) {
821e817cc15SEmil Constantinescu         ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
822e817cc15SEmil Constantinescu       }
82364b5d2f7SDebojyoti Ghosh       /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
82464b5d2f7SDebojyoti Ghosh       if (ark->init_guess_extrp) {
82564b5d2f7SDebojyoti Ghosh         for (i = 0; i<s; i++) {
82664b5d2f7SDebojyoti Ghosh           ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr);
82764b5d2f7SDebojyoti Ghosh           ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr);
82864b5d2f7SDebojyoti Ghosh           ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr);
82964b5d2f7SDebojyoti Ghosh         }
830*9eef816dSJed Brown         ark->prev_step_valid = PETSC_TRUE;
83164b5d2f7SDebojyoti Ghosh       }
832108c343cSJed Brown       break;
833108c343cSJed Brown     } else {                    /* Roll back the current step */
8342c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*bt[j];
835108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotI);CHKERRQ(ierr);
8362c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
837108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotRHS);CHKERRQ(ierr);
838108c343cSJed Brown       ts->time_step = next_time_step;
839108c343cSJed Brown       ark->status   = TS_STEP_INCOMPLETE;
840108c343cSJed Brown     }
841476b6736SJed Brown reject_step: continue;
842108c343cSJed Brown   }
843b2ce242eSJed Brown   if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
8448a381b04SJed Brown   PetscFunctionReturn(0);
8458a381b04SJed Brown }
8468a381b04SJed Brown 
847cd652676SJed Brown #undef __FUNCT__
848cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
849cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
850cd652676SJed Brown {
851cd652676SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8524f385281SJed Brown   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
853108c343cSJed Brown   PetscReal       h;
854108c343cSJed Brown   PetscReal       tt,t;
855cd652676SJed Brown   PetscScalar     *bt,*b;
856cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
857cd652676SJed Brown   PetscErrorCode  ierr;
858cd652676SJed Brown 
859cd652676SJed Brown   PetscFunctionBegin;
860ce94432eSBarry Smith   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
861108c343cSJed Brown   switch (ark->status) {
862108c343cSJed Brown   case TS_STEP_INCOMPLETE:
863108c343cSJed Brown   case TS_STEP_PENDING:
864108c343cSJed Brown     h = ts->time_step;
865108c343cSJed Brown     t = (itime - ts->ptime)/h;
866108c343cSJed Brown     break;
867108c343cSJed Brown   case TS_STEP_COMPLETE:
868108c343cSJed Brown     h = ts->time_step_prev;
869108c343cSJed Brown     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
870108c343cSJed Brown     break;
871ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
872108c343cSJed Brown   }
873cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
874cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
8754f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
876cd652676SJed Brown     for (i=0; i<s; i++) {
877c1758d98SDebojyoti Ghosh       bt[i] += h * Bt[i*pinterp+j] * tt;
878108c343cSJed Brown       b[i]  += h * B[i*pinterp+j] * tt;
879cd652676SJed Brown     }
880cd652676SJed Brown   }
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 
88856dcabbaSDebojyoti Ghosh #undef __FUNCT__
88956dcabbaSDebojyoti Ghosh #define __FUNCT__ "TSExtrapolate_ARKIMEX"
89056dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X)
89156dcabbaSDebojyoti Ghosh {
89256dcabbaSDebojyoti Ghosh   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
89356dcabbaSDebojyoti Ghosh   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
89456dcabbaSDebojyoti Ghosh   PetscReal       h;
89556dcabbaSDebojyoti Ghosh   PetscReal       tt,t;
89656dcabbaSDebojyoti Ghosh   PetscScalar     *bt,*b;
89756dcabbaSDebojyoti Ghosh   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
89856dcabbaSDebojyoti Ghosh   PetscErrorCode  ierr;
89956dcabbaSDebojyoti Ghosh 
90056dcabbaSDebojyoti Ghosh   PetscFunctionBegin;
90156dcabbaSDebojyoti Ghosh   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
90256dcabbaSDebojyoti Ghosh   t = 1.0 + (ts->time_step/ts->time_step_prev)*c;
90381d12688SDebojyoti Ghosh   h = ts->time_step;
90456dcabbaSDebojyoti Ghosh   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
90556dcabbaSDebojyoti Ghosh   for (i=0; i<s; i++) bt[i] = b[i] = 0;
90656dcabbaSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
90756dcabbaSDebojyoti Ghosh     for (i=0; i<s; i++) {
90881d12688SDebojyoti Ghosh       bt[i] += h * Bt[i*pinterp+j] * tt;
90956dcabbaSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
91056dcabbaSDebojyoti Ghosh     }
91156dcabbaSDebojyoti Ghosh   }
912*9eef816dSJed Brown   if (!ark->prev_step_valid) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored");
91356dcabbaSDebojyoti Ghosh   ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr);
91456dcabbaSDebojyoti Ghosh   ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr);
91556dcabbaSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr);
91656dcabbaSDebojyoti Ghosh   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
91756dcabbaSDebojyoti Ghosh   PetscFunctionReturn(0);
91856dcabbaSDebojyoti Ghosh }
91956dcabbaSDebojyoti Ghosh 
9208a381b04SJed Brown /*------------------------------------------------------------*/
9218a381b04SJed Brown #undef __FUNCT__
9228a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
9238a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
9248a381b04SJed Brown {
9258a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
9268a381b04SJed Brown   PetscInt       s;
9278a381b04SJed Brown   PetscErrorCode ierr;
9288a381b04SJed Brown 
9298a381b04SJed Brown   PetscFunctionBegin;
9308a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
9318a381b04SJed Brown   s    = ark->tableau->s;
9328a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
9338a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
9348a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
93556dcabbaSDebojyoti Ghosh   if (&ark->init_guess_extrp) {
93656dcabbaSDebojyoti Ghosh     ierr = VecDestroyVecs(s,&ark->Y_prev);CHKERRQ(ierr);
93756dcabbaSDebojyoti Ghosh     ierr = VecDestroyVecs(s,&ark->YdotI_prev);CHKERRQ(ierr);
93856dcabbaSDebojyoti Ghosh     ierr = VecDestroyVecs(s,&ark->YdotRHS_prev);CHKERRQ(ierr);
93956dcabbaSDebojyoti Ghosh   }
9408a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
9418a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
942e817cc15SEmil Constantinescu   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
9438a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
9448a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
9458a381b04SJed Brown   PetscFunctionReturn(0);
9468a381b04SJed Brown }
9478a381b04SJed Brown 
9488a381b04SJed Brown #undef __FUNCT__
9498a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
9508a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
9518a381b04SJed Brown {
9528a381b04SJed Brown   PetscErrorCode ierr;
9538a381b04SJed Brown 
9548a381b04SJed Brown   PetscFunctionBegin;
9558a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
9568a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
957bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr);
958bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr);
959bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr);
9608a381b04SJed Brown   PetscFunctionReturn(0);
9618a381b04SJed Brown }
9628a381b04SJed Brown 
963d5e6173cSPeter Brune 
964d5e6173cSPeter Brune #undef __FUNCT__
965d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXGetVecs"
966d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
967d5e6173cSPeter Brune {
968d5e6173cSPeter Brune   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
969d5e6173cSPeter Brune   PetscErrorCode ierr;
970d5e6173cSPeter Brune 
971d5e6173cSPeter Brune   PetscFunctionBegin;
972d5e6173cSPeter Brune   if (Z) {
973d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
974d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
975d5e6173cSPeter Brune     } else *Z = ax->Z;
976d5e6173cSPeter Brune   }
977d5e6173cSPeter Brune   if (Ydot) {
978d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
979d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
980d5e6173cSPeter Brune     } else *Ydot = ax->Ydot;
981d5e6173cSPeter Brune   }
982d5e6173cSPeter Brune   PetscFunctionReturn(0);
983d5e6173cSPeter Brune }
984d5e6173cSPeter Brune 
985d5e6173cSPeter Brune 
986d5e6173cSPeter Brune #undef __FUNCT__
987d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXRestoreVecs"
988d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
989d5e6173cSPeter Brune {
990d5e6173cSPeter Brune   PetscErrorCode ierr;
991d5e6173cSPeter Brune 
992d5e6173cSPeter Brune   PetscFunctionBegin;
993d5e6173cSPeter Brune   if (Z) {
994d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
995d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
996d5e6173cSPeter Brune     }
997d5e6173cSPeter Brune   }
998d5e6173cSPeter Brune   if (Ydot) {
999d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
1000d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
1001d5e6173cSPeter Brune     }
1002d5e6173cSPeter Brune   }
1003d5e6173cSPeter Brune   PetscFunctionReturn(0);
1004d5e6173cSPeter Brune }
1005d5e6173cSPeter Brune 
10068a381b04SJed Brown /*
10078a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
10088a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
10098a381b04SJed Brown */
10108a381b04SJed Brown #undef __FUNCT__
10118a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
10128a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
10138a381b04SJed Brown {
10148a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1015d5e6173cSPeter Brune   DM             dm,dmsave;
1016d5e6173cSPeter Brune   Vec            Z,Ydot;
1017b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10188a381b04SJed Brown   PetscErrorCode ierr;
10198a381b04SJed Brown 
10208a381b04SJed Brown   PetscFunctionBegin;
1021d5e6173cSPeter Brune   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1022d5e6173cSPeter Brune   ierr   = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
1023b296d7d5SJed Brown   ierr   = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
1024d5e6173cSPeter Brune   dmsave = ts->dm;
1025d5e6173cSPeter Brune   ts->dm = dm;
1026740132f1SEmil Constantinescu 
1027d5e6173cSPeter Brune   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
1028e817cc15SEmil Constantinescu 
1029d5e6173cSPeter Brune   ts->dm = dmsave;
1030d5e6173cSPeter Brune   ierr   = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
10318a381b04SJed Brown   PetscFunctionReturn(0);
10328a381b04SJed Brown }
10338a381b04SJed Brown 
10348a381b04SJed Brown #undef __FUNCT__
10358a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
10368a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
10378a381b04SJed Brown {
10388a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1039d5e6173cSPeter Brune   DM             dm,dmsave;
1040d5e6173cSPeter Brune   Vec            Ydot;
1041b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10428a381b04SJed Brown   PetscErrorCode ierr;
10438a381b04SJed Brown 
10448a381b04SJed Brown   PetscFunctionBegin;
1045d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
10460298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
10478a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1048d5e6173cSPeter Brune   dmsave = ts->dm;
1049d5e6173cSPeter Brune   ts->dm = dm;
1050740132f1SEmil Constantinescu 
1051b296d7d5SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,str,ark->imex);CHKERRQ(ierr);
1052740132f1SEmil Constantinescu 
1053d5e6173cSPeter Brune   ts->dm = dmsave;
10540298fd71SBarry Smith   ierr   = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
1055d5e6173cSPeter Brune   PetscFunctionReturn(0);
1056d5e6173cSPeter Brune }
1057d5e6173cSPeter Brune 
1058d5e6173cSPeter Brune #undef __FUNCT__
1059d5e6173cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSARKIMEX"
1060d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1061d5e6173cSPeter Brune {
1062d5e6173cSPeter Brune   PetscFunctionBegin;
1063d5e6173cSPeter Brune   PetscFunctionReturn(0);
1064d5e6173cSPeter Brune }
1065d5e6173cSPeter Brune 
1066d5e6173cSPeter Brune #undef __FUNCT__
1067d5e6173cSPeter Brune #define __FUNCT__ "DMRestrictHook_TSARKIMEX"
1068d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1069d5e6173cSPeter Brune {
1070d5e6173cSPeter Brune   TS             ts = (TS)ctx;
1071d5e6173cSPeter Brune   PetscErrorCode ierr;
1072d5e6173cSPeter Brune   Vec            Z,Z_c;
1073d5e6173cSPeter Brune 
1074d5e6173cSPeter Brune   PetscFunctionBegin;
10750298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10760298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
1077d5e6173cSPeter Brune   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1078d5e6173cSPeter Brune   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
10790298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10800298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
10818a381b04SJed Brown   PetscFunctionReturn(0);
10828a381b04SJed Brown }
10838a381b04SJed Brown 
1084cdb298fcSPeter Brune 
1085cdb298fcSPeter Brune #undef __FUNCT__
1086cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainHook_TSARKIMEX"
1087cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1088cdb298fcSPeter Brune {
1089cdb298fcSPeter Brune   PetscFunctionBegin;
1090cdb298fcSPeter Brune   PetscFunctionReturn(0);
1091cdb298fcSPeter Brune }
1092cdb298fcSPeter Brune 
1093cdb298fcSPeter Brune #undef __FUNCT__
1094cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX"
1095cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1096cdb298fcSPeter Brune {
1097cdb298fcSPeter Brune   TS             ts = (TS)ctx;
1098cdb298fcSPeter Brune   PetscErrorCode ierr;
1099cdb298fcSPeter Brune   Vec            Z,Z_c;
1100cdb298fcSPeter Brune 
1101cdb298fcSPeter Brune   PetscFunctionBegin;
11020298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
11030298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1104cdb298fcSPeter Brune 
1105cdb298fcSPeter Brune   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1106cdb298fcSPeter Brune   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1107cdb298fcSPeter Brune 
11080298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
11090298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1110cdb298fcSPeter Brune   PetscFunctionReturn(0);
1111cdb298fcSPeter Brune }
1112cdb298fcSPeter Brune 
11138a381b04SJed Brown #undef __FUNCT__
11148a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
11158a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
11168a381b04SJed Brown {
11178a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1118f2c2a1b9SBarry Smith   ARKTableau     tab;
1119f2c2a1b9SBarry Smith   PetscInt       s;
11208a381b04SJed Brown   PetscErrorCode ierr;
1121d5e6173cSPeter Brune   DM             dm;
1122f9c1d6abSBarry Smith 
11238a381b04SJed Brown   PetscFunctionBegin;
11248a381b04SJed Brown   if (!ark->tableau) {
1125e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
11268a381b04SJed Brown   }
1127f2c2a1b9SBarry Smith   tab  = ark->tableau;
1128f2c2a1b9SBarry Smith   s    = tab->s;
11298a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
11308a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
11318a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
113256dcabbaSDebojyoti Ghosh   if (ark->init_guess_extrp) {
113356dcabbaSDebojyoti Ghosh     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y_prev);CHKERRQ(ierr);
113456dcabbaSDebojyoti Ghosh     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI_prev);CHKERRQ(ierr);
113556dcabbaSDebojyoti Ghosh     ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS_prev);CHKERRQ(ierr);
113656dcabbaSDebojyoti Ghosh   }
11378a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
11388a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
1139e817cc15SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
11408a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
11418a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
1142d5e6173cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1143d5e6173cSPeter Brune   if (dm) {
1144d5e6173cSPeter Brune     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1145cdb298fcSPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1146d5e6173cSPeter Brune   }
11478a381b04SJed Brown   PetscFunctionReturn(0);
11488a381b04SJed Brown }
11498a381b04SJed Brown /*------------------------------------------------------------*/
11508a381b04SJed Brown 
11518a381b04SJed Brown #undef __FUNCT__
11528a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
11538a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
11548a381b04SJed Brown {
11554cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11568a381b04SJed Brown   PetscErrorCode ierr;
11578a381b04SJed Brown   char           arktype[256];
11588a381b04SJed Brown 
11598a381b04SJed Brown   PetscFunctionBegin;
11608a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
11618a381b04SJed Brown   {
11628a381b04SJed Brown     ARKTableauLink link;
11638a381b04SJed Brown     PetscInt       count,choice;
11648a381b04SJed Brown     PetscBool      flg;
11658a381b04SJed Brown     const char     **namelist;
11668caf3d72SBarry Smith     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr);
11678a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
11688a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
11698a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11708a381b04SJed Brown     ierr      = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
11718a381b04SJed Brown     ierr      = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
11728a381b04SJed Brown     ierr      = PetscFree(namelist);CHKERRQ(ierr);
11734cc180ffSJed Brown     flg       = (PetscBool) !ark->imex;
11740298fd71SBarry Smith     ierr      = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr);
11754cc180ffSJed Brown     ark->imex = (PetscBool) !flg;
117656dcabbaSDebojyoti Ghosh     ark->init_guess_extrp = PETSC_FALSE;
117756dcabbaSDebojyoti Ghosh     ierr      = PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->init_guess_extrp,&ark->init_guess_extrp,NULL);CHKERRQ(ierr);
1178d52bd9f3SBarry Smith     ierr      = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
11798a381b04SJed Brown   }
11808a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
11818a381b04SJed Brown   PetscFunctionReturn(0);
11828a381b04SJed Brown }
11838a381b04SJed Brown 
11848a381b04SJed Brown #undef __FUNCT__
11858a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
11868a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
11878a381b04SJed Brown {
1188257d2499SJed Brown   PetscErrorCode ierr;
1189f1d86077SJed Brown   PetscInt       i;
1190f1d86077SJed Brown   size_t         left,count;
11918a381b04SJed Brown   char           *p;
11928a381b04SJed Brown 
11938a381b04SJed Brown   PetscFunctionBegin;
1194f1d86077SJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
1195f1d86077SJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
11968a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
11978a381b04SJed Brown     left -= count;
11988a381b04SJed Brown     p    += count;
11998a381b04SJed Brown     *p++  = ' ';
12008a381b04SJed Brown   }
12018a381b04SJed Brown   p[i ? 0 : -1] = 0;
12028a381b04SJed Brown   PetscFunctionReturn(0);
12038a381b04SJed Brown }
12048a381b04SJed Brown 
12058a381b04SJed Brown #undef __FUNCT__
12068a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
12078a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
12088a381b04SJed Brown {
12098a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
12108a381b04SJed Brown   ARKTableau     tab  = ark->tableau;
12118a381b04SJed Brown   PetscBool      iascii;
12128a381b04SJed Brown   PetscErrorCode ierr;
1213559eea31SJed Brown   TSAdapt        adapt;
12148a381b04SJed Brown 
12158a381b04SJed Brown   PetscFunctionBegin;
1216251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12178a381b04SJed Brown   if (iascii) {
121819fd82e9SBarry Smith     TSARKIMEXType arktype;
12198a381b04SJed Brown     char          buf[512];
12208a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
12218a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
12228caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
122331f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
12248caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1225e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
1226e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr);
1227e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr);
122831f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
12298a381b04SJed Brown   }
1230552698daSJed Brown   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
1231559eea31SJed Brown   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
1232d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
12338a381b04SJed Brown   PetscFunctionReturn(0);
12348a381b04SJed Brown }
12358a381b04SJed Brown 
12368a381b04SJed Brown #undef __FUNCT__
1237f2c2a1b9SBarry Smith #define __FUNCT__ "TSLoad_ARKIMEX"
1238f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1239f2c2a1b9SBarry Smith {
1240f2c2a1b9SBarry Smith   PetscErrorCode ierr;
1241f2c2a1b9SBarry Smith   SNES           snes;
1242ad6bc421SBarry Smith   TSAdapt        tsadapt;
1243f2c2a1b9SBarry Smith 
1244f2c2a1b9SBarry Smith   PetscFunctionBegin;
1245552698daSJed Brown   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
1246ad6bc421SBarry Smith   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
1247f2c2a1b9SBarry Smith   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1248f2c2a1b9SBarry Smith   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1249ad6bc421SBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
12500298fd71SBarry Smith   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
12510298fd71SBarry Smith   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
1252f2c2a1b9SBarry Smith   PetscFunctionReturn(0);
1253f2c2a1b9SBarry Smith }
1254f2c2a1b9SBarry Smith 
1255f2c2a1b9SBarry Smith #undef __FUNCT__
12568a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
12578a381b04SJed Brown /*@C
12588a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
12598a381b04SJed Brown 
12608a381b04SJed Brown   Logically collective
12618a381b04SJed Brown 
12628a381b04SJed Brown   Input Parameter:
12638a381b04SJed Brown +  ts - timestepping context
12648a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
12658a381b04SJed Brown 
12668a381b04SJed Brown   Level: intermediate
12678a381b04SJed Brown 
1268020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
12698a381b04SJed Brown @*/
127019fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
12718a381b04SJed Brown {
12728a381b04SJed Brown   PetscErrorCode ierr;
12738a381b04SJed Brown 
12748a381b04SJed Brown   PetscFunctionBegin;
12758a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
127619fd82e9SBarry Smith   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
12778a381b04SJed Brown   PetscFunctionReturn(0);
12788a381b04SJed Brown }
12798a381b04SJed Brown 
12808a381b04SJed Brown #undef __FUNCT__
12818a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
12828a381b04SJed Brown /*@C
12838a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
12848a381b04SJed Brown 
12858a381b04SJed Brown   Logically collective
12868a381b04SJed Brown 
12878a381b04SJed Brown   Input Parameter:
12888a381b04SJed Brown .  ts - timestepping context
12898a381b04SJed Brown 
12908a381b04SJed Brown   Output Parameter:
12918a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
12928a381b04SJed Brown 
12938a381b04SJed Brown   Level: intermediate
12948a381b04SJed Brown 
12958a381b04SJed Brown .seealso: TSARKIMEXGetType()
12968a381b04SJed Brown @*/
129719fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
12988a381b04SJed Brown {
12998a381b04SJed Brown   PetscErrorCode ierr;
13008a381b04SJed Brown 
13018a381b04SJed Brown   PetscFunctionBegin;
13028a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
130319fd82e9SBarry Smith   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
13048a381b04SJed Brown   PetscFunctionReturn(0);
13058a381b04SJed Brown }
13068a381b04SJed Brown 
13074cc180ffSJed Brown #undef __FUNCT__
13084cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
13094cc180ffSJed Brown /*@C
13104cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
13114cc180ffSJed Brown 
13124cc180ffSJed Brown   Logically collective
13134cc180ffSJed Brown 
13144cc180ffSJed Brown   Input Parameter:
13154cc180ffSJed Brown +  ts - timestepping context
13164cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
13174cc180ffSJed Brown 
13184cc180ffSJed Brown   Level: intermediate
13194cc180ffSJed Brown 
13204cc180ffSJed Brown .seealso: TSARKIMEXGetType()
13214cc180ffSJed Brown @*/
13224cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
13234cc180ffSJed Brown {
13244cc180ffSJed Brown   PetscErrorCode ierr;
13254cc180ffSJed Brown 
13264cc180ffSJed Brown   PetscFunctionBegin;
13274cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13284cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
13294cc180ffSJed Brown   PetscFunctionReturn(0);
13304cc180ffSJed Brown }
13314cc180ffSJed Brown 
13328a381b04SJed Brown #undef __FUNCT__
13338a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
133419fd82e9SBarry Smith PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
13358a381b04SJed Brown {
13368a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
13378a381b04SJed Brown   PetscErrorCode ierr;
13388a381b04SJed Brown 
13398a381b04SJed Brown   PetscFunctionBegin;
1340f2c2a1b9SBarry Smith   if (!ark->tableau) {
1341f2c2a1b9SBarry Smith     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1342f2c2a1b9SBarry Smith   }
13438a381b04SJed Brown   *arktype = ark->tableau->name;
13448a381b04SJed Brown   PetscFunctionReturn(0);
13458a381b04SJed Brown }
13468a381b04SJed Brown #undef __FUNCT__
13478a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
134819fd82e9SBarry Smith PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
13498a381b04SJed Brown {
13508a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
13518a381b04SJed Brown   PetscErrorCode ierr;
13528a381b04SJed Brown   PetscBool      match;
13538a381b04SJed Brown   ARKTableauLink link;
13548a381b04SJed Brown 
13558a381b04SJed Brown   PetscFunctionBegin;
13568a381b04SJed Brown   if (ark->tableau) {
13578a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
13588a381b04SJed Brown     if (match) PetscFunctionReturn(0);
13598a381b04SJed Brown   }
13608a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
13618a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
13628a381b04SJed Brown     if (match) {
13638a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
13648a381b04SJed Brown       ark->tableau = &link->tab;
13658a381b04SJed Brown       PetscFunctionReturn(0);
13668a381b04SJed Brown     }
13678a381b04SJed Brown   }
1368ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
13698a381b04SJed Brown   PetscFunctionReturn(0);
13708a381b04SJed Brown }
13714cc180ffSJed Brown #undef __FUNCT__
13724cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
13734cc180ffSJed Brown PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
13744cc180ffSJed Brown {
13754cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13764cc180ffSJed Brown 
13774cc180ffSJed Brown   PetscFunctionBegin;
13784cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
13794cc180ffSJed Brown   PetscFunctionReturn(0);
13804cc180ffSJed Brown }
13818a381b04SJed Brown 
13828a381b04SJed Brown /* ------------------------------------------------------------ */
13838a381b04SJed Brown /*MC
1384a4386c9eSJed Brown       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
13858a381b04SJed Brown 
1386fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1387fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1388fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1389fca742c7SJed Brown 
1390fca742c7SJed Brown   Notes:
1391a4386c9eSJed Brown   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1392c8058688SBarry Smith 
1393a4386c9eSJed 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).
1394fca742c7SJed Brown 
13958a381b04SJed Brown   Level: beginner
13968a381b04SJed Brown 
1397c8058688SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
1398a4386c9eSJed Brown            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
13998a381b04SJed Brown 
14008a381b04SJed Brown M*/
14018a381b04SJed Brown #undef __FUNCT__
14028a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
14038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
14048a381b04SJed Brown {
14058a381b04SJed Brown   TS_ARKIMEX     *th;
14068a381b04SJed Brown   PetscErrorCode ierr;
14078a381b04SJed Brown 
14088a381b04SJed Brown   PetscFunctionBegin;
14098a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1410607a6623SBarry Smith   ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr);
14118a381b04SJed Brown #endif
14128a381b04SJed Brown 
14138a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
14148a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
14158a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
1416f2c2a1b9SBarry Smith   ts->ops->load           = TSLoad_ARKIMEX;
14178a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
14188a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
1419cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1420108c343cSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
14218a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
14228a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
14238a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
14248a381b04SJed Brown 
14258a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
14268a381b04SJed Brown   ts->data = (void*)th;
14274cc180ffSJed Brown   th->imex = PETSC_TRUE;
14288a381b04SJed Brown 
1429bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
1430bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
1431bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
14328a381b04SJed Brown   PetscFunctionReturn(0);
14338a381b04SJed Brown }
1434