xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision 2c0c504e7053a2ac2ad951e72b13ff570bb496ea)
18a381b04SJed Brown /*
28a381b04SJed Brown   Code for timestepping with additive Runge-Kutta IMEX method
38a381b04SJed Brown 
48a381b04SJed Brown   Notes:
58a381b04SJed Brown   The general system is written as
68a381b04SJed Brown 
7f9c1d6abSBarry Smith   F(t,U,Udot) = G(t,U)
88a381b04SJed Brown 
98a381b04SJed Brown   where F represents the stiff part of the physics and G represents the non-stiff part.
108a381b04SJed Brown 
118a381b04SJed Brown */
12b45d2f2cSJed Brown #include <petsc-private/tsimpl.h>                /*I   "petscts.h"   I*/
138a381b04SJed Brown 
1419fd82e9SBarry Smith static TSARKIMEXType TSARKIMEXDefault = TSARKIMEX3;
158a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled;
168a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized;
17e817cc15SEmil Constantinescu static PetscInt explicit_stage_time_id;
188a381b04SJed Brown 
198a381b04SJed Brown typedef struct _ARKTableau *ARKTableau;
208a381b04SJed Brown struct _ARKTableau {
218a381b04SJed Brown   char      *name;
224f385281SJed Brown   PetscInt  order;               /* Classical approximation order of the method */
234f385281SJed Brown   PetscInt  s;                   /* Number of stages */
24e817cc15SEmil Constantinescu   PetscBool stiffly_accurate;    /* The implicit part is stiffly accurate*/
25e817cc15SEmil Constantinescu   PetscBool FSAL_implicit;       /* The implicit part is FSAL*/
26e817cc15SEmil Constantinescu   PetscBool explicit_first_stage;/* The implicit part has an explicit first stage*/
274f385281SJed Brown   PetscInt  pinterp;             /* Interpolation order */
284f385281SJed Brown   PetscReal *At,*bt,*ct;        /* Stiff tableau */
298a381b04SJed Brown   PetscReal *A,*b,*c;           /* Non-stiff tableau */
30108c343cSJed Brown   PetscReal *bembedt,*bembed;   /* Embedded formula of order one less (order-1) */
31cd652676SJed Brown   PetscReal *binterpt,*binterp; /* Dense output formula */
32108c343cSJed Brown   PetscReal ccfl;               /* Placeholder for CFL coefficient relative to forward Euler */
338a381b04SJed Brown };
348a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink;
358a381b04SJed Brown struct _ARKTableauLink {
368a381b04SJed Brown   struct _ARKTableau tab;
378a381b04SJed Brown   ARKTableauLink next;
388a381b04SJed Brown };
398a381b04SJed Brown static ARKTableauLink ARKTableauList;
408a381b04SJed Brown 
418a381b04SJed Brown typedef struct {
428a381b04SJed Brown   ARKTableau   tableau;
438a381b04SJed Brown   Vec          *Y;               /* States computed during the step */
448a381b04SJed Brown   Vec          *YdotI;           /* Time derivatives for the stiff part */
458a381b04SJed Brown   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
46e817cc15SEmil Constantinescu   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
478a381b04SJed Brown   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
488a381b04SJed Brown   Vec          Work;             /* Generic work vector */
498a381b04SJed Brown   Vec          Z;                /* Ydot = shift(Y-Z) */
508a381b04SJed Brown   PetscScalar  *work;            /* Scalar work */
51b296d7d5SJed Brown   PetscReal    scoeff;           /* shift = scoeff/dt */
528a381b04SJed Brown   PetscReal    stage_time;
534cc180ffSJed Brown   PetscBool    imex;
54e817cc15SEmil Constantinescu   /*PetscBool    init_slope;*/
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};
234e817cc15SEmil 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*/
235*2c0c504eSEmil Constantinescu           ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr);
236e817cc15SEmil Constantinescu   }
2378a381b04SJed Brown   {
2388a381b04SJed Brown     const PetscReal
2391f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2401f80e275SEmil Constantinescu                  {0.5,0.0}},
2411f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2421f80e275SEmil Constantinescu                   {0.0,0.5}},
2431f80e275SEmil Constantinescu         b[2] = {0.0,1.0},
2441f80e275SEmil Constantinescu           bembedt[2] = {0.5,0.5};
2451f80e275SEmil 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*/
2461f80e275SEmil Constantinescu           ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr);
2471f80e275SEmil Constantinescu   }
2481f80e275SEmil Constantinescu   {
2491f80e275SEmil Constantinescu     const PetscReal
2501f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2511f80e275SEmil Constantinescu                  {1.0,0.0}},
2521f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2531f80e275SEmil Constantinescu                   {0.5,0.5}},
2541f80e275SEmil Constantinescu         b[2] = {0.5,0.5},
2551f80e275SEmil Constantinescu           bembedt[2] = {0.0,1.0};
2561f80e275SEmil 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*/
2571f80e275SEmil Constantinescu           ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,1,b,PETSC_NULL);CHKERRQ(ierr);
2581f80e275SEmil Constantinescu   }
2591f80e275SEmil Constantinescu   {
2601f80e275SEmil Constantinescu     const PetscReal us2 = 1.0-1.0/PetscSqrtReal((PetscReal)2.0);
2611f80e275SEmil Constantinescu     const PetscReal
2621f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2631f80e275SEmil Constantinescu                  {1.0,0.0}},
2641f80e275SEmil Constantinescu       At[2][2] = {{us2,0.0},
2651f80e275SEmil Constantinescu                   {1.0-2.0*us2,us2}},
2661f80e275SEmil Constantinescu         b[2] = {0.5,0.5},
2671f80e275SEmil Constantinescu           bembedt[2] = {0.0,1.0},
2681f80e275SEmil Constantinescu             binterpt[2][2] = {{(us2-1.0)/(2.0*us2-1.0),-1/(2.0*(1.0-2.0*us2))},{1-(us2-1.0)/(2.0*us2-1.0),-1/(2.0*(1.0-2.0*us2))}},
2691f80e275SEmil Constantinescu               binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
2701f80e275SEmil Constantinescu               ierr = TSARKIMEXRegister(TSARKIMEXL2,2,2,&At[0][0],b,PETSC_NULL,&A[0][0],b,PETSC_NULL,bembedt,bembedt,2,binterpt[0],binterp[0]);CHKERRQ(ierr);
2711f80e275SEmil Constantinescu   }
2721f80e275SEmil Constantinescu   {
2731f80e275SEmil Constantinescu     const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
2748a381b04SJed Brown       A[3][3] = {{0,0,0},
2751f80e275SEmil Constantinescu                  {2-s2,0,0},
276617a39beSEmil Constantinescu                  {0.5,0.5,0}},
2771f80e275SEmil Constantinescu       At[3][3] = {{0,0,0},
2781f80e275SEmil Constantinescu                   {1-1/s2,1-1/s2,0},
2791f80e275SEmil Constantinescu                   {1/(2*s2),1/(2*s2),1-1/s2}},
280e93ac1c2SEmil Constantinescu         bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)},
281ce4a059fSEmil Constantinescu         binterpt[3][2] = {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}};
2821f80e275SEmil Constantinescu     ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
2831f80e275SEmil Constantinescu   }
2841f80e275SEmil Constantinescu   {
2851f80e275SEmil Constantinescu     const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
2861f80e275SEmil Constantinescu       A[3][3] = {{0,0,0},
2871f80e275SEmil Constantinescu                  {2-s2,0,0},
2888a381b04SJed Brown                  {0.75,0.25,0}},
2898a381b04SJed Brown       At[3][3] = {{0,0,0},
2901f80e275SEmil Constantinescu                   {1-1/s2,1-1/s2,0},
2911f80e275SEmil Constantinescu                   {1/(2*s2),1/(2*s2),1-1/s2}},
292e93ac1c2SEmil Constantinescu       bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)},
293ce4a059fSEmil Constantinescu       binterpt[3][2] =  {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}};
294108c343cSJed Brown       ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
2958a381b04SJed Brown   }
29606db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
29703403c7fSJed Brown     const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),
298a3a57f36SJed Brown       A[3][3] = {{0,0,0},
299a3a57f36SJed Brown                  {2-s2,0,0},
300a3a57f36SJed Brown                  {(3-2*s2)/6,(3+2*s2)/6,0}},
301a3a57f36SJed Brown       At[3][3] = {{0,0,0},
302a3a57f36SJed Brown                   {1-1/s2,1-1/s2,0},
303cd652676SJed Brown                   {1/(2*s2),1/(2*s2),1-1/s2}},
304e93ac1c2SEmil Constantinescu       bembedt[3] = {(4.-s2)/8.,(4.-s2)/8.,1/(2.*s2)},
305ce4a059fSEmil Constantinescu       binterpt[3][2] =  {{1.0/s2,-1.0/(2.0*s2)},{1.0/s2,-1.0/(2.0*s2)},{1.0-s2,1.0/s2}};
3061f80e275SEmil Constantinescu     ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
307a3a57f36SJed Brown   }
3086cf0794eSJed Brown   {                             /* Optimal for linear implicit part */
3096cf0794eSJed Brown     const PetscReal
3106cf0794eSJed Brown       A[3][3] = {{0,0,0},
3116cf0794eSJed Brown                  {0.5,0,0},
3126cf0794eSJed Brown                  {0.5,0.5,0}},
3136cf0794eSJed Brown       At[3][3] = {{0.25,0,0},
3146cf0794eSJed Brown                   {0,0.25,0},
3156cf0794eSJed Brown                   {1./3,1./3,1./3}};
316108c343cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
3176cf0794eSJed Brown   }
318a3a57f36SJed Brown   {
319a3a57f36SJed Brown     const PetscReal
320a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
3214040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
3224040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
3234040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
324a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
3254040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
3264040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
3274040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
328cc46b9d1SJed Brown       bembedt[4] = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
3294040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
3304040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
3314040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
3324040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
333108c343cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
334a3a57f36SJed Brown   }
335a3a57f36SJed Brown   {
336a3a57f36SJed Brown     const PetscReal
337e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3386cf0794eSJed Brown                  {1./2,0,0,0,0},
3396cf0794eSJed Brown                  {11./18,1./18,0,0,0},
3406cf0794eSJed Brown                  {5./6,-5./6,.5,0,0},
3416cf0794eSJed Brown                  {1./4,7./4,3./4,-7./4,0}},
3426cf0794eSJed Brown       At[5][5] = {{0,0,0,0,0},
3436cf0794eSJed Brown                   {0,1./2,0,0,0},
3446cf0794eSJed Brown                   {0,1./6,1./2,0,0},
3456cf0794eSJed Brown                   {0,-1./2,1./2,1./2,0},
346108c343cSJed Brown                   {0,3./2,-3./2,1./2,1./2}},
347108c343cSJed Brown       *bembedt = PETSC_NULL;
348108c343cSJed Brown       ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
3496cf0794eSJed Brown   }
3506cf0794eSJed Brown   {
3516cf0794eSJed Brown     const PetscReal
352e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3536cf0794eSJed Brown                  {1,0,0,0,0},
3546cf0794eSJed Brown                  {4./9,2./9,0,0,0},
3556cf0794eSJed Brown                  {1./4,0,3./4,0,0},
3566cf0794eSJed Brown                  {1./4,0,3./5,0,0}},
357e74514c0SSatish Balay       At[5][5] = {{0,0,0,0,0},
3586cf0794eSJed Brown                   {.5,.5,0,0,0},
3596cf0794eSJed Brown                   {5./18,-1./9,.5,0,0},
3606cf0794eSJed Brown                   {.5,0,0,.5,0},
361108c343cSJed Brown                   {.25,0,.75,-.5,.5}},
362108c343cSJed Brown       *bembedt = PETSC_NULL;
363108c343cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,0,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
3646cf0794eSJed Brown   }
3656cf0794eSJed Brown   {
3666cf0794eSJed Brown     const PetscReal
367a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
368a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
3694040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
3704040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
3714040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
3724040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
373a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
374a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
3754040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
3764040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
3774040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
3784040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
379cc46b9d1SJed Brown       bembedt[6] = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
3804040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
381cd652676SJed Brown                         {0,0,0},
3824040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
3834040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
3844040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
3854040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
386108c343cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
387a3a57f36SJed Brown   }
388a3a57f36SJed Brown   {
389a3a57f36SJed Brown     const PetscReal
390a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
391a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
3924040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
3934040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
3944040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
3954040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
3964040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
3974040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
398a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
3994040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
4004040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
4014040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
4024040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
4034040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
4044040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
4054040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
406cc46b9d1SJed Brown       bembedt[8] = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
4074040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614. ,  43486358583215./12773830924787. , -9257016797708./5021505065439.},
408cd652676SJed Brown                         {0                               ,  0                              , 0                            },
409cd652676SJed Brown                         {0                               ,  0                              , 0                            },
4104040e9f2SJed Brown                         {65168852399939./7868540260826.  ,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
4114040e9f2SJed Brown                         {15494834004392./5936557850923.  ,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
4124040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473. , 30029262896817./10175596800299.},
4134040e9f2SJed Brown                         {-19024464361622./5461577185407. ,  115839755401235./10719374521269., -26136350496073./3983972220547.},
4144040e9f2SJed Brown                         {-6511271360970./6095937251113.  ,  5843115559534./2180450260947.   , -5289405421727./3760307252460. }};
415108c343cSJed Brown     ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,bembedt,bembedt,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr);
416a3a57f36SJed Brown   }
417a3a57f36SJed Brown 
4188a381b04SJed Brown   PetscFunctionReturn(0);
4198a381b04SJed Brown }
4208a381b04SJed Brown 
4218a381b04SJed Brown #undef __FUNCT__
4228a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy"
4238a381b04SJed Brown /*@C
4248a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
4258a381b04SJed Brown 
4268a381b04SJed Brown    Not Collective
4278a381b04SJed Brown 
4288a381b04SJed Brown    Level: advanced
4298a381b04SJed Brown 
4308a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
4318a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic()
4328a381b04SJed Brown @*/
4338a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
4348a381b04SJed Brown {
4358a381b04SJed Brown   PetscErrorCode ierr;
4368a381b04SJed Brown   ARKTableauLink link;
4378a381b04SJed Brown 
4388a381b04SJed Brown   PetscFunctionBegin;
4398a381b04SJed Brown   while ((link = ARKTableauList)) {
4408a381b04SJed Brown     ARKTableau t = &link->tab;
4418a381b04SJed Brown     ARKTableauList = link->next;
4428a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
443108c343cSJed Brown     ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr);
444cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
4458a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
4468a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
4478a381b04SJed Brown   }
4488a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
4498a381b04SJed Brown   PetscFunctionReturn(0);
4508a381b04SJed Brown }
4518a381b04SJed Brown 
4528a381b04SJed Brown #undef __FUNCT__
4538a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage"
4548a381b04SJed Brown /*@C
4558a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
4568a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
4578a381b04SJed Brown   when using static libraries.
4588a381b04SJed Brown 
4598a381b04SJed Brown   Input Parameter:
4608a381b04SJed Brown   path - The dynamic library path, or PETSC_NULL
4618a381b04SJed Brown 
4628a381b04SJed Brown   Level: developer
4638a381b04SJed Brown 
4648a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
4658a381b04SJed Brown .seealso: PetscInitialize()
4668a381b04SJed Brown @*/
4678a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[])
4688a381b04SJed Brown {
4698a381b04SJed Brown   PetscErrorCode ierr;
4708a381b04SJed Brown 
4718a381b04SJed Brown   PetscFunctionBegin;
4728a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
4738a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
4748a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
475e817cc15SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
4768a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
4778a381b04SJed Brown   PetscFunctionReturn(0);
4788a381b04SJed Brown }
4798a381b04SJed Brown 
4808a381b04SJed Brown #undef __FUNCT__
4818a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage"
4828a381b04SJed Brown /*@C
4838a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
4848a381b04SJed Brown   called from PetscFinalize().
4858a381b04SJed Brown 
4868a381b04SJed Brown   Level: developer
4878a381b04SJed Brown 
4888a381b04SJed Brown .keywords: Petsc, destroy, package
4898a381b04SJed Brown .seealso: PetscFinalize()
4908a381b04SJed Brown @*/
4918a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
4928a381b04SJed Brown {
4938a381b04SJed Brown   PetscErrorCode ierr;
4948a381b04SJed Brown 
4958a381b04SJed Brown   PetscFunctionBegin;
4968a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
4978a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
4988a381b04SJed Brown   PetscFunctionReturn(0);
4998a381b04SJed Brown }
5008a381b04SJed Brown 
5018a381b04SJed Brown #undef __FUNCT__
5028a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister"
503cd652676SJed Brown /*@C
504cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
505cd652676SJed Brown 
506cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
507cd652676SJed Brown 
508cd652676SJed Brown    Input Parameters:
509cd652676SJed Brown +  name - identifier for method
510cd652676SJed Brown .  order - approximation order of method
511cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
512cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
513cd652676SJed Brown .  bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At)
514cd652676SJed Brown .  ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At)
515cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
516cd652676SJed Brown .  b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At)
517cd652676SJed Brown .  c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A)
518108c343cSJed Brown .  bembedt - Stiff part of completion table for embedded method (dimension s; PETSC_NULL if not available)
519108c343cSJed Brown .  bembed - Non-stiff part of completion table for embedded method (dimension s; PETSC_NULL to use bembedt if provided)
520cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
521cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
522cd652676SJed Brown -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt)
523cd652676SJed Brown 
524cd652676SJed Brown    Notes:
525cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
526cd652676SJed Brown 
527cd652676SJed Brown    Level: advanced
528cd652676SJed Brown 
529cd652676SJed Brown .keywords: TS, register
530cd652676SJed Brown 
531cd652676SJed Brown .seealso: TSARKIMEX
532cd652676SJed Brown @*/
53319fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
5348a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
535cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
536108c343cSJed Brown                                  const PetscReal bembedt[],const PetscReal bembed[],
537cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
5388a381b04SJed Brown {
5398a381b04SJed Brown   PetscErrorCode ierr;
5408a381b04SJed Brown   ARKTableauLink link;
5418a381b04SJed Brown   ARKTableau     t;
5428a381b04SJed Brown   PetscInt       i,j;
5438a381b04SJed Brown 
5448a381b04SJed Brown   PetscFunctionBegin;
5458a381b04SJed Brown   ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
546cd652676SJed Brown   ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
5478a381b04SJed Brown   t = &link->tab;
5488a381b04SJed Brown   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5498a381b04SJed Brown   t->order = order;
5508a381b04SJed Brown   t->s = s;
5518a381b04SJed 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);
5528a381b04SJed Brown   ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
5538a381b04SJed Brown   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5548a381b04SJed Brown   if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);}
5558a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
5568a381b04SJed Brown   if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);}
5578a381b04SJed Brown   else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i];
5588a381b04SJed Brown   if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);}
5598a381b04SJed 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];
5608a381b04SJed Brown   if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);}
5618a381b04SJed 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];
562e817cc15SEmil Constantinescu   t->stiffly_accurate = PETSC_TRUE;
563e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
564e817cc15SEmil Constantinescu   t->explicit_first_stage = PETSC_TRUE;
565e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
566e817cc15SEmil Constantinescu   /*def of FSAL can be made more precise*/
5674e9d4bf5SJed Brown   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
568108c343cSJed Brown   if (bembedt) {
569108c343cSJed Brown     ierr = PetscMalloc2(s,PetscReal,&t->bembedt,s,PetscReal,&t->bembed);CHKERRQ(ierr);
570108c343cSJed Brown     ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr);
571108c343cSJed Brown     ierr = PetscMemcpy(t->bembed,bembed?bembed:bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr);
572108c343cSJed Brown   }
573108c343cSJed Brown 
5744f385281SJed Brown   t->pinterp = pinterp;
575cd652676SJed Brown   ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr);
576cd652676SJed Brown   ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
577cd652676SJed Brown   ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
5788a381b04SJed Brown   link->next = ARKTableauList;
5798a381b04SJed Brown   ARKTableauList = link;
5808a381b04SJed Brown   PetscFunctionReturn(0);
5818a381b04SJed Brown }
5828a381b04SJed Brown 
5838a381b04SJed Brown #undef __FUNCT__
584108c343cSJed Brown #define __FUNCT__ "TSEvaluateStep_ARKIMEX"
585108c343cSJed Brown /*
586108c343cSJed Brown  The step completion formula is
587108c343cSJed Brown 
588108c343cSJed Brown  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
589108c343cSJed Brown 
590108c343cSJed Brown  This function can be called before or after ts->vec_sol has been updated.
591108c343cSJed Brown  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
592108c343cSJed Brown  We can write
593108c343cSJed Brown 
594108c343cSJed Brown  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
595108c343cSJed Brown      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
596108c343cSJed Brown      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
597108c343cSJed Brown 
598108c343cSJed Brown  so we can evaluate the method with different order even after the step has been optimistically completed.
599108c343cSJed Brown */
600108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
601108c343cSJed Brown {
602108c343cSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
603108c343cSJed Brown   ARKTableau     tab  = ark->tableau;
604108c343cSJed Brown   PetscScalar    *w = ark->work;
605108c343cSJed Brown   PetscReal      h;
606108c343cSJed Brown   PetscInt       s = tab->s,j;
607108c343cSJed Brown   PetscErrorCode ierr;
608108c343cSJed Brown 
609108c343cSJed Brown   PetscFunctionBegin;
610108c343cSJed Brown   switch (ark->status) {
611108c343cSJed Brown   case TS_STEP_INCOMPLETE:
612108c343cSJed Brown   case TS_STEP_PENDING:
613108c343cSJed Brown     h = ts->time_step; break;
614108c343cSJed Brown   case TS_STEP_COMPLETE:
615108c343cSJed Brown     h = ts->time_step_prev; break;
616b9ce6d65SJed Brown   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
617108c343cSJed Brown   }
618108c343cSJed Brown   if (order == tab->order) {
619e817cc15SEmil Constantinescu     if (ark->status == TS_STEP_INCOMPLETE) {
620e817cc15SEmil Constantinescu       if(!ark->imex && tab->FSAL_implicit) {/* Only the stiffly accurate implicit formula is used */
621e817cc15SEmil Constantinescu         ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr);
622e817cc15SEmil Constantinescu       } else { /* Use the standard completion formula (bt,b) */
623108c343cSJed Brown         ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
624e817cc15SEmil Constantinescu         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
625108c343cSJed Brown         ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
626e817cc15SEmil Constantinescu         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
627108c343cSJed Brown           for (j=0; j<s; j++) w[j] = h*tab->b[j];
628108c343cSJed Brown           ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
629e817cc15SEmil Constantinescu         }
630e817cc15SEmil Constantinescu       }
631108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
632108c343cSJed Brown     if (done) *done = PETSC_TRUE;
633108c343cSJed Brown     PetscFunctionReturn(0);
634108c343cSJed Brown   } else if (order == tab->order-1) {
635108c343cSJed Brown     if (!tab->bembedt) goto unavailable;
636108c343cSJed Brown     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
637108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
638e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
639108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
640108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
641108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
642108c343cSJed Brown     } else {                    /* Rollback and re-complete using (bet-be,be-b) */
643108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
644e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
645108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr);
646108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
647108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
648108c343cSJed Brown     }
649108c343cSJed Brown     if (done) *done = PETSC_TRUE;
650108c343cSJed Brown     PetscFunctionReturn(0);
651108c343cSJed Brown   }
652108c343cSJed Brown   unavailable:
653108c343cSJed Brown   if (done) *done = PETSC_FALSE;
654108c343cSJed Brown   else SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
655108c343cSJed Brown   PetscFunctionReturn(0);
656108c343cSJed Brown }
657108c343cSJed Brown 
658108c343cSJed Brown #undef __FUNCT__
6598a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX"
6608a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
6618a381b04SJed Brown {
6628a381b04SJed Brown   TS_ARKIMEX          *ark = (TS_ARKIMEX*)ts->data;
6638a381b04SJed Brown   ARKTableau          tab  = ark->tableau;
6648a381b04SJed Brown   const PetscInt      s    = tab->s;
6658a381b04SJed Brown   const PetscReal     *At  = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c;
666406d0ec2SJed Brown   PetscScalar         *w   = ark->work;
667e817cc15SEmil Constantinescu   Vec                 *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,W = ark->Work,Z = ark->Z;
668108c343cSJed Brown   TSAdapt             adapt;
6698a381b04SJed Brown   SNES                snes;
670108c343cSJed Brown   PetscInt            i,j,its,lits,reject,next_scheme;
671cdbf8f93SLisandro Dalcin   PetscReal           next_time_step;
672108c343cSJed Brown   PetscReal           t;
673108c343cSJed Brown   PetscBool           accept;
6748a381b04SJed Brown   PetscErrorCode      ierr;
6758a381b04SJed Brown 
6768a381b04SJed Brown   PetscFunctionBegin;
677e817cc15SEmil Constantinescu 
678e817cc15SEmil Constantinescu   /*ark->init_slope=PETSC_FALSE;*/
679e817cc15SEmil Constantinescu 
680e817cc15SEmil Constantinescu 
681e817cc15SEmil Constantinescu   if (ts->equation_type>=TS_EQ_IMPLICIT && tab->explicit_first_stage) {
682e817cc15SEmil Constantinescu     PetscReal valid_time;
683e817cc15SEmil Constantinescu     PetscBool isvalid;
684e817cc15SEmil Constantinescu     ierr = PetscObjectComposedDataGetReal((PetscObject)ts->vec_sol,
685e817cc15SEmil Constantinescu                                           explicit_stage_time_id,
686e817cc15SEmil Constantinescu                                           valid_time,
687e817cc15SEmil Constantinescu                                           isvalid);
688e817cc15SEmil Constantinescu     CHKERRQ(ierr);
689e817cc15SEmil Constantinescu     if (!isvalid || valid_time != ts->ptime) {
690e817cc15SEmil Constantinescu       TS             ts_start;
691e817cc15SEmil Constantinescu       SNES           snes_start;
692e817cc15SEmil Constantinescu       TSRHSFunction  rhsfunction;
693e817cc15SEmil Constantinescu       TSIFunction    ifunction;
694e817cc15SEmil Constantinescu       TSIJacobian    ijacobian;
695e817cc15SEmil Constantinescu       void           *ctxrhs,*ctxif,*ctxij;
696e817cc15SEmil Constantinescu       DM             dm,dm_start;
69719436ca2SJed Brown 
69819436ca2SJed Brown       ierr = TSCreate(PETSC_COMM_WORLD,&ts_start);CHKERRQ(ierr);
69919436ca2SJed Brown       ierr = TSGetSNES(ts,&snes_start);CHKERRQ(ierr);
70019436ca2SJed Brown       ierr = TSSetSNES(ts_start,snes_start);CHKERRQ(ierr);
701e817cc15SEmil Constantinescu       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
702e817cc15SEmil Constantinescu       ierr = TSGetDM(ts_start,&dm_start);CHKERRQ(ierr);
703e817cc15SEmil Constantinescu 
704e817cc15SEmil Constantinescu       ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctxrhs);CHKERRQ(ierr);
705e817cc15SEmil Constantinescu       ierr = DMTSSetRHSFunction(dm_start,rhsfunction,ctxrhs);CHKERRQ(ierr);
706e817cc15SEmil Constantinescu 
707e817cc15SEmil Constantinescu       ierr = DMTSGetIFunction(dm,&ifunction,&ctxif);CHKERRQ(ierr);
708e817cc15SEmil Constantinescu       ierr = DMTSSetIFunction(dm_start,ifunction,ctxif);CHKERRQ(ierr);
709e817cc15SEmil Constantinescu 
710e817cc15SEmil Constantinescu       ierr = DMTSGetIJacobian(dm,&ijacobian,&ctxij);CHKERRQ(ierr);
711e817cc15SEmil Constantinescu       ierr = DMTSSetIJacobian(dm_start,ijacobian,ctxij);CHKERRQ(ierr);
712e817cc15SEmil Constantinescu       ts_start->adapt=ts->adapt;
713e817cc15SEmil Constantinescu       ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
714e817cc15SEmil Constantinescu       ierr = TSSetTime(ts_start,ts->ptime); CHKERRQ(ierr);
715e817cc15SEmil Constantinescu       ierr = TSSetDuration(ts_start,1,ts->time_step);CHKERRQ(ierr);
716e817cc15SEmil Constantinescu       ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
717e817cc15SEmil Constantinescu       ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
718e817cc15SEmil Constantinescu       ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
719e817cc15SEmil Constantinescu       ierr = TSSetEquationType(ts_start,ts->equation_type);CHKERRQ(ierr);
720e817cc15SEmil Constantinescu       ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
721e817cc15SEmil Constantinescu       PetscReal h=-ts->ptime;
722e817cc15SEmil Constantinescu       ierr = TSGetTime(ts_start,&ts->ptime); CHKERRQ(ierr);
723e817cc15SEmil Constantinescu       ts->time_step = h + ts->ptime;
724e817cc15SEmil Constantinescu       ts->steps = 1;
725e817cc15SEmil Constantinescu       ierr = VecCopy(((TS_ARKIMEX *)ts_start->data)->Ydot0,Ydot0);CHKERRQ(ierr);
726e817cc15SEmil Constantinescu       /*ierr = TSDestroy(&ts_start);CHKERRQ(ierr); will this destroy snes as well?*/
727e817cc15SEmil Constantinescu     }
728e817cc15SEmil Constantinescu   }
729e817cc15SEmil Constantinescu 
7308a381b04SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
731cdbf8f93SLisandro Dalcin   next_time_step = ts->time_step;
7328a381b04SJed Brown   t = ts->ptime;
733108c343cSJed Brown   accept = PETSC_TRUE;
734108c343cSJed Brown   ark->status = TS_STEP_INCOMPLETE;
7358a381b04SJed Brown 
736e817cc15SEmil Constantinescu 
73797335746SJed Brown   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
738108c343cSJed Brown     PetscReal h = ts->time_step;
739b8123daeSJed Brown     ierr = TSPreStep(ts);CHKERRQ(ierr);
7408a381b04SJed Brown     for (i=0; i<s; i++) {
7418a381b04SJed Brown       if (At[i*s+i] == 0) {           /* This stage is explicit */
742e817cc15SEmil Constantinescu         /*if(ts->equation_type>=TS_EQ_IMPLICIT){
743e817cc15SEmil Constantinescu           if(i!=0){
744e817cc15SEmil Constantinescu             printf("Throw an error: we cannot have explicit stages for DAEs other than the first stage when used in FSAL\n");
745e817cc15SEmil Constantinescu           }
746e817cc15SEmil Constantinescu           if(ts->steps==0){ //initialize the slope - needs to be moved outside
747e817cc15SEmil Constantinescu 
748e817cc15SEmil Constantinescu             ark->init_slope=PETSC_TRUE;
749e817cc15SEmil Constantinescu             ierr = VecZeroEntries(Ydot0);CHKERRQ(ierr);
750e817cc15SEmil Constantinescu             ierr = SNESSolve(snes,PETSC_NULL,Ydot0);CHKERRQ(ierr);
751e817cc15SEmil Constantinescu 
752e817cc15SEmil Constantinescu             ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
753e817cc15SEmil Constantinescu             ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
754e817cc15SEmil Constantinescu             ts->snes_its += its; ts->ksp_its += lits;
755e817cc15SEmil Constantinescu             ark->init_slope=PETSC_FALSE;
756e817cc15SEmil Constantinescu           }
757e817cc15SEmil Constantinescu         }*/
7588a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
759e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7608a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
7618a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7628a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
7638a381b04SJed Brown       } else {
7648a381b04SJed Brown         ark->stage_time = t + h*ct[i];
765b296d7d5SJed Brown         ark->scoeff = 1./At[i*s+i];
766b8123daeSJed Brown         ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
7678a381b04SJed Brown         /* Affine part */
7688a381b04SJed Brown         ierr = VecZeroEntries(W);CHKERRQ(ierr);
7698a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
7708a381b04SJed Brown         ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr);
771b296d7d5SJed Brown         ierr = VecScale(W, ark->scoeff/h);CHKERRQ(ierr);
772f16577ceSEmil Constantinescu 
7738a381b04SJed Brown         /* Ydot = shift*(Y-Z) */
7748a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
775e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7764f385281SJed Brown         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
777f16577ceSEmil Constantinescu 
7788a381b04SJed Brown         /* Initial guess taken from last stage */
7798a381b04SJed Brown         ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr);
7808a381b04SJed Brown         ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr);
781e817cc15SEmil Constantinescu         ierr = (ts->ops->snesfunction)(snes,Y[i],W,ts);CHKERRQ(ierr);
7828a381b04SJed Brown         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
7838a381b04SJed Brown         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
7845ef26d82SJed Brown         ts->snes_its += its; ts->ksp_its += lits;
785ad6bc421SBarry Smith         ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
78697335746SJed Brown         ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
78797335746SJed Brown         if (!accept) goto reject_step;
7888a381b04SJed Brown       }
789e817cc15SEmil Constantinescu       if(ts->equation_type>=TS_EQ_IMPLICIT){
790e817cc15SEmil Constantinescu         if(i==0 && tab->explicit_first_stage){
791e817cc15SEmil Constantinescu           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);
792e817cc15SEmil Constantinescu         } else {
793e817cc15SEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
794e817cc15SEmil Constantinescu         }
795e817cc15SEmil Constantinescu       }else{
7968a381b04SJed Brown         ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
7974cc180ffSJed Brown         ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);
798e817cc15SEmil Constantinescu         ierr = VecScale(YdotI[i], -1.0);CHKERRQ(ierr);
7994cc180ffSJed Brown         if (ark->imex) {
8008a381b04SJed Brown           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
8014cc180ffSJed Brown         } else {
8024cc180ffSJed Brown           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
8034cc180ffSJed Brown         }
8048a381b04SJed Brown       }
805e817cc15SEmil Constantinescu     }
806108c343cSJed Brown     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
807108c343cSJed Brown     ark->status = TS_STEP_PENDING;
8088a381b04SJed Brown 
809108c343cSJed Brown     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
810ad6bc421SBarry Smith     ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
811108c343cSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
812108c343cSJed Brown     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
813108c343cSJed Brown     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
814108c343cSJed Brown     if (accept) {
815108c343cSJed Brown       /* ignore next_scheme for now */
8168a381b04SJed Brown       ts->ptime += ts->time_step;
817cdbf8f93SLisandro Dalcin       ts->time_step = next_time_step;
8188a381b04SJed Brown       ts->steps++;
819e817cc15SEmil Constantinescu       if(ts->equation_type>=TS_EQ_IMPLICIT){/* save the initial slope for the next step*/
820e817cc15SEmil Constantinescu         ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
821e817cc15SEmil Constantinescu       }
822108c343cSJed Brown       ark->status = TS_STEP_COMPLETE;
823e817cc15SEmil Constantinescu       if (tab->explicit_first_stage) {
824e817cc15SEmil Constantinescu         ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
825e817cc15SEmil Constantinescu       }
826e817cc15SEmil Constantinescu 
827108c343cSJed Brown       break;
828108c343cSJed Brown     } else {                    /* Roll back the current step */
829*2c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*bt[j];
830108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotI);CHKERRQ(ierr);
831*2c0c504eSEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
832108c343cSJed Brown       ierr = VecMAXPY(ts->vec_sol,s,w,ark->YdotRHS);CHKERRQ(ierr);
833108c343cSJed Brown       ts->time_step = next_time_step;
834108c343cSJed Brown       ark->status = TS_STEP_INCOMPLETE;
835108c343cSJed Brown     }
836476b6736SJed Brown     reject_step: continue;
837108c343cSJed Brown   }
838b2ce242eSJed Brown   if (ark->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
8398a381b04SJed Brown   PetscFunctionReturn(0);
8408a381b04SJed Brown }
8418a381b04SJed Brown 
842cd652676SJed Brown #undef __FUNCT__
843cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX"
844cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
845cd652676SJed Brown {
846cd652676SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8474f385281SJed Brown   PetscInt        s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
848108c343cSJed Brown   PetscReal       h;
849108c343cSJed Brown   PetscReal       tt,t;
850cd652676SJed Brown   PetscScalar     *bt,*b;
851cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
852cd652676SJed Brown   PetscErrorCode  ierr;
853cd652676SJed Brown 
854cd652676SJed Brown   PetscFunctionBegin;
855cd652676SJed Brown   if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
856108c343cSJed Brown   switch (ark->status) {
857108c343cSJed Brown   case TS_STEP_INCOMPLETE:
858108c343cSJed Brown   case TS_STEP_PENDING:
859108c343cSJed Brown     h = ts->time_step;
860108c343cSJed Brown     t = (itime - ts->ptime)/h;
861108c343cSJed Brown     break;
862108c343cSJed Brown   case TS_STEP_COMPLETE:
863108c343cSJed Brown     h = ts->time_step_prev;
864108c343cSJed Brown     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
865108c343cSJed Brown     break;
866b9ce6d65SJed Brown   default: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Invalid TSStepStatus");
867108c343cSJed Brown   }
868cd652676SJed Brown   ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr);
869cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
8704f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
871cd652676SJed Brown     for (i=0; i<s; i++) {
872108c343cSJed Brown       bt[i] += h * Bt[i*pinterp+j] * tt * -1.0;
873108c343cSJed Brown       b[i]  += h * B[i*pinterp+j] * tt;
874cd652676SJed Brown     }
875cd652676SJed Brown   }
876cd652676SJed Brown   if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved");
877cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
878cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
879cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
880cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
881cd652676SJed Brown   PetscFunctionReturn(0);
882cd652676SJed Brown }
883cd652676SJed Brown 
8848a381b04SJed Brown /*------------------------------------------------------------*/
8858a381b04SJed Brown #undef __FUNCT__
8868a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX"
8878a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
8888a381b04SJed Brown {
8898a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8908a381b04SJed Brown   PetscInt        s;
8918a381b04SJed Brown   PetscErrorCode  ierr;
8928a381b04SJed Brown 
8938a381b04SJed Brown   PetscFunctionBegin;
8948a381b04SJed Brown   if (!ark->tableau) PetscFunctionReturn(0);
8958a381b04SJed Brown   s = ark->tableau->s;
8968a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr);
8978a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr);
8988a381b04SJed Brown   ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr);
8998a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
9008a381b04SJed Brown   ierr = VecDestroy(&ark->Work);CHKERRQ(ierr);
901e817cc15SEmil Constantinescu   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
9028a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
9038a381b04SJed Brown   ierr = PetscFree(ark->work);CHKERRQ(ierr);
9048a381b04SJed Brown   PetscFunctionReturn(0);
9058a381b04SJed Brown }
9068a381b04SJed Brown 
9078a381b04SJed Brown #undef __FUNCT__
9088a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX"
9098a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
9108a381b04SJed Brown {
9118a381b04SJed Brown   PetscErrorCode  ierr;
9128a381b04SJed Brown 
9138a381b04SJed Brown   PetscFunctionBegin;
9148a381b04SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
9158a381b04SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
9168a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr);
9178a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr);
918995b3938SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","",PETSC_NULL);CHKERRQ(ierr);
9198a381b04SJed Brown   PetscFunctionReturn(0);
9208a381b04SJed Brown }
9218a381b04SJed Brown 
922d5e6173cSPeter Brune 
923d5e6173cSPeter Brune #undef __FUNCT__
924d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXGetVecs"
925d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
926d5e6173cSPeter Brune {
927d5e6173cSPeter Brune   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
928d5e6173cSPeter Brune   PetscErrorCode ierr;
929d5e6173cSPeter Brune 
930d5e6173cSPeter Brune   PetscFunctionBegin;
931d5e6173cSPeter Brune   if (Z) {
932d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
933d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
934d5e6173cSPeter Brune     } else *Z = ax->Z;
935d5e6173cSPeter Brune   }
936d5e6173cSPeter Brune   if (Ydot) {
937d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
938d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
939d5e6173cSPeter Brune     } else *Ydot = ax->Ydot;
940d5e6173cSPeter Brune   }
941d5e6173cSPeter Brune   PetscFunctionReturn(0);
942d5e6173cSPeter Brune }
943d5e6173cSPeter Brune 
944d5e6173cSPeter Brune 
945d5e6173cSPeter Brune #undef __FUNCT__
946d5e6173cSPeter Brune #define __FUNCT__ "TSARKIMEXRestoreVecs"
947d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
948d5e6173cSPeter Brune {
949d5e6173cSPeter Brune   PetscErrorCode ierr;
950d5e6173cSPeter Brune 
951d5e6173cSPeter Brune   PetscFunctionBegin;
952d5e6173cSPeter Brune   if (Z) {
953d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
954d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
955d5e6173cSPeter Brune     }
956d5e6173cSPeter Brune   }
957d5e6173cSPeter Brune   if (Ydot) {
958d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
959d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
960d5e6173cSPeter Brune     }
961d5e6173cSPeter Brune   }
962d5e6173cSPeter Brune   PetscFunctionReturn(0);
963d5e6173cSPeter Brune }
964d5e6173cSPeter Brune 
9658a381b04SJed Brown /*
9668a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
9678a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
9688a381b04SJed Brown */
9698a381b04SJed Brown #undef __FUNCT__
9708a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX"
9718a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
9728a381b04SJed Brown {
9738a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
974d5e6173cSPeter Brune   DM             dm,dmsave;
975d5e6173cSPeter Brune   Vec            Z,Ydot;
976b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
9778a381b04SJed Brown   PetscErrorCode ierr;
9788a381b04SJed Brown 
9798a381b04SJed Brown   PetscFunctionBegin;
980d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
981d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
982b296d7d5SJed Brown   ierr = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
983d5e6173cSPeter Brune   dmsave = ts->dm;
984d5e6173cSPeter Brune   ts->dm = dm;
985e817cc15SEmil Constantinescu   /*  if(!ark->init_slope){*/
986d5e6173cSPeter Brune   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
987e817cc15SEmil Constantinescu     /*  }else{
988e817cc15SEmil Constantinescu     ierr = TSComputeIFunction(ts,ark->stage_time,ts->vec_sol,X,F,ark->imex);CHKERRQ(ierr);
989e817cc15SEmil Constantinescu      }*/
990e817cc15SEmil Constantinescu 
991d5e6173cSPeter Brune   ts->dm = dmsave;
992d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
9938a381b04SJed Brown   PetscFunctionReturn(0);
9948a381b04SJed Brown }
9958a381b04SJed Brown 
9968a381b04SJed Brown #undef __FUNCT__
9978a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX"
9988a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts)
9998a381b04SJed Brown {
10008a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1001d5e6173cSPeter Brune   DM             dm,dmsave;
1002d5e6173cSPeter Brune   Vec            Ydot;
1003b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10048a381b04SJed Brown   PetscErrorCode ierr;
10058a381b04SJed Brown 
10068a381b04SJed Brown   PetscFunctionBegin;
1007d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1008d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,PETSC_NULL,&Ydot);CHKERRQ(ierr);
10098a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1010d5e6173cSPeter Brune   dmsave = ts->dm;
1011d5e6173cSPeter Brune   ts->dm = dm;
1012e817cc15SEmil Constantinescu   /*if(!ark->init_slope){*/
1013b296d7d5SJed Brown   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,str,ark->imex);CHKERRQ(ierr);
1014e817cc15SEmil Constantinescu     /*  }else{
1015e817cc15SEmil Constantinescu     ierr = VecZeroEntries(ark->Work);CHKERRQ(ierr);
1016e817cc15SEmil Constantinescu     ierr = TSComputeIJacobian(ts,ark->stage_time,ark->Work,X,1.0,A,B,str,ark->imex);CHKERRQ(ierr);
1017e817cc15SEmil Constantinescu      }*/
1018d5e6173cSPeter Brune   ts->dm = dmsave;
1019d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,PETSC_NULL,&Ydot);CHKERRQ(ierr);
1020d5e6173cSPeter Brune   PetscFunctionReturn(0);
1021d5e6173cSPeter Brune }
1022d5e6173cSPeter Brune 
1023d5e6173cSPeter Brune #undef __FUNCT__
1024d5e6173cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSARKIMEX"
1025d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1026d5e6173cSPeter Brune {
1027d5e6173cSPeter Brune 
1028d5e6173cSPeter Brune   PetscFunctionBegin;
1029d5e6173cSPeter Brune   PetscFunctionReturn(0);
1030d5e6173cSPeter Brune }
1031d5e6173cSPeter Brune 
1032d5e6173cSPeter Brune #undef __FUNCT__
1033d5e6173cSPeter Brune #define __FUNCT__ "DMRestrictHook_TSARKIMEX"
1034d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1035d5e6173cSPeter Brune {
1036d5e6173cSPeter Brune   TS             ts = (TS)ctx;
1037d5e6173cSPeter Brune   PetscErrorCode ierr;
1038d5e6173cSPeter Brune   Vec            Z,Z_c;
1039d5e6173cSPeter Brune 
1040d5e6173cSPeter Brune   PetscFunctionBegin;
1041d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,fine,&Z,PETSC_NULL);CHKERRQ(ierr);
1042d5e6173cSPeter Brune   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1043d5e6173cSPeter Brune   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1044d5e6173cSPeter Brune   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
1045d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,PETSC_NULL);CHKERRQ(ierr);
1046d5e6173cSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,PETSC_NULL);CHKERRQ(ierr);
10478a381b04SJed Brown   PetscFunctionReturn(0);
10488a381b04SJed Brown }
10498a381b04SJed Brown 
1050cdb298fcSPeter Brune 
1051cdb298fcSPeter Brune #undef __FUNCT__
1052cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainHook_TSARKIMEX"
1053cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1054cdb298fcSPeter Brune {
1055cdb298fcSPeter Brune 
1056cdb298fcSPeter Brune   PetscFunctionBegin;
1057cdb298fcSPeter Brune   PetscFunctionReturn(0);
1058cdb298fcSPeter Brune }
1059cdb298fcSPeter Brune 
1060cdb298fcSPeter Brune #undef __FUNCT__
1061cdb298fcSPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_TSARKIMEX"
1062cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1063cdb298fcSPeter Brune {
1064cdb298fcSPeter Brune   TS             ts = (TS)ctx;
1065cdb298fcSPeter Brune   PetscErrorCode ierr;
1066cdb298fcSPeter Brune   Vec            Z,Z_c;
1067cdb298fcSPeter Brune 
1068cdb298fcSPeter Brune   PetscFunctionBegin;
1069cdb298fcSPeter Brune   ierr = TSARKIMEXGetVecs(ts,dm,&Z,PETSC_NULL);CHKERRQ(ierr);
1070cdb298fcSPeter Brune   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1071cdb298fcSPeter Brune 
1072cdb298fcSPeter Brune   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1073cdb298fcSPeter Brune   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1074cdb298fcSPeter Brune 
1075cdb298fcSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,PETSC_NULL);CHKERRQ(ierr);
1076cdb298fcSPeter Brune   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,PETSC_NULL);CHKERRQ(ierr);
1077cdb298fcSPeter Brune   PetscFunctionReturn(0);
1078cdb298fcSPeter Brune }
1079cdb298fcSPeter Brune 
10808a381b04SJed Brown #undef __FUNCT__
10818a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX"
10828a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
10838a381b04SJed Brown {
10848a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1085f2c2a1b9SBarry Smith   ARKTableau     tab;
1086f2c2a1b9SBarry Smith   PetscInt       s;
10878a381b04SJed Brown   PetscErrorCode ierr;
1088d5e6173cSPeter Brune   DM             dm;
1089f9c1d6abSBarry Smith 
10908a381b04SJed Brown   PetscFunctionBegin;
10918a381b04SJed Brown   if (!ark->tableau) {
1092e24355feSJed Brown     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
10938a381b04SJed Brown   }
1094f2c2a1b9SBarry Smith   tab  = ark->tableau;
1095f2c2a1b9SBarry Smith   s    = tab->s;
10968a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr);
10978a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr);
10988a381b04SJed Brown   ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr);
10998a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
11008a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr);
1101e817cc15SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
11028a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
11038a381b04SJed Brown   ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr);
1104d5e6173cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1105d5e6173cSPeter Brune   if (dm) {
1106d5e6173cSPeter Brune     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1107cdb298fcSPeter Brune     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1108d5e6173cSPeter Brune   }
11098a381b04SJed Brown   PetscFunctionReturn(0);
11108a381b04SJed Brown }
11118a381b04SJed Brown /*------------------------------------------------------------*/
11128a381b04SJed Brown 
11138a381b04SJed Brown #undef __FUNCT__
11148a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX"
11158a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts)
11168a381b04SJed Brown {
11174cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11188a381b04SJed Brown   PetscErrorCode ierr;
11198a381b04SJed Brown   char           arktype[256];
11208a381b04SJed Brown 
11218a381b04SJed Brown   PetscFunctionBegin;
11228a381b04SJed Brown   ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr);
11238a381b04SJed Brown   {
11248a381b04SJed Brown     ARKTableauLink link;
11258a381b04SJed Brown     PetscInt       count,choice;
11268a381b04SJed Brown     PetscBool      flg;
11278a381b04SJed Brown     const char     **namelist;
11288caf3d72SBarry Smith     ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof(arktype));CHKERRQ(ierr);
11298a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
11308a381b04SJed Brown     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
11318a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
11328a381b04SJed Brown     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr);
11338a381b04SJed Brown     ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr);
11348a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
11354cc180ffSJed Brown     flg = (PetscBool)!ark->imex;
11364cc180ffSJed Brown     ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
11374cc180ffSJed Brown     ark->imex = (PetscBool)!flg;
1138d52bd9f3SBarry Smith     ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr);
11398a381b04SJed Brown   }
11408a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
11418a381b04SJed Brown   PetscFunctionReturn(0);
11428a381b04SJed Brown }
11438a381b04SJed Brown 
11448a381b04SJed Brown #undef __FUNCT__
11458a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray"
11468a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
11478a381b04SJed Brown {
1148257d2499SJed Brown   PetscErrorCode ierr;
1149f1d86077SJed Brown   PetscInt       i;
1150f1d86077SJed Brown   size_t         left,count;
11518a381b04SJed Brown   char           *p;
11528a381b04SJed Brown 
11538a381b04SJed Brown   PetscFunctionBegin;
1154f1d86077SJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
1155f1d86077SJed Brown     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
11568a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
11578a381b04SJed Brown     left -= count;
11588a381b04SJed Brown     p += count;
11598a381b04SJed Brown     *p++ = ' ';
11608a381b04SJed Brown   }
11618a381b04SJed Brown   p[i ? 0 : -1] = 0;
11628a381b04SJed Brown   PetscFunctionReturn(0);
11638a381b04SJed Brown }
11648a381b04SJed Brown 
11658a381b04SJed Brown #undef __FUNCT__
11668a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX"
11678a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
11688a381b04SJed Brown {
11698a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11708a381b04SJed Brown   ARKTableau     tab = ark->tableau;
11718a381b04SJed Brown   PetscBool      iascii;
11728a381b04SJed Brown   PetscErrorCode ierr;
1173559eea31SJed Brown   TSAdapt        adapt;
11748a381b04SJed Brown 
11758a381b04SJed Brown   PetscFunctionBegin;
1176251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11778a381b04SJed Brown   if (iascii) {
117819fd82e9SBarry Smith     TSARKIMEXType arktype;
11798a381b04SJed Brown     char buf[512];
11808a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
11818a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
11828caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
118331f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
11848caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1185e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate?"yes":"no");CHKERRQ(ierr);
1186e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage?"yes":"no");CHKERRQ(ierr);
1187e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit?"yes":"no");CHKERRQ(ierr);
118831f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
11898a381b04SJed Brown   }
1190ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&adapt);CHKERRQ(ierr);
1191559eea31SJed Brown   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
1192d52bd9f3SBarry Smith   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
11938a381b04SJed Brown   PetscFunctionReturn(0);
11948a381b04SJed Brown }
11958a381b04SJed Brown 
11968a381b04SJed Brown #undef __FUNCT__
1197f2c2a1b9SBarry Smith #define __FUNCT__ "TSLoad_ARKIMEX"
1198f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1199f2c2a1b9SBarry Smith {
1200f2c2a1b9SBarry Smith   PetscErrorCode ierr;
1201f2c2a1b9SBarry Smith   SNES           snes;
1202ad6bc421SBarry Smith   TSAdapt        tsadapt;
1203f2c2a1b9SBarry Smith 
1204f2c2a1b9SBarry Smith   PetscFunctionBegin;
1205ad6bc421SBarry Smith   ierr = TSGetTSAdapt(ts,&tsadapt);CHKERRQ(ierr);
1206ad6bc421SBarry Smith   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
1207f2c2a1b9SBarry Smith   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1208f2c2a1b9SBarry Smith   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1209ad6bc421SBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
1210ad6bc421SBarry Smith   ierr = SNESSetFunction(snes,PETSC_NULL,PETSC_NULL,ts);CHKERRQ(ierr);
1211ad6bc421SBarry Smith   ierr = SNESSetJacobian(snes,PETSC_NULL,PETSC_NULL,PETSC_NULL,ts);CHKERRQ(ierr);
1212f2c2a1b9SBarry Smith   PetscFunctionReturn(0);
1213f2c2a1b9SBarry Smith }
1214f2c2a1b9SBarry Smith 
1215f2c2a1b9SBarry Smith #undef __FUNCT__
12168a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType"
12178a381b04SJed Brown /*@C
12188a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
12198a381b04SJed Brown 
12208a381b04SJed Brown   Logically collective
12218a381b04SJed Brown 
12228a381b04SJed Brown   Input Parameter:
12238a381b04SJed Brown +  ts - timestepping context
12248a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
12258a381b04SJed Brown 
12268a381b04SJed Brown   Level: intermediate
12278a381b04SJed Brown 
1228020d8f30SJed Brown .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2, TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
12298a381b04SJed Brown @*/
123019fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
12318a381b04SJed Brown {
12328a381b04SJed Brown   PetscErrorCode ierr;
12338a381b04SJed Brown 
12348a381b04SJed Brown   PetscFunctionBegin;
12358a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
123619fd82e9SBarry Smith   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
12378a381b04SJed Brown   PetscFunctionReturn(0);
12388a381b04SJed Brown }
12398a381b04SJed Brown 
12408a381b04SJed Brown #undef __FUNCT__
12418a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType"
12428a381b04SJed Brown /*@C
12438a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
12448a381b04SJed Brown 
12458a381b04SJed Brown   Logically collective
12468a381b04SJed Brown 
12478a381b04SJed Brown   Input Parameter:
12488a381b04SJed Brown .  ts - timestepping context
12498a381b04SJed Brown 
12508a381b04SJed Brown   Output Parameter:
12518a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
12528a381b04SJed Brown 
12538a381b04SJed Brown   Level: intermediate
12548a381b04SJed Brown 
12558a381b04SJed Brown .seealso: TSARKIMEXGetType()
12568a381b04SJed Brown @*/
125719fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
12588a381b04SJed Brown {
12598a381b04SJed Brown   PetscErrorCode ierr;
12608a381b04SJed Brown 
12618a381b04SJed Brown   PetscFunctionBegin;
12628a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
126319fd82e9SBarry Smith   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
12648a381b04SJed Brown   PetscFunctionReturn(0);
12658a381b04SJed Brown }
12668a381b04SJed Brown 
12674cc180ffSJed Brown #undef __FUNCT__
12684cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit"
12694cc180ffSJed Brown /*@C
12704cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
12714cc180ffSJed Brown 
12724cc180ffSJed Brown   Logically collective
12734cc180ffSJed Brown 
12744cc180ffSJed Brown   Input Parameter:
12754cc180ffSJed Brown +  ts - timestepping context
12764cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
12774cc180ffSJed Brown 
12784cc180ffSJed Brown   Level: intermediate
12794cc180ffSJed Brown 
12804cc180ffSJed Brown .seealso: TSARKIMEXGetType()
12814cc180ffSJed Brown @*/
12824cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
12834cc180ffSJed Brown {
12844cc180ffSJed Brown   PetscErrorCode ierr;
12854cc180ffSJed Brown 
12864cc180ffSJed Brown   PetscFunctionBegin;
12874cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
12884cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
12894cc180ffSJed Brown   PetscFunctionReturn(0);
12904cc180ffSJed Brown }
12914cc180ffSJed Brown 
12928a381b04SJed Brown EXTERN_C_BEGIN
12938a381b04SJed Brown #undef __FUNCT__
12948a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX"
129519fd82e9SBarry Smith PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
12968a381b04SJed Brown {
12978a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
12988a381b04SJed Brown   PetscErrorCode ierr;
12998a381b04SJed Brown 
13008a381b04SJed Brown   PetscFunctionBegin;
1301f2c2a1b9SBarry Smith   if (!ark->tableau) {
1302f2c2a1b9SBarry Smith     ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
1303f2c2a1b9SBarry Smith   }
13048a381b04SJed Brown   *arktype = ark->tableau->name;
13058a381b04SJed Brown   PetscFunctionReturn(0);
13068a381b04SJed Brown }
13078a381b04SJed Brown #undef __FUNCT__
13088a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX"
130919fd82e9SBarry Smith PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
13108a381b04SJed Brown {
13118a381b04SJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13128a381b04SJed Brown   PetscErrorCode ierr;
13138a381b04SJed Brown   PetscBool match;
13148a381b04SJed Brown   ARKTableauLink link;
13158a381b04SJed Brown 
13168a381b04SJed Brown   PetscFunctionBegin;
13178a381b04SJed Brown   if (ark->tableau) {
13188a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
13198a381b04SJed Brown     if (match) PetscFunctionReturn(0);
13208a381b04SJed Brown   }
13218a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
13228a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
13238a381b04SJed Brown     if (match) {
13248a381b04SJed Brown       ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
13258a381b04SJed Brown       ark->tableau = &link->tab;
13268a381b04SJed Brown       PetscFunctionReturn(0);
13278a381b04SJed Brown     }
13288a381b04SJed Brown   }
13298a381b04SJed Brown   SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
13308a381b04SJed Brown   PetscFunctionReturn(0);
13318a381b04SJed Brown }
13324cc180ffSJed Brown #undef __FUNCT__
13334cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX"
13344cc180ffSJed Brown PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
13354cc180ffSJed Brown {
13364cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13374cc180ffSJed Brown 
13384cc180ffSJed Brown   PetscFunctionBegin;
13394cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
13404cc180ffSJed Brown   PetscFunctionReturn(0);
13414cc180ffSJed Brown }
13428a381b04SJed Brown EXTERN_C_END
13438a381b04SJed Brown 
13448a381b04SJed Brown /* ------------------------------------------------------------ */
13458a381b04SJed Brown /*MC
1346a4386c9eSJed Brown       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
13478a381b04SJed Brown 
1348fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1349fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1350fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1351fca742c7SJed Brown 
1352fca742c7SJed Brown   Notes:
1353a4386c9eSJed Brown   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1354c8058688SBarry Smith 
1355a4386c9eSJed 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).
1356fca742c7SJed Brown 
13578a381b04SJed Brown   Level: beginner
13588a381b04SJed Brown 
1359c8058688SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3,
1360a4386c9eSJed Brown            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
13618a381b04SJed Brown 
13628a381b04SJed Brown M*/
13638a381b04SJed Brown EXTERN_C_BEGIN
13648a381b04SJed Brown #undef __FUNCT__
13658a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX"
13668a381b04SJed Brown PetscErrorCode  TSCreate_ARKIMEX(TS ts)
13678a381b04SJed Brown {
13688a381b04SJed Brown   TS_ARKIMEX     *th;
13698a381b04SJed Brown   PetscErrorCode ierr;
13708a381b04SJed Brown 
13718a381b04SJed Brown   PetscFunctionBegin;
13728a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
13738a381b04SJed Brown   ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr);
13748a381b04SJed Brown #endif
13758a381b04SJed Brown 
13768a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
13778a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
13788a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
1379f2c2a1b9SBarry Smith   ts->ops->load           = TSLoad_ARKIMEX;
13808a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
13818a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
1382cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1383108c343cSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
13848a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
13858a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
13868a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
13878a381b04SJed Brown 
13888a381b04SJed Brown   ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr);
13898a381b04SJed Brown   ts->data = (void*)th;
13904cc180ffSJed Brown   th->imex = PETSC_TRUE;
13918a381b04SJed Brown 
13928a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
13938a381b04SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
13944cc180ffSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
13958a381b04SJed Brown   PetscFunctionReturn(0);
13968a381b04SJed Brown }
13978a381b04SJed Brown EXTERN_C_END
1398