xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 166a683482ae06e3a60148a2fc00bd06bb4e50b0)
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;
198a381b04SJed Brown 
208a381b04SJed Brown typedef struct _ARKTableau *ARKTableau;
218a381b04SJed Brown struct _ARKTableau {
228a381b04SJed Brown   char      *name;
234f385281SJed Brown   PetscInt  order;                /* Classical approximation order of the method */
244f385281SJed Brown   PetscInt  s;                    /* Number of stages */
25e817cc15SEmil Constantinescu   PetscBool stiffly_accurate;     /* The implicit part is stiffly accurate*/
26e817cc15SEmil Constantinescu   PetscBool FSAL_implicit;        /* The implicit part is FSAL*/
27e817cc15SEmil Constantinescu   PetscBool explicit_first_stage; /* The implicit part has an explicit first stage*/
284f385281SJed Brown   PetscInt  pinterp;              /* Interpolation order */
294f385281SJed Brown   PetscReal *At,*bt,*ct;          /* Stiff tableau */
308a381b04SJed Brown   PetscReal *A,*b,*c;             /* Non-stiff tableau */
31108c343cSJed Brown   PetscReal *bembedt,*bembed;     /* Embedded formula of order one less (order-1) */
32cd652676SJed Brown   PetscReal *binterpt,*binterp;   /* Dense output formula */
33108c343cSJed Brown   PetscReal ccfl;                 /* Placeholder for CFL coefficient relative to forward Euler */
348a381b04SJed Brown };
358a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
368a381b04SJed Brown struct _ARKTableauLink {
378a381b04SJed Brown   struct _ARKTableau tab;
388a381b04SJed Brown   ARKTableauLink     next;
398a381b04SJed Brown };
408a381b04SJed Brown static ARKTableauLink ARKTableauList;
418a381b04SJed Brown 
428a381b04SJed Brown typedef struct {
438a381b04SJed Brown   ARKTableau   tableau;
448a381b04SJed Brown   Vec          *Y;               /* States computed during the step */
458a381b04SJed Brown   Vec          *YdotI;           /* Time derivatives for the stiff part */
468a381b04SJed Brown   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
47e817cc15SEmil Constantinescu   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
488a381b04SJed Brown   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
498a381b04SJed Brown   Vec          Work;             /* Generic work vector */
508a381b04SJed Brown   Vec          Z;                /* Ydot = shift(Y-Z) */
518a381b04SJed Brown   PetscScalar  *work;            /* Scalar work */
52b296d7d5SJed Brown   PetscReal    scoeff;           /* shift = scoeff/dt */
538a381b04SJed Brown   PetscReal    stage_time;
544cc180ffSJed Brown   PetscBool    imex;
55108c343cSJed Brown   TSStepStatus status;
568a381b04SJed Brown } TS_ARKIMEX;
571f80e275SEmil Constantinescu /*MC
581f80e275SEmil Constantinescu      TSARKIMEXARS122 - Second order ARK IMEX scheme.
598a381b04SJed Brown 
601f80e275SEmil Constantinescu      This method has one explicit stage and one implicit stage.
611f80e275SEmil Constantinescu 
621f80e275SEmil Constantinescu      References:
631f80e275SEmil 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.
641f80e275SEmil Constantinescu 
651f80e275SEmil Constantinescu      Level: advanced
661f80e275SEmil Constantinescu 
671f80e275SEmil Constantinescu .seealso: TSARKIMEX
681f80e275SEmil Constantinescu M*/
691f80e275SEmil Constantinescu /*MC
701f80e275SEmil Constantinescu      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.
711f80e275SEmil Constantinescu 
721f80e275SEmil 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.
731f80e275SEmil Constantinescu 
741f80e275SEmil Constantinescu      Level: advanced
751f80e275SEmil Constantinescu 
761f80e275SEmil Constantinescu .seealso: TSARKIMEX
771f80e275SEmil Constantinescu M*/
781f80e275SEmil Constantinescu /*MC
791f80e275SEmil Constantinescu      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.
801f80e275SEmil Constantinescu 
811f80e275SEmil Constantinescu      This method has two implicit stages, and L-stable implicit scheme.
821f80e275SEmil Constantinescu 
831f80e275SEmil Constantinescu     References:
841f80e275SEmil 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
851f80e275SEmil Constantinescu 
861f80e275SEmil Constantinescu      Level: advanced
871f80e275SEmil Constantinescu 
881f80e275SEmil Constantinescu .seealso: TSARKIMEX
891f80e275SEmil Constantinescu M*/
901f80e275SEmil Constantinescu /*MC
91e817cc15SEmil Constantinescu      TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.
92e817cc15SEmil Constantinescu 
93e817cc15SEmil Constantinescu      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.
94e817cc15SEmil Constantinescu 
95e817cc15SEmil Constantinescu      Level: advanced
96e817cc15SEmil Constantinescu 
97e817cc15SEmil Constantinescu .seealso: TSARKIMEX
98e817cc15SEmil Constantinescu M*/
99e817cc15SEmil Constantinescu /*MC
1001f80e275SEmil Constantinescu      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.
1011f80e275SEmil Constantinescu 
1021f80e275SEmil 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.
1031f80e275SEmil Constantinescu 
1041f80e275SEmil Constantinescu      Level: advanced
1051f80e275SEmil Constantinescu 
1061f80e275SEmil Constantinescu .seealso: TSARKIMEX
1071f80e275SEmil Constantinescu M*/
10864f491ddSJed Brown /*MC
10964f491ddSJed Brown      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
11064f491ddSJed Brown 
111617a39beSEmil 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.
11264f491ddSJed Brown 
113b330ce4dSSatish Balay      Level: advanced
114b330ce4dSSatish Balay 
11564f491ddSJed Brown .seealso: TSARKIMEX
11664f491ddSJed Brown M*/
11764f491ddSJed Brown /*MC
11864f491ddSJed Brown      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
11964f491ddSJed Brown 
12064f491ddSJed Brown      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
12164f491ddSJed Brown 
122b330ce4dSSatish Balay      Level: advanced
123b330ce4dSSatish Balay 
12464f491ddSJed Brown .seealso: TSARKIMEX
12564f491ddSJed Brown M*/
12664f491ddSJed Brown /*MC
1276cf0794eSJed Brown      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.
1286cf0794eSJed Brown 
1296cf0794eSJed Brown      This method has three implicit stages.
1306cf0794eSJed Brown 
1316cf0794eSJed Brown      References:
1326cf0794eSJed 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
1336cf0794eSJed Brown 
1346cf0794eSJed Brown      This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375
1356cf0794eSJed Brown 
1366cf0794eSJed Brown      Level: advanced
1376cf0794eSJed Brown 
1386cf0794eSJed Brown .seealso: TSARKIMEX
1396cf0794eSJed Brown M*/
1406cf0794eSJed Brown /*MC
14164f491ddSJed Brown      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
14264f491ddSJed Brown 
14364f491ddSJed Brown      This method has one explicit stage and three implicit stages.
14464f491ddSJed Brown 
14564f491ddSJed Brown      References:
14664f491ddSJed Brown      Kennedy and Carpenter 2003.
14764f491ddSJed Brown 
148b330ce4dSSatish Balay      Level: advanced
149b330ce4dSSatish Balay 
15064f491ddSJed Brown .seealso: TSARKIMEX
15164f491ddSJed Brown M*/
15264f491ddSJed Brown /*MC
1536cf0794eSJed Brown      TSARKIMEXARS443 - Third order ARK IMEX scheme.
1546cf0794eSJed Brown 
1556cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1566cf0794eSJed Brown 
1576cf0794eSJed Brown      References:
1586cf0794eSJed 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.
1596cf0794eSJed Brown 
1606cf0794eSJed Brown      This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375
1616cf0794eSJed Brown 
1626cf0794eSJed Brown      Level: advanced
1636cf0794eSJed Brown 
1646cf0794eSJed Brown .seealso: TSARKIMEX
1656cf0794eSJed Brown M*/
1666cf0794eSJed Brown /*MC
1676cf0794eSJed Brown      TSARKIMEXBPR3 - Third order ARK IMEX scheme.
1686cf0794eSJed Brown 
1696cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1706cf0794eSJed Brown 
1716cf0794eSJed Brown      References:
1726cf0794eSJed Brown      This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375
1736cf0794eSJed Brown 
1746cf0794eSJed Brown      Level: advanced
1756cf0794eSJed Brown 
1766cf0794eSJed Brown .seealso: TSARKIMEX
1776cf0794eSJed Brown M*/
1786cf0794eSJed Brown /*MC
17964f491ddSJed Brown      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
18064f491ddSJed Brown 
18164f491ddSJed Brown      This method has one explicit stage and four implicit stages.
18264f491ddSJed Brown 
18364f491ddSJed Brown      References:
18464f491ddSJed Brown      Kennedy and Carpenter 2003.
18564f491ddSJed Brown 
186b330ce4dSSatish Balay      Level: advanced
187b330ce4dSSatish Balay 
18864f491ddSJed Brown .seealso: TSARKIMEX
18964f491ddSJed Brown M*/
19064f491ddSJed Brown /*MC
19164f491ddSJed Brown      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
19264f491ddSJed Brown 
19364f491ddSJed Brown      This method has one explicit stage and five implicit stages.
19464f491ddSJed Brown 
19564f491ddSJed Brown      References:
19664f491ddSJed Brown      Kennedy and Carpenter 2003.
19764f491ddSJed Brown 
198b330ce4dSSatish Balay      Level: advanced
199b330ce4dSSatish Balay 
20064f491ddSJed Brown .seealso: TSARKIMEX
20164f491ddSJed Brown M*/
20264f491ddSJed Brown 
2038a381b04SJed Brown #undef __FUNCT__
2048a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll"
2058a381b04SJed Brown /*@C
2068a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
2078a381b04SJed Brown 
208fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
2098a381b04SJed Brown 
2108a381b04SJed Brown   Level: advanced
2118a381b04SJed Brown 
2128a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
2138a381b04SJed Brown 
2148a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
2158a381b04SJed Brown @*/
2168a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
2178a381b04SJed Brown {
2188a381b04SJed Brown   PetscErrorCode ierr;
2198a381b04SJed Brown 
2208a381b04SJed Brown   PetscFunctionBegin;
2218a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
2228a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
223e817cc15SEmil Constantinescu 
224e817cc15SEmil Constantinescu   {
225e817cc15SEmil Constantinescu     const PetscReal
226e817cc15SEmil Constantinescu       A[3][3] = {{0.0,0.0,0.0},
227e817cc15SEmil Constantinescu                  {0.0,0.0,0.0},
228748ad121SEmil Constantinescu                  {0.0,0.5,0.0}},
229e817cc15SEmil Constantinescu       At[3][3] = {{1.0,0.0,0.0},
230e817cc15SEmil Constantinescu                   {0.0,0.5,0.0},
231e817cc15SEmil Constantinescu                   {0.0,0.5,0.5}},
232e817cc15SEmil Constantinescu       b[3]       = {0.0,0.5,0.5},
233e817cc15SEmil Constantinescu       bembedt[3] = {1.0,0.0,0.0};
2340298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
235e817cc15SEmil Constantinescu   }
2368a381b04SJed Brown   {
2378a381b04SJed Brown     const PetscReal
2381f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2391f80e275SEmil Constantinescu                  {0.5,0.0}},
2401f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2411f80e275SEmil Constantinescu                   {0.0,0.5}},
2421f80e275SEmil Constantinescu       b[2]       = {0.0,1.0},
2431f80e275SEmil Constantinescu       bembedt[2] = {0.5,0.5};
2441f80e275SEmil 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*/
2450298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2461f80e275SEmil Constantinescu   }
2471f80e275SEmil Constantinescu   {
2481f80e275SEmil Constantinescu     const PetscReal
2491f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2501f80e275SEmil Constantinescu                  {1.0,0.0}},
2511f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2521f80e275SEmil Constantinescu                   {0.5,0.5}},
2531f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
2541f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0};
2551f80e275SEmil 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*/
2560298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2571f80e275SEmil Constantinescu   }
2581f80e275SEmil Constantinescu   {
259da80777bSKarl 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   */
2601f80e275SEmil Constantinescu     const PetscReal
2611f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2621f80e275SEmil Constantinescu                  {1.0,0.0}},
263da80777bSKarl Rupp       At[2][2] = {{0.2928932188134524755992,0.0},
264da80777bSKarl Rupp                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
2651f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
2661f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0},
267da80777bSKarl Rupp       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
268da80777bSKarl Rupp                         {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}},
2691f80e275SEmil Constantinescu       binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
2700298fd71SBarry 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);
2711f80e275SEmil Constantinescu   }
2721f80e275SEmil Constantinescu   {
273da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
274da80777bSKarl Rupp     const PetscReal
2758a381b04SJed Brown       A[3][3] = {{0,0,0},
276da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
277617a39beSEmil Constantinescu                  {0.5,0.5,0}},
278da80777bSKarl Rupp       At[3][3] = {{0,0,0},
279da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
280da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
281da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
282da80777bSKarl Rupp       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
283da80777bSKarl Rupp                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
284da80777bSKarl Rupp                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
2850298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
2861f80e275SEmil Constantinescu   }
2871f80e275SEmil Constantinescu   {
288da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
289da80777bSKarl Rupp     const PetscReal
2901f80e275SEmil Constantinescu       A[3][3] = {{0,0,0},
291da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
2928a381b04SJed Brown                  {0.75,0.25,0}},
293da80777bSKarl Rupp       At[3][3] = {{0,0,0},
294da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
295da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
296da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
297da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
298da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
299da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3000298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
3018a381b04SJed Brown   }
30206db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
303da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
304da80777bSKarl Rupp     const PetscReal
305da80777bSKarl Rupp       A[3][3] = {{0,0,0},
306da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
307da80777bSKarl Rupp                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
308da80777bSKarl Rupp       At[3][3] = {{0,0,0},
309da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
310da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
311da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
312da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
313da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
314da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3150298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
316a3a57f36SJed Brown   }
3176cf0794eSJed Brown   {                             /* Optimal for linear implicit part */
3186cf0794eSJed Brown     const PetscReal
3196cf0794eSJed Brown       A[3][3] = {{0,0,0},
3206cf0794eSJed Brown                  {0.5,0,0},
3216cf0794eSJed Brown                  {0.5,0.5,0}},
3226cf0794eSJed Brown       At[3][3] = {{0.25,0,0},
3236cf0794eSJed Brown                   {0,0.25,0},
3246cf0794eSJed Brown                   {1./3,1./3,1./3}};
3250298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr);
3266cf0794eSJed Brown   }
327a3a57f36SJed Brown   {
328a3a57f36SJed Brown     const PetscReal
329a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
3304040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
3314040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
3324040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
333a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
3344040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
3354040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
3364040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
337cc46b9d1SJed Brown       bembedt[4]     = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
3384040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
3394040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
3404040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
3414040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
3420298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
343a3a57f36SJed Brown   }
344a3a57f36SJed Brown   {
345a3a57f36SJed Brown     const PetscReal
346e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3476cf0794eSJed Brown                  {1./2,0,0,0,0},
3486cf0794eSJed Brown                  {11./18,1./18,0,0,0},
3496cf0794eSJed Brown                  {5./6,-5./6,.5,0,0},
3506cf0794eSJed Brown                  {1./4,7./4,3./4,-7./4,0}},
3516cf0794eSJed Brown       At[5][5] = {{0,0,0,0,0},
3526cf0794eSJed Brown                   {0,1./2,0,0,0},
3536cf0794eSJed Brown                   {0,1./6,1./2,0,0},
3546cf0794eSJed Brown                   {0,-1./2,1./2,1./2,0},
355108c343cSJed Brown                   {0,3./2,-3./2,1./2,1./2}},
3560298fd71SBarry Smith     *bembedt = NULL;
3570298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
3586cf0794eSJed Brown   }
3596cf0794eSJed Brown   {
3606cf0794eSJed Brown     const PetscReal
361e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3626cf0794eSJed Brown                  {1,0,0,0,0},
3636cf0794eSJed Brown                  {4./9,2./9,0,0,0},
3646cf0794eSJed Brown                  {1./4,0,3./4,0,0},
3656cf0794eSJed Brown                  {1./4,0,3./5,0,0}},
366e74514c0SSatish Balay       At[5][5] = {{0,0,0,0,0},
3676cf0794eSJed Brown                   {.5,.5,0,0,0},
3686cf0794eSJed Brown                   {5./18,-1./9,.5,0,0},
3696cf0794eSJed Brown                   {.5,0,0,.5,0},
370108c343cSJed Brown                   {.25,0,.75,-.5,.5}},
3710298fd71SBarry Smith     *bembedt = NULL;
3720298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
3736cf0794eSJed Brown   }
3746cf0794eSJed Brown   {
3756cf0794eSJed Brown     const PetscReal
376a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
377a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
3784040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
3794040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
3804040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
3814040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
382a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
383a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
3844040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
3854040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
3864040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
3874040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
388cc46b9d1SJed Brown       bembedt[6]     = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
3894040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
390cd652676SJed Brown                         {0,0,0},
3914040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
3924040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
3934040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
3944040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
3950298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
396a3a57f36SJed Brown   }
397a3a57f36SJed Brown   {
398a3a57f36SJed Brown     const PetscReal
399a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
400a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
4014040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
4024040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
4034040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
4044040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
4054040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
4064040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
407a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
4084040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
4094040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
4104040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
4114040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
4124040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
4134040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
4144040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
415cc46b9d1SJed Brown       bembedt[8]     = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
4164040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614.,  43486358583215./12773830924787., -9257016797708./5021505065439.},
417cd652676SJed Brown                         {0,  0, 0                            },
418cd652676SJed Brown                         {0,  0, 0                            },
4194040e9f2SJed Brown                         {65168852399939./7868540260826.,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
4204040e9f2SJed Brown                         {15494834004392./5936557850923.,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
4214040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473., 30029262896817./10175596800299.},
4224040e9f2SJed Brown                         {-19024464361622./5461577185407.,  115839755401235./10719374521269., -26136350496073./3983972220547.},
4234040e9f2SJed Brown                         {-6511271360970./6095937251113.,  5843115559534./2180450260947., -5289405421727./3760307252460. }};
4240298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
425a3a57f36SJed Brown   }
4268a381b04SJed Brown   PetscFunctionReturn(0);
4278a381b04SJed Brown }
4288a381b04SJed Brown 
4298a381b04SJed Brown #undef __FUNCT__
4308a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
4318a381b04SJed Brown /*@C
4328a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
4338a381b04SJed Brown 
4348a381b04SJed Brown    Not Collective
4358a381b04SJed Brown 
4368a381b04SJed Brown    Level: advanced
4378a381b04SJed Brown 
4388a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
4398a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
4408a381b04SJed Brown @*/
4418a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
4428a381b04SJed Brown {
4438a381b04SJed Brown   PetscErrorCode ierr;
4448a381b04SJed Brown   ARKTableauLink link;
4458a381b04SJed Brown 
4468a381b04SJed Brown   PetscFunctionBegin;
4478a381b04SJed Brown   while ((link = ARKTableauList)) {
4488a381b04SJed Brown     ARKTableau t = &link->tab;
4498a381b04SJed Brown     ARKTableauList = link->next;
4508a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
451108c343cSJed Brown     ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr);
452cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
4538a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
4548a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
4558a381b04SJed Brown   }
4568a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
4578a381b04SJed Brown   PetscFunctionReturn(0);
4588a381b04SJed Brown }
4598a381b04SJed Brown 
4608a381b04SJed Brown #undef __FUNCT__
4618a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
4628a381b04SJed Brown /*@C
4638a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
4648a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
4658a381b04SJed Brown   when using static libraries.
4668a381b04SJed Brown 
4678a381b04SJed Brown   Input Parameter:
4680298fd71SBarry Smith   path - The dynamic library path, or NULL
4698a381b04SJed Brown 
4708a381b04SJed Brown   Level: developer
4718a381b04SJed Brown 
4728a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
4738a381b04SJed Brown .seealso: PetscInitialize()
4748a381b04SJed Brown @*/
4758a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
4768a381b04SJed Brown {
4778a381b04SJed Brown   PetscErrorCode ierr;
4788a381b04SJed Brown 
4798a381b04SJed Brown   PetscFunctionBegin;
4808a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
4818a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
4828a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
483e817cc15SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
4848a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
4858a381b04SJed Brown   PetscFunctionReturn(0);
4868a381b04SJed Brown }
4878a381b04SJed Brown 
4888a381b04SJed Brown #undef __FUNCT__
4898a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
4908a381b04SJed Brown /*@C
4918a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
4928a381b04SJed Brown   called from PetscFinalize().
4938a381b04SJed Brown 
4948a381b04SJed Brown   Level: developer
4958a381b04SJed Brown 
4968a381b04SJed Brown .keywords: Petsc, destroy, package
4978a381b04SJed Brown .seealso: PetscFinalize()
4988a381b04SJed Brown @*/
4998a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
5008a381b04SJed Brown {
5018a381b04SJed Brown   PetscErrorCode ierr;
5028a381b04SJed Brown 
5038a381b04SJed Brown   PetscFunctionBegin;
5048a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
5058a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
5068a381b04SJed Brown   PetscFunctionReturn(0);
5078a381b04SJed Brown }
5088a381b04SJed Brown 
5098a381b04SJed Brown #undef __FUNCT__
5108a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
511cd652676SJed Brown /*@C
512cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
513cd652676SJed Brown 
514cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
515cd652676SJed Brown 
516cd652676SJed Brown    Input Parameters:
517cd652676SJed Brown +  name - identifier for method
518cd652676SJed Brown .  order - approximation order of method
519cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
520cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
5210298fd71SBarry Smith .  bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At)
5220298fd71SBarry Smith .  ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At)
523cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
5240298fd71SBarry Smith .  b - Non-stiff step completion table (dimension s; NULL to use last row of At)
5250298fd71SBarry Smith .  c - Non-stiff abscissa (dimension s; NULL to use row sums of A)
5260298fd71SBarry Smith .  bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available)
5270298fd71SBarry Smith .  bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided)
528cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
529cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
5300298fd71SBarry Smith -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt)
531cd652676SJed Brown 
532cd652676SJed Brown    Notes:
533cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
534cd652676SJed Brown 
535cd652676SJed Brown    Level: advanced
536cd652676SJed Brown 
537cd652676SJed Brown .keywords: TS, register
538cd652676SJed Brown 
539cd652676SJed Brown .seealso: TSARKIMEX
540cd652676SJed Brown @*/
54119fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
5428a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
543cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
544108c343cSJed Brown                                  const PetscReal bembedt[],const PetscReal bembed[],
545cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
5468a381b04SJed Brown {
5478a381b04SJed Brown   PetscErrorCode ierr;
5488a381b04SJed Brown   ARKTableauLink link;
5498a381b04SJed Brown   ARKTableau     t;
5508a381b04SJed Brown   PetscInt       i,j;
5518a381b04SJed Brown 
5528a381b04SJed Brown   PetscFunctionBegin;
5538a381b04SJed Brown   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
554cd652676SJed Brown   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
5558a381b04SJed Brown   t        = &link->tab;
5568a381b04SJed Brown   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5578a381b04SJed Brown   t->order = order;
5588a381b04SJed Brown   t->s     = s;
5598a381b04SJed 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);
5608a381b04SJed Brown   ierr     = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
5618a381b04SJed Brown   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5628a381b04SJed Brown   if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); }
5638a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
5648a381b04SJed Brown   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
5658a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
5668a381b04SJed Brown   if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); }
5678a381b04SJed 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];
5688a381b04SJed Brown   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
5698a381b04SJed 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];
570e817cc15SEmil Constantinescu   t->stiffly_accurate = PETSC_TRUE;
571e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
572e817cc15SEmil Constantinescu   t->explicit_first_stage = PETSC_TRUE;
573e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
574e817cc15SEmil Constantinescu   /*def of FSAL can be made more precise*/
5754e9d4bf5SJed Brown   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
576108c343cSJed Brown   if (bembedt) {
577108c343cSJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembedt,s,PetscReal,&t->bembed);CHKERRQ(ierr);
578108c343cSJed Brown     ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr);
579108c343cSJed Brown     ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr);
580108c343cSJed Brown   }
581108c343cSJed Brown 
5824f385281SJed Brown   t->pinterp     = pinterp;
583cd652676SJed Brown   ierr           = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
584cd652676SJed Brown   ierr           = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
585cd652676SJed Brown   ierr           = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
5868a381b04SJed Brown   link->next     = ARKTableauList;
5878a381b04SJed Brown   ARKTableauList = link;
5888a381b04SJed Brown   PetscFunctionReturn(0);
5898a381b04SJed Brown }
5908a381b04SJed Brown 
5918a381b04SJed Brown #undef __FUNCT__
592108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX"
593108c343cSJed Brown /*
594108c343cSJed Brown  The step completion formula is
595108c343cSJed Brown 
596108c343cSJed Brown  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
597108c343cSJed Brown 
598108c343cSJed Brown  This function can be called before or after ts->vec_sol has been updated.
599108c343cSJed Brown  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
600108c343cSJed Brown  We can write
601108c343cSJed Brown 
602108c343cSJed Brown  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
603108c343cSJed Brown      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
604108c343cSJed Brown      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
605108c343cSJed Brown 
606108c343cSJed Brown  so we can evaluate the method with different order even after the step has been optimistically completed.
607108c343cSJed Brown */
608108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
609108c343cSJed Brown {
610108c343cSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
611108c343cSJed Brown   ARKTableau     tab  = ark->tableau;
612108c343cSJed Brown   PetscScalar    *w   = ark->work;
613108c343cSJed Brown   PetscReal      h;
614108c343cSJed Brown   PetscInt       s = tab->s,j;
615108c343cSJed Brown   PetscErrorCode ierr;
616108c343cSJed Brown 
617108c343cSJed Brown   PetscFunctionBegin;
618108c343cSJed Brown   switch (ark->status) {
619108c343cSJed Brown   case TS_STEP_INCOMPLETE:
620108c343cSJed Brown   case TS_STEP_PENDING:
621108c343cSJed Brown     h = ts->time_step; break;
622108c343cSJed Brown   case TS_STEP_COMPLETE:
623108c343cSJed Brown     h = ts->time_step_prev; break;
624ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
625108c343cSJed Brown   }
626108c343cSJed Brown   if (order == tab->order) {
627e817cc15SEmil Constantinescu     if (ark->status == TS_STEP_INCOMPLETE) {
628740132f1SEmil Constantinescu       if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */
629e817cc15SEmil Constantinescu         ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr);
630e817cc15SEmil Constantinescu       } else { /* Use the standard completion formula (bt,b) */
631108c343cSJed Brown         ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
632e817cc15SEmil Constantinescu         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
633108c343cSJed Brown         ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
634e817cc15SEmil Constantinescu         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
635108c343cSJed Brown           for (j=0; j<s; j++) w[j] = h*tab->b[j];
636108c343cSJed Brown           ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
637e817cc15SEmil Constantinescu         }
638e817cc15SEmil Constantinescu       }
639108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
640108c343cSJed Brown     if (done) *done = PETSC_TRUE;
641108c343cSJed Brown     PetscFunctionReturn(0);
642108c343cSJed Brown   } else if (order == tab->order-1) {
643108c343cSJed Brown     if (!tab->bembedt) goto unavailable;
644108c343cSJed Brown     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
645108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
646e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
647108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
648108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
649108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
650108c343cSJed Brown     } else {                    /* Rollback and re-complete using (bet-be,be-b) */
651108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
652e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
653108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr);
654108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
655108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
656108c343cSJed Brown     }
657108c343cSJed Brown     if (done) *done = PETSC_TRUE;
658108c343cSJed Brown     PetscFunctionReturn(0);
659108c343cSJed Brown   }
660108c343cSJed Brown unavailable:
661108c343cSJed Brown   if (done) *done = PETSC_FALSE;
662ce94432eSBarry 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);
663108c343cSJed Brown   PetscFunctionReturn(0);
664108c343cSJed Brown }
665108c343cSJed Brown 
666108c343cSJed Brown #undef __FUNCT__
6678a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
6688a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
6698a381b04SJed Brown {
6708a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
6718a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
6728a381b04SJed Brown   const PetscInt  s    = tab->s;
6738a381b04SJed Brown   const PetscReal *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
674406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
675e817cc15SEmil Constantinescu   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,W = ark->Work,Z = ark->Z;
676108c343cSJed Brown   TSAdapt         adapt;
6778a381b04SJed Brown   SNES            snes;
678108c343cSJed Brown   PetscInt        i,j,its,lits,reject,next_scheme;
679cdbf8f93SLisandro Dalcin   PetscReal       next_time_step;
680108c343cSJed Brown   PetscReal       t;
681108c343cSJed Brown   PetscBool       accept;
6828a381b04SJed Brown   PetscErrorCode  ierr;
6838a381b04SJed Brown 
6848a381b04SJed Brown   PetscFunctionBegin;
685e817cc15SEmil Constantinescu   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage) {
686e817cc15SEmil Constantinescu     PetscReal valid_time;
687e817cc15SEmil Constantinescu     PetscBool isvalid;
688e817cc15SEmil Constantinescu     ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,
689e817cc15SEmil Constantinescu                                           explicit_stage_time_id,
690e817cc15SEmil Constantinescu                                           valid_time,
691e817cc15SEmil Constantinescu                                           isvalid);
692e817cc15SEmil Constantinescu     CHKERRQ(ierr);
693e817cc15SEmil Constantinescu     if (!isvalid || valid_time != ts->ptime) {
694e817cc15SEmil Constantinescu       TS        ts_start;
695e817cc15SEmil Constantinescu       SNES      snes_start;
696740132f1SEmil Constantinescu       DM        dm;
697740132f1SEmil Constantinescu       PetscReal atol;
698740132f1SEmil Constantinescu       Vec       vatol;
699740132f1SEmil Constantinescu       PetscReal rtol;
700740132f1SEmil Constantinescu       Vec       vrtol;
70119436ca2SJed Brown 
70219436ca2SJed Brown       ierr = TSCreate(PETSC_COMM_WORLD,&ts_start);CHKERRQ(ierr);
70319436ca2SJed Brown       ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr);
70419436ca2SJed Brown       ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr);
705e817cc15SEmil Constantinescu       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
706740132f1SEmil Constantinescu       ierr = TSSetDM(ts_start,dm);CHKERRQ(ierr);
707bbd56ea5SKarl Rupp 
708e817cc15SEmil Constantinescu       ts_start->adapt=ts->adapt;
709740132f1SEmil Constantinescu       PetscObjectReference((PetscObject)ts_start->adapt);
710bbd56ea5SKarl Rupp 
711e817cc15SEmil Constantinescu       ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
712e817cc15SEmil Constantinescu       ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr);
713eb082435SEmil Constantinescu       ierr = TSSetDuration(ts_start,1,ts->ptime+ts->time_step);CHKERRQ(ierr);
714740132f1SEmil Constantinescu       ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr);
715e817cc15SEmil Constantinescu       ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
716740132f1SEmil Constantinescu       ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
717e817cc15SEmil Constantinescu       ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
718e817cc15SEmil Constantinescu       ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr);
719740132f1SEmil Constantinescu       ierr = TSGetTolerances(ts,&atol,&vatol,&rtol,&vrtol);CHKERRQ(ierr);
720740132f1SEmil Constantinescu       ierr = TSSetTolerances(ts_start,atol,vatol,rtol,vrtol);CHKERRQ(ierr);
721e817cc15SEmil Constantinescu       ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
722e817cc15SEmil Constantinescu       ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr);
723bbd56ea5SKarl Rupp 
724740132f1SEmil Constantinescu       ts->time_step = ts_start->time_step;
725740132f1SEmil Constantinescu       ts->steps++;
726e817cc15SEmil Constantinescu       ierr = VecCopy(((TS_ARKIMEX*)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr);
727*166a6834SEmil Constantinescu       ts_start->snes=NULL;
728740132f1SEmil Constantinescu       ierr = TSSetSNES(ts,snes_start);CHKERRQ(ierr);
729*166a6834SEmil Constantinescu       ierr = SNESDestroy(&snes_start);CHKERRQ(ierr);
730*166a6834SEmil Constantinescu       ierr = TSDestroy(&ts_start);CHKERRQ(ierr);
731e817cc15SEmil Constantinescu     }
732e817cc15SEmil Constantinescu   }
733e817cc15SEmil Constantinescu 
7348a381b04SJed Brown   ierr           = TSGetSNES(ts,&snes);CHKERRQ(ierr);
735cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
7368a381b04SJed Brown   t              = ts->ptime;
737108c343cSJed Brown   accept         = PETSC_TRUE;
738108c343cSJed Brown   ark->status    = TS_STEP_INCOMPLETE;
7398a381b04SJed Brown 
740e817cc15SEmil Constantinescu 
74197335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
742108c343cSJed Brown     PetscReal h = ts->time_step;
743b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
7448a381b04SJed Brown     for (i=0; i<s; i++) {
7458a381b04SJed Brown       if (At[i*s+i] == 0) {           /* This stage is explicit */
7468a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
747e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7488a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
7498a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7508a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7518a381b04SJed Brown       } else {
7528a381b04SJed Brown         ark->stage_time = t + h*ct[i];
753b296d7d5SJed Brown         ark->scoeff     = 1./At[i*s+i];
754b8123daeSJed Brown         ierr            = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
7558a381b04SJed Brown         /* Affine part */
7568a381b04SJed Brown         ierr = VecZeroEntries(W);CHKERRQ(ierr);
7578a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7588a381b04SJed Brown         ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
759b296d7d5SJed Brown         ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr);
760f16577ceSEmil Constantinescu 
7618a381b04SJed Brown         /* Ydot = shift*(Y-Z) */
7628a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
763e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7644f385281SJed Brown         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
765f16577ceSEmil Constantinescu 
7668a381b04SJed Brown         /* Initial guess taken from last stage */
7678a381b04SJed Brown         ierr          = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr);
7688a381b04SJed Brown         ierr          = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
769e817cc15SEmil Constantinescu         ierr          = (ts->ops->snesfunction)(snes,Y[i],W,ts);CHKERRQ(ierr);
7708a381b04SJed Brown         ierr          = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
7718a381b04SJed Brown         ierr          = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
7725ef26d82SJed Brown         ts->snes_its += its; ts->ksp_its += lits;
773ad6bc421SBarry Smith         ierr          = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
77497335746SJed Brown         ierr          = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
77597335746SJed Brown         if (!accept) goto reject_step;
7768a381b04SJed Brown       }
777e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) {
778e817cc15SEmil Constantinescu         if (i==0 && tab->explicit_first_stage) {
779e817cc15SEmil Constantinescu           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);
780e817cc15SEmil Constantinescu         } else {
781e817cc15SEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
782e817cc15SEmil Constantinescu         }
783e817cc15SEmil Constantinescu       } else {
7848a381b04SJed Brown         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
7854cc180ffSJed Brown         ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
786e817cc15SEmil Constantinescu         ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr);
7874cc180ffSJed Brown         if (ark->imex) {
7888a381b04SJed Brown           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
7894cc180ffSJed Brown         } else {
7904cc180ffSJed Brown           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
7914cc180ffSJed Brown         }
7928a381b04SJed Brown       }
793e817cc15SEmil Constantinescu     }
7940298fd71SBarry Smith     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
795108c343cSJed Brown     ark->status = TS_STEP_PENDING;
7968a381b04SJed Brown 
797108c343cSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
798ad6bc421SBarry Smith     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
799108c343cSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
800108c343cSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
801108c343cSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
802108c343cSJed Brown     if (accept) {
803108c343cSJed Brown       /* ignore next_scheme for now */
8048a381b04SJed Brown       ts->ptime    += ts->time_step;
805cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
8068a381b04SJed Brown       ts->steps++;
807e817cc15SEmil Constantinescu       if (ts->equation_type>=TS_EQ_IMPLICIT) { /* save the initial slope for the next step*/
808e817cc15SEmil Constantinescu         ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
809e817cc15SEmil Constantinescu       }
810108c343cSJed Brown       ark->status = TS_STEP_COMPLETE;
811e817cc15SEmil Constantinescu       if (tab->explicit_first_stage) {
812e817cc15SEmil Constantinescu         ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
813e817cc15SEmil Constantinescu       }
814e817cc15SEmil Constantinescu 
815108c343cSJed Brown       break;
816108c343cSJed Brown     } else {                    /* Roll back the current step */
8172c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*bt[j];
818108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotI);CHKERRQ(ierr);
8192c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
820108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotRHS);CHKERRQ(ierr);
821108c343cSJed Brown       ts->time_step = next_time_step;
822108c343cSJed Brown       ark->status   = TS_STEP_INCOMPLETE;
823108c343cSJed Brown     }
824476b6736SJed Brown reject_step: continue;
825108c343cSJed Brown   }
826b2ce242eSJed Brown   if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
8278a381b04SJed Brown   PetscFunctionReturn(0);
8288a381b04SJed Brown }
8298a381b04SJed Brown 
830cd652676SJed Brown #undef __FUNCT__
831cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
832cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
833cd652676SJed Brown {
834cd652676SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8354f385281SJed Brown   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
836108c343cSJed Brown   PetscReal       h;
837108c343cSJed Brown   PetscReal       tt,t;
838cd652676SJed Brown   PetscScalar     *bt,*b;
839cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
840cd652676SJed Brown   PetscErrorCode  ierr;
841cd652676SJed Brown 
842cd652676SJed Brown   PetscFunctionBegin;
843ce94432eSBarry Smith   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
844108c343cSJed Brown   switch (ark->status) {
845108c343cSJed Brown   case TS_STEP_INCOMPLETE:
846108c343cSJed Brown   case TS_STEP_PENDING:
847108c343cSJed Brown     h = ts->time_step;
848108c343cSJed Brown     t = (itime - ts->ptime)/h;
849108c343cSJed Brown     break;
850108c343cSJed Brown   case TS_STEP_COMPLETE:
851108c343cSJed Brown     h = ts->time_step_prev;
852108c343cSJed Brown     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
853108c343cSJed Brown     break;
854ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
855108c343cSJed Brown   }
856cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
857cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
8584f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
859cd652676SJed Brown     for (i=0; i<s; i++) {
860108c343cSJed Brown       bt[i] += h * Bt[i*pinterp+j] * tt * -1.0;
861108c343cSJed Brown       b[i]  += h * B[i*pinterp+j] * tt;
862cd652676SJed Brown     }
863cd652676SJed Brown   }
864ce94432eSBarry Smith   if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"First stage not explicit so starting stage not saved");
865cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
866cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
867cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
868cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
869cd652676SJed Brown   PetscFunctionReturn(0);
870cd652676SJed Brown }
871cd652676SJed Brown 
8728a381b04SJed Brown /*------------------------------------------------------------*/
8738a381b04SJed Brown #undef __FUNCT__
8748a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
8758a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
8768a381b04SJed Brown {
8778a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
8788a381b04SJed Brown   PetscInt       s;
8798a381b04SJed Brown   PetscErrorCode ierr;
8808a381b04SJed Brown 
8818a381b04SJed Brown   PetscFunctionBegin;
8828a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
8838a381b04SJed Brown   s    = ark->tableau->s;
8848a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
8858a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
8868a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
8878a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
8888a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
889e817cc15SEmil Constantinescu   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
8908a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
8918a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
8928a381b04SJed Brown   PetscFunctionReturn(0);
8938a381b04SJed Brown }
8948a381b04SJed Brown 
8958a381b04SJed Brown #undef __FUNCT__
8968a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
8978a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
8988a381b04SJed Brown {
8998a381b04SJed Brown   PetscErrorCode ierr;
9008a381b04SJed Brown 
9018a381b04SJed Brown   PetscFunctionBegin;
9028a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
9038a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
90400de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C","",NULL);CHKERRQ(ierr);
90500de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C","",NULL);CHKERRQ(ierr);
90600de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",NULL);CHKERRQ(ierr);
9078a381b04SJed Brown   PetscFunctionReturn(0);
9088a381b04SJed Brown }
9098a381b04SJed Brown 
910d5e6173cSPeter Brune 
911d5e6173cSPeter Brune #undef __FUNCT__
912d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXGetVecs"
913d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
914d5e6173cSPeter Brune {
915d5e6173cSPeter Brune   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
916d5e6173cSPeter Brune   PetscErrorCode ierr;
917d5e6173cSPeter Brune 
918d5e6173cSPeter Brune   PetscFunctionBegin;
919d5e6173cSPeter Brune   if (Z) {
920d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
921d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
922d5e6173cSPeter Brune     } else *Z = ax->Z;
923d5e6173cSPeter Brune   }
924d5e6173cSPeter Brune   if (Ydot) {
925d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
926d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
927d5e6173cSPeter Brune     } else *Ydot = ax->Ydot;
928d5e6173cSPeter Brune   }
929d5e6173cSPeter Brune   PetscFunctionReturn(0);
930d5e6173cSPeter Brune }
931d5e6173cSPeter Brune 
932d5e6173cSPeter Brune 
933d5e6173cSPeter Brune #undef __FUNCT__
934d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXRestoreVecs"
935d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
936d5e6173cSPeter Brune {
937d5e6173cSPeter Brune   PetscErrorCode ierr;
938d5e6173cSPeter Brune 
939d5e6173cSPeter Brune   PetscFunctionBegin;
940d5e6173cSPeter Brune   if (Z) {
941d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
942d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
943d5e6173cSPeter Brune     }
944d5e6173cSPeter Brune   }
945d5e6173cSPeter Brune   if (Ydot) {
946d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
947d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
948d5e6173cSPeter Brune     }
949d5e6173cSPeter Brune   }
950d5e6173cSPeter Brune   PetscFunctionReturn(0);
951d5e6173cSPeter Brune }
952d5e6173cSPeter Brune 
9538a381b04SJed Brown /*
9548a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9558a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9568a381b04SJed Brown */
9578a381b04SJed Brown #undef __FUNCT__
9588a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
9598a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
9608a381b04SJed Brown {
9618a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
962d5e6173cSPeter Brune   DM             dm,dmsave;
963d5e6173cSPeter Brune   Vec            Z,Ydot;
964b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
9658a381b04SJed Brown   PetscErrorCode ierr;
9668a381b04SJed Brown 
9678a381b04SJed Brown   PetscFunctionBegin;
968d5e6173cSPeter Brune   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
969d5e6173cSPeter Brune   ierr   = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
970b296d7d5SJed Brown   ierr   = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
971d5e6173cSPeter Brune   dmsave = ts->dm;
972d5e6173cSPeter Brune   ts->dm = dm;
973740132f1SEmil Constantinescu 
974d5e6173cSPeter Brune   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
975e817cc15SEmil Constantinescu 
976d5e6173cSPeter Brune   ts->dm = dmsave;
977d5e6173cSPeter Brune   ierr   = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
9788a381b04SJed Brown   PetscFunctionReturn(0);
9798a381b04SJed Brown }
9808a381b04SJed Brown 
9818a381b04SJed Brown #undef __FUNCT__
9828a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
9838a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
9848a381b04SJed Brown {
9858a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
986d5e6173cSPeter Brune   DM             dm,dmsave;
987d5e6173cSPeter Brune   Vec            Ydot;
988b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
9898a381b04SJed Brown   PetscErrorCode ierr;
9908a381b04SJed Brown 
9918a381b04SJed Brown   PetscFunctionBegin;
992d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
9930298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
9948a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
995d5e6173cSPeter Brune   dmsave = ts->dm;
996d5e6173cSPeter Brune   ts->dm = dm;
997740132f1SEmil Constantinescu 
998b296d7d5SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,str,ark->imex);CHKERRQ(ierr);
999740132f1SEmil Constantinescu 
1000d5e6173cSPeter Brune   ts->dm = dmsave;
10010298fd71SBarry Smith   ierr   = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
1002d5e6173cSPeter Brune   PetscFunctionReturn(0);
1003d5e6173cSPeter Brune }
1004d5e6173cSPeter Brune 
1005d5e6173cSPeter Brune #undef __FUNCT__
1006d5e6173cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSARKIMEX"
1007d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1008d5e6173cSPeter Brune {
1009d5e6173cSPeter Brune   PetscFunctionBegin;
1010d5e6173cSPeter Brune   PetscFunctionReturn(0);
1011d5e6173cSPeter Brune }
1012d5e6173cSPeter Brune 
1013d5e6173cSPeter Brune #undef __FUNCT__
1014d5e6173cSPeter Brune #define __FUNCT__ "DMRestrictHook_TSARKIMEX"
1015d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1016d5e6173cSPeter Brune {
1017d5e6173cSPeter Brune   TS             ts = (TS)ctx;
1018d5e6173cSPeter Brune   PetscErrorCode ierr;
1019d5e6173cSPeter Brune   Vec            Z,Z_c;
1020d5e6173cSPeter Brune 
1021d5e6173cSPeter Brune   PetscFunctionBegin;
10220298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10230298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
1024d5e6173cSPeter Brune   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1025d5e6173cSPeter Brune   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
10260298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10270298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
10288a381b04SJed Brown   PetscFunctionReturn(0);
10298a381b04SJed Brown }
10308a381b04SJed Brown 
1031cdb298fcSPeter Brune 
1032cdb298fcSPeter Brune #undef __FUNCT__
1033cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainHook_TSARKIMEX"
1034cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1035cdb298fcSPeter Brune {
1036cdb298fcSPeter Brune   PetscFunctionBegin;
1037cdb298fcSPeter Brune   PetscFunctionReturn(0);
1038cdb298fcSPeter Brune }
1039cdb298fcSPeter Brune 
1040cdb298fcSPeter Brune #undef __FUNCT__
1041cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX"
1042cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1043cdb298fcSPeter Brune {
1044cdb298fcSPeter Brune   TS             ts = (TS)ctx;
1045cdb298fcSPeter Brune   PetscErrorCode ierr;
1046cdb298fcSPeter Brune   Vec            Z,Z_c;
1047cdb298fcSPeter Brune 
1048cdb298fcSPeter Brune   PetscFunctionBegin;
10490298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
10500298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1051cdb298fcSPeter Brune 
1052cdb298fcSPeter Brune   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1053cdb298fcSPeter Brune   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1054cdb298fcSPeter Brune 
10550298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
10560298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1057cdb298fcSPeter Brune   PetscFunctionReturn(0);
1058cdb298fcSPeter Brune }
1059cdb298fcSPeter Brune 
10608a381b04SJed Brown #undef __FUNCT__
10618a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
10628a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
10638a381b04SJed Brown {
10648a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1065f2c2a1b9SBarry Smith   ARKTableau     tab;
1066f2c2a1b9SBarry Smith   PetscInt       s;
10678a381b04SJed Brown   PetscErrorCode ierr;
1068d5e6173cSPeter Brune   DM             dm;
1069f9c1d6abSBarry Smith 
10708a381b04SJed Brown   PetscFunctionBegin;
10718a381b04SJed Brown   if (!ark->tableau) {
1072e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
10738a381b04SJed Brown   }
1074f2c2a1b9SBarry Smith   tab  = ark->tableau;
1075f2c2a1b9SBarry Smith   s    = tab->s;
10768a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
10778a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
10788a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
10798a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
10808a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
1081e817cc15SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
10828a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
10838a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
1084d5e6173cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1085d5e6173cSPeter Brune   if (dm) {
1086d5e6173cSPeter Brune     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1087cdb298fcSPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1088d5e6173cSPeter Brune   }
10898a381b04SJed Brown   PetscFunctionReturn(0);
10908a381b04SJed Brown }
10918a381b04SJed Brown /*------------------------------------------------------------*/
10928a381b04SJed Brown 
10938a381b04SJed Brown #undef __FUNCT__
10948a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
10958a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
10968a381b04SJed Brown {
10974cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
10988a381b04SJed Brown   PetscErrorCode ierr;
10998a381b04SJed Brown   char           arktype[256];
11008a381b04SJed Brown 
11018a381b04SJed Brown   PetscFunctionBegin;
11028a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
11038a381b04SJed Brown   {
11048a381b04SJed Brown     ARKTableauLink link;
11058a381b04SJed Brown     PetscInt       count,choice;
11068a381b04SJed Brown     PetscBool      flg;
11078a381b04SJed Brown     const char     **namelist;
11088caf3d72SBarry Smith     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr);
11098a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
11108a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
11118a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11128a381b04SJed Brown     ierr      = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
11138a381b04SJed Brown     ierr      = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
11148a381b04SJed Brown     ierr      = PetscFree(namelist);CHKERRQ(ierr);
11154cc180ffSJed Brown     flg       = (PetscBool) !ark->imex;
11160298fd71SBarry Smith     ierr      = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr);
11174cc180ffSJed Brown     ark->imex = (PetscBool) !flg;
1118d52bd9f3SBarry Smith     ierr      = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
11198a381b04SJed Brown   }
11208a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
11218a381b04SJed Brown   PetscFunctionReturn(0);
11228a381b04SJed Brown }
11238a381b04SJed Brown 
11248a381b04SJed Brown #undef __FUNCT__
11258a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
11268a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
11278a381b04SJed Brown {
1128257d2499SJed Brown   PetscErrorCode ierr;
1129f1d86077SJed Brown   PetscInt       i;
1130f1d86077SJed Brown   size_t         left,count;
11318a381b04SJed Brown   char           *p;
11328a381b04SJed Brown 
11338a381b04SJed Brown   PetscFunctionBegin;
1134f1d86077SJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
1135f1d86077SJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
11368a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
11378a381b04SJed Brown     left -= count;
11388a381b04SJed Brown     p    += count;
11398a381b04SJed Brown     *p++  = ' ';
11408a381b04SJed Brown   }
11418a381b04SJed Brown   p[i ? 0 : -1] = 0;
11428a381b04SJed Brown   PetscFunctionReturn(0);
11438a381b04SJed Brown }
11448a381b04SJed Brown 
11458a381b04SJed Brown #undef __FUNCT__
11468a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
11478a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
11488a381b04SJed Brown {
11498a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11508a381b04SJed Brown   ARKTableau     tab  = ark->tableau;
11518a381b04SJed Brown   PetscBool      iascii;
11528a381b04SJed Brown   PetscErrorCode ierr;
1153559eea31SJed Brown   TSAdapt        adapt;
11548a381b04SJed Brown 
11558a381b04SJed Brown   PetscFunctionBegin;
1156251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11578a381b04SJed Brown   if (iascii) {
115819fd82e9SBarry Smith     TSARKIMEXType arktype;
11598a381b04SJed Brown     char          buf[512];
11608a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
11618a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
11628caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
116331f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
11648caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1165e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
1166e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr);
1167e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr);
116831f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
11698a381b04SJed Brown   }
1170ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
1171559eea31SJed Brown   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
1172d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
11738a381b04SJed Brown   PetscFunctionReturn(0);
11748a381b04SJed Brown }
11758a381b04SJed Brown 
11768a381b04SJed Brown #undef __FUNCT__
1177f2c2a1b9SBarry Smith #define __FUNCT__ "TSLoad_ARKIMEX"
1178f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1179f2c2a1b9SBarry Smith {
1180f2c2a1b9SBarry Smith   PetscErrorCode ierr;
1181f2c2a1b9SBarry Smith   SNES           snes;
1182ad6bc421SBarry Smith   TSAdapt        tsadapt;
1183f2c2a1b9SBarry Smith 
1184f2c2a1b9SBarry Smith   PetscFunctionBegin;
1185ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&tsadapt);CHKERRQ(ierr);
1186ad6bc421SBarry Smith   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
1187f2c2a1b9SBarry Smith   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1188f2c2a1b9SBarry Smith   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1189ad6bc421SBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
11900298fd71SBarry Smith   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
11910298fd71SBarry Smith   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
1192f2c2a1b9SBarry Smith   PetscFunctionReturn(0);
1193f2c2a1b9SBarry Smith }
1194f2c2a1b9SBarry Smith 
1195f2c2a1b9SBarry Smith #undef __FUNCT__
11968a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
11978a381b04SJed Brown /*@C
11988a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
11998a381b04SJed Brown 
12008a381b04SJed Brown   Logically collective
12018a381b04SJed Brown 
12028a381b04SJed Brown   Input Parameter:
12038a381b04SJed Brown +  ts - timestepping context
12048a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
12058a381b04SJed Brown 
12068a381b04SJed Brown   Level: intermediate
12078a381b04SJed Brown 
1208020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
12098a381b04SJed Brown @*/
121019fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
12118a381b04SJed Brown {
12128a381b04SJed Brown   PetscErrorCode ierr;
12138a381b04SJed Brown 
12148a381b04SJed Brown   PetscFunctionBegin;
12158a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
121619fd82e9SBarry Smith   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
12178a381b04SJed Brown   PetscFunctionReturn(0);
12188a381b04SJed Brown }
12198a381b04SJed Brown 
12208a381b04SJed Brown #undef __FUNCT__
12218a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
12228a381b04SJed Brown /*@C
12238a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
12248a381b04SJed Brown 
12258a381b04SJed Brown   Logically collective
12268a381b04SJed Brown 
12278a381b04SJed Brown   Input Parameter:
12288a381b04SJed Brown .  ts - timestepping context
12298a381b04SJed Brown 
12308a381b04SJed Brown   Output Parameter:
12318a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
12328a381b04SJed Brown 
12338a381b04SJed Brown   Level: intermediate
12348a381b04SJed Brown 
12358a381b04SJed Brown .seealso: TSARKIMEXGetType()
12368a381b04SJed Brown @*/
123719fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
12388a381b04SJed Brown {
12398a381b04SJed Brown   PetscErrorCode ierr;
12408a381b04SJed Brown 
12418a381b04SJed Brown   PetscFunctionBegin;
12428a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
124319fd82e9SBarry Smith   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
12448a381b04SJed Brown   PetscFunctionReturn(0);
12458a381b04SJed Brown }
12468a381b04SJed Brown 
12474cc180ffSJed Brown #undef __FUNCT__
12484cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
12494cc180ffSJed Brown /*@C
12504cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
12514cc180ffSJed Brown 
12524cc180ffSJed Brown   Logically collective
12534cc180ffSJed Brown 
12544cc180ffSJed Brown   Input Parameter:
12554cc180ffSJed Brown +  ts - timestepping context
12564cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
12574cc180ffSJed Brown 
12584cc180ffSJed Brown   Level: intermediate
12594cc180ffSJed Brown 
12604cc180ffSJed Brown .seealso: TSARKIMEXGetType()
12614cc180ffSJed Brown @*/
12624cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
12634cc180ffSJed Brown {
12644cc180ffSJed Brown   PetscErrorCode ierr;
12654cc180ffSJed Brown 
12664cc180ffSJed Brown   PetscFunctionBegin;
12674cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12684cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
12694cc180ffSJed Brown   PetscFunctionReturn(0);
12704cc180ffSJed Brown }
12714cc180ffSJed Brown 
12728a381b04SJed Brown #undef __FUNCT__
12738a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
127419fd82e9SBarry Smith PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
12758a381b04SJed Brown {
12768a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
12778a381b04SJed Brown   PetscErrorCode ierr;
12788a381b04SJed Brown 
12798a381b04SJed Brown   PetscFunctionBegin;
1280f2c2a1b9SBarry Smith   if (!ark->tableau) {
1281f2c2a1b9SBarry Smith     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1282f2c2a1b9SBarry Smith   }
12838a381b04SJed Brown   *arktype = ark->tableau->name;
12848a381b04SJed Brown   PetscFunctionReturn(0);
12858a381b04SJed Brown }
12868a381b04SJed Brown #undef __FUNCT__
12878a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
128819fd82e9SBarry Smith PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
12898a381b04SJed Brown {
12908a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
12918a381b04SJed Brown   PetscErrorCode ierr;
12928a381b04SJed Brown   PetscBool      match;
12938a381b04SJed Brown   ARKTableauLink link;
12948a381b04SJed Brown 
12958a381b04SJed Brown   PetscFunctionBegin;
12968a381b04SJed Brown   if (ark->tableau) {
12978a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
12988a381b04SJed Brown     if (match) PetscFunctionReturn(0);
12998a381b04SJed Brown   }
13008a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
13018a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
13028a381b04SJed Brown     if (match) {
13038a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
13048a381b04SJed Brown       ark->tableau = &link->tab;
13058a381b04SJed Brown       PetscFunctionReturn(0);
13068a381b04SJed Brown     }
13078a381b04SJed Brown   }
1308ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
13098a381b04SJed Brown   PetscFunctionReturn(0);
13108a381b04SJed Brown }
13114cc180ffSJed Brown #undef __FUNCT__
13124cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
13134cc180ffSJed Brown PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
13144cc180ffSJed Brown {
13154cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13164cc180ffSJed Brown 
13174cc180ffSJed Brown   PetscFunctionBegin;
13184cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
13194cc180ffSJed Brown   PetscFunctionReturn(0);
13204cc180ffSJed Brown }
13218a381b04SJed Brown 
13228a381b04SJed Brown /* ------------------------------------------------------------ */
13238a381b04SJed Brown /*MC
1324a4386c9eSJed Brown       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
13258a381b04SJed Brown 
1326fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1327fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1328fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1329fca742c7SJed Brown 
1330fca742c7SJed Brown   Notes:
1331a4386c9eSJed Brown   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1332c8058688SBarry Smith 
1333a4386c9eSJed 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).
1334fca742c7SJed Brown 
13358a381b04SJed Brown   Level: beginner
13368a381b04SJed Brown 
1337c8058688SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
1338a4386c9eSJed Brown            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
13398a381b04SJed Brown 
13408a381b04SJed Brown M*/
13418a381b04SJed Brown #undef __FUNCT__
13428a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
13438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
13448a381b04SJed Brown {
13458a381b04SJed Brown   TS_ARKIMEX     *th;
13468a381b04SJed Brown   PetscErrorCode ierr;
13478a381b04SJed Brown 
13488a381b04SJed Brown   PetscFunctionBegin;
13498a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
13500298fd71SBarry Smith   ierr = TSARKIMEXInitializePackage(NULL);CHKERRQ(ierr);
13518a381b04SJed Brown #endif
13528a381b04SJed Brown 
13538a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
13548a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
13558a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
1356f2c2a1b9SBarry Smith   ts->ops->load           = TSLoad_ARKIMEX;
13578a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
13588a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
1359cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1360108c343cSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
13618a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
13628a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
13638a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
13648a381b04SJed Brown 
13658a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
13668a381b04SJed Brown   ts->data = (void*)th;
13674cc180ffSJed Brown   th->imex = PETSC_TRUE;
13688a381b04SJed Brown 
136900de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
137000de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
137100de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
13728a381b04SJed Brown   PetscFunctionReturn(0);
13738a381b04SJed Brown }
1374