xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
1e4dd521cSBarry Smith /*
22e698f8bSDebojyoti Ghosh   Code for time stepping with the Runge-Kutta method
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
5f68a32c8SEmil Constantinescu   The general system is written as
6f68a32c8SEmil Constantinescu 
72e698f8bSDebojyoti Ghosh   Udot = F(t,U)
8f68a32c8SEmil Constantinescu 
9e4dd521cSBarry Smith */
10af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
11f68a32c8SEmil Constantinescu #include <petscdm.h>
12f68a32c8SEmil Constantinescu 
13484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
14f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
15f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
16f68a32c8SEmil Constantinescu 
17f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau;
18f68a32c8SEmil Constantinescu struct _RKTableau {
19f68a32c8SEmil Constantinescu   char      *name;
20d760c35bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i              */
21f68a32c8SEmil Constantinescu   PetscInt   s;                   /* Number of stages                                           */
223a8a9803SLisandro Dalcin   PetscInt   p;                   /* Interpolation order                                        */
23d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
24f68a32c8SEmil Constantinescu   PetscReal *A,*b,*c;             /* Tableau                                                    */
25f68a32c8SEmil Constantinescu   PetscReal *bembed;              /* Embedded formula of order one less (order-1)               */
26f68a32c8SEmil Constantinescu   PetscReal *binterp;             /* Dense output formula                                       */
27f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
28f68a32c8SEmil Constantinescu };
29f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
30f68a32c8SEmil Constantinescu struct _RKTableauLink {
31f68a32c8SEmil Constantinescu   struct _RKTableau tab;
32f68a32c8SEmil Constantinescu   RKTableauLink     next;
33f68a32c8SEmil Constantinescu };
34f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
35e4dd521cSBarry Smith 
36e4dd521cSBarry Smith typedef struct {
37f68a32c8SEmil Constantinescu   RKTableau    tableau;
38f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step */
39f68a32c8SEmil Constantinescu   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
40ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage */
41ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage */
42b18ea86cSHong Zhang   Vec          *VecSensiTemp;    /* Vector to be timed with Jacobian transpose */
430f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events */
44f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
45f68a32c8SEmil Constantinescu   PetscReal    stage_time;
46f68a32c8SEmil Constantinescu   TSStepStatus status;
470f7a1166SHong Zhang   PetscReal    ptime;
480f7a1166SHong Zhang   PetscReal    time_step;
495f70b478SJed Brown } TS_RK;
50e4dd521cSBarry Smith 
51f68a32c8SEmil Constantinescu /*MC
52f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
53262ff7bbSBarry Smith 
54f68a32c8SEmil Constantinescu      This method has one stage.
55f68a32c8SEmil Constantinescu 
56f68a32c8SEmil Constantinescu      Level: advanced
57f68a32c8SEmil Constantinescu 
58f68a32c8SEmil Constantinescu .seealso: TSRK
59f68a32c8SEmil Constantinescu M*/
60f68a32c8SEmil Constantinescu /*MC
612109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
62f68a32c8SEmil Constantinescu 
63f68a32c8SEmil Constantinescu      This method has two stages.
64f68a32c8SEmil Constantinescu 
65f68a32c8SEmil Constantinescu      Level: advanced
66f68a32c8SEmil Constantinescu 
67f68a32c8SEmil Constantinescu .seealso: TSRK
68f68a32c8SEmil Constantinescu M*/
69f68a32c8SEmil Constantinescu /*MC
70f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
71f68a32c8SEmil Constantinescu 
72f68a32c8SEmil Constantinescu      This method has three stages.
73f68a32c8SEmil Constantinescu 
74f68a32c8SEmil Constantinescu      Level: advanced
75f68a32c8SEmil Constantinescu 
76f68a32c8SEmil Constantinescu .seealso: TSRK
77f68a32c8SEmil Constantinescu M*/
78f68a32c8SEmil Constantinescu /*MC
792109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
802109b73fSDebojyoti Ghosh 
81d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
822109b73fSDebojyoti Ghosh 
832109b73fSDebojyoti Ghosh      Level: advanced
842109b73fSDebojyoti Ghosh 
8598164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
86d1905564SLisandro Dalcin 
872109b73fSDebojyoti Ghosh .seealso: TSRK
882109b73fSDebojyoti Ghosh M*/
892109b73fSDebojyoti Ghosh /*MC
90f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
91f68a32c8SEmil Constantinescu 
922e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
93f68a32c8SEmil Constantinescu 
94f68a32c8SEmil Constantinescu      Level: advanced
95f68a32c8SEmil Constantinescu 
96f68a32c8SEmil Constantinescu .seealso: TSRK
97f68a32c8SEmil Constantinescu M*/
98f68a32c8SEmil Constantinescu /*MC
992e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
100f68a32c8SEmil Constantinescu 
101f68a32c8SEmil Constantinescu      This method has six stages.
102f68a32c8SEmil Constantinescu 
103f68a32c8SEmil Constantinescu      Level: advanced
104f68a32c8SEmil Constantinescu 
105f68a32c8SEmil Constantinescu .seealso: TSRK
106f68a32c8SEmil Constantinescu M*/
1072109b73fSDebojyoti Ghosh /*MC
1082e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1092109b73fSDebojyoti Ghosh 
110d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1112109b73fSDebojyoti Ghosh 
1122109b73fSDebojyoti Ghosh      Level: advanced
1132109b73fSDebojyoti Ghosh 
11498164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
115d1905564SLisandro Dalcin 
1162109b73fSDebojyoti Ghosh .seealso: TSRK
1172109b73fSDebojyoti Ghosh M*/
11805e23783SLisandro Dalcin /*MC
11905e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
12005e23783SLisandro Dalcin 
12105e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
12205e23783SLisandro Dalcin 
12305e23783SLisandro Dalcin      Level: advanced
12405e23783SLisandro Dalcin 
12598164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
12605e23783SLisandro Dalcin 
12705e23783SLisandro Dalcin .seealso: TSRK
12805e23783SLisandro Dalcin M*/
129262ff7bbSBarry Smith 
130f68a32c8SEmil Constantinescu /*@C
131f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
132262ff7bbSBarry Smith 
133f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
134262ff7bbSBarry Smith 
135f68a32c8SEmil Constantinescu   Level: advanced
136262ff7bbSBarry Smith 
137f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
138262ff7bbSBarry Smith 
139f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
140262ff7bbSBarry Smith @*/
141f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
142262ff7bbSBarry Smith {
1434ac538c5SBarry Smith   PetscErrorCode ierr;
144262ff7bbSBarry Smith 
145262ff7bbSBarry Smith   PetscFunctionBegin;
146f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
147f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
148e4dd521cSBarry Smith 
1494e82626cSLisandro Dalcin #define RC PetscRealConstant
150e4dd521cSBarry Smith   {
151f68a32c8SEmil Constantinescu     const PetscReal
1524e82626cSLisandro Dalcin       A[1][1] = {{0}},
1534e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1543a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
155e4dd521cSBarry Smith   }
156db046809SBarry Smith   {
157f68a32c8SEmil Constantinescu     const PetscReal
1584e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1594e82626cSLisandro Dalcin                    {RC(1.0),0}},
1604e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1614e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1623a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
163db046809SBarry Smith   }
164f68a32c8SEmil Constantinescu   {
165f68a32c8SEmil Constantinescu     const PetscReal
16617f6384fSBarry Smith       A[3][3] = {{0,0,0},
1674e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
1684e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
1694e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
1703a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
171fdefd5e5SDebojyoti Ghosh   }
172fdefd5e5SDebojyoti Ghosh   {
173fdefd5e5SDebojyoti Ghosh     const PetscReal
17417f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
1754e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
1764e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
1774e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
1784e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
1794e82626cSLisandro Dalcin       bembed[4] =  {RC(7.0)/RC(24.0),RC(1.0)/RC(4.0),RC(1.0)/RC(3.0),RC(1.0)/RC(8.0)};
1803a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
181db046809SBarry Smith   }
182f68a32c8SEmil Constantinescu   {
183f68a32c8SEmil Constantinescu     const PetscReal
184f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
1854e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
1864e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
1874e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
1884e82626cSLisandro Dalcin       b[4]    =  {RC(1.0)/RC(6.0),RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),RC(1.0)/RC(6.0)};
1893a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
190db046809SBarry Smith   }
191f68a32c8SEmil Constantinescu   {
192f68a32c8SEmil Constantinescu     const PetscReal
19317f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
1944e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
1954e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
1964e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
1974e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
1984e82626cSLisandro Dalcin                    {RC(-8.0)/RC(27.0),RC(2.0),RC(-3544.0)/RC(2565.0),RC(1859.0)/RC(4104.0),RC(-11.0)/RC(40.0),0}},
1994e82626cSLisandro Dalcin       b[6]      =  {RC(16.0)/RC(135.0),0,RC(6656.0)/RC(12825.0),RC(28561.0)/RC(56430.0),RC(-9.0)/RC(50.0),RC(2.0)/RC(55.0)},
2004e82626cSLisandro Dalcin       bembed[6] =  {RC(25.0)/RC(216.0),0,RC(1408.0)/RC(2565.0),RC(2197.0)/RC(4104.0),RC(-1.0)/RC(5.0),0};
2013a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
202fdefd5e5SDebojyoti Ghosh   }
203fdefd5e5SDebojyoti Ghosh   {
204fdefd5e5SDebojyoti Ghosh     const PetscReal
20517f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2064e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2074e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2084e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2094e82626cSLisandro Dalcin                    {RC(19372.0)/RC(6561.0),RC(-25360.0)/RC(2187.0),RC(64448.0)/RC(6561.0),RC(-212.0)/RC(729.0),0,0,0},
2104e82626cSLisandro Dalcin                    {RC(9017.0)/RC(3168.0),RC(-355.0)/RC(33.0),RC(46732.0)/RC(5247.0),RC(49.0)/RC(176.0),RC(-5103.0)/RC(18656.0),0,0},
2114e82626cSLisandro Dalcin                    {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0}},
2124e82626cSLisandro Dalcin       b[7]      =  {RC(35.0)/RC(384.0),0,RC(500.0)/RC(1113.0),RC(125.0)/RC(192.0),RC(-2187.0)/RC(6784.0),RC(11.0)/RC(84.0),0},
2134e82626cSLisandro Dalcin       bembed[7] =  {RC(5179.0)/RC(57600.0),0,RC(7571.0)/RC(16695.0),RC(393.0)/RC(640.0),RC(-92097.0)/RC(339200.0),RC(187.0)/RC(2100.0),RC(1.0)/RC(40.0)};
2143a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
215f68a32c8SEmil Constantinescu   }
21605e23783SLisandro Dalcin   {
21705e23783SLisandro Dalcin     const PetscReal
21817f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2194e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2204e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2214e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2224e82626cSLisandro Dalcin                    {RC(68.0)/RC(297.0),RC(-4.0)/RC(11.0),RC(42.0)/RC(143.0),RC(1960.0)/RC(3861.0),0,0,0,0},
2234e82626cSLisandro Dalcin                    {RC(597.0)/RC(22528.0),RC(81.0)/RC(352.0),RC(63099.0)/RC(585728.0),RC(58653.0)/RC(366080.0),RC(4617.0)/RC(20480.0),0,0,0},
2244e82626cSLisandro Dalcin                    {RC(174197.0)/RC(959244.0),RC(-30942.0)/RC(79937.0),RC(8152137.0)/RC(19744439.0),RC(666106.0)/RC(1039181.0),RC(-29421.0)/RC(29068.0),RC(482048.0)/RC(414219.0),0,0},
2254e82626cSLisandro Dalcin                    {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0}},
2264e82626cSLisandro Dalcin       b[8]      =  {RC(587.0)/RC(8064.0),0,RC(4440339.0)/RC(15491840.0),RC(24353.0)/RC(124800.0),RC(387.0)/RC(44800.0),RC(2152.0)/RC(5985.0),RC(7267.0)/RC(94080.0),0},
2274e82626cSLisandro Dalcin       bembed[8] =  {RC(2479.0)/RC(34992.0),0,RC(123.0)/RC(416.0),RC(612941.0)/RC(3411720.0),RC(43.0)/RC(1440.0),RC(2272.0)/RC(6561.0),RC(79937.0)/RC(1113912.0),RC(3293.0)/RC(556956.0)};
22805e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
22905e23783SLisandro Dalcin   }
2304e82626cSLisandro Dalcin #undef RC
231db046809SBarry Smith   PetscFunctionReturn(0);
232db046809SBarry Smith }
233db046809SBarry Smith 
234f68a32c8SEmil Constantinescu /*@C
235f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
236f68a32c8SEmil Constantinescu 
237f68a32c8SEmil Constantinescu    Not Collective
238f68a32c8SEmil Constantinescu 
239f68a32c8SEmil Constantinescu    Level: advanced
240f68a32c8SEmil Constantinescu 
241f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
242f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
243f68a32c8SEmil Constantinescu @*/
244f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
245e4dd521cSBarry Smith {
246f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
247f68a32c8SEmil Constantinescu   RKTableauLink link;
248f68a32c8SEmil Constantinescu 
249f68a32c8SEmil Constantinescu   PetscFunctionBegin;
250f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
251f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
252f68a32c8SEmil Constantinescu     RKTableauList = link->next;
253f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
254f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
255f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
256f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
258f68a32c8SEmil Constantinescu   }
259f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
260f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
261f68a32c8SEmil Constantinescu }
262f68a32c8SEmil Constantinescu 
263f68a32c8SEmil Constantinescu /*@C
264f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
265f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
266f68a32c8SEmil Constantinescu   when using static libraries.
267f68a32c8SEmil Constantinescu 
268f68a32c8SEmil Constantinescu   Level: developer
269f68a32c8SEmil Constantinescu 
270f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
271f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
272f68a32c8SEmil Constantinescu @*/
273f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
274f68a32c8SEmil Constantinescu {
275186e87acSLisandro Dalcin   PetscErrorCode ierr;
276e4dd521cSBarry Smith 
277e4dd521cSBarry Smith   PetscFunctionBegin;
278f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
279f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
280f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
281f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
282f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
283f68a32c8SEmil Constantinescu }
284186e87acSLisandro Dalcin 
285f68a32c8SEmil Constantinescu /*@C
286f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
287f68a32c8SEmil Constantinescu   called from PetscFinalize().
288186e87acSLisandro Dalcin 
289f68a32c8SEmil Constantinescu   Level: developer
290f68a32c8SEmil Constantinescu 
291f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
292f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
293f68a32c8SEmil Constantinescu @*/
294f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
295f68a32c8SEmil Constantinescu {
296f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu   PetscFunctionBegin;
299f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
300f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
301f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
302f68a32c8SEmil Constantinescu }
303f68a32c8SEmil Constantinescu 
304f68a32c8SEmil Constantinescu /*@C
305f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
306f68a32c8SEmil Constantinescu 
307f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
308f68a32c8SEmil Constantinescu 
309f68a32c8SEmil Constantinescu    Input Parameters:
310f68a32c8SEmil Constantinescu +  name - identifier for method
311f68a32c8SEmil Constantinescu .  order - approximation order of method
312f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
313f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
314f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
315f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
316f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3173a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3183a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
319f68a32c8SEmil Constantinescu 
320f68a32c8SEmil Constantinescu    Notes:
321f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
322f68a32c8SEmil Constantinescu 
323f68a32c8SEmil Constantinescu    Level: advanced
324f68a32c8SEmil Constantinescu 
325f68a32c8SEmil Constantinescu .keywords: TS, register
326f68a32c8SEmil Constantinescu 
327f68a32c8SEmil Constantinescu .seealso: TSRK
328f68a32c8SEmil Constantinescu @*/
329f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
330f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3313a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
332f68a32c8SEmil Constantinescu {
333f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
334f68a32c8SEmil Constantinescu   RKTableauLink   link;
335f68a32c8SEmil Constantinescu   RKTableau       t;
336f68a32c8SEmil Constantinescu   PetscInt        i,j;
337f68a32c8SEmil Constantinescu 
338f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3393a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3403a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3413a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3423a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3433a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3443a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3453a8a9803SLisandro Dalcin 
34695dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
347f68a32c8SEmil Constantinescu   t = &link->tab;
3483a8a9803SLisandro Dalcin 
349f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
350f68a32c8SEmil Constantinescu   t->order = order;
351f68a32c8SEmil Constantinescu   t->s = s;
352dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
353f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
354f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
355f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
356f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
357f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
358d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
359d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3603a8a9803SLisandro Dalcin 
361f68a32c8SEmil Constantinescu   if (bembed) {
362785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
363f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
364f68a32c8SEmil Constantinescu   }
365f68a32c8SEmil Constantinescu 
3663a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3673a8a9803SLisandro Dalcin   t->p = p;
3683a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3693a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3703a8a9803SLisandro Dalcin 
371f68a32c8SEmil Constantinescu   link->next = RKTableauList;
372f68a32c8SEmil Constantinescu   RKTableauList = link;
373f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
374f68a32c8SEmil Constantinescu }
375f68a32c8SEmil Constantinescu 
376e4dd521cSBarry Smith /*
377f68a32c8SEmil Constantinescu  The step completion formula is
378e4dd521cSBarry Smith 
379f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
380f68a32c8SEmil Constantinescu 
381f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
382f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
383f68a32c8SEmil Constantinescu  We can write
384f68a32c8SEmil Constantinescu 
385f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
386f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
387f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
388f68a32c8SEmil Constantinescu 
389f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
390e4dd521cSBarry Smith */
391f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
392f68a32c8SEmil Constantinescu {
393f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
394f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
395f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
396f68a32c8SEmil Constantinescu   PetscReal      h;
397f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
398f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
399f68a32c8SEmil Constantinescu 
400f68a32c8SEmil Constantinescu   PetscFunctionBegin;
401f68a32c8SEmil Constantinescu   switch (rk->status) {
402f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
403f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
404f68a32c8SEmil Constantinescu     h = ts->time_step; break;
405f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
406be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
407f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
408f68a32c8SEmil Constantinescu   }
409f68a32c8SEmil Constantinescu   if (order == tab->order) {
410f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
411f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
412f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
413f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
414f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
415f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
416f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
417f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
418f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
419f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
420f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
421f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
422f68a32c8SEmil Constantinescu     } else { /* Rollback and re-complete using (be-b) */
423f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
424f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
425f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
4260f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
4270f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4280f7a1166SHong Zhang       }
429f68a32c8SEmil Constantinescu     }
430f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
431f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
432f68a32c8SEmil Constantinescu   }
433f68a32c8SEmil Constantinescu unavailable:
434f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
435a7fac7c2SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%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);
436f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
437f68a32c8SEmil Constantinescu }
438f68a32c8SEmil Constantinescu 
4390f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4400f7a1166SHong Zhang {
4410f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4420f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4430f7a1166SHong Zhang   const PetscInt  s = tab->s;
4440f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4450f7a1166SHong Zhang   Vec             *Y = rk->Y;
4460f7a1166SHong Zhang   PetscInt        i;
4470f7a1166SHong Zhang   PetscErrorCode  ierr;
4480f7a1166SHong Zhang 
4490f7a1166SHong Zhang   PetscFunctionBegin;
4500f7a1166SHong Zhang   /* backup cost integral */
4510f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4520f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4530f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4540f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4550f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4560f7a1166SHong Zhang   }
4570f7a1166SHong Zhang   PetscFunctionReturn(0);
4580f7a1166SHong Zhang }
4590f7a1166SHong Zhang 
4600f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4610f7a1166SHong Zhang {
4620f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4630f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4640f7a1166SHong Zhang   const PetscInt  s = tab->s;
4650f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4660f7a1166SHong Zhang   Vec             *Y = rk->Y;
4670f7a1166SHong Zhang   PetscInt        i;
4680f7a1166SHong Zhang   PetscErrorCode  ierr;
4690f7a1166SHong Zhang 
4700f7a1166SHong Zhang   PetscFunctionBegin;
4710f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4720f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4730f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4740f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4750f7a1166SHong Zhang   }
4760f7a1166SHong Zhang   PetscFunctionReturn(0);
4770f7a1166SHong Zhang }
4780f7a1166SHong Zhang 
479fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
480fecfb714SLisandro Dalcin {
481fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
482fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
483fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
484fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
485fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
486fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
487fecfb714SLisandro Dalcin   PetscInt        j;
488fecfb714SLisandro Dalcin   PetscReal       h;
489fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
490fecfb714SLisandro Dalcin 
491fecfb714SLisandro Dalcin   PetscFunctionBegin;
492fecfb714SLisandro Dalcin   switch (rk->status) {
493fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
494fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
495fecfb714SLisandro Dalcin     h = ts->time_step; break;
496fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
497fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
498fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
499fecfb714SLisandro Dalcin   }
500fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
501fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
502fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
503fecfb714SLisandro Dalcin }
504fecfb714SLisandro Dalcin 
505f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
506f68a32c8SEmil Constantinescu {
507f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
508f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
509f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
510fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
511f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
512f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
513d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
514f68a32c8SEmil Constantinescu   TSAdapt          adapt;
515fecfb714SLisandro Dalcin   PetscInt         i,j;
516be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
517be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
518be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
519f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
520f68a32c8SEmil Constantinescu 
521f68a32c8SEmil Constantinescu   PetscFunctionBegin;
522d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
523d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
524d1905564SLisandro Dalcin 
525f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
526be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
527be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
528f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
529f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5309be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5319be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
532f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
533f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
534f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5359be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
536f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
537be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
538fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5398f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
540f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
541f68a32c8SEmil Constantinescu     }
542be5899b3SLisandro Dalcin 
543be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
544f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
545f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
546f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
547f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5481917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
549fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
550be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
551be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
552fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
553be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
554be5899b3SLisandro Dalcin       goto reject_step;
555be5899b3SLisandro Dalcin     }
556be5899b3SLisandro Dalcin 
5570f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
5580f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5590f7a1166SHong Zhang       rk->time_step = ts->time_step;
560493ed95fSHong Zhang     }
561be5899b3SLisandro Dalcin 
562f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
563f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
564f68a32c8SEmil Constantinescu     break;
565be5899b3SLisandro Dalcin 
566be5899b3SLisandro Dalcin   reject_step:
567fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
568be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
569be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
570be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
571f68a32c8SEmil Constantinescu     }
572f68a32c8SEmil Constantinescu   }
573f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
574e4dd521cSBarry Smith }
575e4dd521cSBarry Smith 
576f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
577f6a906c0SBarry Smith {
578f6a906c0SBarry Smith   TS_RK         *rk  = (TS_RK*)ts->data;
579be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
580be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
581f6a906c0SBarry Smith   PetscErrorCode ierr;
582f6a906c0SBarry Smith 
583f6a906c0SBarry Smith   PetscFunctionBegin;
584f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
585abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
586f6a906c0SBarry Smith   if(ts->vecs_sensip) {
587abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
588f6a906c0SBarry Smith   }
589abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
590f6a906c0SBarry Smith   PetscFunctionReturn(0);
591f6a906c0SBarry Smith }
592f6a906c0SBarry Smith 
59342f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
594d2daff3dSHong Zhang {
595c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
596c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
597c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
598c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
599c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
600ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
601b18ea86cSHong Zhang   PetscInt         i,j,nadj;
6023d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
603d2daff3dSHong Zhang   PetscErrorCode   ierr;
604cef76868SBarry Smith   PetscReal        h = ts->time_step;
605cef76868SBarry Smith   Mat              J,Jp;
606c235aa8dSHong Zhang 
607d2daff3dSHong Zhang   PetscFunctionBegin;
608c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
609c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
610c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
611abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
612b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
613302440fdSBarry Smith       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
614c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
615302440fdSBarry Smith         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
616b18ea86cSHong Zhang       }
617c235aa8dSHong Zhang     }
618ad8e2604SHong Zhang     /* Stage values of lambda */
619c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
620c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
621493ed95fSHong Zhang     if (ts->vec_costintegral) {
622493ed95fSHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
623493ed95fSHong Zhang     }
624abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
625b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
626493ed95fSHong Zhang       if (ts->vec_costintegral) {
627493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
628493ed95fSHong Zhang       }
629b18ea86cSHong Zhang     }
6306497ce14SHong Zhang 
631ad8e2604SHong Zhang     /* Stage values of mu */
6326497ce14SHong Zhang     if(ts->vecs_sensip) {
6335bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
634493ed95fSHong Zhang       if (ts->vec_costintegral) {
635493ed95fSHong Zhang         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
636493ed95fSHong Zhang       }
637493ed95fSHong Zhang 
638abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
639ad8e2604SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
640493ed95fSHong Zhang         if (ts->vec_costintegral) {
641493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
642493ed95fSHong Zhang         }
643ad8e2604SHong Zhang       }
644c235aa8dSHong Zhang     }
6456497ce14SHong Zhang   }
646c235aa8dSHong Zhang 
647c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
648abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
649b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6506497ce14SHong Zhang     if(ts->vecs_sensip) {
651ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
652b18ea86cSHong Zhang     }
6536497ce14SHong Zhang   }
654c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
655d2daff3dSHong Zhang   PetscFunctionReturn(0);
656d2daff3dSHong Zhang }
657d2daff3dSHong Zhang 
658f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
659f68a32c8SEmil Constantinescu {
660f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
6613a8a9803SLisandro Dalcin   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
662f68a32c8SEmil Constantinescu   PetscReal        h;
663f68a32c8SEmil Constantinescu   PetscReal        tt,t;
664f68a32c8SEmil Constantinescu   PetscScalar     *b;
665f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
666f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
667e4dd521cSBarry Smith 
668f68a32c8SEmil Constantinescu   PetscFunctionBegin;
669f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
670be5899b3SLisandro Dalcin 
671f68a32c8SEmil Constantinescu   switch (rk->status) {
672f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
673f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
674f68a32c8SEmil Constantinescu     h = ts->time_step;
675f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
676f68a32c8SEmil Constantinescu     break;
677f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
678be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
679f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
680f68a32c8SEmil Constantinescu     break;
681f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
682e4dd521cSBarry Smith   }
683785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
684f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
6853a8a9803SLisandro Dalcin   for (j=0,tt=t; j<p; j++,tt*=t) {
686f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
6873a8a9803SLisandro Dalcin       b[i]  += h * B[i*p+j] * tt;
688e4dd521cSBarry Smith     }
689f68a32c8SEmil Constantinescu   }
690be5899b3SLisandro Dalcin 
691f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
692f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
693be5899b3SLisandro Dalcin 
694f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
695e4dd521cSBarry Smith   PetscFunctionReturn(0);
696e4dd521cSBarry Smith }
697e4dd521cSBarry Smith 
698e4dd521cSBarry Smith /*------------------------------------------------------------*/
699be5899b3SLisandro Dalcin 
700be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
701be5899b3SLisandro Dalcin {
702be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
703be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
704be5899b3SLisandro Dalcin   PetscErrorCode ierr;
705be5899b3SLisandro Dalcin 
706be5899b3SLisandro Dalcin   PetscFunctionBegin;
707be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
708be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
709be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
710be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
711be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
712be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
713be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
714be5899b3SLisandro Dalcin }
715be5899b3SLisandro Dalcin 
716277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
717e4dd521cSBarry Smith {
7185f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7196849ba73SBarry Smith   PetscErrorCode ierr;
720e4dd521cSBarry Smith 
721e4dd521cSBarry Smith   PetscFunctionBegin;
722be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7230f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
724abc2977eSBarry Smith   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
725277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
726e4dd521cSBarry Smith }
727277b19d0SLisandro Dalcin 
728277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
729277b19d0SLisandro Dalcin {
730277b19d0SLisandro Dalcin   PetscErrorCode ierr;
731277b19d0SLisandro Dalcin 
732277b19d0SLisandro Dalcin   PetscFunctionBegin;
733277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
734277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
735f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
736f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
737f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
738f68a32c8SEmil Constantinescu }
739f68a32c8SEmil Constantinescu 
740f68a32c8SEmil Constantinescu 
741f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
742f68a32c8SEmil Constantinescu {
743f68a32c8SEmil Constantinescu   PetscFunctionBegin;
744f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
745f68a32c8SEmil Constantinescu }
746f68a32c8SEmil Constantinescu 
747f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
748f68a32c8SEmil Constantinescu {
749f68a32c8SEmil Constantinescu   PetscFunctionBegin;
750f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
751f68a32c8SEmil Constantinescu }
752f68a32c8SEmil Constantinescu 
753f68a32c8SEmil Constantinescu 
754f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
755f68a32c8SEmil Constantinescu {
756f68a32c8SEmil Constantinescu   PetscFunctionBegin;
757f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
758f68a32c8SEmil Constantinescu }
759f68a32c8SEmil Constantinescu 
760f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
761f68a32c8SEmil Constantinescu {
762f68a32c8SEmil Constantinescu 
763f68a32c8SEmil Constantinescu   PetscFunctionBegin;
764f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
765f68a32c8SEmil Constantinescu }
766c235aa8dSHong Zhang /*
767d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
768d2daff3dSHong Zhang {
769d2daff3dSHong Zhang   PetscReal *A,*b;
770d2daff3dSHong Zhang   PetscInt        s,i,j;
771d2daff3dSHong Zhang   PetscErrorCode  ierr;
772d2daff3dSHong Zhang 
773d2daff3dSHong Zhang   PetscFunctionBegin;
774d2daff3dSHong Zhang   s    = tab->s;
775d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
776d2daff3dSHong Zhang 
777d2daff3dSHong Zhang   for (i=0; i<s; i++)
778d2daff3dSHong Zhang     for (j=0; j<s; j++) {
779d2daff3dSHong Zhang       A[i*s+j] = (tab->b[s-1-i]==0) ? 0: -tab->A[s-1-i+(s-1-j)*s] * tab->b[s-1-j] / tab->b[s-1-i];
780d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
781d2daff3dSHong Zhang     }
782d2daff3dSHong Zhang 
783d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
784d2daff3dSHong Zhang 
785d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
786d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
787d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
788d2daff3dSHong Zhang   PetscFunctionReturn(0);
789d2daff3dSHong Zhang }
790c235aa8dSHong Zhang */
791be5899b3SLisandro Dalcin 
792be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
793be5899b3SLisandro Dalcin {
794be5899b3SLisandro Dalcin   TS_RK         *rk  = (TS_RK*)ts->data;
795be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
796be5899b3SLisandro Dalcin   PetscErrorCode ierr;
797be5899b3SLisandro Dalcin 
798be5899b3SLisandro Dalcin   PetscFunctionBegin;
799be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
800be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
801be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
802be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
803be5899b3SLisandro Dalcin }
804be5899b3SLisandro Dalcin 
805be5899b3SLisandro Dalcin 
806f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
807f68a32c8SEmil Constantinescu {
808f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
809f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
810f68a32c8SEmil Constantinescu   DM             dm;
811f68a32c8SEmil Constantinescu 
812f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8133dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
814be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8150f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8160f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8170f7a1166SHong Zhang   }
818f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
819f68a32c8SEmil Constantinescu   if (dm) {
820f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
821f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
822f68a32c8SEmil Constantinescu   }
823e4dd521cSBarry Smith   PetscFunctionReturn(0);
824e4dd521cSBarry Smith }
825d2daff3dSHong Zhang 
826f6a906c0SBarry Smith 
827e4dd521cSBarry Smith /*------------------------------------------------------------*/
828e4dd521cSBarry Smith 
8294416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
830e4dd521cSBarry Smith {
831be5899b3SLisandro Dalcin   TS_RK         *rk = (TS_RK*)ts->data;
832dfbe8321SBarry Smith   PetscErrorCode ierr;
833262ff7bbSBarry Smith 
834e4dd521cSBarry Smith   PetscFunctionBegin;
835e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
836f68a32c8SEmil Constantinescu   {
837f68a32c8SEmil Constantinescu     RKTableauLink  link;
838f68a32c8SEmil Constantinescu     PetscInt       count,choice;
839f68a32c8SEmil Constantinescu     PetscBool      flg;
840f68a32c8SEmil Constantinescu     const char   **namelist;
841f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
842785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
843f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
844be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
845be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
846f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
847f68a32c8SEmil Constantinescu   }
848262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
849e4dd521cSBarry Smith   PetscFunctionReturn(0);
850e4dd521cSBarry Smith }
851e4dd521cSBarry Smith 
852f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
853f68a32c8SEmil Constantinescu {
854f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
855f68a32c8SEmil Constantinescu   PetscInt       i;
856f68a32c8SEmil Constantinescu   size_t         left,count;
857f68a32c8SEmil Constantinescu   char           *p;
858f68a32c8SEmil Constantinescu 
859f68a32c8SEmil Constantinescu   PetscFunctionBegin;
860f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
861f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
862f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
863f68a32c8SEmil Constantinescu     left -= count;
864f68a32c8SEmil Constantinescu     p    += count;
865f68a32c8SEmil Constantinescu     *p++  = ' ';
866f68a32c8SEmil Constantinescu   }
867f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
868f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
869f68a32c8SEmil Constantinescu }
870f68a32c8SEmil Constantinescu 
8715f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
872e4dd521cSBarry Smith {
8735f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8748370ee3bSLisandro Dalcin   PetscBool      iascii;
875dfbe8321SBarry Smith   PetscErrorCode ierr;
8762cdf8120Sasbjorn 
877e4dd521cSBarry Smith   PetscFunctionBegin;
878251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8798370ee3bSLisandro Dalcin   if (iascii) {
8809c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
881f68a32c8SEmil Constantinescu     TSRKType  rktype;
882f68a32c8SEmil Constantinescu     char      buf[512];
883f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
884*efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
8850c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
8860c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
887f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
888f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
8898370ee3bSLisandro Dalcin   }
890f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
891f68a32c8SEmil Constantinescu }
892f68a32c8SEmil Constantinescu 
893f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
894f68a32c8SEmil Constantinescu {
895f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
8969c334d8fSLisandro Dalcin   TSAdapt        adapt;
897f68a32c8SEmil Constantinescu 
898f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8999c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9009c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
901f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
902f68a32c8SEmil Constantinescu }
903f68a32c8SEmil Constantinescu 
904f68a32c8SEmil Constantinescu /*@C
905f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
906f68a32c8SEmil Constantinescu 
907f68a32c8SEmil Constantinescu   Logically collective
908f68a32c8SEmil Constantinescu 
909f68a32c8SEmil Constantinescu   Input Parameter:
910f68a32c8SEmil Constantinescu +  ts - timestepping context
911f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
912f68a32c8SEmil Constantinescu 
913f68a32c8SEmil Constantinescu   Level: intermediate
914f68a32c8SEmil Constantinescu 
915f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
916f68a32c8SEmil Constantinescu @*/
917f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
918f68a32c8SEmil Constantinescu {
919f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
920f68a32c8SEmil Constantinescu 
921f68a32c8SEmil Constantinescu   PetscFunctionBegin;
922f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
923b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
924f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
925f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
926f68a32c8SEmil Constantinescu }
927f68a32c8SEmil Constantinescu 
928f68a32c8SEmil Constantinescu /*@C
929f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
930f68a32c8SEmil Constantinescu 
931f68a32c8SEmil Constantinescu   Logically collective
932f68a32c8SEmil Constantinescu 
933f68a32c8SEmil Constantinescu   Input Parameter:
934f68a32c8SEmil Constantinescu .  ts - timestepping context
935f68a32c8SEmil Constantinescu 
936f68a32c8SEmil Constantinescu   Output Parameter:
937f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
938f68a32c8SEmil Constantinescu 
939f68a32c8SEmil Constantinescu   Level: intermediate
940f68a32c8SEmil Constantinescu 
941f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
942f68a32c8SEmil Constantinescu @*/
943f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
944f68a32c8SEmil Constantinescu {
945f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
946f68a32c8SEmil Constantinescu 
947f68a32c8SEmil Constantinescu   PetscFunctionBegin;
948f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
949f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
950f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
951f68a32c8SEmil Constantinescu }
952f68a32c8SEmil Constantinescu 
953560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
954f68a32c8SEmil Constantinescu {
955f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
956f68a32c8SEmil Constantinescu 
957f68a32c8SEmil Constantinescu   PetscFunctionBegin;
958f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
959f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
960f68a32c8SEmil Constantinescu }
961560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
962f68a32c8SEmil Constantinescu {
963f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
964f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
965f68a32c8SEmil Constantinescu   PetscBool      match;
966f68a32c8SEmil Constantinescu   RKTableauLink  link;
967f68a32c8SEmil Constantinescu 
968f68a32c8SEmil Constantinescu   PetscFunctionBegin;
969f68a32c8SEmil Constantinescu   if (rk->tableau) {
970f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
971f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
972f68a32c8SEmil Constantinescu   }
973f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
974f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
975f68a32c8SEmil Constantinescu     if (match) {
976be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
977f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
978be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
9792ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
980f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
981f68a32c8SEmil Constantinescu     }
982f68a32c8SEmil Constantinescu   }
983f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
984e4dd521cSBarry Smith   PetscFunctionReturn(0);
985e4dd521cSBarry Smith }
986e4dd521cSBarry Smith 
987ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
988ff22ae23SHong Zhang {
989ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
990ff22ae23SHong Zhang 
991ff22ae23SHong Zhang   PetscFunctionBegin;
992ff22ae23SHong Zhang   *ns = rk->tableau->s;
993d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
994ff22ae23SHong Zhang   PetscFunctionReturn(0);
995ff22ae23SHong Zhang }
996ff22ae23SHong Zhang 
997ff22ae23SHong Zhang 
998e4dd521cSBarry Smith /* ------------------------------------------------------------ */
999ebe3b25bSBarry Smith /*MC
1000f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1001ebe3b25bSBarry Smith 
1002f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1003f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1004ebe3b25bSBarry Smith 
1005f68a32c8SEmil Constantinescu   Notes:
100698164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1007ebe3b25bSBarry Smith 
1008d41222bbSBarry Smith   Level: beginner
1009d41222bbSBarry Smith 
1010f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1011f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
1012ebe3b25bSBarry Smith 
1013ebe3b25bSBarry Smith M*/
10148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1015e4dd521cSBarry Smith {
10161566a47fSLisandro Dalcin   TS_RK          *rk;
1017dfbe8321SBarry Smith   PetscErrorCode ierr;
1018e4dd521cSBarry Smith 
1019e4dd521cSBarry Smith   PetscFunctionBegin;
1020f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1021f68a32c8SEmil Constantinescu 
1022f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10235f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10245f70b478SJed Brown   ts->ops->view           = TSView_RK;
1025f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1026f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
102742f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1028f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1029f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1030f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1031fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1032f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1033ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
103442f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1035e4dd521cSBarry Smith 
10360f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10370f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10380f7a1166SHong Zhang 
10391566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10401566a47fSLisandro Dalcin   ts->data = (void*)rk;
1041e4dd521cSBarry Smith 
1042f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1043f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1044be5899b3SLisandro Dalcin 
1045be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1046e4dd521cSBarry Smith   PetscFunctionReturn(0);
1047e4dd521cSBarry Smith }
1048