xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 3389c7d920c03cb50d8c32cd94278d87a56215bc)
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 */
420f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events */
43f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
44f68a32c8SEmil Constantinescu   PetscReal    stage_time;
45f68a32c8SEmil Constantinescu   TSStepStatus status;
460f7a1166SHong Zhang   PetscReal    ptime;
470f7a1166SHong Zhang   PetscReal    time_step;
485f70b478SJed Brown } TS_RK;
49e4dd521cSBarry Smith 
50f68a32c8SEmil Constantinescu /*MC
51287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
52262ff7bbSBarry Smith 
53f68a32c8SEmil Constantinescu      This method has one stage.
54f68a32c8SEmil Constantinescu 
55287c2655SBarry Smith      Options database:
569bd3e852SBarry Smith .     -ts_rk_type 1fe
57287c2655SBarry Smith 
58f68a32c8SEmil Constantinescu      Level: advanced
59f68a32c8SEmil Constantinescu 
60287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
61f68a32c8SEmil Constantinescu M*/
62f68a32c8SEmil Constantinescu /*MC
632109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
64f68a32c8SEmil Constantinescu 
65f68a32c8SEmil Constantinescu      This method has two stages.
66f68a32c8SEmil Constantinescu 
67287c2655SBarry Smith      Options database:
689bd3e852SBarry Smith .     -ts_rk_type 2a
69287c2655SBarry Smith 
70f68a32c8SEmil Constantinescu      Level: advanced
71f68a32c8SEmil Constantinescu 
72287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
73f68a32c8SEmil Constantinescu M*/
74f68a32c8SEmil Constantinescu /*MC
75f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
76f68a32c8SEmil Constantinescu 
77f68a32c8SEmil Constantinescu      This method has three stages.
78f68a32c8SEmil Constantinescu 
79287c2655SBarry Smith      Options database:
809bd3e852SBarry Smith .     -ts_rk_type 3
81287c2655SBarry Smith 
82f68a32c8SEmil Constantinescu      Level: advanced
83f68a32c8SEmil Constantinescu 
84287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
85f68a32c8SEmil Constantinescu M*/
86f68a32c8SEmil Constantinescu /*MC
872109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
882109b73fSDebojyoti Ghosh 
89d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
902109b73fSDebojyoti Ghosh 
91287c2655SBarry Smith      Options database:
929bd3e852SBarry Smith .     -ts_rk_type 3bs
93287c2655SBarry Smith 
942109b73fSDebojyoti Ghosh      Level: advanced
952109b73fSDebojyoti Ghosh 
9698164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
97d1905564SLisandro Dalcin 
98287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
992109b73fSDebojyoti Ghosh M*/
1002109b73fSDebojyoti Ghosh /*MC
101f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
102f68a32c8SEmil Constantinescu 
1032e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
104f68a32c8SEmil Constantinescu 
105287c2655SBarry Smith      Options database:
1069bd3e852SBarry Smith .     -ts_rk_type 4
107287c2655SBarry Smith 
108f68a32c8SEmil Constantinescu      Level: advanced
109f68a32c8SEmil Constantinescu 
110287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
111f68a32c8SEmil Constantinescu M*/
112f68a32c8SEmil Constantinescu /*MC
1132e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
114f68a32c8SEmil Constantinescu 
115f68a32c8SEmil Constantinescu      This method has six stages.
116f68a32c8SEmil Constantinescu 
117287c2655SBarry Smith      Options database:
1189bd3e852SBarry Smith .     -ts_rk_type 5f
119287c2655SBarry Smith 
120f68a32c8SEmil Constantinescu      Level: advanced
121f68a32c8SEmil Constantinescu 
122287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
123f68a32c8SEmil Constantinescu M*/
1242109b73fSDebojyoti Ghosh /*MC
1252e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1262109b73fSDebojyoti Ghosh 
127d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1282109b73fSDebojyoti Ghosh 
129287c2655SBarry Smith      Options database:
1309bd3e852SBarry Smith .     -ts_rk_type 5dp
131287c2655SBarry Smith 
1322109b73fSDebojyoti Ghosh      Level: advanced
1332109b73fSDebojyoti Ghosh 
13498164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
135d1905564SLisandro Dalcin 
136287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1372109b73fSDebojyoti Ghosh M*/
13805e23783SLisandro Dalcin /*MC
13905e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
14005e23783SLisandro Dalcin 
14105e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
14205e23783SLisandro Dalcin 
143287c2655SBarry Smith      Options database:
1449bd3e852SBarry Smith .     -ts_rk_type 5bs
145287c2655SBarry Smith 
14605e23783SLisandro Dalcin      Level: advanced
14705e23783SLisandro Dalcin 
14898164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
14905e23783SLisandro Dalcin 
150287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
15105e23783SLisandro Dalcin M*/
152262ff7bbSBarry Smith 
153f68a32c8SEmil Constantinescu /*@C
154f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
155262ff7bbSBarry Smith 
156f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
157262ff7bbSBarry Smith 
158f68a32c8SEmil Constantinescu   Level: advanced
159262ff7bbSBarry Smith 
160f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
161262ff7bbSBarry Smith 
162f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
163262ff7bbSBarry Smith @*/
164f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
165262ff7bbSBarry Smith {
1664ac538c5SBarry Smith   PetscErrorCode ierr;
167262ff7bbSBarry Smith 
168262ff7bbSBarry Smith   PetscFunctionBegin;
169f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
170f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
171e4dd521cSBarry Smith 
1724e82626cSLisandro Dalcin #define RC PetscRealConstant
173e4dd521cSBarry Smith   {
174f68a32c8SEmil Constantinescu     const PetscReal
1754e82626cSLisandro Dalcin       A[1][1] = {{0}},
1764e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1773a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
178e4dd521cSBarry Smith   }
179db046809SBarry Smith   {
180f68a32c8SEmil Constantinescu     const PetscReal
1814e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1824e82626cSLisandro Dalcin                    {RC(1.0),0}},
1834e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1844e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1853a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
186db046809SBarry Smith   }
187f68a32c8SEmil Constantinescu   {
188f68a32c8SEmil Constantinescu     const PetscReal
18917f6384fSBarry Smith       A[3][3] = {{0,0,0},
1904e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
1914e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
1924e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
1933a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
194fdefd5e5SDebojyoti Ghosh   }
195fdefd5e5SDebojyoti Ghosh   {
196fdefd5e5SDebojyoti Ghosh     const PetscReal
19717f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
1984e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
1994e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2004e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2014e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2024e82626cSLisandro 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)};
2033a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
204db046809SBarry Smith   }
205f68a32c8SEmil Constantinescu   {
206f68a32c8SEmil Constantinescu     const PetscReal
207f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2084e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2094e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2104e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2114e82626cSLisandro 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)};
2123a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
213db046809SBarry Smith   }
214f68a32c8SEmil Constantinescu   {
215f68a32c8SEmil Constantinescu     const PetscReal
21617f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2174e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2184e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2194e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2204e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2214e82626cSLisandro 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}},
2224e82626cSLisandro 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)},
2234e82626cSLisandro 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};
2243a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
225fdefd5e5SDebojyoti Ghosh   }
226fdefd5e5SDebojyoti Ghosh   {
227fdefd5e5SDebojyoti Ghosh     const PetscReal
22817f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2294e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2304e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2314e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2324e82626cSLisandro 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},
2334e82626cSLisandro 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},
2344e82626cSLisandro 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}},
2354e82626cSLisandro 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},
2364e82626cSLisandro 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)};
2373a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
238f68a32c8SEmil Constantinescu   }
23905e23783SLisandro Dalcin   {
24005e23783SLisandro Dalcin     const PetscReal
24117f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2424e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2434e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2444e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2454e82626cSLisandro 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},
2464e82626cSLisandro 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},
2474e82626cSLisandro 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},
2484e82626cSLisandro 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}},
2494e82626cSLisandro 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},
2504e82626cSLisandro 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)};
25105e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
25205e23783SLisandro Dalcin   }
2534e82626cSLisandro Dalcin #undef RC
254db046809SBarry Smith   PetscFunctionReturn(0);
255db046809SBarry Smith }
256db046809SBarry Smith 
257f68a32c8SEmil Constantinescu /*@C
258f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
259f68a32c8SEmil Constantinescu 
260f68a32c8SEmil Constantinescu    Not Collective
261f68a32c8SEmil Constantinescu 
262f68a32c8SEmil Constantinescu    Level: advanced
263f68a32c8SEmil Constantinescu 
264f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
265f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
266f68a32c8SEmil Constantinescu @*/
267f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
268e4dd521cSBarry Smith {
269f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
270f68a32c8SEmil Constantinescu   RKTableauLink link;
271f68a32c8SEmil Constantinescu 
272f68a32c8SEmil Constantinescu   PetscFunctionBegin;
273f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
274f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
275f68a32c8SEmil Constantinescu     RKTableauList = link->next;
276f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
277f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
278f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
279f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
280f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
281f68a32c8SEmil Constantinescu   }
282f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
283f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
284f68a32c8SEmil Constantinescu }
285f68a32c8SEmil Constantinescu 
286f68a32c8SEmil Constantinescu /*@C
287f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
288f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
289f68a32c8SEmil Constantinescu   when using static libraries.
290f68a32c8SEmil Constantinescu 
291f68a32c8SEmil Constantinescu   Level: developer
292f68a32c8SEmil Constantinescu 
293f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
294f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
295f68a32c8SEmil Constantinescu @*/
296f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
297f68a32c8SEmil Constantinescu {
298186e87acSLisandro Dalcin   PetscErrorCode ierr;
299e4dd521cSBarry Smith 
300e4dd521cSBarry Smith   PetscFunctionBegin;
301f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
302f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
303f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
304f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
305f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
306f68a32c8SEmil Constantinescu }
307186e87acSLisandro Dalcin 
308f68a32c8SEmil Constantinescu /*@C
309f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
310f68a32c8SEmil Constantinescu   called from PetscFinalize().
311186e87acSLisandro Dalcin 
312f68a32c8SEmil Constantinescu   Level: developer
313f68a32c8SEmil Constantinescu 
314f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
315f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
316f68a32c8SEmil Constantinescu @*/
317f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
318f68a32c8SEmil Constantinescu {
319f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
320f68a32c8SEmil Constantinescu 
321f68a32c8SEmil Constantinescu   PetscFunctionBegin;
322f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
323f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
324f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
325f68a32c8SEmil Constantinescu }
326f68a32c8SEmil Constantinescu 
327f68a32c8SEmil Constantinescu /*@C
328f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
329f68a32c8SEmil Constantinescu 
330f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
331f68a32c8SEmil Constantinescu 
332f68a32c8SEmil Constantinescu    Input Parameters:
333f68a32c8SEmil Constantinescu +  name - identifier for method
334f68a32c8SEmil Constantinescu .  order - approximation order of method
335f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
336f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
337f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
338f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
339f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3403a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3413a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
342f68a32c8SEmil Constantinescu 
343f68a32c8SEmil Constantinescu    Notes:
344f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
345f68a32c8SEmil Constantinescu 
346f68a32c8SEmil Constantinescu    Level: advanced
347f68a32c8SEmil Constantinescu 
348f68a32c8SEmil Constantinescu .keywords: TS, register
349f68a32c8SEmil Constantinescu 
350f68a32c8SEmil Constantinescu .seealso: TSRK
351f68a32c8SEmil Constantinescu @*/
352f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
353f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3543a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
355f68a32c8SEmil Constantinescu {
356f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
357f68a32c8SEmil Constantinescu   RKTableauLink   link;
358f68a32c8SEmil Constantinescu   RKTableau       t;
359f68a32c8SEmil Constantinescu   PetscInt        i,j;
360f68a32c8SEmil Constantinescu 
361f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3623a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3633a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3643a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3653a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3663a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3673a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3683a8a9803SLisandro Dalcin 
36995dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
370f68a32c8SEmil Constantinescu   t = &link->tab;
3713a8a9803SLisandro Dalcin 
372f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
373f68a32c8SEmil Constantinescu   t->order = order;
374f68a32c8SEmil Constantinescu   t->s = s;
375dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
376f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
377f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
378f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
379f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
380f68a32c8SEmil 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];
381d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
382d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3833a8a9803SLisandro Dalcin 
384f68a32c8SEmil Constantinescu   if (bembed) {
385785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
386f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
387f68a32c8SEmil Constantinescu   }
388f68a32c8SEmil Constantinescu 
3893a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3903a8a9803SLisandro Dalcin   t->p = p;
3913a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3923a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3933a8a9803SLisandro Dalcin 
394f68a32c8SEmil Constantinescu   link->next = RKTableauList;
395f68a32c8SEmil Constantinescu   RKTableauList = link;
396f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
397f68a32c8SEmil Constantinescu }
398f68a32c8SEmil Constantinescu 
399e4dd521cSBarry Smith /*
400f68a32c8SEmil Constantinescu  The step completion formula is
401e4dd521cSBarry Smith 
402f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
403f68a32c8SEmil Constantinescu 
404f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
405f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
406f68a32c8SEmil Constantinescu  We can write
407f68a32c8SEmil Constantinescu 
408f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
409f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
410f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
411f68a32c8SEmil Constantinescu 
412f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
413e4dd521cSBarry Smith */
414f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
415f68a32c8SEmil Constantinescu {
416f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
417f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
418f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
419f68a32c8SEmil Constantinescu   PetscReal      h;
420f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
421f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
422f68a32c8SEmil Constantinescu 
423f68a32c8SEmil Constantinescu   PetscFunctionBegin;
424f68a32c8SEmil Constantinescu   switch (rk->status) {
425f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
426f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
427f68a32c8SEmil Constantinescu     h = ts->time_step; break;
428f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
429be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
430f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
431f68a32c8SEmil Constantinescu   }
432f68a32c8SEmil Constantinescu   if (order == tab->order) {
433f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
434f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
435f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
436f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
437f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
438f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
439f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
440f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
441f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
442f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
443f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
444f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
445f68a32c8SEmil Constantinescu     } else { /* Rollback and re-complete using (be-b) */
446f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
447f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
448f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
4490f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
4500f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4510f7a1166SHong Zhang       }
452f68a32c8SEmil Constantinescu     }
453f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
454f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
455f68a32c8SEmil Constantinescu   }
456f68a32c8SEmil Constantinescu unavailable:
457f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
458a7fac7c2SEmil 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);
459f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
460f68a32c8SEmil Constantinescu }
461f68a32c8SEmil Constantinescu 
4620f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4630f7a1166SHong Zhang {
4640f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4650f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4660f7a1166SHong Zhang   const PetscInt  s = tab->s;
4670f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4680f7a1166SHong Zhang   Vec             *Y = rk->Y;
4690f7a1166SHong Zhang   PetscInt        i;
4700f7a1166SHong Zhang   PetscErrorCode  ierr;
4710f7a1166SHong Zhang 
4720f7a1166SHong Zhang   PetscFunctionBegin;
4730f7a1166SHong Zhang   /* backup cost integral */
4740f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4750f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4760f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
477022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4780f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4790f7a1166SHong Zhang   }
4800f7a1166SHong Zhang   PetscFunctionReturn(0);
4810f7a1166SHong Zhang }
4820f7a1166SHong Zhang 
4830f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4840f7a1166SHong Zhang {
4850f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4860f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4870f7a1166SHong Zhang   const PetscInt  s = tab->s;
4880f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4890f7a1166SHong Zhang   Vec             *Y = rk->Y;
4900f7a1166SHong Zhang   PetscInt        i;
4910f7a1166SHong Zhang   PetscErrorCode  ierr;
4920f7a1166SHong Zhang 
4930f7a1166SHong Zhang   PetscFunctionBegin;
4940f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4950f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
496022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4970f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4980f7a1166SHong Zhang   }
4990f7a1166SHong Zhang   PetscFunctionReturn(0);
5000f7a1166SHong Zhang }
5010f7a1166SHong Zhang 
502fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
503fecfb714SLisandro Dalcin {
504fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
505fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
506fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
507fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
508fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
509fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
510fecfb714SLisandro Dalcin   PetscInt        j;
511fecfb714SLisandro Dalcin   PetscReal       h;
512fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
513fecfb714SLisandro Dalcin 
514fecfb714SLisandro Dalcin   PetscFunctionBegin;
515fecfb714SLisandro Dalcin   switch (rk->status) {
516fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
517fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
518fecfb714SLisandro Dalcin     h = ts->time_step; break;
519fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
520fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
521fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
522fecfb714SLisandro Dalcin   }
523fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
524fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
525fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
526fecfb714SLisandro Dalcin }
527fecfb714SLisandro Dalcin 
528f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
529f68a32c8SEmil Constantinescu {
530f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
531f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
532f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
533fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
534f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
535f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
536d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
537f68a32c8SEmil Constantinescu   TSAdapt          adapt;
538fecfb714SLisandro Dalcin   PetscInt         i,j;
539be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
540be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
541be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
542f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
543f68a32c8SEmil Constantinescu 
544f68a32c8SEmil Constantinescu   PetscFunctionBegin;
545d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
546d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
547d1905564SLisandro Dalcin 
548f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
549be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
550be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
551f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
552f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5539be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5549be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
555f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
556f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
557f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5589be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
559f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
560be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
561fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5628f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
563f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
564f68a32c8SEmil Constantinescu     }
565be5899b3SLisandro Dalcin 
566be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
567f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
568f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
569f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
570f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5711917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
572fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
573be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
574be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
575fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
576be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
577be5899b3SLisandro Dalcin       goto reject_step;
578be5899b3SLisandro Dalcin     }
579be5899b3SLisandro Dalcin 
5800f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
5810f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5820f7a1166SHong Zhang       rk->time_step = ts->time_step;
583493ed95fSHong Zhang     }
584be5899b3SLisandro Dalcin 
585f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
586f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
587f68a32c8SEmil Constantinescu     break;
588be5899b3SLisandro Dalcin 
589be5899b3SLisandro Dalcin   reject_step:
590fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
591be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
592be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
593be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
594f68a32c8SEmil Constantinescu     }
595f68a32c8SEmil Constantinescu   }
596f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
597e4dd521cSBarry Smith }
598e4dd521cSBarry Smith 
599f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
600f6a906c0SBarry Smith {
601f6a906c0SBarry Smith   TS_RK         *rk  = (TS_RK*)ts->data;
602be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
603be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
604f6a906c0SBarry Smith   PetscErrorCode ierr;
605f6a906c0SBarry Smith 
606f6a906c0SBarry Smith   PetscFunctionBegin;
607f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
608abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
609f6a906c0SBarry Smith   if(ts->vecs_sensip) {
610abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
611f6a906c0SBarry Smith   }
612f6a906c0SBarry Smith   PetscFunctionReturn(0);
613f6a906c0SBarry Smith }
614f6a906c0SBarry Smith 
61542f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
616d2daff3dSHong Zhang {
617c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
618c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
619c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
620c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
621c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
622*3389c7d9SStefano Zampini   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu;
623b18ea86cSHong Zhang   PetscInt         i,j,nadj;
6243d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
625d2daff3dSHong Zhang   PetscErrorCode   ierr;
626cef76868SBarry Smith   PetscReal        h = ts->time_step;
627c235aa8dSHong Zhang 
628d2daff3dSHong Zhang   PetscFunctionBegin;
629c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
630c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
631*3389c7d9SStefano Zampini     Mat       J;
632*3389c7d9SStefano Zampini     PetscReal stage_time = t + h*(1.0-c[i]);
633*3389c7d9SStefano Zampini     PetscBool zero = PETSC_FALSE;
634*3389c7d9SStefano Zampini 
635*3389c7d9SStefano Zampini     ierr = TSGetRHSJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr);
636*3389c7d9SStefano Zampini     ierr = TSComputeRHSJacobian(ts,stage_time,Y[i],J,J);CHKERRQ(ierr);
637493ed95fSHong Zhang     if (ts->vec_costintegral) {
638*3389c7d9SStefano Zampini       ierr = TSAdjointComputeDRDYFunction(ts,stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
639493ed95fSHong Zhang     }
640*3389c7d9SStefano Zampini     /* Stage values of mu */
641*3389c7d9SStefano Zampini     if (ts->vecs_sensip) {
642*3389c7d9SStefano Zampini       ierr = TSAdjointComputeRHSJacobian(ts,stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
643*3389c7d9SStefano Zampini       if (ts->vec_costintegral) {
644*3389c7d9SStefano Zampini         ierr = TSAdjointComputeDRDPFunction(ts,stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
645*3389c7d9SStefano Zampini       }
646*3389c7d9SStefano Zampini     }
647*3389c7d9SStefano Zampini 
648*3389c7d9SStefano Zampini     if (b[i] == 0 && i == s-1) zero = PETSC_TRUE;
649abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
650*3389c7d9SStefano Zampini       DM  dm;
651*3389c7d9SStefano Zampini       Vec VecSensiTemp;
652*3389c7d9SStefano Zampini 
653*3389c7d9SStefano Zampini       ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
654*3389c7d9SStefano Zampini       ierr = DMGetGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
655*3389c7d9SStefano Zampini       /* Stage values of lambda */
656*3389c7d9SStefano Zampini       if (!zero) {
657*3389c7d9SStefano Zampini         ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp);CHKERRQ(ierr);
658*3389c7d9SStefano Zampini         ierr = VecScale(VecSensiTemp,-h*b[i]);CHKERRQ(ierr);
659*3389c7d9SStefano Zampini         for (j=i+1; j<s; j++) w[j-i-1]  = -h*A[j*s+i];
660*3389c7d9SStefano Zampini         ierr = VecMAXPY(VecSensiTemp,s-i-1,w,&VecDeltaLam[nadj*s+i+1]);CHKERRQ(ierr);
661*3389c7d9SStefano Zampini         ierr = MatMultTranspose(J,VecSensiTemp,VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
662*3389c7d9SStefano Zampini       } else {
663*3389c7d9SStefano Zampini         ierr = VecSet(VecDeltaLam[nadj*s+i],0);CHKERRQ(ierr);
664*3389c7d9SStefano Zampini       }
665493ed95fSHong Zhang       if (ts->vec_costintegral) {
666493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
667493ed95fSHong Zhang       }
6686497ce14SHong Zhang 
669ad8e2604SHong Zhang       /* Stage values of mu */
6706497ce14SHong Zhang       if (ts->vecs_sensip) {
671*3389c7d9SStefano Zampini         if (!zero) {
672*3389c7d9SStefano Zampini           ierr = MatMultTranspose(ts->Jacp,VecSensiTemp,VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
673*3389c7d9SStefano Zampini         } else {
674*3389c7d9SStefano Zampini           ierr = VecSet(VecDeltaMu[nadj*s+i],0);CHKERRQ(ierr);
675493ed95fSHong Zhang         }
676493ed95fSHong Zhang         if (ts->vec_costintegral) {
677493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
678493ed95fSHong Zhang         }
679ad8e2604SHong Zhang       }
680*3389c7d9SStefano Zampini       ierr = DMRestoreGlobalVector(dm,&VecSensiTemp);CHKERRQ(ierr);
681c235aa8dSHong Zhang     }
6826497ce14SHong Zhang   }
683c235aa8dSHong Zhang 
684c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
685abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
686b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6876497ce14SHong Zhang     if (ts->vecs_sensip) {
688ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
689b18ea86cSHong Zhang     }
6906497ce14SHong Zhang   }
691c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
692d2daff3dSHong Zhang   PetscFunctionReturn(0);
693d2daff3dSHong Zhang }
694d2daff3dSHong Zhang 
695f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
696f68a32c8SEmil Constantinescu {
697f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
6983a8a9803SLisandro Dalcin   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
699f68a32c8SEmil Constantinescu   PetscReal        h;
700f68a32c8SEmil Constantinescu   PetscReal        tt,t;
701f68a32c8SEmil Constantinescu   PetscScalar     *b;
702f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
703f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
704e4dd521cSBarry Smith 
705f68a32c8SEmil Constantinescu   PetscFunctionBegin;
706f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
707be5899b3SLisandro Dalcin 
708f68a32c8SEmil Constantinescu   switch (rk->status) {
709f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
710f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
711f68a32c8SEmil Constantinescu     h = ts->time_step;
712f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
713f68a32c8SEmil Constantinescu     break;
714f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
715be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
716f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
717f68a32c8SEmil Constantinescu     break;
718f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
719e4dd521cSBarry Smith   }
720785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
721f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
7223a8a9803SLisandro Dalcin   for (j=0,tt=t; j<p; j++,tt*=t) {
723f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
7243a8a9803SLisandro Dalcin       b[i]  += h * B[i*p+j] * tt;
725e4dd521cSBarry Smith     }
726f68a32c8SEmil Constantinescu   }
727be5899b3SLisandro Dalcin 
728f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
729f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
730be5899b3SLisandro Dalcin 
731f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
732e4dd521cSBarry Smith   PetscFunctionReturn(0);
733e4dd521cSBarry Smith }
734e4dd521cSBarry Smith 
735e4dd521cSBarry Smith /*------------------------------------------------------------*/
736be5899b3SLisandro Dalcin 
737be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
738be5899b3SLisandro Dalcin {
739be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
740be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
741be5899b3SLisandro Dalcin   PetscErrorCode ierr;
742be5899b3SLisandro Dalcin 
743be5899b3SLisandro Dalcin   PetscFunctionBegin;
744be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
745be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
746be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
747be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
748be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
749be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
750be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
751be5899b3SLisandro Dalcin }
752be5899b3SLisandro Dalcin 
753277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
754e4dd521cSBarry Smith {
7555f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7566849ba73SBarry Smith   PetscErrorCode ierr;
757e4dd521cSBarry Smith 
758e4dd521cSBarry Smith   PetscFunctionBegin;
759be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7600f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
761277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
762e4dd521cSBarry Smith }
763277b19d0SLisandro Dalcin 
764f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
765f68a32c8SEmil Constantinescu {
766f68a32c8SEmil Constantinescu   PetscFunctionBegin;
767f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
768f68a32c8SEmil Constantinescu }
769f68a32c8SEmil Constantinescu 
770f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
771f68a32c8SEmil Constantinescu {
772f68a32c8SEmil Constantinescu   PetscFunctionBegin;
773f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
774f68a32c8SEmil Constantinescu }
775f68a32c8SEmil Constantinescu 
776f68a32c8SEmil Constantinescu 
777f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
778f68a32c8SEmil Constantinescu {
779f68a32c8SEmil Constantinescu   PetscFunctionBegin;
780f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
781f68a32c8SEmil Constantinescu }
782f68a32c8SEmil Constantinescu 
783f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
784f68a32c8SEmil Constantinescu {
785f68a32c8SEmil Constantinescu 
786f68a32c8SEmil Constantinescu   PetscFunctionBegin;
787f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
788f68a32c8SEmil Constantinescu }
789c235aa8dSHong Zhang /*
790d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
791d2daff3dSHong Zhang {
792d2daff3dSHong Zhang   PetscReal *A,*b;
793d2daff3dSHong Zhang   PetscInt        s,i,j;
794d2daff3dSHong Zhang   PetscErrorCode  ierr;
795d2daff3dSHong Zhang 
796d2daff3dSHong Zhang   PetscFunctionBegin;
797d2daff3dSHong Zhang   s    = tab->s;
798d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
799d2daff3dSHong Zhang 
800d2daff3dSHong Zhang   for (i=0; i<s; i++)
801d2daff3dSHong Zhang     for (j=0; j<s; j++) {
802d2daff3dSHong 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];
803d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
804d2daff3dSHong Zhang     }
805d2daff3dSHong Zhang 
806d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
807d2daff3dSHong Zhang 
808d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
809d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
810d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
811d2daff3dSHong Zhang   PetscFunctionReturn(0);
812d2daff3dSHong Zhang }
813c235aa8dSHong Zhang */
814be5899b3SLisandro Dalcin 
815be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
816be5899b3SLisandro Dalcin {
817be5899b3SLisandro Dalcin   TS_RK         *rk  = (TS_RK*)ts->data;
818be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
819be5899b3SLisandro Dalcin   PetscErrorCode ierr;
820be5899b3SLisandro Dalcin 
821be5899b3SLisandro Dalcin   PetscFunctionBegin;
822be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
823be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
824be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
825be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
826be5899b3SLisandro Dalcin }
827be5899b3SLisandro Dalcin 
828be5899b3SLisandro Dalcin 
829f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
830f68a32c8SEmil Constantinescu {
831f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
832f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
833f68a32c8SEmil Constantinescu   DM             dm;
834f68a32c8SEmil Constantinescu 
835f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8363dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
837be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8380f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8390f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8400f7a1166SHong Zhang   }
841f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
842f68a32c8SEmil Constantinescu   ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
843f68a32c8SEmil Constantinescu   ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
844e4dd521cSBarry Smith   PetscFunctionReturn(0);
845e4dd521cSBarry Smith }
846d2daff3dSHong Zhang 
847f6a906c0SBarry Smith 
848e4dd521cSBarry Smith /*------------------------------------------------------------*/
849e4dd521cSBarry Smith 
8504416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
851e4dd521cSBarry Smith {
852be5899b3SLisandro Dalcin   TS_RK         *rk = (TS_RK*)ts->data;
853dfbe8321SBarry Smith   PetscErrorCode ierr;
854262ff7bbSBarry Smith 
855e4dd521cSBarry Smith   PetscFunctionBegin;
856e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
857f68a32c8SEmil Constantinescu   {
858f68a32c8SEmil Constantinescu     RKTableauLink  link;
859f68a32c8SEmil Constantinescu     PetscInt       count,choice;
860f68a32c8SEmil Constantinescu     PetscBool      flg;
861f68a32c8SEmil Constantinescu     const char   **namelist;
862f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
863f489ac74SBarry Smith     ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr);
864f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
865be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
866be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
867f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
868f68a32c8SEmil Constantinescu   }
869262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
870e4dd521cSBarry Smith   PetscFunctionReturn(0);
871e4dd521cSBarry Smith }
872e4dd521cSBarry Smith 
873f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
874f68a32c8SEmil Constantinescu {
875f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
876f68a32c8SEmil Constantinescu   PetscInt       i;
877f68a32c8SEmil Constantinescu   size_t         left,count;
878f68a32c8SEmil Constantinescu   char           *p;
879f68a32c8SEmil Constantinescu 
880f68a32c8SEmil Constantinescu   PetscFunctionBegin;
881f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
882f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
883f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
884f68a32c8SEmil Constantinescu     left -= count;
885f68a32c8SEmil Constantinescu     p    += count;
886f68a32c8SEmil Constantinescu     *p++  = ' ';
887f68a32c8SEmil Constantinescu   }
888f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
889f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
890f68a32c8SEmil Constantinescu }
891f68a32c8SEmil Constantinescu 
8925f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
893e4dd521cSBarry Smith {
8945f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8958370ee3bSLisandro Dalcin   PetscBool      iascii;
896dfbe8321SBarry Smith   PetscErrorCode ierr;
8972cdf8120Sasbjorn 
898e4dd521cSBarry Smith   PetscFunctionBegin;
899251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9008370ee3bSLisandro Dalcin   if (iascii) {
9019c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
902f68a32c8SEmil Constantinescu     TSRKType  rktype;
903f68a32c8SEmil Constantinescu     char      buf[512];
904f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
905efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
9060c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
9070c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
908f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
909f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
9108370ee3bSLisandro Dalcin   }
911f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
912f68a32c8SEmil Constantinescu }
913f68a32c8SEmil Constantinescu 
914f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
915f68a32c8SEmil Constantinescu {
916f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
9179c334d8fSLisandro Dalcin   TSAdapt        adapt;
918f68a32c8SEmil Constantinescu 
919f68a32c8SEmil Constantinescu   PetscFunctionBegin;
9209c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9219c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
922f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
923f68a32c8SEmil Constantinescu }
924f68a32c8SEmil Constantinescu 
925f68a32c8SEmil Constantinescu /*@C
926f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
927f68a32c8SEmil Constantinescu 
928f68a32c8SEmil Constantinescu   Logically collective
929f68a32c8SEmil Constantinescu 
930f68a32c8SEmil Constantinescu   Input Parameter:
931f68a32c8SEmil Constantinescu +  ts - timestepping context
932f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
933f68a32c8SEmil Constantinescu 
934287c2655SBarry Smith   Options Database:
9359bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
936287c2655SBarry Smith 
937f68a32c8SEmil Constantinescu   Level: intermediate
938f68a32c8SEmil Constantinescu 
939287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
940f68a32c8SEmil Constantinescu @*/
941f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
942f68a32c8SEmil Constantinescu {
943f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
944f68a32c8SEmil Constantinescu 
945f68a32c8SEmil Constantinescu   PetscFunctionBegin;
946f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
947b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
948f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
949f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
950f68a32c8SEmil Constantinescu }
951f68a32c8SEmil Constantinescu 
952f68a32c8SEmil Constantinescu /*@C
953f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
954f68a32c8SEmil Constantinescu 
955f68a32c8SEmil Constantinescu   Logically collective
956f68a32c8SEmil Constantinescu 
957f68a32c8SEmil Constantinescu   Input Parameter:
958f68a32c8SEmil Constantinescu .  ts - timestepping context
959f68a32c8SEmil Constantinescu 
960f68a32c8SEmil Constantinescu   Output Parameter:
961f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
962f68a32c8SEmil Constantinescu 
963f68a32c8SEmil Constantinescu   Level: intermediate
964f68a32c8SEmil Constantinescu 
965f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
966f68a32c8SEmil Constantinescu @*/
967f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
968f68a32c8SEmil Constantinescu {
969f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
970f68a32c8SEmil Constantinescu 
971f68a32c8SEmil Constantinescu   PetscFunctionBegin;
972f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
973f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
974f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
975f68a32c8SEmil Constantinescu }
976f68a32c8SEmil Constantinescu 
977560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
978f68a32c8SEmil Constantinescu {
979f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
980f68a32c8SEmil Constantinescu 
981f68a32c8SEmil Constantinescu   PetscFunctionBegin;
982f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
983f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
984f68a32c8SEmil Constantinescu }
985560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
986f68a32c8SEmil Constantinescu {
987f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
988f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
989f68a32c8SEmil Constantinescu   PetscBool      match;
990f68a32c8SEmil Constantinescu   RKTableauLink  link;
991f68a32c8SEmil Constantinescu 
992f68a32c8SEmil Constantinescu   PetscFunctionBegin;
993f68a32c8SEmil Constantinescu   if (rk->tableau) {
994f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
995f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
996f68a32c8SEmil Constantinescu   }
997f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
998f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
999f68a32c8SEmil Constantinescu     if (match) {
1000be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1001f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1002be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
10032ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1004f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1005f68a32c8SEmil Constantinescu     }
1006f68a32c8SEmil Constantinescu   }
1007f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1008e4dd521cSBarry Smith   PetscFunctionReturn(0);
1009e4dd521cSBarry Smith }
1010e4dd521cSBarry Smith 
1011ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1012ff22ae23SHong Zhang {
1013ff22ae23SHong Zhang   TS_RK *rk = (TS_RK*)ts->data;
1014ff22ae23SHong Zhang 
1015ff22ae23SHong Zhang   PetscFunctionBegin;
1016ff22ae23SHong Zhang   *ns = rk->tableau->s;
1017d2daff3dSHong Zhang   if (Y) *Y = rk->Y;
1018ff22ae23SHong Zhang   PetscFunctionReturn(0);
1019ff22ae23SHong Zhang }
1020ff22ae23SHong Zhang 
1021b3a6b972SJed Brown static PetscErrorCode TSDestroy_RK(TS ts)
1022b3a6b972SJed Brown {
1023b3a6b972SJed Brown   PetscErrorCode ierr;
1024b3a6b972SJed Brown 
1025b3a6b972SJed Brown   PetscFunctionBegin;
1026b3a6b972SJed Brown   ierr = TSReset_RK(ts);CHKERRQ(ierr);
1027b3a6b972SJed Brown   if (ts->dm) {
1028b3a6b972SJed Brown     ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
1029b3a6b972SJed Brown     ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
1030b3a6b972SJed Brown   }
1031b3a6b972SJed Brown   ierr = PetscFree(ts->data);CHKERRQ(ierr);
1032b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
1033b3a6b972SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
1034b3a6b972SJed Brown   PetscFunctionReturn(0);
1035b3a6b972SJed Brown }
1036ff22ae23SHong Zhang 
1037e4dd521cSBarry Smith /* ------------------------------------------------------------ */
1038ebe3b25bSBarry Smith /*MC
1039f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1040ebe3b25bSBarry Smith 
1041f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1042f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1043ebe3b25bSBarry Smith 
1044f68a32c8SEmil Constantinescu   Notes:
104598164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1046ebe3b25bSBarry Smith 
1047d41222bbSBarry Smith   Level: beginner
1048d41222bbSBarry Smith 
1049f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1050f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
1051ebe3b25bSBarry Smith 
1052ebe3b25bSBarry Smith M*/
10538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1054e4dd521cSBarry Smith {
10551566a47fSLisandro Dalcin   TS_RK          *rk;
1056dfbe8321SBarry Smith   PetscErrorCode ierr;
1057e4dd521cSBarry Smith 
1058e4dd521cSBarry Smith   PetscFunctionBegin;
1059f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1060f68a32c8SEmil Constantinescu 
1061f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10625f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10635f70b478SJed Brown   ts->ops->view           = TSView_RK;
1064f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1065f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
106642f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1067f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1068f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1069f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1070fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1071f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1072ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
107342f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1074e4dd521cSBarry Smith 
10750f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10760f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10770f7a1166SHong Zhang 
10781566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10791566a47fSLisandro Dalcin   ts->data = (void*)rk;
1080e4dd521cSBarry Smith 
1081f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1082f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1083be5899b3SLisandro Dalcin 
1084be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1085e4dd521cSBarry Smith   PetscFunctionReturn(0);
1086e4dd521cSBarry Smith }
1087