xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 022c081a6d8b02a670b429ae64c0e7826e2eb534)
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
52287c2655SBarry Smith      TSRK1FE - First order forward Euler scheme.
53262ff7bbSBarry Smith 
54f68a32c8SEmil Constantinescu      This method has one stage.
55f68a32c8SEmil Constantinescu 
56287c2655SBarry Smith      Options database:
579bd3e852SBarry Smith .     -ts_rk_type 1fe
58287c2655SBarry Smith 
59f68a32c8SEmil Constantinescu      Level: advanced
60f68a32c8SEmil Constantinescu 
61287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
62f68a32c8SEmil Constantinescu M*/
63f68a32c8SEmil Constantinescu /*MC
642109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
65f68a32c8SEmil Constantinescu 
66f68a32c8SEmil Constantinescu      This method has two stages.
67f68a32c8SEmil Constantinescu 
68287c2655SBarry Smith      Options database:
699bd3e852SBarry Smith .     -ts_rk_type 2a
70287c2655SBarry Smith 
71f68a32c8SEmil Constantinescu      Level: advanced
72f68a32c8SEmil Constantinescu 
73287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
74f68a32c8SEmil Constantinescu M*/
75f68a32c8SEmil Constantinescu /*MC
76f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
77f68a32c8SEmil Constantinescu 
78f68a32c8SEmil Constantinescu      This method has three stages.
79f68a32c8SEmil Constantinescu 
80287c2655SBarry Smith      Options database:
819bd3e852SBarry Smith .     -ts_rk_type 3
82287c2655SBarry Smith 
83f68a32c8SEmil Constantinescu      Level: advanced
84f68a32c8SEmil Constantinescu 
85287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
86f68a32c8SEmil Constantinescu M*/
87f68a32c8SEmil Constantinescu /*MC
882109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
892109b73fSDebojyoti Ghosh 
90d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
912109b73fSDebojyoti Ghosh 
92287c2655SBarry Smith      Options database:
939bd3e852SBarry Smith .     -ts_rk_type 3bs
94287c2655SBarry Smith 
952109b73fSDebojyoti Ghosh      Level: advanced
962109b73fSDebojyoti Ghosh 
9798164e13SLisandro Dalcin      References: https://doi.org/10.1016/0893-9659(89)90079-7
98d1905564SLisandro Dalcin 
99287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1002109b73fSDebojyoti Ghosh M*/
1012109b73fSDebojyoti Ghosh /*MC
102f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
103f68a32c8SEmil Constantinescu 
1042e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
105f68a32c8SEmil Constantinescu 
106287c2655SBarry Smith      Options database:
1079bd3e852SBarry Smith .     -ts_rk_type 4
108287c2655SBarry Smith 
109f68a32c8SEmil Constantinescu      Level: advanced
110f68a32c8SEmil Constantinescu 
111287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
112f68a32c8SEmil Constantinescu M*/
113f68a32c8SEmil Constantinescu /*MC
1142e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
115f68a32c8SEmil Constantinescu 
116f68a32c8SEmil Constantinescu      This method has six stages.
117f68a32c8SEmil Constantinescu 
118287c2655SBarry Smith      Options database:
1199bd3e852SBarry Smith .     -ts_rk_type 5f
120287c2655SBarry Smith 
121f68a32c8SEmil Constantinescu      Level: advanced
122f68a32c8SEmil Constantinescu 
123287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
124f68a32c8SEmil Constantinescu M*/
1252109b73fSDebojyoti Ghosh /*MC
1262e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1272109b73fSDebojyoti Ghosh 
128d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1292109b73fSDebojyoti Ghosh 
130287c2655SBarry Smith      Options database:
1319bd3e852SBarry Smith .     -ts_rk_type 5dp
132287c2655SBarry Smith 
1332109b73fSDebojyoti Ghosh      Level: advanced
1342109b73fSDebojyoti Ghosh 
13598164e13SLisandro Dalcin      References: https://doi.org/10.1016/0771-050X(80)90013-3
136d1905564SLisandro Dalcin 
137287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
1382109b73fSDebojyoti Ghosh M*/
13905e23783SLisandro Dalcin /*MC
14005e23783SLisandro Dalcin      TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method.
14105e23783SLisandro Dalcin 
14205e23783SLisandro Dalcin      This method has eight stages with the First Same As Last (FSAL) property.
14305e23783SLisandro Dalcin 
144287c2655SBarry Smith      Options database:
1459bd3e852SBarry Smith .     -ts_rk_type 5bs
146287c2655SBarry Smith 
14705e23783SLisandro Dalcin      Level: advanced
14805e23783SLisandro Dalcin 
14998164e13SLisandro Dalcin      References: https://doi.org/10.1016/0898-1221(96)00141-1
15005e23783SLisandro Dalcin 
151287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType()
15205e23783SLisandro Dalcin M*/
153262ff7bbSBarry Smith 
154f68a32c8SEmil Constantinescu /*@C
155f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
156262ff7bbSBarry Smith 
157f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
158262ff7bbSBarry Smith 
159f68a32c8SEmil Constantinescu   Level: advanced
160262ff7bbSBarry Smith 
161f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
162262ff7bbSBarry Smith 
163f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
164262ff7bbSBarry Smith @*/
165f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
166262ff7bbSBarry Smith {
1674ac538c5SBarry Smith   PetscErrorCode ierr;
168262ff7bbSBarry Smith 
169262ff7bbSBarry Smith   PetscFunctionBegin;
170f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
171f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
172e4dd521cSBarry Smith 
1734e82626cSLisandro Dalcin #define RC PetscRealConstant
174e4dd521cSBarry Smith   {
175f68a32c8SEmil Constantinescu     const PetscReal
1764e82626cSLisandro Dalcin       A[1][1] = {{0}},
1774e82626cSLisandro Dalcin       b[1]    = {RC(1.0)};
1783a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
179e4dd521cSBarry Smith   }
180db046809SBarry Smith   {
181f68a32c8SEmil Constantinescu     const PetscReal
1824e82626cSLisandro Dalcin       A[2][2]   = {{0,0},
1834e82626cSLisandro Dalcin                    {RC(1.0),0}},
1844e82626cSLisandro Dalcin       b[2]      =  {RC(0.5),RC(0.5)},
1854e82626cSLisandro Dalcin       bembed[2] =  {RC(1.0),0};
1863a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
187db046809SBarry Smith   }
188f68a32c8SEmil Constantinescu   {
189f68a32c8SEmil Constantinescu     const PetscReal
19017f6384fSBarry Smith       A[3][3] = {{0,0,0},
1914e82626cSLisandro Dalcin                  {RC(2.0)/RC(3.0),0,0},
1924e82626cSLisandro Dalcin                  {RC(-1.0)/RC(3.0),RC(1.0),0}},
1934e82626cSLisandro Dalcin       b[3]    =  {RC(0.25),RC(0.5),RC(0.25)};
1943a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
195fdefd5e5SDebojyoti Ghosh   }
196fdefd5e5SDebojyoti Ghosh   {
197fdefd5e5SDebojyoti Ghosh     const PetscReal
19817f6384fSBarry Smith       A[4][4]   = {{0,0,0,0},
1994e82626cSLisandro Dalcin                    {RC(1.0)/RC(2.0),0,0,0},
2004e82626cSLisandro Dalcin                    {0,RC(3.0)/RC(4.0),0,0},
2014e82626cSLisandro Dalcin                    {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0}},
2024e82626cSLisandro Dalcin       b[4]      =  {RC(2.0)/RC(9.0),RC(1.0)/RC(3.0),RC(4.0)/RC(9.0),0},
2034e82626cSLisandro 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)};
2043a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
205db046809SBarry Smith   }
206f68a32c8SEmil Constantinescu   {
207f68a32c8SEmil Constantinescu     const PetscReal
208f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
2094e82626cSLisandro Dalcin                  {RC(0.5),0,0,0},
2104e82626cSLisandro Dalcin                  {0,RC(0.5),0,0},
2114e82626cSLisandro Dalcin                  {0,0,RC(1.0),0}},
2124e82626cSLisandro 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)};
2133a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
214db046809SBarry Smith   }
215f68a32c8SEmil Constantinescu   {
216f68a32c8SEmil Constantinescu     const PetscReal
21717f6384fSBarry Smith       A[6][6]   = {{0,0,0,0,0,0},
2184e82626cSLisandro Dalcin                    {RC(0.25),0,0,0,0,0},
2194e82626cSLisandro Dalcin                    {RC(3.0)/RC(32.0),RC(9.0)/RC(32.0),0,0,0,0},
2204e82626cSLisandro Dalcin                    {RC(1932.0)/RC(2197.0),RC(-7200.0)/RC(2197.0),RC(7296.0)/RC(2197.0),0,0,0},
2214e82626cSLisandro Dalcin                    {RC(439.0)/RC(216.0),RC(-8.0),RC(3680.0)/RC(513.0),RC(-845.0)/RC(4104.0),0,0},
2224e82626cSLisandro 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}},
2234e82626cSLisandro 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)},
2244e82626cSLisandro 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};
2253a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
226fdefd5e5SDebojyoti Ghosh   }
227fdefd5e5SDebojyoti Ghosh   {
228fdefd5e5SDebojyoti Ghosh     const PetscReal
22917f6384fSBarry Smith       A[7][7]   = {{0,0,0,0,0,0,0},
2304e82626cSLisandro Dalcin                    {RC(1.0)/RC(5.0),0,0,0,0,0,0},
2314e82626cSLisandro Dalcin                    {RC(3.0)/RC(40.0),RC(9.0)/RC(40.0),0,0,0,0,0},
2324e82626cSLisandro Dalcin                    {RC(44.0)/RC(45.0),RC(-56.0)/RC(15.0),RC(32.0)/RC(9.0),0,0,0,0},
2334e82626cSLisandro 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},
2344e82626cSLisandro 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},
2354e82626cSLisandro 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}},
2364e82626cSLisandro 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},
2374e82626cSLisandro 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)};
2383a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
239f68a32c8SEmil Constantinescu   }
24005e23783SLisandro Dalcin   {
24105e23783SLisandro Dalcin     const PetscReal
24217f6384fSBarry Smith       A[8][8]   = {{0,0,0,0,0,0,0,0},
2434e82626cSLisandro Dalcin                    {RC(1.0)/RC(6.0),0,0,0,0,0,0,0},
2444e82626cSLisandro Dalcin                    {RC(2.0)/RC(27.0),RC(4.0)/RC(27.0),0,0,0,0,0,0},
2454e82626cSLisandro Dalcin                    {RC(183.0)/RC(1372.0),RC(-162.0)/RC(343.0),RC(1053.0)/RC(1372.0),0,0,0,0,0},
2464e82626cSLisandro 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},
2474e82626cSLisandro 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},
2484e82626cSLisandro 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},
2494e82626cSLisandro 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}},
2504e82626cSLisandro 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},
2514e82626cSLisandro 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)};
25205e23783SLisandro Dalcin     ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
25305e23783SLisandro Dalcin   }
2544e82626cSLisandro Dalcin #undef RC
255db046809SBarry Smith   PetscFunctionReturn(0);
256db046809SBarry Smith }
257db046809SBarry Smith 
258f68a32c8SEmil Constantinescu /*@C
259f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
260f68a32c8SEmil Constantinescu 
261f68a32c8SEmil Constantinescu    Not Collective
262f68a32c8SEmil Constantinescu 
263f68a32c8SEmil Constantinescu    Level: advanced
264f68a32c8SEmil Constantinescu 
265f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
266f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
267f68a32c8SEmil Constantinescu @*/
268f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
269e4dd521cSBarry Smith {
270f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
271f68a32c8SEmil Constantinescu   RKTableauLink link;
272f68a32c8SEmil Constantinescu 
273f68a32c8SEmil Constantinescu   PetscFunctionBegin;
274f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
275f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
276f68a32c8SEmil Constantinescu     RKTableauList = link->next;
277f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
278f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
279f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
280f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
281f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
282f68a32c8SEmil Constantinescu   }
283f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
284f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
285f68a32c8SEmil Constantinescu }
286f68a32c8SEmil Constantinescu 
287f68a32c8SEmil Constantinescu /*@C
288f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
289f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
290f68a32c8SEmil Constantinescu   when using static libraries.
291f68a32c8SEmil Constantinescu 
292f68a32c8SEmil Constantinescu   Level: developer
293f68a32c8SEmil Constantinescu 
294f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
295f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
296f68a32c8SEmil Constantinescu @*/
297f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
298f68a32c8SEmil Constantinescu {
299186e87acSLisandro Dalcin   PetscErrorCode ierr;
300e4dd521cSBarry Smith 
301e4dd521cSBarry Smith   PetscFunctionBegin;
302f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
303f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
304f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
305f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
306f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
307f68a32c8SEmil Constantinescu }
308186e87acSLisandro Dalcin 
309f68a32c8SEmil Constantinescu /*@C
310f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
311f68a32c8SEmil Constantinescu   called from PetscFinalize().
312186e87acSLisandro Dalcin 
313f68a32c8SEmil Constantinescu   Level: developer
314f68a32c8SEmil Constantinescu 
315f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
316f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
317f68a32c8SEmil Constantinescu @*/
318f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
319f68a32c8SEmil Constantinescu {
320f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
321f68a32c8SEmil Constantinescu 
322f68a32c8SEmil Constantinescu   PetscFunctionBegin;
323f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
324f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
325f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
326f68a32c8SEmil Constantinescu }
327f68a32c8SEmil Constantinescu 
328f68a32c8SEmil Constantinescu /*@C
329f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
330f68a32c8SEmil Constantinescu 
331f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
332f68a32c8SEmil Constantinescu 
333f68a32c8SEmil Constantinescu    Input Parameters:
334f68a32c8SEmil Constantinescu +  name - identifier for method
335f68a32c8SEmil Constantinescu .  order - approximation order of method
336f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
337f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
338f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
339f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
340f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
3413a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
3423a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
343f68a32c8SEmil Constantinescu 
344f68a32c8SEmil Constantinescu    Notes:
345f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
346f68a32c8SEmil Constantinescu 
347f68a32c8SEmil Constantinescu    Level: advanced
348f68a32c8SEmil Constantinescu 
349f68a32c8SEmil Constantinescu .keywords: TS, register
350f68a32c8SEmil Constantinescu 
351f68a32c8SEmil Constantinescu .seealso: TSRK
352f68a32c8SEmil Constantinescu @*/
353f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
354f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3553a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
356f68a32c8SEmil Constantinescu {
357f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
358f68a32c8SEmil Constantinescu   RKTableauLink   link;
359f68a32c8SEmil Constantinescu   RKTableau       t;
360f68a32c8SEmil Constantinescu   PetscInt        i,j;
361f68a32c8SEmil Constantinescu 
362f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3633a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3643a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3653a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3663a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3673a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3683a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3693a8a9803SLisandro Dalcin 
37095dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
371f68a32c8SEmil Constantinescu   t = &link->tab;
3723a8a9803SLisandro Dalcin 
373f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
374f68a32c8SEmil Constantinescu   t->order = order;
375f68a32c8SEmil Constantinescu   t->s = s;
376dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
377f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
378f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
379f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
380f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
381f68a32c8SEmil 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];
382d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
383d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3843a8a9803SLisandro Dalcin 
385f68a32c8SEmil Constantinescu   if (bembed) {
386785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
387f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
388f68a32c8SEmil Constantinescu   }
389f68a32c8SEmil Constantinescu 
3903a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3913a8a9803SLisandro Dalcin   t->p = p;
3923a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3933a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3943a8a9803SLisandro Dalcin 
395f68a32c8SEmil Constantinescu   link->next = RKTableauList;
396f68a32c8SEmil Constantinescu   RKTableauList = link;
397f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
398f68a32c8SEmil Constantinescu }
399f68a32c8SEmil Constantinescu 
400e4dd521cSBarry Smith /*
401f68a32c8SEmil Constantinescu  The step completion formula is
402e4dd521cSBarry Smith 
403f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
404f68a32c8SEmil Constantinescu 
405f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
406f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
407f68a32c8SEmil Constantinescu  We can write
408f68a32c8SEmil Constantinescu 
409f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
410f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
411f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
412f68a32c8SEmil Constantinescu 
413f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
414e4dd521cSBarry Smith */
415f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
416f68a32c8SEmil Constantinescu {
417f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
418f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
419f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
420f68a32c8SEmil Constantinescu   PetscReal      h;
421f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
422f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu   PetscFunctionBegin;
425f68a32c8SEmil Constantinescu   switch (rk->status) {
426f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
427f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
428f68a32c8SEmil Constantinescu     h = ts->time_step; break;
429f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
430be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
431f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
432f68a32c8SEmil Constantinescu   }
433f68a32c8SEmil Constantinescu   if (order == tab->order) {
434f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
435f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
436f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
437f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
438f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
439f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
440f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
441f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
442f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
443f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
444f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
445f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
446f68a32c8SEmil Constantinescu     } else { /* Rollback and re-complete using (be-b) */
447f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
448f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
449f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
4500f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
4510f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4520f7a1166SHong Zhang       }
453f68a32c8SEmil Constantinescu     }
454f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
455f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
456f68a32c8SEmil Constantinescu   }
457f68a32c8SEmil Constantinescu unavailable:
458f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
459a7fac7c2SEmil 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);
460f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
461f68a32c8SEmil Constantinescu }
462f68a32c8SEmil Constantinescu 
4630f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4640f7a1166SHong Zhang {
4650f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4660f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4670f7a1166SHong Zhang   const PetscInt  s = tab->s;
4680f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4690f7a1166SHong Zhang   Vec             *Y = rk->Y;
4700f7a1166SHong Zhang   PetscInt        i;
4710f7a1166SHong Zhang   PetscErrorCode  ierr;
4720f7a1166SHong Zhang 
4730f7a1166SHong Zhang   PetscFunctionBegin;
4740f7a1166SHong Zhang   /* backup cost integral */
4750f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4760f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4770f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
478*022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4790f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4800f7a1166SHong Zhang   }
4810f7a1166SHong Zhang   PetscFunctionReturn(0);
4820f7a1166SHong Zhang }
4830f7a1166SHong Zhang 
4840f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4850f7a1166SHong Zhang {
4860f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4870f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4880f7a1166SHong Zhang   const PetscInt  s = tab->s;
4890f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4900f7a1166SHong Zhang   Vec             *Y = rk->Y;
4910f7a1166SHong Zhang   PetscInt        i;
4920f7a1166SHong Zhang   PetscErrorCode  ierr;
4930f7a1166SHong Zhang 
4940f7a1166SHong Zhang   PetscFunctionBegin;
4950f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4960f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
497*022c081aSHong Zhang     ierr = TSComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4980f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4990f7a1166SHong Zhang   }
5000f7a1166SHong Zhang   PetscFunctionReturn(0);
5010f7a1166SHong Zhang }
5020f7a1166SHong Zhang 
503fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
504fecfb714SLisandro Dalcin {
505fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
506fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
507fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
508fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
509fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
510fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
511fecfb714SLisandro Dalcin   PetscInt        j;
512fecfb714SLisandro Dalcin   PetscReal       h;
513fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
514fecfb714SLisandro Dalcin 
515fecfb714SLisandro Dalcin   PetscFunctionBegin;
516fecfb714SLisandro Dalcin   switch (rk->status) {
517fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
518fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
519fecfb714SLisandro Dalcin     h = ts->time_step; break;
520fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
521fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
522fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
523fecfb714SLisandro Dalcin   }
524fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
525fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
526fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
527fecfb714SLisandro Dalcin }
528fecfb714SLisandro Dalcin 
529f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
530f68a32c8SEmil Constantinescu {
531f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
532f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
533f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
534fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
535f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
536f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
537d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
538f68a32c8SEmil Constantinescu   TSAdapt          adapt;
539fecfb714SLisandro Dalcin   PetscInt         i,j;
540be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
541be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
542be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
543f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
544f68a32c8SEmil Constantinescu 
545f68a32c8SEmil Constantinescu   PetscFunctionBegin;
546d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
547d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
548d1905564SLisandro Dalcin 
549f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
550be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
551be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
552f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
553f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5549be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5559be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
556f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
557f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
558f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5599be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
560f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
561be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
562fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
5638f738a7cSLisandro Dalcin       if (FSAL && !i) continue;
564f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
565f68a32c8SEmil Constantinescu     }
566be5899b3SLisandro Dalcin 
567be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
568f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
569f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
570f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
571f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
5721917a363SLisandro Dalcin     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr);
573fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
574be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
575be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
576fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
577be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
578be5899b3SLisandro Dalcin       goto reject_step;
579be5899b3SLisandro Dalcin     }
580be5899b3SLisandro Dalcin 
5810f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
5820f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5830f7a1166SHong Zhang       rk->time_step = ts->time_step;
584493ed95fSHong Zhang     }
585be5899b3SLisandro Dalcin 
586f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
587f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
588f68a32c8SEmil Constantinescu     break;
589be5899b3SLisandro Dalcin 
590be5899b3SLisandro Dalcin   reject_step:
591fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
592be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
593be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
594be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
595f68a32c8SEmil Constantinescu     }
596f68a32c8SEmil Constantinescu   }
597f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
598e4dd521cSBarry Smith }
599e4dd521cSBarry Smith 
600f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
601f6a906c0SBarry Smith {
602f6a906c0SBarry Smith   TS_RK         *rk  = (TS_RK*)ts->data;
603be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
604be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
605f6a906c0SBarry Smith   PetscErrorCode ierr;
606f6a906c0SBarry Smith 
607f6a906c0SBarry Smith   PetscFunctionBegin;
608f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
609abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
610f6a906c0SBarry Smith   if(ts->vecs_sensip) {
611abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
612f6a906c0SBarry Smith   }
613abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
614f6a906c0SBarry Smith   PetscFunctionReturn(0);
615f6a906c0SBarry Smith }
616f6a906c0SBarry Smith 
61742f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
618d2daff3dSHong Zhang {
619c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
620c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
621c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
622c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
623c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
624ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
625b18ea86cSHong Zhang   PetscInt         i,j,nadj;
6263d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
627d2daff3dSHong Zhang   PetscErrorCode   ierr;
628cef76868SBarry Smith   PetscReal        h = ts->time_step;
629cef76868SBarry Smith   Mat              J,Jp;
630c235aa8dSHong Zhang 
631d2daff3dSHong Zhang   PetscFunctionBegin;
632c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
633c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
634c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
635abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
636b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
637302440fdSBarry Smith       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
638c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
639302440fdSBarry Smith         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
640b18ea86cSHong Zhang       }
641c235aa8dSHong Zhang     }
642ad8e2604SHong Zhang     /* Stage values of lambda */
643c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
644c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
645493ed95fSHong Zhang     if (ts->vec_costintegral) {
646493ed95fSHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
647493ed95fSHong Zhang     }
648abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
649b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
650493ed95fSHong Zhang       if (ts->vec_costintegral) {
651493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
652493ed95fSHong Zhang       }
653b18ea86cSHong Zhang     }
6546497ce14SHong Zhang 
655ad8e2604SHong Zhang     /* Stage values of mu */
6566497ce14SHong Zhang     if(ts->vecs_sensip) {
6575bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
658493ed95fSHong Zhang       if (ts->vec_costintegral) {
659493ed95fSHong Zhang         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
660493ed95fSHong Zhang       }
661493ed95fSHong Zhang 
662abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
663ad8e2604SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
664493ed95fSHong Zhang         if (ts->vec_costintegral) {
665493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
666493ed95fSHong Zhang         }
667ad8e2604SHong Zhang       }
668c235aa8dSHong Zhang     }
6696497ce14SHong Zhang   }
670c235aa8dSHong Zhang 
671c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
672abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
673b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6746497ce14SHong Zhang     if(ts->vecs_sensip) {
675ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
676b18ea86cSHong Zhang     }
6776497ce14SHong Zhang   }
678c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
679d2daff3dSHong Zhang   PetscFunctionReturn(0);
680d2daff3dSHong Zhang }
681d2daff3dSHong Zhang 
682f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
683f68a32c8SEmil Constantinescu {
684f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
6853a8a9803SLisandro Dalcin   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
686f68a32c8SEmil Constantinescu   PetscReal        h;
687f68a32c8SEmil Constantinescu   PetscReal        tt,t;
688f68a32c8SEmil Constantinescu   PetscScalar     *b;
689f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
690f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
691e4dd521cSBarry Smith 
692f68a32c8SEmil Constantinescu   PetscFunctionBegin;
693f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
694be5899b3SLisandro Dalcin 
695f68a32c8SEmil Constantinescu   switch (rk->status) {
696f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
697f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
698f68a32c8SEmil Constantinescu     h = ts->time_step;
699f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
700f68a32c8SEmil Constantinescu     break;
701f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
702be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
703f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
704f68a32c8SEmil Constantinescu     break;
705f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
706e4dd521cSBarry Smith   }
707785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
708f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
7093a8a9803SLisandro Dalcin   for (j=0,tt=t; j<p; j++,tt*=t) {
710f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
7113a8a9803SLisandro Dalcin       b[i]  += h * B[i*p+j] * tt;
712e4dd521cSBarry Smith     }
713f68a32c8SEmil Constantinescu   }
714be5899b3SLisandro Dalcin 
715f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
716f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
717be5899b3SLisandro Dalcin 
718f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
719e4dd521cSBarry Smith   PetscFunctionReturn(0);
720e4dd521cSBarry Smith }
721e4dd521cSBarry Smith 
722e4dd521cSBarry Smith /*------------------------------------------------------------*/
723be5899b3SLisandro Dalcin 
724be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
725be5899b3SLisandro Dalcin {
726be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
727be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
728be5899b3SLisandro Dalcin   PetscErrorCode ierr;
729be5899b3SLisandro Dalcin 
730be5899b3SLisandro Dalcin   PetscFunctionBegin;
731be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
732be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
733be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
734be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
735be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
736be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
737be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
738be5899b3SLisandro Dalcin }
739be5899b3SLisandro Dalcin 
740277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
741e4dd521cSBarry Smith {
7425f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7436849ba73SBarry Smith   PetscErrorCode ierr;
744e4dd521cSBarry Smith 
745e4dd521cSBarry Smith   PetscFunctionBegin;
746be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7470f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
748abc2977eSBarry Smith   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
749277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
750e4dd521cSBarry Smith }
751277b19d0SLisandro Dalcin 
752277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
753277b19d0SLisandro Dalcin {
754277b19d0SLisandro Dalcin   PetscErrorCode ierr;
755277b19d0SLisandro Dalcin 
756277b19d0SLisandro Dalcin   PetscFunctionBegin;
757277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
758277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
759f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
760f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
761f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
762f68a32c8SEmil Constantinescu }
763f68a32c8SEmil Constantinescu 
764f68a32c8SEmil Constantinescu 
765f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
766f68a32c8SEmil Constantinescu {
767f68a32c8SEmil Constantinescu   PetscFunctionBegin;
768f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
769f68a32c8SEmil Constantinescu }
770f68a32c8SEmil Constantinescu 
771f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
772f68a32c8SEmil Constantinescu {
773f68a32c8SEmil Constantinescu   PetscFunctionBegin;
774f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
775f68a32c8SEmil Constantinescu }
776f68a32c8SEmil Constantinescu 
777f68a32c8SEmil Constantinescu 
778f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
779f68a32c8SEmil Constantinescu {
780f68a32c8SEmil Constantinescu   PetscFunctionBegin;
781f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
782f68a32c8SEmil Constantinescu }
783f68a32c8SEmil Constantinescu 
784f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
785f68a32c8SEmil Constantinescu {
786f68a32c8SEmil Constantinescu 
787f68a32c8SEmil Constantinescu   PetscFunctionBegin;
788f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
789f68a32c8SEmil Constantinescu }
790c235aa8dSHong Zhang /*
791d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
792d2daff3dSHong Zhang {
793d2daff3dSHong Zhang   PetscReal *A,*b;
794d2daff3dSHong Zhang   PetscInt        s,i,j;
795d2daff3dSHong Zhang   PetscErrorCode  ierr;
796d2daff3dSHong Zhang 
797d2daff3dSHong Zhang   PetscFunctionBegin;
798d2daff3dSHong Zhang   s    = tab->s;
799d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
800d2daff3dSHong Zhang 
801d2daff3dSHong Zhang   for (i=0; i<s; i++)
802d2daff3dSHong Zhang     for (j=0; j<s; j++) {
803d2daff3dSHong 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];
804d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
805d2daff3dSHong Zhang     }
806d2daff3dSHong Zhang 
807d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
808d2daff3dSHong Zhang 
809d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
810d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
811d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
812d2daff3dSHong Zhang   PetscFunctionReturn(0);
813d2daff3dSHong Zhang }
814c235aa8dSHong Zhang */
815be5899b3SLisandro Dalcin 
816be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
817be5899b3SLisandro Dalcin {
818be5899b3SLisandro Dalcin   TS_RK         *rk  = (TS_RK*)ts->data;
819be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
820be5899b3SLisandro Dalcin   PetscErrorCode ierr;
821be5899b3SLisandro Dalcin 
822be5899b3SLisandro Dalcin   PetscFunctionBegin;
823be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
824be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
825be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
826be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
827be5899b3SLisandro Dalcin }
828be5899b3SLisandro Dalcin 
829be5899b3SLisandro Dalcin 
830f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
831f68a32c8SEmil Constantinescu {
832f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
833f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
834f68a32c8SEmil Constantinescu   DM             dm;
835f68a32c8SEmil Constantinescu 
836f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8373dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
838be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8390f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8400f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8410f7a1166SHong Zhang   }
842f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
843f68a32c8SEmil Constantinescu   if (dm) {
844f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
845f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
846f68a32c8SEmil Constantinescu   }
847e4dd521cSBarry Smith   PetscFunctionReturn(0);
848e4dd521cSBarry Smith }
849d2daff3dSHong Zhang 
850f6a906c0SBarry Smith 
851e4dd521cSBarry Smith /*------------------------------------------------------------*/
852e4dd521cSBarry Smith 
8534416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
854e4dd521cSBarry Smith {
855be5899b3SLisandro Dalcin   TS_RK         *rk = (TS_RK*)ts->data;
856dfbe8321SBarry Smith   PetscErrorCode ierr;
857262ff7bbSBarry Smith 
858e4dd521cSBarry Smith   PetscFunctionBegin;
859e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
860f68a32c8SEmil Constantinescu   {
861f68a32c8SEmil Constantinescu     RKTableauLink  link;
862f68a32c8SEmil Constantinescu     PetscInt       count,choice;
863f68a32c8SEmil Constantinescu     PetscBool      flg;
864f68a32c8SEmil Constantinescu     const char   **namelist;
865f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
866785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
867f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
868be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
869be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
870f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
871f68a32c8SEmil Constantinescu   }
872262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
873e4dd521cSBarry Smith   PetscFunctionReturn(0);
874e4dd521cSBarry Smith }
875e4dd521cSBarry Smith 
876f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
877f68a32c8SEmil Constantinescu {
878f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
879f68a32c8SEmil Constantinescu   PetscInt       i;
880f68a32c8SEmil Constantinescu   size_t         left,count;
881f68a32c8SEmil Constantinescu   char           *p;
882f68a32c8SEmil Constantinescu 
883f68a32c8SEmil Constantinescu   PetscFunctionBegin;
884f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
885f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
886f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
887f68a32c8SEmil Constantinescu     left -= count;
888f68a32c8SEmil Constantinescu     p    += count;
889f68a32c8SEmil Constantinescu     *p++  = ' ';
890f68a32c8SEmil Constantinescu   }
891f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
892f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
893f68a32c8SEmil Constantinescu }
894f68a32c8SEmil Constantinescu 
8955f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
896e4dd521cSBarry Smith {
8975f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8988370ee3bSLisandro Dalcin   PetscBool      iascii;
899dfbe8321SBarry Smith   PetscErrorCode ierr;
9002cdf8120Sasbjorn 
901e4dd521cSBarry Smith   PetscFunctionBegin;
902251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
9038370ee3bSLisandro Dalcin   if (iascii) {
9049c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
905f68a32c8SEmil Constantinescu     TSRKType  rktype;
906f68a32c8SEmil Constantinescu     char      buf[512];
907f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
908efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  RK type %s\n",rktype);CHKERRQ(ierr);
9090c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
9100c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
911f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
912f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
9138370ee3bSLisandro Dalcin   }
914f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
915f68a32c8SEmil Constantinescu }
916f68a32c8SEmil Constantinescu 
917f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
918f68a32c8SEmil Constantinescu {
919f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
9209c334d8fSLisandro Dalcin   TSAdapt        adapt;
921f68a32c8SEmil Constantinescu 
922f68a32c8SEmil Constantinescu   PetscFunctionBegin;
9239c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9249c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
925f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
926f68a32c8SEmil Constantinescu }
927f68a32c8SEmil Constantinescu 
928f68a32c8SEmil Constantinescu /*@C
929f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
930f68a32c8SEmil Constantinescu 
931f68a32c8SEmil Constantinescu   Logically collective
932f68a32c8SEmil Constantinescu 
933f68a32c8SEmil Constantinescu   Input Parameter:
934f68a32c8SEmil Constantinescu +  ts - timestepping context
935f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
936f68a32c8SEmil Constantinescu 
937287c2655SBarry Smith   Options Database:
9389bd3e852SBarry Smith .   -ts_rk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs>
939287c2655SBarry Smith 
940f68a32c8SEmil Constantinescu   Level: intermediate
941f68a32c8SEmil Constantinescu 
942287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS
943f68a32c8SEmil Constantinescu @*/
944f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
945f68a32c8SEmil Constantinescu {
946f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
947f68a32c8SEmil Constantinescu 
948f68a32c8SEmil Constantinescu   PetscFunctionBegin;
949f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
950b92453a8SLisandro Dalcin   PetscValidCharPointer(rktype,2);
951f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
952f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
953f68a32c8SEmil Constantinescu }
954f68a32c8SEmil Constantinescu 
955f68a32c8SEmil Constantinescu /*@C
956f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
957f68a32c8SEmil Constantinescu 
958f68a32c8SEmil Constantinescu   Logically collective
959f68a32c8SEmil Constantinescu 
960f68a32c8SEmil Constantinescu   Input Parameter:
961f68a32c8SEmil Constantinescu .  ts - timestepping context
962f68a32c8SEmil Constantinescu 
963f68a32c8SEmil Constantinescu   Output Parameter:
964f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
965f68a32c8SEmil Constantinescu 
966f68a32c8SEmil Constantinescu   Level: intermediate
967f68a32c8SEmil Constantinescu 
968f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
969f68a32c8SEmil Constantinescu @*/
970f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
971f68a32c8SEmil Constantinescu {
972f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
973f68a32c8SEmil Constantinescu 
974f68a32c8SEmil Constantinescu   PetscFunctionBegin;
975f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
976f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
977f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
978f68a32c8SEmil Constantinescu }
979f68a32c8SEmil Constantinescu 
980560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
981f68a32c8SEmil Constantinescu {
982f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
983f68a32c8SEmil Constantinescu 
984f68a32c8SEmil Constantinescu   PetscFunctionBegin;
985f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
986f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
987f68a32c8SEmil Constantinescu }
988560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
989f68a32c8SEmil Constantinescu {
990f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
991f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
992f68a32c8SEmil Constantinescu   PetscBool      match;
993f68a32c8SEmil Constantinescu   RKTableauLink  link;
994f68a32c8SEmil Constantinescu 
995f68a32c8SEmil Constantinescu   PetscFunctionBegin;
996f68a32c8SEmil Constantinescu   if (rk->tableau) {
997f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
998f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
999f68a32c8SEmil Constantinescu   }
1000f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
1001f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
1002f68a32c8SEmil Constantinescu     if (match) {
1003be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
1004f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
1005be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
10062ffb9264SLisandro Dalcin       ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1007f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
1008f68a32c8SEmil Constantinescu     }
1009f68a32c8SEmil Constantinescu   }
1010f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
1011e4dd521cSBarry Smith   PetscFunctionReturn(0);
1012e4dd521cSBarry Smith }
1013e4dd521cSBarry Smith 
1014ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1015ff22ae23SHong Zhang {
1016ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1017ff22ae23SHong Zhang 
1018ff22ae23SHong Zhang   PetscFunctionBegin;
1019ff22ae23SHong Zhang   *ns = rk->tableau->s;
1020d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
1021ff22ae23SHong Zhang   PetscFunctionReturn(0);
1022ff22ae23SHong Zhang }
1023ff22ae23SHong Zhang 
1024ff22ae23SHong Zhang 
1025e4dd521cSBarry Smith /* ------------------------------------------------------------ */
1026ebe3b25bSBarry Smith /*MC
1027f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1028ebe3b25bSBarry Smith 
1029f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1030f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1031ebe3b25bSBarry Smith 
1032f68a32c8SEmil Constantinescu   Notes:
103398164e13SLisandro Dalcin   The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type
1034ebe3b25bSBarry Smith 
1035d41222bbSBarry Smith   Level: beginner
1036d41222bbSBarry Smith 
1037f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1038f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
1039ebe3b25bSBarry Smith 
1040ebe3b25bSBarry Smith M*/
10418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1042e4dd521cSBarry Smith {
10431566a47fSLisandro Dalcin   TS_RK          *rk;
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
1045e4dd521cSBarry Smith 
1046e4dd521cSBarry Smith   PetscFunctionBegin;
1047f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1048f68a32c8SEmil Constantinescu 
1049f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10505f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10515f70b478SJed Brown   ts->ops->view           = TSView_RK;
1052f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1053f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
105442f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1055f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1056f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1057f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1058fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1059f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1060ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
106142f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1062e4dd521cSBarry Smith 
10630f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10640f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10650f7a1166SHong Zhang 
10661566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10671566a47fSLisandro Dalcin   ts->data = (void*)rk;
1068e4dd521cSBarry Smith 
1069f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1070f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1071be5899b3SLisandro Dalcin 
1072be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1073e4dd521cSBarry Smith   PetscFunctionReturn(0);
1074e4dd521cSBarry Smith }
1075