xref: /petsc/src/ts/impls/arkimex/arkimex.c (revision b3a6b9728fcfbb58cc066441aaceaf50a6f816d2)
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 */
12af0996ceSBarry Smith #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;
1856dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS,PetscReal,Vec);
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 */
4756dcabbaSDebojyoti Ghosh   Vec          *Y_prev;          /* States computed during the previous time step */
4856dcabbaSDebojyoti Ghosh   Vec          *YdotI_prev;      /* Time derivatives for the stiff part for the previous time step*/
4956dcabbaSDebojyoti Ghosh   Vec          *YdotRHS_prev;    /* Function evaluations for the non-stiff part for the previous time step*/
50e817cc15SEmil Constantinescu   Vec          Ydot0;            /* Holds the slope from the previous step in FSAL case */
518a381b04SJed Brown   Vec          Ydot;             /* Work vector holding Ydot during residual evaluation */
528a381b04SJed Brown   Vec          Z;                /* Ydot = shift(Y-Z) */
538a381b04SJed Brown   PetscScalar  *work;            /* Scalar work */
54b296d7d5SJed Brown   PetscReal    scoeff;           /* shift = scoeff/dt */
558a381b04SJed Brown   PetscReal    stage_time;
564cc180ffSJed Brown   PetscBool    imex;
5796400bd6SLisandro Dalcin   PetscBool    extrapolate;      /* Extrapolate initial guess from previous time-step stage values */
58108c343cSJed Brown   TSStepStatus status;
598a381b04SJed Brown } TS_ARKIMEX;
601f80e275SEmil Constantinescu /*MC
611f80e275SEmil Constantinescu      TSARKIMEXARS122 - Second order ARK IMEX scheme.
628a381b04SJed Brown 
631f80e275SEmil Constantinescu      This method has one explicit stage and one implicit stage.
641f80e275SEmil Constantinescu 
659bd3e852SBarry Smith      Options Database:
669bd3e852SBarry Smith .      -ts_arkimex_type ars122
679bd3e852SBarry Smith 
681f80e275SEmil Constantinescu      References:
6996a0c994SBarry Smith .   1. -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).
701f80e275SEmil Constantinescu 
711f80e275SEmil Constantinescu      Level: advanced
721f80e275SEmil Constantinescu 
739bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
741f80e275SEmil Constantinescu M*/
751f80e275SEmil Constantinescu /*MC
761f80e275SEmil Constantinescu      TSARKIMEXA2 - Second order ARK IMEX scheme with A-stable implicit part.
771f80e275SEmil Constantinescu 
781f80e275SEmil Constantinescu      This method has an explicit stage and one implicit stage, and has an A-stable implicit scheme. This method was provided by Emil Constantinescu.
791f80e275SEmil Constantinescu 
809bd3e852SBarry Smith      Options Database:
819bd3e852SBarry Smith .      -ts_arkimex_type a2
829bd3e852SBarry Smith 
831f80e275SEmil Constantinescu      Level: advanced
841f80e275SEmil Constantinescu 
859bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
861f80e275SEmil Constantinescu M*/
871f80e275SEmil Constantinescu /*MC
881f80e275SEmil Constantinescu      TSARKIMEXL2 - Second order ARK IMEX scheme with L-stable implicit part.
891f80e275SEmil Constantinescu 
901f80e275SEmil Constantinescu      This method has two implicit stages, and L-stable implicit scheme.
911f80e275SEmil Constantinescu 
929bd3e852SBarry Smith      Options Database:
939bd3e852SBarry Smith .      -ts_arkimex_type l2
949bd3e852SBarry Smith 
951f80e275SEmil Constantinescu     References:
9696a0c994SBarry Smith .   1. -  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.
971f80e275SEmil Constantinescu 
981f80e275SEmil Constantinescu      Level: advanced
991f80e275SEmil Constantinescu 
1009bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
1011f80e275SEmil Constantinescu M*/
1021f80e275SEmil Constantinescu /*MC
103e817cc15SEmil Constantinescu      TSARKIMEX1BEE - First order Backward Euler represented as an ARK IMEX scheme with extrapolation as error estimator. This is a 3-stage method.
104e817cc15SEmil Constantinescu 
105e817cc15SEmil Constantinescu      This method is aimed at starting the integration of implicit DAEs when explicit first-stage ARK methods are used.
106e817cc15SEmil Constantinescu 
1079bd3e852SBarry Smith      Options Database:
1089bd3e852SBarry Smith .      -ts_arkimex_type 1bee
1099bd3e852SBarry Smith 
110e817cc15SEmil Constantinescu      Level: advanced
111e817cc15SEmil Constantinescu 
1129bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
113e817cc15SEmil Constantinescu M*/
114e817cc15SEmil Constantinescu /*MC
1151f80e275SEmil Constantinescu      TSARKIMEX2C - Second order ARK IMEX scheme with L-stable implicit part.
1161f80e275SEmil Constantinescu 
1171f80e275SEmil 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.
1181f80e275SEmil Constantinescu 
1199bd3e852SBarry Smith      Options Database:
1209bd3e852SBarry Smith .      -ts_arkimex_type 2c
1219bd3e852SBarry Smith 
1221f80e275SEmil Constantinescu      Level: advanced
1231f80e275SEmil Constantinescu 
1249bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
1251f80e275SEmil Constantinescu M*/
12664f491ddSJed Brown /*MC
12764f491ddSJed Brown      TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part.
12864f491ddSJed Brown 
129617a39beSEmil 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.
13064f491ddSJed Brown 
1319bd3e852SBarry Smith      Options Database:
1329bd3e852SBarry Smith .      -ts_arkimex_type 2d
1339bd3e852SBarry Smith 
134b330ce4dSSatish Balay      Level: advanced
135b330ce4dSSatish Balay 
1369bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
13764f491ddSJed Brown M*/
13864f491ddSJed Brown /*MC
13964f491ddSJed Brown      TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part.
14064f491ddSJed Brown 
14164f491ddSJed Brown      This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu.
14264f491ddSJed Brown 
1439bd3e852SBarry Smith      Options Database:
1449bd3e852SBarry Smith .      -ts_arkimex_type 2e
1459bd3e852SBarry Smith 
146b330ce4dSSatish Balay     Level: advanced
147b330ce4dSSatish Balay 
1489bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
14964f491ddSJed Brown M*/
15064f491ddSJed Brown /*MC
1516cf0794eSJed Brown      TSARKIMEXPRSSP2 - Second order SSP ARK IMEX scheme.
1526cf0794eSJed Brown 
1536cf0794eSJed Brown      This method has three implicit stages.
1546cf0794eSJed Brown 
1556cf0794eSJed Brown      References:
15696a0c994SBarry Smith .   1. -  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.
1576cf0794eSJed Brown 
1586cf0794eSJed Brown      This method is referred to as SSP2-(3,3,2) in http://arxiv.org/abs/1110.4375
1596cf0794eSJed Brown 
1609bd3e852SBarry Smith      Options Database:
1619bd3e852SBarry Smith .      -ts_arkimex_type prssp2
1629bd3e852SBarry Smith 
1636cf0794eSJed Brown      Level: advanced
1646cf0794eSJed Brown 
1659bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
1666cf0794eSJed Brown M*/
1676cf0794eSJed Brown /*MC
16864f491ddSJed Brown      TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part.
16964f491ddSJed Brown 
17064f491ddSJed Brown      This method has one explicit stage and three implicit stages.
17164f491ddSJed Brown 
1729bd3e852SBarry Smith      Options Database:
1739bd3e852SBarry Smith .      -ts_arkimex_type 3
1749bd3e852SBarry Smith 
17564f491ddSJed Brown      References:
17696a0c994SBarry Smith .   1. -  Kennedy and Carpenter 2003.
17764f491ddSJed Brown 
178b330ce4dSSatish Balay      Level: advanced
179b330ce4dSSatish Balay 
1809bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
18164f491ddSJed Brown M*/
18264f491ddSJed Brown /*MC
1836cf0794eSJed Brown      TSARKIMEXARS443 - Third order ARK IMEX scheme.
1846cf0794eSJed Brown 
1856cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
1866cf0794eSJed Brown 
1879bd3e852SBarry Smith      Options Database:
1889bd3e852SBarry Smith .      -ts_arkimex_type ars443
1899bd3e852SBarry Smith 
1906cf0794eSJed Brown      References:
19196a0c994SBarry Smith +   1. -  U. Ascher, S. Ruuth, R. J. Spiteri, Implicit explicit Runge Kutta methods for time dependent Partial Differential Equations. Appl. Numer. Math. 25, (1997).
19296a0c994SBarry Smith -   2. -  This method is referred to as ARS(4,4,3) in http://arxiv.org/abs/1110.4375
1936cf0794eSJed Brown 
1946cf0794eSJed Brown      Level: advanced
1956cf0794eSJed Brown 
1969bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
1976cf0794eSJed Brown M*/
1986cf0794eSJed Brown /*MC
1996cf0794eSJed Brown      TSARKIMEXBPR3 - Third order ARK IMEX scheme.
2006cf0794eSJed Brown 
2016cf0794eSJed Brown      This method has one explicit stage and four implicit stages.
2026cf0794eSJed Brown 
2039bd3e852SBarry Smith      Options Database:
2049bd3e852SBarry Smith .      -ts_arkimex_type bpr3
2059bd3e852SBarry Smith 
2066cf0794eSJed Brown      References:
20796a0c994SBarry Smith  .    This method is referred to as ARK3 in http://arxiv.org/abs/1110.4375
2086cf0794eSJed Brown 
2096cf0794eSJed Brown      Level: advanced
2106cf0794eSJed Brown 
2119bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
2126cf0794eSJed Brown M*/
2136cf0794eSJed Brown /*MC
21464f491ddSJed Brown      TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part.
21564f491ddSJed Brown 
21664f491ddSJed Brown      This method has one explicit stage and four implicit stages.
21764f491ddSJed Brown 
2189bd3e852SBarry Smith      Options Database:
2199bd3e852SBarry Smith .      -ts_arkimex_type 4
2209bd3e852SBarry Smith 
22164f491ddSJed Brown      References:
22296a0c994SBarry Smith .     Kennedy and Carpenter 2003.
22364f491ddSJed Brown 
224b330ce4dSSatish Balay      Level: advanced
225b330ce4dSSatish Balay 
2269bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
22764f491ddSJed Brown M*/
22864f491ddSJed Brown /*MC
22964f491ddSJed Brown      TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part.
23064f491ddSJed Brown 
23164f491ddSJed Brown      This method has one explicit stage and five implicit stages.
23264f491ddSJed Brown 
2339bd3e852SBarry Smith      Options Database:
2349bd3e852SBarry Smith .      -ts_arkimex_type 5
2359bd3e852SBarry Smith 
23664f491ddSJed Brown      References:
23796a0c994SBarry Smith .     Kennedy and Carpenter 2003.
23864f491ddSJed Brown 
239b330ce4dSSatish Balay      Level: advanced
240b330ce4dSSatish Balay 
2419bd3e852SBarry Smith .seealso: TSARKIMEX, TSARKIMEXType, TSARKIMEXSetType()
24264f491ddSJed Brown M*/
24364f491ddSJed Brown 
2448a381b04SJed Brown /*@C
2458a381b04SJed Brown   TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX
2468a381b04SJed Brown 
247fca742c7SJed Brown   Not Collective, but should be called by all processes which will need the schemes to be registered
2488a381b04SJed Brown 
2498a381b04SJed Brown   Level: advanced
2508a381b04SJed Brown 
2518a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all
2528a381b04SJed Brown 
2538a381b04SJed Brown .seealso:  TSARKIMEXRegisterDestroy()
2548a381b04SJed Brown @*/
2558a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void)
2568a381b04SJed Brown {
2578a381b04SJed Brown   PetscErrorCode ierr;
2588a381b04SJed Brown 
2598a381b04SJed Brown   PetscFunctionBegin;
2608a381b04SJed Brown   if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0);
2618a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_TRUE;
262e817cc15SEmil Constantinescu 
263e817cc15SEmil Constantinescu   {
264e817cc15SEmil Constantinescu     const PetscReal
265e817cc15SEmil Constantinescu       A[3][3] = {{0.0,0.0,0.0},
266e817cc15SEmil Constantinescu                  {0.0,0.0,0.0},
267748ad121SEmil Constantinescu                  {0.0,0.5,0.0}},
268e817cc15SEmil Constantinescu       At[3][3] = {{1.0,0.0,0.0},
269e817cc15SEmil Constantinescu                   {0.0,0.5,0.0},
270e817cc15SEmil Constantinescu                   {0.0,0.5,0.5}},
271e817cc15SEmil Constantinescu       b[3]       = {0.0,0.5,0.5},
272e817cc15SEmil Constantinescu       bembedt[3] = {1.0,0.0,0.0};
2730298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX1BEE,2,3,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
274e817cc15SEmil Constantinescu   }
2758a381b04SJed Brown   {
2768a381b04SJed Brown     const PetscReal
2771f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2781f80e275SEmil Constantinescu                  {0.5,0.0}},
2791f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2801f80e275SEmil Constantinescu                   {0.0,0.5}},
2811f80e275SEmil Constantinescu       b[2]       = {0.0,1.0},
2821f80e275SEmil Constantinescu       bembedt[2] = {0.5,0.5};
2831f80e275SEmil 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*/
2840298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS122,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2851f80e275SEmil Constantinescu   }
2861f80e275SEmil Constantinescu   {
2871f80e275SEmil Constantinescu     const PetscReal
2881f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
2891f80e275SEmil Constantinescu                  {1.0,0.0}},
2901f80e275SEmil Constantinescu       At[2][2] = {{0.0,0.0},
2911f80e275SEmil Constantinescu                   {0.5,0.5}},
2921f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
2931f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0};
2941f80e275SEmil 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*/
2950298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXA2,2,2,&At[0][0],b,NULL,&A[0][0],b,NULL,bembedt,bembedt,1,b,NULL);CHKERRQ(ierr);
2961f80e275SEmil Constantinescu   }
2971f80e275SEmil Constantinescu   {
298da80777bSKarl 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   */
2991f80e275SEmil Constantinescu     const PetscReal
3001f80e275SEmil Constantinescu       A[2][2] = {{0.0,0.0},
3011f80e275SEmil Constantinescu                  {1.0,0.0}},
302da80777bSKarl Rupp       At[2][2] = {{0.2928932188134524755992,0.0},
303da80777bSKarl Rupp                   {1.0-2.0*0.2928932188134524755992,0.2928932188134524755992}},
3041f80e275SEmil Constantinescu       b[2]       = {0.5,0.5},
3051f80e275SEmil Constantinescu       bembedt[2] = {0.0,1.0},
306da80777bSKarl Rupp       binterpt[2][2] = {{  (0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))},
307da80777bSKarl Rupp                         {1-(0.2928932188134524755992-1.0)/(2.0*0.2928932188134524755992-1.0),-1/(2.0*(1.0-2.0*0.2928932188134524755992))}},
3081f80e275SEmil Constantinescu       binterp[2][2] = {{1.0,-0.5},{0.0,0.5}};
3090298fd71SBarry 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);
3101f80e275SEmil Constantinescu   }
3111f80e275SEmil Constantinescu   {
312da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
313da80777bSKarl Rupp     const PetscReal
3148a381b04SJed Brown       A[3][3] = {{0,0,0},
315da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
316617a39beSEmil Constantinescu                  {0.5,0.5,0}},
317da80777bSKarl Rupp       At[3][3] = {{0,0,0},
318da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
319da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
320da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
321da80777bSKarl Rupp       binterpt[3][2] = {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
322da80777bSKarl Rupp                         {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
323da80777bSKarl Rupp                         {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3240298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2C,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
3251f80e275SEmil Constantinescu   }
3261f80e275SEmil Constantinescu   {
327da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
328da80777bSKarl Rupp     const PetscReal
3291f80e275SEmil Constantinescu       A[3][3] = {{0,0,0},
330da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
3318a381b04SJed Brown                  {0.75,0.25,0}},
332da80777bSKarl Rupp       At[3][3] = {{0,0,0},
333da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
334da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
335da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
336da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
337da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
338da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3390298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
3408a381b04SJed Brown   }
34106db7b1cSJed Brown   {                             /* Optimal for linear implicit part */
342da80777bSKarl Rupp     /* const PetscReal s2 = PetscSqrtReal((PetscReal)2.0),  Direct evaluation: 1.414213562373095048802. Used below to ensure all values are available at compile time   */
343da80777bSKarl Rupp     const PetscReal
344da80777bSKarl Rupp       A[3][3] = {{0,0,0},
345da80777bSKarl Rupp                  {2-1.414213562373095048802,0,0},
346da80777bSKarl Rupp                  {(3-2*1.414213562373095048802)/6,(3+2*1.414213562373095048802)/6,0}},
347da80777bSKarl Rupp       At[3][3] = {{0,0,0},
348da80777bSKarl Rupp                   {1-1/1.414213562373095048802,1-1/1.414213562373095048802,0},
349da80777bSKarl Rupp                   {1/(2*1.414213562373095048802),1/(2*1.414213562373095048802),1-1/1.414213562373095048802}},
350da80777bSKarl Rupp       bembedt[3] = {(4.-1.414213562373095048802)/8.,(4.-1.414213562373095048802)/8.,1/(2.*1.414213562373095048802)},
351da80777bSKarl Rupp       binterpt[3][2] =  {{1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
352da80777bSKarl Rupp                          {1.0/1.414213562373095048802,-1.0/(2.0*1.414213562373095048802)},
353da80777bSKarl Rupp                          {1.0-1.414213562373095048802,1.0/1.414213562373095048802}};
3540298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
355a3a57f36SJed Brown   }
3566cf0794eSJed Brown   {                             /* Optimal for linear implicit part */
3576cf0794eSJed Brown     const PetscReal
3586cf0794eSJed Brown       A[3][3] = {{0,0,0},
3596cf0794eSJed Brown                  {0.5,0,0},
3606cf0794eSJed Brown                  {0.5,0.5,0}},
3616cf0794eSJed Brown       At[3][3] = {{0.25,0,0},
3626cf0794eSJed Brown                   {0,0.25,0},
3636cf0794eSJed Brown                   {1./3,1./3,1./3}};
3640298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXPRSSP2,2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);CHKERRQ(ierr);
3656cf0794eSJed Brown   }
366a3a57f36SJed Brown   {
367a3a57f36SJed Brown     const PetscReal
368a3a57f36SJed Brown       A[4][4] = {{0,0,0,0},
3694040e9f2SJed Brown                  {1767732205903./2027836641118.,0,0,0},
3704040e9f2SJed Brown                  {5535828885825./10492691773637.,788022342437./10882634858940.,0,0},
3714040e9f2SJed Brown                  {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}},
372a3a57f36SJed Brown       At[4][4] = {{0,0,0,0},
3734040e9f2SJed Brown                   {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0},
3744040e9f2SJed Brown                   {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0},
3754040e9f2SJed Brown                   {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}},
376cc46b9d1SJed Brown       bembedt[4]     = {2756255671327./12835298489170.,-10771552573575./22201958757719.,9247589265047./10645013368117.,2193209047091./5459859503100.},
3774040e9f2SJed Brown       binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.},
3784040e9f2SJed Brown                         {-18682724506714./9892148508045.,17870216137069./13817060693119.},
3794040e9f2SJed Brown                         {34259539580243./13192909600954.,-28141676662227./17317692491321.},
3804040e9f2SJed Brown                         {584795268549./6622622206610.,   2508943948391./7218656332882.}};
3810298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,2,binterpt[0],NULL);CHKERRQ(ierr);
382a3a57f36SJed Brown   }
383a3a57f36SJed Brown   {
384a3a57f36SJed Brown     const PetscReal
385e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
3866cf0794eSJed Brown                  {1./2,0,0,0,0},
3876cf0794eSJed Brown                  {11./18,1./18,0,0,0},
3886cf0794eSJed Brown                  {5./6,-5./6,.5,0,0},
3896cf0794eSJed Brown                  {1./4,7./4,3./4,-7./4,0}},
3906cf0794eSJed Brown       At[5][5] = {{0,0,0,0,0},
3916cf0794eSJed Brown                   {0,1./2,0,0,0},
3926cf0794eSJed Brown                   {0,1./6,1./2,0,0},
3936cf0794eSJed Brown                   {0,-1./2,1./2,1./2,0},
394108c343cSJed Brown                   {0,3./2,-3./2,1./2,1./2}},
3950298fd71SBarry Smith     *bembedt = NULL;
3960298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXARS443,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
3976cf0794eSJed Brown   }
3986cf0794eSJed Brown   {
3996cf0794eSJed Brown     const PetscReal
400e74514c0SSatish Balay       A[5][5] = {{0,0,0,0,0},
4016cf0794eSJed Brown                  {1,0,0,0,0},
4026cf0794eSJed Brown                  {4./9,2./9,0,0,0},
4036cf0794eSJed Brown                  {1./4,0,3./4,0,0},
4046cf0794eSJed Brown                  {1./4,0,3./5,0,0}},
405e74514c0SSatish Balay       At[5][5] = {{0,0,0,0,0},
4066cf0794eSJed Brown                   {.5,.5,0,0,0},
4076cf0794eSJed Brown                   {5./18,-1./9,.5,0,0},
4086cf0794eSJed Brown                   {.5,0,0,.5,0},
409108c343cSJed Brown                   {.25,0,.75,-.5,.5}},
4100298fd71SBarry Smith     *bembedt = NULL;
4110298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEXBPR3,3,5,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,0,NULL,NULL);CHKERRQ(ierr);
4126cf0794eSJed Brown   }
4136cf0794eSJed Brown   {
4146cf0794eSJed Brown     const PetscReal
415a3a57f36SJed Brown       A[6][6] = {{0,0,0,0,0,0},
416a3a57f36SJed Brown                  {1./2,0,0,0,0,0},
4174040e9f2SJed Brown                  {13861./62500.,6889./62500.,0,0,0,0},
4184040e9f2SJed Brown                  {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0},
4194040e9f2SJed Brown                  {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0},
4204040e9f2SJed Brown                  {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}},
421a3a57f36SJed Brown       At[6][6] = {{0,0,0,0,0,0},
422a3a57f36SJed Brown                   {1./4,1./4,0,0,0,0},
4234040e9f2SJed Brown                   {8611./62500.,-1743./31250.,1./4,0,0,0},
4244040e9f2SJed Brown                   {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0},
4254040e9f2SJed Brown                   {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0},
4264040e9f2SJed Brown                   {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}},
427cc46b9d1SJed Brown       bembedt[6]     = {4586570599./29645900160.,0,178811875./945068544.,814220225./1159782912.,-3700637./11593932.,61727./225920.},
4284040e9f2SJed Brown       binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.},
429cd652676SJed Brown                         {0,0,0},
4304040e9f2SJed Brown                         {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.},
4314040e9f2SJed Brown                         {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.},
4324040e9f2SJed Brown                         {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.},
4334040e9f2SJed Brown                         {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}};
4340298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
435a3a57f36SJed Brown   }
436a3a57f36SJed Brown   {
437a3a57f36SJed Brown     const PetscReal
438a3a57f36SJed Brown       A[8][8] = {{0,0,0,0,0,0,0,0},
439a3a57f36SJed Brown                  {41./100,0,0,0,0,0,0,0},
4404040e9f2SJed Brown                  {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0},
4414040e9f2SJed Brown                  {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0},
4424040e9f2SJed Brown                  {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0},
4434040e9f2SJed Brown                  {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0},
4444040e9f2SJed Brown                  {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0},
4454040e9f2SJed Brown                  {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}},
446a3a57f36SJed Brown       At[8][8] = {{0,0,0,0,0,0,0,0},
4474040e9f2SJed Brown                   {41./200.,41./200.,0,0,0,0,0,0},
4484040e9f2SJed Brown                   {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0},
4494040e9f2SJed Brown                   {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0},
4504040e9f2SJed Brown                   {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0},
4514040e9f2SJed Brown                   {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0},
4524040e9f2SJed Brown                   {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0},
4534040e9f2SJed Brown                   {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}},
454cc46b9d1SJed Brown       bembedt[8]     = {-975461918565./9796059967033.,0,0,78070527104295./32432590147079.,-548382580838./3424219808633.,-33438840321285./15594753105479.,3629800801594./4656183773603.,4035322873751./18575991585200.},
4554040e9f2SJed Brown       binterpt[8][3] = {{-17674230611817./10670229744614.,  43486358583215./12773830924787., -9257016797708./5021505065439.},
456cd652676SJed Brown                         {0,  0, 0                            },
457cd652676SJed Brown                         {0,  0, 0                            },
4584040e9f2SJed Brown                         {65168852399939./7868540260826.,  -91478233927265./11067650958493., 26096422576131./11239449250142.},
4594040e9f2SJed Brown                         {15494834004392./5936557850923.,  -79368583304911./10890268929626., 92396832856987./20362823103730.},
4604040e9f2SJed Brown                         {-99329723586156./26959484932159.,  -12239297817655./9152339842473., 30029262896817./10175596800299.},
4614040e9f2SJed Brown                         {-19024464361622./5461577185407.,  115839755401235./10719374521269., -26136350496073./3983972220547.},
4624040e9f2SJed Brown                         {-6511271360970./6095937251113.,  5843115559534./2180450260947., -5289405421727./3760307252460. }};
4630298fd71SBarry Smith     ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembedt,3,binterpt[0],NULL);CHKERRQ(ierr);
464a3a57f36SJed Brown   }
4658a381b04SJed Brown   PetscFunctionReturn(0);
4668a381b04SJed Brown }
4678a381b04SJed Brown 
4688a381b04SJed Brown /*@C
4698a381b04SJed Brown    TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister().
4708a381b04SJed Brown 
4718a381b04SJed Brown    Not Collective
4728a381b04SJed Brown 
4738a381b04SJed Brown    Level: advanced
4748a381b04SJed Brown 
4758a381b04SJed Brown .keywords: TSARKIMEX, register, destroy
476607a6623SBarry Smith .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll()
4778a381b04SJed Brown @*/
4788a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void)
4798a381b04SJed Brown {
4808a381b04SJed Brown   PetscErrorCode ierr;
4818a381b04SJed Brown   ARKTableauLink link;
4828a381b04SJed Brown 
4838a381b04SJed Brown   PetscFunctionBegin;
4848a381b04SJed Brown   while ((link = ARKTableauList)) {
4858a381b04SJed Brown     ARKTableau t = &link->tab;
4868a381b04SJed Brown     ARKTableauList = link->next;
4878a381b04SJed Brown     ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr);
488108c343cSJed Brown     ierr = PetscFree2(t->bembedt,t->bembed);CHKERRQ(ierr);
489cd652676SJed Brown     ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr);
4908a381b04SJed Brown     ierr = PetscFree(t->name);CHKERRQ(ierr);
4918a381b04SJed Brown     ierr = PetscFree(link);CHKERRQ(ierr);
4928a381b04SJed Brown   }
4938a381b04SJed Brown   TSARKIMEXRegisterAllCalled = PETSC_FALSE;
4948a381b04SJed Brown   PetscFunctionReturn(0);
4958a381b04SJed Brown }
4968a381b04SJed Brown 
4978a381b04SJed Brown /*@C
4988a381b04SJed Brown   TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called
4998a381b04SJed Brown   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX()
5008a381b04SJed Brown   when using static libraries.
5018a381b04SJed Brown 
5028a381b04SJed Brown   Level: developer
5038a381b04SJed Brown 
5048a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package
5058a381b04SJed Brown .seealso: PetscInitialize()
5068a381b04SJed Brown @*/
507607a6623SBarry Smith PetscErrorCode TSARKIMEXInitializePackage(void)
5088a381b04SJed Brown {
5098a381b04SJed Brown   PetscErrorCode ierr;
5108a381b04SJed Brown 
5118a381b04SJed Brown   PetscFunctionBegin;
5128a381b04SJed Brown   if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0);
5138a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_TRUE;
5148a381b04SJed Brown   ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr);
5158a381b04SJed Brown   ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr);
5168a381b04SJed Brown   PetscFunctionReturn(0);
5178a381b04SJed Brown }
5188a381b04SJed Brown 
5198a381b04SJed Brown /*@C
5208a381b04SJed Brown   TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is
5218a381b04SJed Brown   called from PetscFinalize().
5228a381b04SJed Brown 
5238a381b04SJed Brown   Level: developer
5248a381b04SJed Brown 
5258a381b04SJed Brown .keywords: Petsc, destroy, package
5268a381b04SJed Brown .seealso: PetscFinalize()
5278a381b04SJed Brown @*/
5288a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void)
5298a381b04SJed Brown {
5308a381b04SJed Brown   PetscErrorCode ierr;
5318a381b04SJed Brown 
5328a381b04SJed Brown   PetscFunctionBegin;
5338a381b04SJed Brown   TSARKIMEXPackageInitialized = PETSC_FALSE;
5348a381b04SJed Brown   ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr);
5358a381b04SJed Brown   PetscFunctionReturn(0);
5368a381b04SJed Brown }
5378a381b04SJed Brown 
538cd652676SJed Brown /*@C
539cd652676SJed Brown    TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
540cd652676SJed Brown 
541cd652676SJed Brown    Not Collective, but the same schemes should be registered on all processes on which they will be used
542cd652676SJed Brown 
543cd652676SJed Brown    Input Parameters:
544cd652676SJed Brown +  name - identifier for method
545cd652676SJed Brown .  order - approximation order of method
546cd652676SJed Brown .  s - number of stages, this is the dimension of the matrices below
547cd652676SJed Brown .  At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major)
5480298fd71SBarry Smith .  bt - Butcher table for completing the stiff part of the step (dimension s; NULL to use the last row of At)
5490298fd71SBarry Smith .  ct - Abscissa of each stiff stage (dimension s, NULL to use row sums of At)
550cd652676SJed Brown .  A - Non-stiff stage coefficients (dimension s*s, row-major)
5510298fd71SBarry Smith .  b - Non-stiff step completion table (dimension s; NULL to use last row of At)
5520298fd71SBarry Smith .  c - Non-stiff abscissa (dimension s; NULL to use row sums of A)
5530298fd71SBarry Smith .  bembedt - Stiff part of completion table for embedded method (dimension s; NULL if not available)
5540298fd71SBarry Smith .  bembed - Non-stiff part of completion table for embedded method (dimension s; NULL to use bembedt if provided)
555cd652676SJed Brown .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp
556cd652676SJed Brown .  binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp)
5570298fd71SBarry Smith -  binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; NULL to reuse binterpt)
558cd652676SJed Brown 
559cd652676SJed Brown    Notes:
560cd652676SJed Brown    Several ARK IMEX methods are provided, this function is only needed to create new methods.
561cd652676SJed Brown 
562cd652676SJed Brown    Level: advanced
563cd652676SJed Brown 
564cd652676SJed Brown .keywords: TS, register
565cd652676SJed Brown 
566cd652676SJed Brown .seealso: TSARKIMEX
567cd652676SJed Brown @*/
56819fd82e9SBarry Smith PetscErrorCode TSARKIMEXRegister(TSARKIMEXType name,PetscInt order,PetscInt s,
5698a381b04SJed Brown                                  const PetscReal At[],const PetscReal bt[],const PetscReal ct[],
570cd652676SJed Brown                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
571108c343cSJed Brown                                  const PetscReal bembedt[],const PetscReal bembed[],
572cd652676SJed Brown                                  PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[])
5738a381b04SJed Brown {
5748a381b04SJed Brown   PetscErrorCode ierr;
5758a381b04SJed Brown   ARKTableauLink link;
5768a381b04SJed Brown   ARKTableau     t;
5778a381b04SJed Brown   PetscInt       i,j;
5788a381b04SJed Brown 
5798a381b04SJed Brown   PetscFunctionBegin;
5801795a4d1SJed Brown   ierr     = PetscCalloc1(1,&link);CHKERRQ(ierr);
5818a381b04SJed Brown   t        = &link->tab;
5828a381b04SJed Brown   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
5838a381b04SJed Brown   t->order = order;
5848a381b04SJed Brown   t->s     = s;
585dcca6d9dSJed Brown   ierr     = PetscMalloc6(s*s,&t->At,s,&t->bt,s,&t->ct,s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
5868a381b04SJed Brown   ierr     = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr);
5878a381b04SJed Brown   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
5888a381b04SJed Brown   if (bt) { ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr); }
5898a381b04SJed Brown   else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i];
5908a381b04SJed Brown   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
5915dceddf7SDebojyoti Ghosh   else for (i=0; i<s; i++) t->b[i] = t->bt[i];
5928a381b04SJed Brown   if (ct) { ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr); }
5938a381b04SJed 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];
5948a381b04SJed Brown   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
5958a381b04SJed 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];
596e817cc15SEmil Constantinescu   t->stiffly_accurate = PETSC_TRUE;
597e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[(s-1)*s+i] != t->bt[i]) t->stiffly_accurate = PETSC_FALSE;
598e817cc15SEmil Constantinescu   t->explicit_first_stage = PETSC_TRUE;
599e817cc15SEmil Constantinescu   for (i=0; i<s; i++) if (t->At[i] != 0.0) t->explicit_first_stage = PETSC_FALSE;
600e817cc15SEmil Constantinescu   /*def of FSAL can be made more precise*/
6014e9d4bf5SJed Brown   t->FSAL_implicit = (PetscBool)(t->explicit_first_stage && t->stiffly_accurate);
602108c343cSJed Brown   if (bembedt) {
603dcca6d9dSJed Brown     ierr = PetscMalloc2(s,&t->bembedt,s,&t->bembed);CHKERRQ(ierr);
604108c343cSJed Brown     ierr = PetscMemcpy(t->bembedt,bembedt,s*sizeof(bembedt[0]));CHKERRQ(ierr);
605108c343cSJed Brown     ierr = PetscMemcpy(t->bembed,bembed ? bembed : bembedt,s*sizeof(bembed[0]));CHKERRQ(ierr);
606108c343cSJed Brown   }
607108c343cSJed Brown 
6084f385281SJed Brown   t->pinterp     = pinterp;
609dcca6d9dSJed Brown   ierr           = PetscMalloc2(s*pinterp,&t->binterpt,s*pinterp,&t->binterp);CHKERRQ(ierr);
610cd652676SJed Brown   ierr           = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
611cd652676SJed Brown   ierr           = PetscMemcpy(t->binterp,binterp ? binterp : binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr);
6128a381b04SJed Brown   link->next     = ARKTableauList;
6138a381b04SJed Brown   ARKTableauList = link;
6148a381b04SJed Brown   PetscFunctionReturn(0);
6158a381b04SJed Brown }
6168a381b04SJed Brown 
617108c343cSJed Brown /*
618108c343cSJed Brown  The step completion formula is
619108c343cSJed Brown 
620108c343cSJed Brown  x1 = x0 - h bt^T YdotI + h b^T YdotRHS
621108c343cSJed Brown 
622108c343cSJed Brown  This function can be called before or after ts->vec_sol has been updated.
623108c343cSJed Brown  Suppose we have a completion formula (bt,b) and an embedded formula (bet,be) of different order.
624108c343cSJed Brown  We can write
625108c343cSJed Brown 
626108c343cSJed Brown  x1e = x0 - h bet^T YdotI + h be^T YdotRHS
627108c343cSJed Brown      = x1 + h bt^T YdotI - h b^T YdotRHS - h bet^T YdotI + h be^T YdotRHS
628108c343cSJed Brown      = x1 - h (bet - bt)^T YdotI + h (be - b)^T YdotRHS
629108c343cSJed Brown 
630108c343cSJed Brown  so we can evaluate the method with different order even after the step has been optimistically completed.
631108c343cSJed Brown */
632108c343cSJed Brown static PetscErrorCode TSEvaluateStep_ARKIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
633108c343cSJed Brown {
634108c343cSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
635108c343cSJed Brown   ARKTableau     tab  = ark->tableau;
636108c343cSJed Brown   PetscScalar    *w   = ark->work;
637108c343cSJed Brown   PetscReal      h;
638108c343cSJed Brown   PetscInt       s = tab->s,j;
639108c343cSJed Brown   PetscErrorCode ierr;
640108c343cSJed Brown 
641108c343cSJed Brown   PetscFunctionBegin;
642108c343cSJed Brown   switch (ark->status) {
643108c343cSJed Brown   case TS_STEP_INCOMPLETE:
644108c343cSJed Brown   case TS_STEP_PENDING:
645108c343cSJed Brown     h = ts->time_step; break;
646108c343cSJed Brown   case TS_STEP_COMPLETE:
647be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
648ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
649108c343cSJed Brown   }
650108c343cSJed Brown   if (order == tab->order) {
651e817cc15SEmil Constantinescu     if (ark->status == TS_STEP_INCOMPLETE) {
652740132f1SEmil Constantinescu       if (!ark->imex && tab->stiffly_accurate) { /* Only the stiffly accurate implicit formula is used */
653e817cc15SEmil Constantinescu         ierr = VecCopy(ark->Y[s-1],X);CHKERRQ(ierr);
654e817cc15SEmil Constantinescu       } else { /* Use the standard completion formula (bt,b) */
655108c343cSJed Brown         ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
656e817cc15SEmil Constantinescu         for (j=0; j<s; j++) w[j] = h*tab->bt[j];
657108c343cSJed Brown         ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
658e817cc15SEmil Constantinescu         if (ark->imex) { /* Method is IMEX, complete the explicit formula */
659108c343cSJed Brown           for (j=0; j<s; j++) w[j] = h*tab->b[j];
660108c343cSJed Brown           ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
661e817cc15SEmil Constantinescu         }
662e817cc15SEmil Constantinescu       }
663108c343cSJed Brown     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
664108c343cSJed Brown     if (done) *done = PETSC_TRUE;
665108c343cSJed Brown     PetscFunctionReturn(0);
666108c343cSJed Brown   } else if (order == tab->order-1) {
667108c343cSJed Brown     if (!tab->bembedt) goto unavailable;
668108c343cSJed Brown     if (ark->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (bet,be) */
669108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
670e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembedt[j];
671108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotI);CHKERRQ(ierr);
672108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
673108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
674108c343cSJed Brown     } else { /* Rollback and re-complete using (bet-be,be-b) */
675108c343cSJed Brown       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
676e817cc15SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembedt[j] - tab->bt[j]);
677108c343cSJed Brown       ierr = VecMAXPY(X,tab->s,w,ark->YdotI);CHKERRQ(ierr);
678108c343cSJed Brown       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
679108c343cSJed Brown       ierr = VecMAXPY(X,s,w,ark->YdotRHS);CHKERRQ(ierr);
680108c343cSJed Brown     }
681108c343cSJed Brown     if (done) *done = PETSC_TRUE;
682108c343cSJed Brown     PetscFunctionReturn(0);
683108c343cSJed Brown   }
684108c343cSJed Brown unavailable:
685108c343cSJed Brown   if (done) *done = PETSC_FALSE;
686a7fac7c2SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"ARKIMEX '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
687108c343cSJed Brown   PetscFunctionReturn(0);
688108c343cSJed Brown }
689108c343cSJed Brown 
69024655328SShri static PetscErrorCode TSRollBack_ARKIMEX(TS ts)
69124655328SShri {
69224655328SShri   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
69324655328SShri   ARKTableau      tab  = ark->tableau;
69424655328SShri   const PetscInt  s    = tab->s;
69524655328SShri   const PetscReal *bt  = tab->bt,*b = tab->b;
69624655328SShri   PetscScalar     *w   = ark->work;
69724655328SShri   Vec             *YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS;
69824655328SShri   PetscInt        j;
699be5899b3SLisandro Dalcin   PetscReal       h;
70024655328SShri   PetscErrorCode  ierr;
70124655328SShri 
70224655328SShri   PetscFunctionBegin;
703be5899b3SLisandro Dalcin   switch (ark->status) {
704be5899b3SLisandro Dalcin   case TS_STEP_INCOMPLETE:
705be5899b3SLisandro Dalcin   case TS_STEP_PENDING:
706be5899b3SLisandro Dalcin     h = ts->time_step; break;
707be5899b3SLisandro Dalcin   case TS_STEP_COMPLETE:
708be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
709be5899b3SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
710be5899b3SLisandro Dalcin   }
71124655328SShri   for (j=0; j<s; j++) w[j] = -h*bt[j];
71224655328SShri   ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr);
71324655328SShri   for (j=0; j<s; j++) w[j] = -h*b[j];
71424655328SShri   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
71524655328SShri   PetscFunctionReturn(0);
71624655328SShri }
71724655328SShri 
7188a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts)
7198a381b04SJed Brown {
7208a381b04SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
7218a381b04SJed Brown   ARKTableau      tab  = ark->tableau;
7228a381b04SJed Brown   const PetscInt  s    = tab->s;
72324655328SShri   const PetscReal *At  = tab->At,*A = tab->A,*ct = tab->ct,*c = tab->c;
724406d0ec2SJed Brown   PetscScalar     *w   = ark->work;
7251297b224SEmil Constantinescu   Vec             *Y   = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,Ydot0 = ark->Ydot0,Z = ark->Z;
72696400bd6SLisandro Dalcin   PetscBool       extrapolate = ark->extrapolate;
727108c343cSJed Brown   TSAdapt         adapt;
7288a381b04SJed Brown   SNES            snes;
729fecfb714SLisandro Dalcin   PetscInt        i,j,its,lits;
730be5899b3SLisandro Dalcin   PetscInt        rejections = 0;
73196400bd6SLisandro Dalcin   PetscBool       stageok,accept = PETSC_TRUE;
73296400bd6SLisandro Dalcin   PetscReal       next_time_step = ts->time_step;
7338a381b04SJed Brown   PetscErrorCode  ierr;
7348a381b04SJed Brown 
7358a381b04SJed Brown   PetscFunctionBegin;
73696400bd6SLisandro Dalcin   if (ark->extrapolate && !ark->Y_prev) {
73796400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr);
73896400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr);
73996400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr);
74096400bd6SLisandro Dalcin   }
74196400bd6SLisandro Dalcin 
74296400bd6SLisandro Dalcin   if (!ts->steprollback) {
74396400bd6SLisandro Dalcin     if (ts->equation_type >= TS_EQ_IMPLICIT) { /* Save the initial slope for the next step */
74496400bd6SLisandro Dalcin       ierr = VecCopy(YdotI[s-1],Ydot0);CHKERRQ(ierr);
74596400bd6SLisandro Dalcin     }
746fecfb714SLisandro Dalcin     if (ark->extrapolate && !ts->steprestart) { /* Save the Y, YdotI, YdotRHS for extrapolation initial guess */
74796400bd6SLisandro Dalcin       for (i = 0; i<s; i++) {
74896400bd6SLisandro Dalcin         ierr = VecCopy(Y[i],ark->Y_prev[i]);CHKERRQ(ierr);
74996400bd6SLisandro Dalcin         ierr = VecCopy(YdotRHS[i],ark->YdotRHS_prev[i]);CHKERRQ(ierr);
75096400bd6SLisandro Dalcin         ierr = VecCopy(YdotI[i],ark->YdotI_prev[i]);CHKERRQ(ierr);
75196400bd6SLisandro Dalcin       }
75296400bd6SLisandro Dalcin     }
75396400bd6SLisandro Dalcin   }
75496400bd6SLisandro Dalcin 
755fecfb714SLisandro Dalcin   if (ts->equation_type >= TS_EQ_IMPLICIT && tab->explicit_first_stage && ts->steprestart) {
75696400bd6SLisandro Dalcin     TS ts_start;
757baa10174SEmil Constantinescu     ierr = TSClone(ts,&ts_start);CHKERRQ(ierr);
758e817cc15SEmil Constantinescu     ierr = TSSetSolution(ts_start,ts->vec_sol);CHKERRQ(ierr);
759e817cc15SEmil Constantinescu     ierr = TSSetTime(ts_start,ts->ptime);CHKERRQ(ierr);
76019eac22cSLisandro Dalcin     ierr = TSSetMaxSteps(ts_start,ts->steps+1);CHKERRQ(ierr);
76119eac22cSLisandro Dalcin     ierr = TSSetMaxTime(ts_start,ts->ptime+ts->time_step);CHKERRQ(ierr);
762feed9e9dSBarry Smith     ierr = TSSetExactFinalTime(ts_start,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
763740132f1SEmil Constantinescu     ierr = TSSetTimeStep(ts_start,ts->time_step);CHKERRQ(ierr);
764e817cc15SEmil Constantinescu     ierr = TSSetType(ts_start,TSARKIMEX);CHKERRQ(ierr);
765740132f1SEmil Constantinescu     ierr = TSARKIMEXSetFullyImplicit(ts_start,PETSC_TRUE);CHKERRQ(ierr);
766e817cc15SEmil Constantinescu     ierr = TSARKIMEXSetType(ts_start,TSARKIMEX1BEE);CHKERRQ(ierr);
76734561852SEmil Constantinescu 
768dcb233daSLisandro Dalcin     ierr = TSRestartStep(ts_start);CHKERRQ(ierr);
769e817cc15SEmil Constantinescu     ierr = TSSolve(ts_start,ts->vec_sol);CHKERRQ(ierr);
770e817cc15SEmil Constantinescu     ierr = TSGetTime(ts_start,&ts->ptime);CHKERRQ(ierr);
77196400bd6SLisandro Dalcin     ierr = TSGetTimeStep(ts_start,&ts->time_step);CHKERRQ(ierr);
772bbd56ea5SKarl Rupp 
77385fc7851SLisandro Dalcin     { /* Save the initial slope for the next step */
77485fc7851SLisandro Dalcin       TS_ARKIMEX *ark_start = (TS_ARKIMEX*)ts_start->data;
77585fc7851SLisandro Dalcin       ierr = VecCopy(ark_start->YdotI[ark_start->tableau->s-1],Ydot0);CHKERRQ(ierr);
77685fc7851SLisandro Dalcin     }
77796400bd6SLisandro Dalcin     ts->steps++;
77834561852SEmil Constantinescu 
779d15a3a53SEmil Constantinescu     /* Set the correct TS in SNES */
780d15a3a53SEmil Constantinescu     /* We'll try to bypass this by changing the method on the fly */
78196400bd6SLisandro Dalcin     {
78296400bd6SLisandro Dalcin       ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
78396400bd6SLisandro Dalcin       ierr = TSSetSNES(ts,snes);CHKERRQ(ierr);
78496400bd6SLisandro Dalcin     }
785166a6834SEmil Constantinescu     ierr = TSDestroy(&ts_start);CHKERRQ(ierr);
786e817cc15SEmil Constantinescu   }
787e817cc15SEmil Constantinescu 
788108c343cSJed Brown   ark->status = TS_STEP_INCOMPLETE;
78996400bd6SLisandro Dalcin   while (!ts->reason && ark->status != TS_STEP_COMPLETE) {
79096400bd6SLisandro Dalcin     PetscReal t = ts->ptime;
791108c343cSJed Brown     PetscReal h = ts->time_step;
7928a381b04SJed Brown     for (i=0; i<s; i++) {
7939be3e283SDebojyoti Ghosh       ark->stage_time = t + h*ct[i];
79496400bd6SLisandro Dalcin       ierr = TSPreStage(ts,ark->stage_time);CHKERRQ(ierr);
7958a381b04SJed Brown       if (At[i*s+i] == 0) { /* This stage is explicit */
7966c4ed002SBarry Smith         if (i!=0 && ts->equation_type >= TS_EQ_IMPLICIT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Explicit stages other than the first one are not supported for implicit problems");
7978a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
798e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
7998a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr);
8008a381b04SJed Brown         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
8018a381b04SJed Brown         ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
8028a381b04SJed Brown       } else {
803b296d7d5SJed Brown         ark->scoeff = 1./At[i*s+i];
8048a381b04SJed Brown         /* Ydot = shift*(Y-Z) */
8058a381b04SJed Brown         ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr);
806e817cc15SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*At[i*s+j];
8074f385281SJed Brown         ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr);
808c58d1302SEmil Constantinescu         for (j=0; j<i; j++) w[j] = h*A[i*s+j];
809c58d1302SEmil Constantinescu         ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr);
810fecfb714SLisandro Dalcin         if (extrapolate && !ts->steprestart) {
81156dcabbaSDebojyoti Ghosh           /* Initial guess extrapolated from previous time step stage values */
81256dcabbaSDebojyoti Ghosh           ierr = TSExtrapolate_ARKIMEX(ts,c[i],Y[i]);CHKERRQ(ierr);
81356dcabbaSDebojyoti Ghosh         } else {
8148a381b04SJed Brown           /* Initial guess taken from last stage */
8158a381b04SJed Brown           ierr = VecCopy(i>0 ? Y[i-1] : ts->vec_sol,Y[i]);CHKERRQ(ierr);
81656dcabbaSDebojyoti Ghosh         }
81796400bd6SLisandro Dalcin         ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
818baa10174SEmil Constantinescu         ierr = SNESSolve(snes,NULL,Y[i]);CHKERRQ(ierr);
8198a381b04SJed Brown         ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
8208a381b04SJed Brown         ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
8215ef26d82SJed Brown         ts->snes_its += its; ts->ksp_its += lits;
822552698daSJed Brown         ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
82396400bd6SLisandro Dalcin         ierr = TSAdaptCheckStage(adapt,ts,ark->stage_time,Y[i],&stageok);CHKERRQ(ierr);
82496400bd6SLisandro Dalcin         if (!stageok) {
8251be93e3eSJed Brown           /* We are likely rejecting the step because of solver or function domain problems so we should not attempt to
8261be93e3eSJed Brown            * use extrapolation to initialize the solves on the next attempt. */
82796400bd6SLisandro Dalcin           extrapolate = PETSC_FALSE;
8281be93e3eSJed Brown           goto reject_step;
8291be93e3eSJed Brown         }
8308a381b04SJed Brown       }
831e817cc15SEmil Constantinescu       if (ts->equation_type >= TS_EQ_IMPLICIT) {
832e817cc15SEmil Constantinescu         if (i==0 && tab->explicit_first_stage) {
8336c4ed002SBarry Smith           if (!tab->stiffly_accurate ) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s is not stiffly accurate and therefore explicit-first stage methods cannot be used if the equation is implicit because the slope cannot be evaluated",ark->tableau->name);
834df5e1e3dSEmil Constantinescu           ierr = VecCopy(Ydot0,YdotI[0]);CHKERRQ(ierr);                                      /* YdotI = YdotI(tn-1) */
835e817cc15SEmil Constantinescu         } else {
836df5e1e3dSEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr);  /* YdotI = shift*(X-Z) */
837e817cc15SEmil Constantinescu         }
838e817cc15SEmil Constantinescu       } else {
8395eca1a21SEmil Constantinescu         if (i==0 && tab->explicit_first_stage) {
8408a381b04SJed Brown           ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
841df5e1e3dSEmil Constantinescu           ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr);/* YdotI = -G(t,Y,0)   */
842e817cc15SEmil Constantinescu           ierr = VecScale(YdotI[i],-1.0);CHKERRQ(ierr);
8435eca1a21SEmil Constantinescu         } else {
844df5e1e3dSEmil Constantinescu           ierr = VecAXPBYPCZ(YdotI[i],-ark->scoeff/h,ark->scoeff/h,0,Z,Y[i]);CHKERRQ(ierr);  /* YdotI = shift*(X-Z) */
8455eca1a21SEmil Constantinescu         }
8464cc180ffSJed Brown         if (ark->imex) {
8478a381b04SJed Brown           ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
8484cc180ffSJed Brown         } else {
8494cc180ffSJed Brown           ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr);
8504cc180ffSJed Brown         }
8518a381b04SJed Brown       }
85296400bd6SLisandro Dalcin       ierr = TSPostStage(ts,ark->stage_time,i,Y); CHKERRQ(ierr);
853e817cc15SEmil Constantinescu     }
85496400bd6SLisandro Dalcin 
855be5899b3SLisandro Dalcin     ark->status = TS_STEP_INCOMPLETE;
856fecfb714SLisandro Dalcin     ierr = TSEvaluateStep_ARKIMEX(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
857108c343cSJed Brown     ark->status = TS_STEP_PENDING;
858552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
859108c343cSJed Brown     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
860fecfb714SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
861fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
86296400bd6SLisandro Dalcin     ark->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
86396400bd6SLisandro Dalcin     if (!accept) { /* Roll back the current step */
86496400bd6SLisandro Dalcin       ierr = TSRollBack_ARKIMEX(ts);CHKERRQ(ierr);
865be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
866be5899b3SLisandro Dalcin       goto reject_step;
86796400bd6SLisandro Dalcin     }
86896400bd6SLisandro Dalcin 
8698a381b04SJed Brown     ts->ptime += ts->time_step;
870cdbf8f93SLisandro Dalcin     ts->time_step = next_time_step;
871108c343cSJed Brown     break;
87296400bd6SLisandro Dalcin 
87396400bd6SLisandro Dalcin   reject_step:
874fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
875be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
87696400bd6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
877be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
878108c343cSJed Brown     }
879f85781f1SEmil Constantinescu   }
8808a381b04SJed Brown   PetscFunctionReturn(0);
8818a381b04SJed Brown }
8828a381b04SJed Brown 
883cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X)
884cd652676SJed Brown {
885cd652676SJed Brown   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
8864f385281SJed Brown   PetscInt        s    = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
887108c343cSJed Brown   PetscReal       h;
888108c343cSJed Brown   PetscReal       tt,t;
889cd652676SJed Brown   PetscScalar     *bt,*b;
890cd652676SJed Brown   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
891cd652676SJed Brown   PetscErrorCode  ierr;
892cd652676SJed Brown 
893cd652676SJed Brown   PetscFunctionBegin;
894ce94432eSBarry Smith   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
895108c343cSJed Brown   switch (ark->status) {
896108c343cSJed Brown   case TS_STEP_INCOMPLETE:
897108c343cSJed Brown   case TS_STEP_PENDING:
898108c343cSJed Brown     h = ts->time_step;
899108c343cSJed Brown     t = (itime - ts->ptime)/h;
900108c343cSJed Brown     break;
901108c343cSJed Brown   case TS_STEP_COMPLETE:
902be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
903108c343cSJed Brown     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
904108c343cSJed Brown     break;
905ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
906108c343cSJed Brown   }
907dcca6d9dSJed Brown   ierr = PetscMalloc2(s,&bt,s,&b);CHKERRQ(ierr);
908cd652676SJed Brown   for (i=0; i<s; i++) bt[i] = b[i] = 0;
9094f385281SJed Brown   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
910cd652676SJed Brown     for (i=0; i<s; i++) {
911c1758d98SDebojyoti Ghosh       bt[i] += h * Bt[i*pinterp+j] * tt;
912108c343cSJed Brown       b[i]  += h * B[i*pinterp+j] * tt;
913cd652676SJed Brown     }
914cd652676SJed Brown   }
915cd652676SJed Brown   ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr);
916cd652676SJed Brown   ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr);
917cd652676SJed Brown   ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr);
918cd652676SJed Brown   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
919cd652676SJed Brown   PetscFunctionReturn(0);
920cd652676SJed Brown }
921cd652676SJed Brown 
92256dcabbaSDebojyoti Ghosh static PetscErrorCode TSExtrapolate_ARKIMEX(TS ts,PetscReal c,Vec X)
92356dcabbaSDebojyoti Ghosh {
92456dcabbaSDebojyoti Ghosh   TS_ARKIMEX      *ark = (TS_ARKIMEX*)ts->data;
92556dcabbaSDebojyoti Ghosh   PetscInt        s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j;
926be5899b3SLisandro Dalcin   PetscReal       h,h_prev,t,tt;
92756dcabbaSDebojyoti Ghosh   PetscScalar     *bt,*b;
92856dcabbaSDebojyoti Ghosh   const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp;
92956dcabbaSDebojyoti Ghosh   PetscErrorCode  ierr;
93056dcabbaSDebojyoti Ghosh 
93156dcabbaSDebojyoti Ghosh   PetscFunctionBegin;
93256dcabbaSDebojyoti Ghosh   if (!Bt || !B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name);
933be5899b3SLisandro Dalcin   ierr = PetscCalloc2(s,&bt,s,&b);CHKERRQ(ierr);
93481d12688SDebojyoti Ghosh   h = ts->time_step;
935be5899b3SLisandro Dalcin   h_prev = ts->ptime - ts->ptime_prev;
936be5899b3SLisandro Dalcin   t = 1 + h/h_prev*c;
93756dcabbaSDebojyoti Ghosh   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
93856dcabbaSDebojyoti Ghosh     for (i=0; i<s; i++) {
93981d12688SDebojyoti Ghosh       bt[i] += h * Bt[i*pinterp+j] * tt;
94056dcabbaSDebojyoti Ghosh       b[i]  += h * B[i*pinterp+j] * tt;
94156dcabbaSDebojyoti Ghosh     }
94256dcabbaSDebojyoti Ghosh   }
94396400bd6SLisandro Dalcin   if (!ark->Y_prev) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Stages from previous step have not been stored");
94456dcabbaSDebojyoti Ghosh   ierr = VecCopy(ark->Y_prev[0],X);CHKERRQ(ierr);
94556dcabbaSDebojyoti Ghosh   ierr = VecMAXPY(X,s,bt,ark->YdotI_prev);CHKERRQ(ierr);
94656dcabbaSDebojyoti Ghosh   ierr = VecMAXPY(X,s,b,ark->YdotRHS_prev);CHKERRQ(ierr);
94756dcabbaSDebojyoti Ghosh   ierr = PetscFree2(bt,b);CHKERRQ(ierr);
94856dcabbaSDebojyoti Ghosh   PetscFunctionReturn(0);
94956dcabbaSDebojyoti Ghosh }
95056dcabbaSDebojyoti Ghosh 
9518a381b04SJed Brown /*------------------------------------------------------------*/
95296400bd6SLisandro Dalcin 
95396400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauReset(TS ts)
95496400bd6SLisandro Dalcin {
95596400bd6SLisandro Dalcin   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
95696400bd6SLisandro Dalcin   ARKTableau     tab  = ark->tableau;
95796400bd6SLisandro Dalcin   PetscErrorCode ierr;
95896400bd6SLisandro Dalcin 
95996400bd6SLisandro Dalcin   PetscFunctionBegin;
96096400bd6SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
96196400bd6SLisandro Dalcin   ierr = PetscFree(ark->work);CHKERRQ(ierr);
96296400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->Y);CHKERRQ(ierr);
96396400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->YdotI);CHKERRQ(ierr);
96496400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->YdotRHS);CHKERRQ(ierr);
96596400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->Y_prev);CHKERRQ(ierr);
96696400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->YdotI_prev);CHKERRQ(ierr);
96796400bd6SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr);
96896400bd6SLisandro Dalcin   PetscFunctionReturn(0);
96996400bd6SLisandro Dalcin }
97096400bd6SLisandro Dalcin 
9718a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts)
9728a381b04SJed Brown {
9738a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
9748a381b04SJed Brown   PetscErrorCode ierr;
9758a381b04SJed Brown 
9768a381b04SJed Brown   PetscFunctionBegin;
97796400bd6SLisandro Dalcin   ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr);
9788a381b04SJed Brown   ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr);
979e817cc15SEmil Constantinescu   ierr = VecDestroy(&ark->Ydot0);CHKERRQ(ierr);
9808a381b04SJed Brown   ierr = VecDestroy(&ark->Z);CHKERRQ(ierr);
9818a381b04SJed Brown   PetscFunctionReturn(0);
9828a381b04SJed Brown }
9838a381b04SJed Brown 
984d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
985d5e6173cSPeter Brune {
986d5e6173cSPeter Brune   TS_ARKIMEX     *ax = (TS_ARKIMEX*)ts->data;
987d5e6173cSPeter Brune   PetscErrorCode ierr;
988d5e6173cSPeter Brune 
989d5e6173cSPeter Brune   PetscFunctionBegin;
990d5e6173cSPeter Brune   if (Z) {
991d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
992d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
993d5e6173cSPeter Brune     } else *Z = ax->Z;
994d5e6173cSPeter Brune   }
995d5e6173cSPeter Brune   if (Ydot) {
996d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
997d5e6173cSPeter Brune       ierr = DMGetNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
998d5e6173cSPeter Brune     } else *Ydot = ax->Ydot;
999d5e6173cSPeter Brune   }
1000d5e6173cSPeter Brune   PetscFunctionReturn(0);
1001d5e6173cSPeter Brune }
1002d5e6173cSPeter Brune 
1003d5e6173cSPeter Brune 
1004d5e6173cSPeter Brune static PetscErrorCode TSARKIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot)
1005d5e6173cSPeter Brune {
1006d5e6173cSPeter Brune   PetscErrorCode ierr;
1007d5e6173cSPeter Brune 
1008d5e6173cSPeter Brune   PetscFunctionBegin;
1009d5e6173cSPeter Brune   if (Z) {
1010d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
1011d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Z",Z);CHKERRQ(ierr);
1012d5e6173cSPeter Brune     }
1013d5e6173cSPeter Brune   }
1014d5e6173cSPeter Brune   if (Ydot) {
1015d5e6173cSPeter Brune     if (dm && dm != ts->dm) {
1016d5e6173cSPeter Brune       ierr = DMRestoreNamedGlobalVector(dm,"TSARKIMEX_Ydot",Ydot);CHKERRQ(ierr);
1017d5e6173cSPeter Brune     }
1018d5e6173cSPeter Brune   }
1019d5e6173cSPeter Brune   PetscFunctionReturn(0);
1020d5e6173cSPeter Brune }
1021d5e6173cSPeter Brune 
10228a381b04SJed Brown /*
10238a381b04SJed Brown   This defines the nonlinear equation that is to be solved with SNES
10248a381b04SJed Brown   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
10258a381b04SJed Brown */
10268a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts)
10278a381b04SJed Brown {
10288a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1029d5e6173cSPeter Brune   DM             dm,dmsave;
1030d5e6173cSPeter Brune   Vec            Z,Ydot;
1031b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10328a381b04SJed Brown   PetscErrorCode ierr;
10338a381b04SJed Brown 
10348a381b04SJed Brown   PetscFunctionBegin;
1035d5e6173cSPeter Brune   ierr   = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1036d5e6173cSPeter Brune   ierr   = TSARKIMEXGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
1037b296d7d5SJed Brown   ierr   = VecAXPBYPCZ(Ydot,-shift,shift,0,Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */
1038d5e6173cSPeter Brune   dmsave = ts->dm;
1039d5e6173cSPeter Brune   ts->dm = dm;
1040740132f1SEmil Constantinescu 
1041d5e6173cSPeter Brune   ierr = TSComputeIFunction(ts,ark->stage_time,X,Ydot,F,ark->imex);CHKERRQ(ierr);
1042e817cc15SEmil Constantinescu 
1043d5e6173cSPeter Brune   ts->dm = dmsave;
1044d5e6173cSPeter Brune   ierr   = TSARKIMEXRestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr);
10458a381b04SJed Brown   PetscFunctionReturn(0);
10468a381b04SJed Brown }
10478a381b04SJed Brown 
1048d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts)
10498a381b04SJed Brown {
10508a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
1051d5e6173cSPeter Brune   DM             dm,dmsave;
1052d5e6173cSPeter Brune   Vec            Ydot;
1053b296d7d5SJed Brown   PetscReal      shift = ark->scoeff / ts->time_step;
10548a381b04SJed Brown   PetscErrorCode ierr;
10558a381b04SJed Brown 
10568a381b04SJed Brown   PetscFunctionBegin;
1057d5e6173cSPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
10580298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
10598a381b04SJed Brown   /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */
1060d5e6173cSPeter Brune   dmsave = ts->dm;
1061d5e6173cSPeter Brune   ts->dm = dm;
1062740132f1SEmil Constantinescu 
1063d1e9a80fSBarry Smith   ierr = TSComputeIJacobian(ts,ark->stage_time,X,Ydot,shift,A,B,ark->imex);CHKERRQ(ierr);
1064740132f1SEmil Constantinescu 
1065d5e6173cSPeter Brune   ts->dm = dmsave;
10660298fd71SBarry Smith   ierr   = TSARKIMEXRestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr);
1067d5e6173cSPeter Brune   PetscFunctionReturn(0);
1068d5e6173cSPeter Brune }
1069d5e6173cSPeter Brune 
1070d5e6173cSPeter Brune static PetscErrorCode DMCoarsenHook_TSARKIMEX(DM fine,DM coarse,void *ctx)
1071d5e6173cSPeter Brune {
1072d5e6173cSPeter Brune   PetscFunctionBegin;
1073d5e6173cSPeter Brune   PetscFunctionReturn(0);
1074d5e6173cSPeter Brune }
1075d5e6173cSPeter Brune 
1076d5e6173cSPeter Brune static PetscErrorCode DMRestrictHook_TSARKIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
1077d5e6173cSPeter Brune {
1078d5e6173cSPeter Brune   TS             ts = (TS)ctx;
1079d5e6173cSPeter Brune   PetscErrorCode ierr;
1080d5e6173cSPeter Brune   Vec            Z,Z_c;
1081d5e6173cSPeter Brune 
1082d5e6173cSPeter Brune   PetscFunctionBegin;
10830298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10840298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
1085d5e6173cSPeter Brune   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
1086d5e6173cSPeter Brune   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
10870298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,fine,&Z,NULL);CHKERRQ(ierr);
10880298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,coarse,&Z_c,NULL);CHKERRQ(ierr);
10898a381b04SJed Brown   PetscFunctionReturn(0);
10908a381b04SJed Brown }
10918a381b04SJed Brown 
1092cdb298fcSPeter Brune 
1093cdb298fcSPeter Brune static PetscErrorCode DMSubDomainHook_TSARKIMEX(DM dm,DM subdm,void *ctx)
1094cdb298fcSPeter Brune {
1095cdb298fcSPeter Brune   PetscFunctionBegin;
1096cdb298fcSPeter Brune   PetscFunctionReturn(0);
1097cdb298fcSPeter Brune }
1098cdb298fcSPeter Brune 
1099cdb298fcSPeter Brune static PetscErrorCode DMSubDomainRestrictHook_TSARKIMEX(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
1100cdb298fcSPeter Brune {
1101cdb298fcSPeter Brune   TS             ts = (TS)ctx;
1102cdb298fcSPeter Brune   PetscErrorCode ierr;
1103cdb298fcSPeter Brune   Vec            Z,Z_c;
1104cdb298fcSPeter Brune 
1105cdb298fcSPeter Brune   PetscFunctionBegin;
11060298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
11070298fd71SBarry Smith   ierr = TSARKIMEXGetVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1108cdb298fcSPeter Brune 
1109cdb298fcSPeter Brune   ierr = VecScatterBegin(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1110cdb298fcSPeter Brune   ierr = VecScatterEnd(gscat,Z,Z_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1111cdb298fcSPeter Brune 
11120298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,dm,&Z,NULL);CHKERRQ(ierr);
11130298fd71SBarry Smith   ierr = TSARKIMEXRestoreVecs(ts,subdm,&Z_c,NULL);CHKERRQ(ierr);
1114cdb298fcSPeter Brune   PetscFunctionReturn(0);
1115cdb298fcSPeter Brune }
1116cdb298fcSPeter Brune 
111796400bd6SLisandro Dalcin static PetscErrorCode TSARKIMEXTableauSetUp(TS ts)
111896400bd6SLisandro Dalcin {
111996400bd6SLisandro Dalcin   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
112096400bd6SLisandro Dalcin   ARKTableau     tab  = ark->tableau;
112196400bd6SLisandro Dalcin   PetscErrorCode ierr;
112296400bd6SLisandro Dalcin 
112396400bd6SLisandro Dalcin   PetscFunctionBegin;
112496400bd6SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&ark->work);CHKERRQ(ierr);
112596400bd6SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y);CHKERRQ(ierr);
112696400bd6SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI);CHKERRQ(ierr);
112796400bd6SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS);CHKERRQ(ierr);
112896400bd6SLisandro Dalcin   if (ark->extrapolate) {
112996400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->Y_prev);CHKERRQ(ierr);
113096400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotI_prev);CHKERRQ(ierr);
113196400bd6SLisandro Dalcin     ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&ark->YdotRHS_prev);CHKERRQ(ierr);
113296400bd6SLisandro Dalcin   }
113396400bd6SLisandro Dalcin   PetscFunctionReturn(0);
113496400bd6SLisandro Dalcin }
113596400bd6SLisandro Dalcin 
11368a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts)
11378a381b04SJed Brown {
11388a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11398a381b04SJed Brown   PetscErrorCode ierr;
1140d5e6173cSPeter Brune   DM             dm;
114196400bd6SLisandro Dalcin   SNES           snes;
1142f9c1d6abSBarry Smith 
11438a381b04SJed Brown   PetscFunctionBegin;
114496400bd6SLisandro Dalcin   ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr);
11458a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr);
1146e817cc15SEmil Constantinescu   ierr = VecDuplicate(ts->vec_sol,&ark->Ydot0);CHKERRQ(ierr);
11478a381b04SJed Brown   ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr);
1148d5e6173cSPeter Brune   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1149d5e6173cSPeter Brune   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1150cdb298fcSPeter Brune   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
115196400bd6SLisandro Dalcin   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
11528a381b04SJed Brown   PetscFunctionReturn(0);
11538a381b04SJed Brown }
11548a381b04SJed Brown /*------------------------------------------------------------*/
11558a381b04SJed Brown 
11564416b707SBarry Smith static PetscErrorCode TSSetFromOptions_ARKIMEX(PetscOptionItems *PetscOptionsObject,TS ts)
11578a381b04SJed Brown {
11584cc180ffSJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
11598a381b04SJed Brown   PetscErrorCode ierr;
11608a381b04SJed Brown 
11618a381b04SJed Brown   PetscFunctionBegin;
1162e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ARKIMEX ODE solver options");CHKERRQ(ierr);
11638a381b04SJed Brown   {
11648a381b04SJed Brown     ARKTableauLink link;
11658a381b04SJed Brown     PetscInt       count,choice;
11668a381b04SJed Brown     PetscBool      flg;
11678a381b04SJed Brown     const char     **namelist;
11688a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) ;
1169785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
11708a381b04SJed Brown     for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
117196400bd6SLisandro Dalcin     ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,ark->tableau->name,&choice,&flg);CHKERRQ(ierr);
117296400bd6SLisandro Dalcin     if (flg) {ierr = TSARKIMEXSetType(ts,namelist[choice]);CHKERRQ(ierr);}
11738a381b04SJed Brown     ierr = PetscFree(namelist);CHKERRQ(ierr);
117496400bd6SLisandro Dalcin 
11754cc180ffSJed Brown     flg  = (PetscBool) !ark->imex;
11760298fd71SBarry Smith     ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,NULL);CHKERRQ(ierr);
11774cc180ffSJed Brown     ark->imex = (PetscBool) !flg;
117803842d09SLisandro Dalcin     ierr = PetscOptionsBool("-ts_arkimex_initial_guess_extrapolate","Extrapolate the initial guess for the stage solution from stage values of the previous time step","",ark->extrapolate,&ark->extrapolate,NULL);CHKERRQ(ierr);
11798a381b04SJed Brown   }
11808a381b04SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
11818a381b04SJed Brown   PetscFunctionReturn(0);
11828a381b04SJed Brown }
11838a381b04SJed Brown 
11848a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
11858a381b04SJed Brown {
1186257d2499SJed Brown   PetscErrorCode ierr;
1187f1d86077SJed Brown   PetscInt       i;
1188f1d86077SJed Brown   size_t         left,count;
11898a381b04SJed Brown   char           *p;
11908a381b04SJed Brown 
11918a381b04SJed Brown   PetscFunctionBegin;
1192f1d86077SJed Brown   for (i=0,p=buf,left=len; i<n; i++) {
1193da649d3eSBarry Smith     ierr = PetscSNPrintfCount(p,left,fmt,&count,(double)x[i]);CHKERRQ(ierr);
11948a381b04SJed Brown     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
11958a381b04SJed Brown     left -= count;
11968a381b04SJed Brown     p    += count;
11978a381b04SJed Brown     *p++  = ' ';
11988a381b04SJed Brown   }
11998a381b04SJed Brown   p[i ? 0 : -1] = 0;
12008a381b04SJed Brown   PetscFunctionReturn(0);
12018a381b04SJed Brown }
12028a381b04SJed Brown 
12038a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer)
12048a381b04SJed Brown {
12058a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
12068a381b04SJed Brown   PetscBool      iascii;
12078a381b04SJed Brown   PetscErrorCode ierr;
12088a381b04SJed Brown 
12098a381b04SJed Brown   PetscFunctionBegin;
1210251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12118a381b04SJed Brown   if (iascii) {
12129c334d8fSLisandro Dalcin     ARKTableau    tab = ark->tableau;
121319fd82e9SBarry Smith     TSARKIMEXType arktype;
12148a381b04SJed Brown     char          buf[512];
12158a381b04SJed Brown     ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr);
12168a381b04SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  ARK IMEX %s\n",arktype);CHKERRQ(ierr);
12178caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr);
121831f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Stiff abscissa       ct = %s\n",buf);CHKERRQ(ierr);
12198caf3d72SBarry Smith     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
1220e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",tab->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr);
1221e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"Explicit first stage: %s\n",tab->explicit_first_stage ? "yes" : "no");CHKERRQ(ierr);
1222e817cc15SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"FSAL property: %s\n",tab->FSAL_implicit ? "yes" : "no");CHKERRQ(ierr);
122331f6fcc0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  Nonstiff abscissa     c = %s\n",buf);CHKERRQ(ierr);
12248a381b04SJed Brown   }
12258a381b04SJed Brown   PetscFunctionReturn(0);
12268a381b04SJed Brown }
12278a381b04SJed Brown 
1228f2c2a1b9SBarry Smith static PetscErrorCode TSLoad_ARKIMEX(TS ts,PetscViewer viewer)
1229f2c2a1b9SBarry Smith {
1230f2c2a1b9SBarry Smith   PetscErrorCode ierr;
1231f2c2a1b9SBarry Smith   SNES           snes;
12329c334d8fSLisandro Dalcin   TSAdapt        adapt;
1233f2c2a1b9SBarry Smith 
1234f2c2a1b9SBarry Smith   PetscFunctionBegin;
12359c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
12369c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
1237f2c2a1b9SBarry Smith   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
1238f2c2a1b9SBarry Smith   ierr = SNESLoad(snes,viewer);CHKERRQ(ierr);
1239ad6bc421SBarry Smith   /* function and Jacobian context for SNES when used with TS is always ts object */
12400298fd71SBarry Smith   ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr);
12410298fd71SBarry Smith   ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr);
1242f2c2a1b9SBarry Smith   PetscFunctionReturn(0);
1243f2c2a1b9SBarry Smith }
1244f2c2a1b9SBarry Smith 
12458a381b04SJed Brown /*@C
12468a381b04SJed Brown   TSARKIMEXSetType - Set the type of ARK IMEX scheme
12478a381b04SJed Brown 
12488a381b04SJed Brown   Logically collective
12498a381b04SJed Brown 
12508a381b04SJed Brown   Input Parameter:
12518a381b04SJed Brown +  ts - timestepping context
12528a381b04SJed Brown -  arktype - type of ARK-IMEX scheme
12538a381b04SJed Brown 
12549bd3e852SBarry Smith   Options Database:
12559bd3e852SBarry Smith .  -ts_arkimex_type <1bee,a2,l2,ars122,2c,2d,2e,prssp2,3,bpr3,ars443,4,5>
12569bd3e852SBarry Smith 
12578a381b04SJed Brown   Level: intermediate
12588a381b04SJed Brown 
12599bd3e852SBarry Smith .seealso: TSARKIMEXGetType(), TSARKIMEX, TSARKIMEXType, TSARKIMEX1BEE, TSARKIMEXA2, TSARKIMEXL2, TSARKIMEXARS122, TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEXPRSSP2,
12609bd3e852SBarry Smith           TSARKIMEX3, TSARKIMEXBPR3, TSARKIMEXARS443, TSARKIMEX4, TSARKIMEX5
12618a381b04SJed Brown @*/
126219fd82e9SBarry Smith PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType arktype)
12638a381b04SJed Brown {
12648a381b04SJed Brown   PetscErrorCode ierr;
12658a381b04SJed Brown 
12668a381b04SJed Brown   PetscFunctionBegin;
12678a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1268b92453a8SLisandro Dalcin   PetscValidCharPointer(arktype,2);
126919fd82e9SBarry Smith   ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,TSARKIMEXType),(ts,arktype));CHKERRQ(ierr);
12708a381b04SJed Brown   PetscFunctionReturn(0);
12718a381b04SJed Brown }
12728a381b04SJed Brown 
12738a381b04SJed Brown /*@C
12748a381b04SJed Brown   TSARKIMEXGetType - Get the type of ARK IMEX scheme
12758a381b04SJed Brown 
12768a381b04SJed Brown   Logically collective
12778a381b04SJed Brown 
12788a381b04SJed Brown   Input Parameter:
12798a381b04SJed Brown .  ts - timestepping context
12808a381b04SJed Brown 
12818a381b04SJed Brown   Output Parameter:
12828a381b04SJed Brown .  arktype - type of ARK-IMEX scheme
12838a381b04SJed Brown 
12848a381b04SJed Brown   Level: intermediate
12858a381b04SJed Brown 
12868a381b04SJed Brown .seealso: TSARKIMEXGetType()
12878a381b04SJed Brown @*/
128819fd82e9SBarry Smith PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType *arktype)
12898a381b04SJed Brown {
12908a381b04SJed Brown   PetscErrorCode ierr;
12918a381b04SJed Brown 
12928a381b04SJed Brown   PetscFunctionBegin;
12938a381b04SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
129419fd82e9SBarry Smith   ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr);
12958a381b04SJed Brown   PetscFunctionReturn(0);
12968a381b04SJed Brown }
12978a381b04SJed Brown 
129816353aafSBarry Smith /*@
12994cc180ffSJed Brown   TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly
13004cc180ffSJed Brown 
13014cc180ffSJed Brown   Logically collective
13024cc180ffSJed Brown 
13034cc180ffSJed Brown   Input Parameter:
13044cc180ffSJed Brown +  ts - timestepping context
13054cc180ffSJed Brown -  flg - PETSC_TRUE for fully implicit
13064cc180ffSJed Brown 
13074cc180ffSJed Brown   Level: intermediate
13084cc180ffSJed Brown 
13094cc180ffSJed Brown .seealso: TSARKIMEXGetType()
13104cc180ffSJed Brown @*/
13114cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg)
13124cc180ffSJed Brown {
13134cc180ffSJed Brown   PetscErrorCode ierr;
13144cc180ffSJed Brown 
13154cc180ffSJed Brown   PetscFunctionBegin;
13164cc180ffSJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
13174cc180ffSJed Brown   ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
13184cc180ffSJed Brown   PetscFunctionReturn(0);
13194cc180ffSJed Brown }
13204cc180ffSJed Brown 
1321e0877f53SBarry Smith static PetscErrorCode  TSARKIMEXGetType_ARKIMEX(TS ts,TSARKIMEXType *arktype)
13228a381b04SJed Brown {
13238a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
13248a381b04SJed Brown 
13258a381b04SJed Brown   PetscFunctionBegin;
13268a381b04SJed Brown   *arktype = ark->tableau->name;
13278a381b04SJed Brown   PetscFunctionReturn(0);
13288a381b04SJed Brown }
1329e0877f53SBarry Smith static PetscErrorCode  TSARKIMEXSetType_ARKIMEX(TS ts,TSARKIMEXType arktype)
13308a381b04SJed Brown {
13318a381b04SJed Brown   TS_ARKIMEX     *ark = (TS_ARKIMEX*)ts->data;
13328a381b04SJed Brown   PetscErrorCode ierr;
13338a381b04SJed Brown   PetscBool      match;
13348a381b04SJed Brown   ARKTableauLink link;
13358a381b04SJed Brown 
13368a381b04SJed Brown   PetscFunctionBegin;
13378a381b04SJed Brown   if (ark->tableau) {
13388a381b04SJed Brown     ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr);
13398a381b04SJed Brown     if (match) PetscFunctionReturn(0);
13408a381b04SJed Brown   }
13418a381b04SJed Brown   for (link = ARKTableauList; link; link=link->next) {
13428a381b04SJed Brown     ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr);
13438a381b04SJed Brown     if (match) {
134496400bd6SLisandro Dalcin       if (ts->setupcalled) {ierr = TSARKIMEXTableauReset(ts);CHKERRQ(ierr);}
13458a381b04SJed Brown       ark->tableau = &link->tab;
134696400bd6SLisandro Dalcin       if (ts->setupcalled) {ierr = TSARKIMEXTableauSetUp(ts);CHKERRQ(ierr);}
13472ffb9264SLisandro Dalcin       ts->default_adapt_type = ark->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
13488a381b04SJed Brown       PetscFunctionReturn(0);
13498a381b04SJed Brown     }
13508a381b04SJed Brown   }
1351ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype);
13528a381b04SJed Brown   PetscFunctionReturn(0);
13538a381b04SJed Brown }
1354e0877f53SBarry Smith 
1355e0877f53SBarry Smith static PetscErrorCode  TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg)
13564cc180ffSJed Brown {
13574cc180ffSJed Brown   TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data;
13584cc180ffSJed Brown 
13594cc180ffSJed Brown   PetscFunctionBegin;
13604cc180ffSJed Brown   ark->imex = (PetscBool)!flg;
13614cc180ffSJed Brown   PetscFunctionReturn(0);
13624cc180ffSJed Brown }
13638a381b04SJed Brown 
1364*b3a6b972SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts)
1365*b3a6b972SJed Brown {
1366*b3a6b972SJed Brown   PetscErrorCode ierr;
1367*b3a6b972SJed Brown 
1368*b3a6b972SJed Brown   PetscFunctionBegin;
1369*b3a6b972SJed Brown   ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr);
1370*b3a6b972SJed Brown   if (ts->dm) {
1371*b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSARKIMEX,DMRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1372*b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSARKIMEX,DMSubDomainRestrictHook_TSARKIMEX,ts);CHKERRQ(ierr);
1373*b3a6b972SJed Brown   }
1374*b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1375*b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",NULL);CHKERRQ(ierr);
1376*b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",NULL);CHKERRQ(ierr);
1377*b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",NULL);CHKERRQ(ierr);
1378*b3a6b972SJed Brown   PetscFunctionReturn(0);
1379*b3a6b972SJed Brown }
1380*b3a6b972SJed Brown 
13818a381b04SJed Brown /* ------------------------------------------------------------ */
13828a381b04SJed Brown /*MC
1383a4386c9eSJed Brown       TSARKIMEX - ODE and DAE solver using Additive Runge-Kutta IMEX schemes
13848a381b04SJed Brown 
1385fca742c7SJed Brown   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1386fca742c7SJed Brown   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1387fca742c7SJed Brown   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
1388fca742c7SJed Brown 
1389fca742c7SJed Brown   Notes:
1390a4386c9eSJed Brown   The default is TSARKIMEX3, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type
1391c8058688SBarry Smith 
13925eca1a21SEmil Constantinescu   If the equation is implicit or a DAE, then TSSetEquationType() needs to be set accordingly. Refer to the manual for further information.
13935eca1a21SEmil Constantinescu 
1394a4386c9eSJed 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).
1395fca742c7SJed Brown 
1396d0685a90SJed Brown   Consider trying TSROSW if the stiff part is linear or weakly nonlinear.
1397d0685a90SJed Brown 
13988a381b04SJed Brown   Level: beginner
13998a381b04SJed Brown 
1400d0685a90SJed Brown .seealso:  TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX1BEE,
1401d0685a90SJed Brown            TSARKIMEX2C, TSARKIMEX2D, TSARKIMEX2E, TSARKIMEX3, TSARKIMEXL2, TSARKIMEXA2, TSARKIMEXARS122,
1402d0685a90SJed Brown            TSARKIMEX4, TSARKIMEX5, TSARKIMEXPRSSP2, TSARKIMEXARS443, TSARKIMEXBPR3, TSARKIMEXType, TSARKIMEXRegister()
14038a381b04SJed Brown 
14048a381b04SJed Brown M*/
14058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_ARKIMEX(TS ts)
14068a381b04SJed Brown {
14078a381b04SJed Brown   TS_ARKIMEX     *th;
14088a381b04SJed Brown   PetscErrorCode ierr;
14098a381b04SJed Brown 
14108a381b04SJed Brown   PetscFunctionBegin;
1411607a6623SBarry Smith   ierr = TSARKIMEXInitializePackage();CHKERRQ(ierr);
14128a381b04SJed Brown 
14138a381b04SJed Brown   ts->ops->reset          = TSReset_ARKIMEX;
14148a381b04SJed Brown   ts->ops->destroy        = TSDestroy_ARKIMEX;
14158a381b04SJed Brown   ts->ops->view           = TSView_ARKIMEX;
1416f2c2a1b9SBarry Smith   ts->ops->load           = TSLoad_ARKIMEX;
14178a381b04SJed Brown   ts->ops->setup          = TSSetUp_ARKIMEX;
14188a381b04SJed Brown   ts->ops->step           = TSStep_ARKIMEX;
1419cd652676SJed Brown   ts->ops->interpolate    = TSInterpolate_ARKIMEX;
1420108c343cSJed Brown   ts->ops->evaluatestep   = TSEvaluateStep_ARKIMEX;
142124655328SShri   ts->ops->rollback       = TSRollBack_ARKIMEX;
14228a381b04SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX;
14238a381b04SJed Brown   ts->ops->snesfunction   = SNESTSFormFunction_ARKIMEX;
14248a381b04SJed Brown   ts->ops->snesjacobian   = SNESTSFormJacobian_ARKIMEX;
14258a381b04SJed Brown 
1426efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
1427efd4aadfSBarry Smith 
1428b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
14298a381b04SJed Brown   ts->data = (void*)th;
14304cc180ffSJed Brown   th->imex = PETSC_TRUE;
14318a381b04SJed Brown 
1432bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXGetType_C",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr);
1433bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetType_C",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr);
1434bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr);
143596400bd6SLisandro Dalcin 
143696400bd6SLisandro Dalcin   ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);
14378a381b04SJed Brown   PetscFunctionReturn(0);
14388a381b04SJed Brown }
1439