xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 3dd83b3843bef2e44404f8ac40307e9aa6509644)
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                                           */
22d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
23f68a32c8SEmil Constantinescu   PetscInt   pinterp;             /* Interpolation order                                        */
24f68a32c8SEmil Constantinescu   PetscReal *A,*b,*c;             /* Tableau                                                    */
25f68a32c8SEmil Constantinescu   PetscReal *bembed;              /* Embedded formula of order one less (order-1)               */
26f68a32c8SEmil Constantinescu   PetscReal *binterp;             /* Dense output formula                                       */
27f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
28f68a32c8SEmil Constantinescu };
29f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
30f68a32c8SEmil Constantinescu struct _RKTableauLink {
31f68a32c8SEmil Constantinescu   struct _RKTableau tab;
32f68a32c8SEmil Constantinescu   RKTableauLink     next;
33f68a32c8SEmil Constantinescu };
34f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
35e4dd521cSBarry Smith 
36e4dd521cSBarry Smith typedef struct {
37f68a32c8SEmil Constantinescu   RKTableau    tableau;
38f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step */
39f68a32c8SEmil Constantinescu   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
40ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage */
41ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage */
42b18ea86cSHong Zhang   Vec          *VecSensiTemp;    /* Vector to be timed with Jacobian transpose */
430f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events */
44f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
45f68a32c8SEmil Constantinescu   PetscReal    stage_time;
46f68a32c8SEmil Constantinescu   TSStepStatus status;
470f7a1166SHong Zhang   PetscReal    ptime;
480f7a1166SHong Zhang   PetscReal    time_step;
495f70b478SJed Brown } TS_RK;
50e4dd521cSBarry Smith 
51f68a32c8SEmil Constantinescu /*MC
52f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
53262ff7bbSBarry Smith 
54f68a32c8SEmil Constantinescu      This method has one stage.
55f68a32c8SEmil Constantinescu 
56f68a32c8SEmil Constantinescu      Level: advanced
57f68a32c8SEmil Constantinescu 
58f68a32c8SEmil Constantinescu .seealso: TSRK
59f68a32c8SEmil Constantinescu M*/
60f68a32c8SEmil Constantinescu /*MC
612109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
62f68a32c8SEmil Constantinescu 
63f68a32c8SEmil Constantinescu      This method has two stages.
64f68a32c8SEmil Constantinescu 
65f68a32c8SEmil Constantinescu      Level: advanced
66f68a32c8SEmil Constantinescu 
67f68a32c8SEmil Constantinescu .seealso: TSRK
68f68a32c8SEmil Constantinescu M*/
69f68a32c8SEmil Constantinescu /*MC
70f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
71f68a32c8SEmil Constantinescu 
72f68a32c8SEmil Constantinescu      This method has three stages.
73f68a32c8SEmil Constantinescu 
74f68a32c8SEmil Constantinescu      Level: advanced
75f68a32c8SEmil Constantinescu 
76f68a32c8SEmil Constantinescu .seealso: TSRK
77f68a32c8SEmil Constantinescu M*/
78f68a32c8SEmil Constantinescu /*MC
792109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
802109b73fSDebojyoti Ghosh 
812109b73fSDebojyoti Ghosh      This method has four stages.
822109b73fSDebojyoti Ghosh 
832109b73fSDebojyoti Ghosh      Level: advanced
842109b73fSDebojyoti Ghosh 
852109b73fSDebojyoti Ghosh .seealso: TSRK
862109b73fSDebojyoti Ghosh M*/
872109b73fSDebojyoti Ghosh /*MC
88f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
89f68a32c8SEmil Constantinescu 
902e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
91f68a32c8SEmil Constantinescu 
92f68a32c8SEmil Constantinescu      Level: advanced
93f68a32c8SEmil Constantinescu 
94f68a32c8SEmil Constantinescu .seealso: TSRK
95f68a32c8SEmil Constantinescu M*/
96f68a32c8SEmil Constantinescu /*MC
972e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
98f68a32c8SEmil Constantinescu 
99f68a32c8SEmil Constantinescu      This method has six stages.
100f68a32c8SEmil Constantinescu 
101f68a32c8SEmil Constantinescu      Level: advanced
102f68a32c8SEmil Constantinescu 
103f68a32c8SEmil Constantinescu .seealso: TSRK
104f68a32c8SEmil Constantinescu M*/
1052109b73fSDebojyoti Ghosh /*MC
1062e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1072109b73fSDebojyoti Ghosh 
1082109b73fSDebojyoti Ghosh      This method has seven stages.
1092109b73fSDebojyoti Ghosh 
1102109b73fSDebojyoti Ghosh      Level: advanced
1112109b73fSDebojyoti Ghosh 
1122109b73fSDebojyoti Ghosh .seealso: TSRK
1132109b73fSDebojyoti Ghosh M*/
114262ff7bbSBarry Smith 
115262ff7bbSBarry Smith #undef __FUNCT__
116f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll"
117f68a32c8SEmil Constantinescu /*@C
118f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
119262ff7bbSBarry Smith 
120f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
121262ff7bbSBarry Smith 
122f68a32c8SEmil Constantinescu   Level: advanced
123262ff7bbSBarry Smith 
124f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
125262ff7bbSBarry Smith 
126f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
127262ff7bbSBarry Smith @*/
128f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
129262ff7bbSBarry Smith {
1304ac538c5SBarry Smith   PetscErrorCode ierr;
131262ff7bbSBarry Smith 
132262ff7bbSBarry Smith   PetscFunctionBegin;
133f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
134f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
135e4dd521cSBarry Smith 
136e4dd521cSBarry Smith   {
137f68a32c8SEmil Constantinescu     const PetscReal
138f68a32c8SEmil Constantinescu       A[1][1] = {{0.0}},
139f68a32c8SEmil Constantinescu       b[1]    = {1.0};
140f68a32c8SEmil Constantinescu     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
141e4dd521cSBarry Smith   }
142db046809SBarry Smith   {
143f68a32c8SEmil Constantinescu     const PetscReal
144f68a32c8SEmil Constantinescu       A[2][2]     = {{0.0,0.0},
145f68a32c8SEmil Constantinescu                     {1.0,0.0}},
146f68a32c8SEmil Constantinescu       b[2]        = {0.5,0.5},
147f68a32c8SEmil Constantinescu       bembed[2]   = {1.0,0};
1488886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr);
149db046809SBarry Smith   }
150f68a32c8SEmil Constantinescu   {
151f68a32c8SEmil Constantinescu     const PetscReal
152f68a32c8SEmil Constantinescu       A[3][3] = {{0,0,0},
153f68a32c8SEmil Constantinescu                  {2.0/3.0,0,0},
154f68a32c8SEmil Constantinescu                  {-1.0/3.0,1.0,0}},
155f68a32c8SEmil Constantinescu       b[3]    = {0.25,0.5,0.25};
1568886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
157fdefd5e5SDebojyoti Ghosh   }
158fdefd5e5SDebojyoti Ghosh   {
159fdefd5e5SDebojyoti Ghosh     const PetscReal
160fdefd5e5SDebojyoti Ghosh       A[4][4] = {{0,0,0,0},
161fdefd5e5SDebojyoti Ghosh                  {1.0/2.0,0,0,0},
162fdefd5e5SDebojyoti Ghosh                  {0,3.0/4.0,0,0},
163fdefd5e5SDebojyoti Ghosh                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
164fdefd5e5SDebojyoti Ghosh       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
165fdefd5e5SDebojyoti Ghosh       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
1668886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr);
167db046809SBarry Smith   }
168f68a32c8SEmil Constantinescu   {
169f68a32c8SEmil Constantinescu     const PetscReal
170f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
171f68a32c8SEmil Constantinescu                  {0.5,0,0,0},
172f68a32c8SEmil Constantinescu                  {0,0.5,0,0},
173f68a32c8SEmil Constantinescu                  {0,0,1.0,0}},
174f68a32c8SEmil Constantinescu       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
1758886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
176db046809SBarry Smith   }
177f68a32c8SEmil Constantinescu   {
178f68a32c8SEmil Constantinescu     const PetscReal
179f68a32c8SEmil Constantinescu       A[6][6]   = {{0,0,0,0,0,0},
180f68a32c8SEmil Constantinescu                    {0.25,0,0,0,0,0},
181f68a32c8SEmil Constantinescu                    {3.0/32.0,9.0/32.0,0,0,0,0},
182f68a32c8SEmil Constantinescu                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
183f68a32c8SEmil Constantinescu                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
184f68a32c8SEmil Constantinescu                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
185f68a32c8SEmil Constantinescu       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
186f68a32c8SEmil Constantinescu       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
1878886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr);
188fdefd5e5SDebojyoti Ghosh   }
189fdefd5e5SDebojyoti Ghosh   {
190fdefd5e5SDebojyoti Ghosh     const PetscReal
191fdefd5e5SDebojyoti Ghosh       A[7][7]   = {{0,0,0,0,0,0,0},
192fdefd5e5SDebojyoti Ghosh                    {1.0/5.0,0,0,0,0,0,0},
193fdefd5e5SDebojyoti Ghosh                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
194fdefd5e5SDebojyoti Ghosh                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
195fdefd5e5SDebojyoti Ghosh                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
196fdefd5e5SDebojyoti Ghosh                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
197fdefd5e5SDebojyoti Ghosh                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
198fdefd5e5SDebojyoti Ghosh       b[7]      = {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0},
199fdefd5e5SDebojyoti Ghosh       bembed[7] = {5179.0/57600.0,0,7571.0/16695.0,393.0/640.0,-92097.0/339200.0,187.0/2100.0,1.0/40.0};
2008886ffe3SEmil Constantinescu     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,1,b);CHKERRQ(ierr);
201f68a32c8SEmil Constantinescu   }
202db046809SBarry Smith   PetscFunctionReturn(0);
203db046809SBarry Smith }
204db046809SBarry Smith 
205e4dd521cSBarry Smith #undef __FUNCT__
206f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy"
207f68a32c8SEmil Constantinescu /*@C
208f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
209f68a32c8SEmil Constantinescu 
210f68a32c8SEmil Constantinescu    Not Collective
211f68a32c8SEmil Constantinescu 
212f68a32c8SEmil Constantinescu    Level: advanced
213f68a32c8SEmil Constantinescu 
214f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
215f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
216f68a32c8SEmil Constantinescu @*/
217f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
218e4dd521cSBarry Smith {
219f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
220f68a32c8SEmil Constantinescu   RKTableauLink link;
221f68a32c8SEmil Constantinescu 
222f68a32c8SEmil Constantinescu   PetscFunctionBegin;
223f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
224f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
225f68a32c8SEmil Constantinescu     RKTableauList = link->next;
226f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
227f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
228f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
229f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
230f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
231f68a32c8SEmil Constantinescu   }
232f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
233f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
234f68a32c8SEmil Constantinescu }
235f68a32c8SEmil Constantinescu 
236f68a32c8SEmil Constantinescu #undef __FUNCT__
237f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage"
238f68a32c8SEmil Constantinescu /*@C
239f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
240f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
241f68a32c8SEmil Constantinescu   when using static libraries.
242f68a32c8SEmil Constantinescu 
243f68a32c8SEmil Constantinescu   Level: developer
244f68a32c8SEmil Constantinescu 
245f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
246f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
247f68a32c8SEmil Constantinescu @*/
248f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
249f68a32c8SEmil Constantinescu {
250186e87acSLisandro Dalcin   PetscErrorCode ierr;
251e4dd521cSBarry Smith 
252e4dd521cSBarry Smith   PetscFunctionBegin;
253f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
254f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
255f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
256f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
258f68a32c8SEmil Constantinescu }
259186e87acSLisandro Dalcin 
260f68a32c8SEmil Constantinescu #undef __FUNCT__
261f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage"
262f68a32c8SEmil Constantinescu /*@C
263f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
264f68a32c8SEmil Constantinescu   called from PetscFinalize().
265186e87acSLisandro Dalcin 
266f68a32c8SEmil Constantinescu   Level: developer
267f68a32c8SEmil Constantinescu 
268f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
269f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
270f68a32c8SEmil Constantinescu @*/
271f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
272f68a32c8SEmil Constantinescu {
273f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
274f68a32c8SEmil Constantinescu 
275f68a32c8SEmil Constantinescu   PetscFunctionBegin;
276f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
277f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
278f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
279f68a32c8SEmil Constantinescu }
280f68a32c8SEmil Constantinescu 
281f68a32c8SEmil Constantinescu #undef __FUNCT__
282f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister"
283f68a32c8SEmil Constantinescu /*@C
284f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
285f68a32c8SEmil Constantinescu 
286f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
287f68a32c8SEmil Constantinescu 
288f68a32c8SEmil Constantinescu    Input Parameters:
289f68a32c8SEmil Constantinescu +  name - identifier for method
290f68a32c8SEmil Constantinescu .  order - approximation order of method
291f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
292f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
293f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
294f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
295f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
296f68a32c8SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterp
297f68a32c8SEmil Constantinescu -  binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt)
298f68a32c8SEmil Constantinescu 
299f68a32c8SEmil Constantinescu    Notes:
300f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
301f68a32c8SEmil Constantinescu 
302f68a32c8SEmil Constantinescu    Level: advanced
303f68a32c8SEmil Constantinescu 
304f68a32c8SEmil Constantinescu .keywords: TS, register
305f68a32c8SEmil Constantinescu 
306f68a32c8SEmil Constantinescu .seealso: TSRK
307f68a32c8SEmil Constantinescu @*/
308f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
309f68a32c8SEmil Constantinescu                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
310f68a32c8SEmil Constantinescu                                  const PetscReal bembed[],
311f68a32c8SEmil Constantinescu                                  PetscInt pinterp,const PetscReal binterp[])
312f68a32c8SEmil Constantinescu {
313f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
314f68a32c8SEmil Constantinescu   RKTableauLink   link;
315f68a32c8SEmil Constantinescu   RKTableau       t;
316f68a32c8SEmil Constantinescu   PetscInt        i,j;
317f68a32c8SEmil Constantinescu 
318f68a32c8SEmil Constantinescu   PetscFunctionBegin;
319f68a32c8SEmil Constantinescu   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
320f68a32c8SEmil Constantinescu   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
321f68a32c8SEmil Constantinescu   t        = &link->tab;
322f68a32c8SEmil Constantinescu   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
323f68a32c8SEmil Constantinescu   t->order = order;
324f68a32c8SEmil Constantinescu   t->s     = s;
325dcca6d9dSJed Brown   ierr     = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
326f68a32c8SEmil Constantinescu   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
327f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
328f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
329f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
330f68a32c8SEmil 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];
331d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
332d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
333f68a32c8SEmil Constantinescu   if (bembed) {
334785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
335f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
336f68a32c8SEmil Constantinescu   }
337f68a32c8SEmil Constantinescu 
338f68a32c8SEmil Constantinescu   t->pinterp     = pinterp;
339785e854fSJed Brown   ierr           = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
340f68a32c8SEmil Constantinescu   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
341f68a32c8SEmil Constantinescu   link->next     = RKTableauList;
342f68a32c8SEmil Constantinescu   RKTableauList = link;
343f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
344f68a32c8SEmil Constantinescu }
345f68a32c8SEmil Constantinescu 
346f68a32c8SEmil Constantinescu #undef __FUNCT__
347f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK"
348e4dd521cSBarry Smith /*
349f68a32c8SEmil Constantinescu  The step completion formula is
350e4dd521cSBarry Smith 
351f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
352f68a32c8SEmil Constantinescu 
353f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
354f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
355f68a32c8SEmil Constantinescu  We can write
356f68a32c8SEmil Constantinescu 
357f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
358f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
359f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
360f68a32c8SEmil Constantinescu 
361f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
362e4dd521cSBarry Smith */
363f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
364f68a32c8SEmil Constantinescu {
365f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
366f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
367f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
368f68a32c8SEmil Constantinescu   PetscReal      h;
369f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
370f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
371f68a32c8SEmil Constantinescu 
372f68a32c8SEmil Constantinescu   PetscFunctionBegin;
373f68a32c8SEmil Constantinescu   switch (rk->status) {
374f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
375f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
376f68a32c8SEmil Constantinescu     h = ts->time_step; break;
377f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
378be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
379f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
380f68a32c8SEmil Constantinescu   }
381f68a32c8SEmil Constantinescu   if (order == tab->order) {
382f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
383f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
384f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
385f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
386f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
387f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
388f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
389f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
390f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
391f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
392f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
393f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
394f68a32c8SEmil Constantinescu     } else { /* Rollback and re-complete using (be-b) */
395f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
396f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
397f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
3980f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
3990f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4000f7a1166SHong Zhang       }
401f68a32c8SEmil Constantinescu     }
402f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
403f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
404f68a32c8SEmil Constantinescu   }
405f68a32c8SEmil Constantinescu unavailable:
406f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
407a7fac7c2SEmil 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);
408f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
409f68a32c8SEmil Constantinescu }
410f68a32c8SEmil Constantinescu 
411f68a32c8SEmil Constantinescu #undef __FUNCT__
4120f7a1166SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_RK"
4130f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4140f7a1166SHong Zhang {
4150f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4160f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4170f7a1166SHong Zhang   const PetscInt  s = tab->s;
4180f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4190f7a1166SHong Zhang   Vec             *Y = rk->Y;
4200f7a1166SHong Zhang   PetscInt        i;
4210f7a1166SHong Zhang   PetscErrorCode  ierr;
4220f7a1166SHong Zhang 
4230f7a1166SHong Zhang   PetscFunctionBegin;
4240f7a1166SHong Zhang   /* backup cost integral */
4250f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4260f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4270f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4280f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4290f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4300f7a1166SHong Zhang   }
4310f7a1166SHong Zhang   PetscFunctionReturn(0);
4320f7a1166SHong Zhang }
4330f7a1166SHong Zhang 
4340f7a1166SHong Zhang #undef __FUNCT__
4350f7a1166SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_RK"
4360f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4370f7a1166SHong Zhang {
4380f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4390f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4400f7a1166SHong Zhang   const PetscInt  s = tab->s;
4410f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4420f7a1166SHong Zhang   Vec             *Y = rk->Y;
4430f7a1166SHong Zhang   PetscInt        i;
4440f7a1166SHong Zhang   PetscErrorCode  ierr;
4450f7a1166SHong Zhang 
4460f7a1166SHong Zhang   PetscFunctionBegin;
4470f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4480f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4490f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4500f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4510f7a1166SHong Zhang   }
4520f7a1166SHong Zhang   PetscFunctionReturn(0);
4530f7a1166SHong Zhang }
4540f7a1166SHong Zhang 
4550f7a1166SHong Zhang #undef __FUNCT__
456fecfb714SLisandro Dalcin #define __FUNCT__ "TSRollBack_RK"
457fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
458fecfb714SLisandro Dalcin {
459fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
460fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
461fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
462fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
463fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
464fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
465fecfb714SLisandro Dalcin   PetscInt        j;
466fecfb714SLisandro Dalcin   PetscReal       h;
467fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
468fecfb714SLisandro Dalcin 
469fecfb714SLisandro Dalcin   PetscFunctionBegin;
470fecfb714SLisandro Dalcin   switch (rk->status) {
471fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
472fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
473fecfb714SLisandro Dalcin     h = ts->time_step; break;
474fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
475fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
476fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
477fecfb714SLisandro Dalcin   }
478fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
479fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
480fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
481fecfb714SLisandro Dalcin }
482fecfb714SLisandro Dalcin 
483fecfb714SLisandro Dalcin #undef __FUNCT__
484f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK"
485f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
486f68a32c8SEmil Constantinescu {
487f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
488f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
489f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
490fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
491f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
492f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
493f68a32c8SEmil Constantinescu   TSAdapt          adapt;
494fecfb714SLisandro Dalcin   PetscInt         i,j;
495be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
496be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
497be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
498f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
499f68a32c8SEmil Constantinescu 
500f68a32c8SEmil Constantinescu   PetscFunctionBegin;
501f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
502be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
503be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
504f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
505f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
5069be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5079be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
508f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
509f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
510f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5119be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
512f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
513be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
514fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
515f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
516f68a32c8SEmil Constantinescu     }
517be5899b3SLisandro Dalcin 
518be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
519f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
520f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
521f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
522f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
523f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
524fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
525be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
526be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
527fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
528be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
529be5899b3SLisandro Dalcin       goto reject_step;
530be5899b3SLisandro Dalcin     }
531be5899b3SLisandro Dalcin 
5320f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
5330f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5340f7a1166SHong Zhang       rk->time_step = ts->time_step;
535493ed95fSHong Zhang     }
536be5899b3SLisandro Dalcin 
537f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
538f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
539f68a32c8SEmil Constantinescu     break;
540be5899b3SLisandro Dalcin 
541be5899b3SLisandro Dalcin   reject_step:
542fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
543be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
544be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
545be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
546f68a32c8SEmil Constantinescu     }
547f68a32c8SEmil Constantinescu   }
548f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
549e4dd521cSBarry Smith }
550e4dd521cSBarry Smith 
551f68a32c8SEmil Constantinescu #undef __FUNCT__
552f6a906c0SBarry Smith #define __FUNCT__ "TSAdjointSetUp_RK"
553f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
554f6a906c0SBarry Smith {
555f6a906c0SBarry Smith   TS_RK         *rk  = (TS_RK*)ts->data;
556be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
557be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
558f6a906c0SBarry Smith   PetscErrorCode ierr;
559f6a906c0SBarry Smith 
560f6a906c0SBarry Smith   PetscFunctionBegin;
561f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
562abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
563f6a906c0SBarry Smith   if(ts->vecs_sensip) {
564abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
565f6a906c0SBarry Smith   }
566abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
567f6a906c0SBarry Smith   PetscFunctionReturn(0);
568f6a906c0SBarry Smith }
569f6a906c0SBarry Smith 
570f6a906c0SBarry Smith #undef __FUNCT__
57142f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_RK"
57242f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
573d2daff3dSHong Zhang {
574c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
575c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
576c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
577c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
578c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
579ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
580b18ea86cSHong Zhang   PetscInt         i,j,nadj;
5813d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
582d2daff3dSHong Zhang   PetscErrorCode   ierr;
583cef76868SBarry Smith   PetscReal        h = ts->time_step;
584cef76868SBarry Smith   Mat              J,Jp;
585c235aa8dSHong Zhang 
586d2daff3dSHong Zhang   PetscFunctionBegin;
587c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
588c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
589c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
590abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
591b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
592302440fdSBarry Smith       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
593c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
594302440fdSBarry Smith         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
595b18ea86cSHong Zhang       }
596c235aa8dSHong Zhang     }
597ad8e2604SHong Zhang     /* Stage values of lambda */
598c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
599c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
600493ed95fSHong Zhang     if (ts->vec_costintegral) {
601493ed95fSHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
602493ed95fSHong Zhang     }
603abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
604b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
605493ed95fSHong Zhang       if (ts->vec_costintegral) {
606493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
607493ed95fSHong Zhang       }
608b18ea86cSHong Zhang     }
6096497ce14SHong Zhang 
610ad8e2604SHong Zhang     /* Stage values of mu */
6116497ce14SHong Zhang     if(ts->vecs_sensip) {
6125bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
613493ed95fSHong Zhang       if (ts->vec_costintegral) {
614493ed95fSHong Zhang         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
615493ed95fSHong Zhang       }
616493ed95fSHong Zhang 
617abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
618ad8e2604SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
619493ed95fSHong Zhang         if (ts->vec_costintegral) {
620493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
621493ed95fSHong Zhang         }
622ad8e2604SHong Zhang       }
623c235aa8dSHong Zhang     }
6246497ce14SHong Zhang   }
625c235aa8dSHong Zhang 
626c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
627abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
628b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6296497ce14SHong Zhang     if(ts->vecs_sensip) {
630ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
631b18ea86cSHong Zhang     }
6326497ce14SHong Zhang   }
633c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
634d2daff3dSHong Zhang   PetscFunctionReturn(0);
635d2daff3dSHong Zhang }
636d2daff3dSHong Zhang 
637d2daff3dSHong Zhang #undef __FUNCT__
638f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK"
639f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
640f68a32c8SEmil Constantinescu {
641f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
642f68a32c8SEmil Constantinescu   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
643f68a32c8SEmil Constantinescu   PetscReal        h;
644f68a32c8SEmil Constantinescu   PetscReal        tt,t;
645f68a32c8SEmil Constantinescu   PetscScalar     *b;
646f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
647f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
648e4dd521cSBarry Smith 
649f68a32c8SEmil Constantinescu   PetscFunctionBegin;
650f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
651be5899b3SLisandro Dalcin 
652f68a32c8SEmil Constantinescu   switch (rk->status) {
653f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
654f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
655f68a32c8SEmil Constantinescu     h = ts->time_step;
656f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
657f68a32c8SEmil Constantinescu     break;
658f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
659be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
660f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
661f68a32c8SEmil Constantinescu     break;
662f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
663e4dd521cSBarry Smith   }
664785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
665f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
666f68a32c8SEmil Constantinescu   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
667f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
668f68a32c8SEmil Constantinescu       b[i]  += h * B[i*pinterp+j] * tt;
669e4dd521cSBarry Smith     }
670f68a32c8SEmil Constantinescu   }
671be5899b3SLisandro Dalcin 
672f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
673f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
674be5899b3SLisandro Dalcin 
675f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
676e4dd521cSBarry Smith   PetscFunctionReturn(0);
677e4dd521cSBarry Smith }
678e4dd521cSBarry Smith 
679e4dd521cSBarry Smith /*------------------------------------------------------------*/
680be5899b3SLisandro Dalcin 
681be5899b3SLisandro Dalcin #undef __FUNCT__
682be5899b3SLisandro Dalcin #define __FUNCT__ "TSRKTableauReset"
683be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
684be5899b3SLisandro Dalcin {
685be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
686be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
687be5899b3SLisandro Dalcin   PetscErrorCode ierr;
688be5899b3SLisandro Dalcin 
689be5899b3SLisandro Dalcin   PetscFunctionBegin;
690be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
691be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
692be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
693be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
694be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
695be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
696be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
697be5899b3SLisandro Dalcin }
698be5899b3SLisandro Dalcin 
699e4dd521cSBarry Smith #undef __FUNCT__
700277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK"
701277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
702e4dd521cSBarry Smith {
7035f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
7046849ba73SBarry Smith   PetscErrorCode ierr;
705e4dd521cSBarry Smith 
706e4dd521cSBarry Smith   PetscFunctionBegin;
707be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
7080f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
709abc2977eSBarry Smith   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
710277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
711e4dd521cSBarry Smith }
712277b19d0SLisandro Dalcin 
713277b19d0SLisandro Dalcin #undef __FUNCT__
714277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK"
715277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
716277b19d0SLisandro Dalcin {
717277b19d0SLisandro Dalcin   PetscErrorCode ierr;
718277b19d0SLisandro Dalcin 
719277b19d0SLisandro Dalcin   PetscFunctionBegin;
720277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
721277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
722f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
723f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
724f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
725f68a32c8SEmil Constantinescu }
726f68a32c8SEmil Constantinescu 
727f68a32c8SEmil Constantinescu 
728f68a32c8SEmil Constantinescu #undef __FUNCT__
729f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK"
730f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
731f68a32c8SEmil Constantinescu {
732f68a32c8SEmil Constantinescu   PetscFunctionBegin;
733f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
734f68a32c8SEmil Constantinescu }
735f68a32c8SEmil Constantinescu 
736f68a32c8SEmil Constantinescu #undef __FUNCT__
737f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK"
738f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
739f68a32c8SEmil Constantinescu {
740f68a32c8SEmil Constantinescu   PetscFunctionBegin;
741f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
742f68a32c8SEmil Constantinescu }
743f68a32c8SEmil Constantinescu 
744f68a32c8SEmil Constantinescu 
745f68a32c8SEmil Constantinescu #undef __FUNCT__
746f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK"
747f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
748f68a32c8SEmil Constantinescu {
749f68a32c8SEmil Constantinescu   PetscFunctionBegin;
750f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
751f68a32c8SEmil Constantinescu }
752f68a32c8SEmil Constantinescu 
753f68a32c8SEmil Constantinescu #undef __FUNCT__
754f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
755f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
756f68a32c8SEmil Constantinescu {
757f68a32c8SEmil Constantinescu 
758f68a32c8SEmil Constantinescu   PetscFunctionBegin;
759f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
760f68a32c8SEmil Constantinescu }
761c235aa8dSHong Zhang /*
762f68a32c8SEmil Constantinescu #undef __FUNCT__
763d2daff3dSHong Zhang #define __FUNCT__ "RKSetAdjCoe"
764d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
765d2daff3dSHong Zhang {
766d2daff3dSHong Zhang   PetscReal *A,*b;
767d2daff3dSHong Zhang   PetscInt        s,i,j;
768d2daff3dSHong Zhang   PetscErrorCode  ierr;
769d2daff3dSHong Zhang 
770d2daff3dSHong Zhang   PetscFunctionBegin;
771d2daff3dSHong Zhang   s    = tab->s;
772d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
773d2daff3dSHong Zhang 
774d2daff3dSHong Zhang   for (i=0; i<s; i++)
775d2daff3dSHong Zhang     for (j=0; j<s; j++) {
776d2daff3dSHong 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];
777d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
778d2daff3dSHong Zhang     }
779d2daff3dSHong Zhang 
780d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
781d2daff3dSHong Zhang 
782d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
783d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
784d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
785d2daff3dSHong Zhang   PetscFunctionReturn(0);
786d2daff3dSHong Zhang }
787c235aa8dSHong Zhang */
788be5899b3SLisandro Dalcin 
789be5899b3SLisandro Dalcin #undef __FUNCT__
790be5899b3SLisandro Dalcin #define __FUNCT__ "TSRKTableauSetUp"
791be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
792be5899b3SLisandro Dalcin {
793be5899b3SLisandro Dalcin   TS_RK         *rk  = (TS_RK*)ts->data;
794be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
795be5899b3SLisandro Dalcin   PetscErrorCode ierr;
796be5899b3SLisandro Dalcin 
797be5899b3SLisandro Dalcin   PetscFunctionBegin;
798be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
799be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
800be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
801be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
802be5899b3SLisandro Dalcin }
803be5899b3SLisandro Dalcin 
804be5899b3SLisandro Dalcin 
805d2daff3dSHong Zhang #undef __FUNCT__
806f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK"
807f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
808f68a32c8SEmil Constantinescu {
809f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
810f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
811f68a32c8SEmil Constantinescu   DM             dm;
812f68a32c8SEmil Constantinescu 
813f68a32c8SEmil Constantinescu   PetscFunctionBegin;
814*3dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
815be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
8160f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
8170f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
8180f7a1166SHong Zhang   }
819f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
820f68a32c8SEmil Constantinescu   if (dm) {
821f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
822f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
823f68a32c8SEmil Constantinescu   }
824e4dd521cSBarry Smith   PetscFunctionReturn(0);
825e4dd521cSBarry Smith }
826d2daff3dSHong Zhang 
827f6a906c0SBarry Smith 
828e4dd521cSBarry Smith /*------------------------------------------------------------*/
829e4dd521cSBarry Smith 
830e4dd521cSBarry Smith #undef __FUNCT__
8315f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK"
8324416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
833e4dd521cSBarry Smith {
834be5899b3SLisandro Dalcin   TS_RK         *rk = (TS_RK*)ts->data;
835dfbe8321SBarry Smith   PetscErrorCode ierr;
836262ff7bbSBarry Smith 
837e4dd521cSBarry Smith   PetscFunctionBegin;
838e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
839f68a32c8SEmil Constantinescu   {
840f68a32c8SEmil Constantinescu     RKTableauLink  link;
841f68a32c8SEmil Constantinescu     PetscInt       count,choice;
842f68a32c8SEmil Constantinescu     PetscBool      flg;
843f68a32c8SEmil Constantinescu     const char   **namelist;
844f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
845785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
846f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
847be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
848be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
849f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
850f68a32c8SEmil Constantinescu   }
851262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
852e4dd521cSBarry Smith   PetscFunctionReturn(0);
853e4dd521cSBarry Smith }
854e4dd521cSBarry Smith 
855e4dd521cSBarry Smith #undef __FUNCT__
856f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray"
857f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
858f68a32c8SEmil Constantinescu {
859f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
860f68a32c8SEmil Constantinescu   PetscInt       i;
861f68a32c8SEmil Constantinescu   size_t         left,count;
862f68a32c8SEmil Constantinescu   char           *p;
863f68a32c8SEmil Constantinescu 
864f68a32c8SEmil Constantinescu   PetscFunctionBegin;
865f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
866f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
867f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
868f68a32c8SEmil Constantinescu     left -= count;
869f68a32c8SEmil Constantinescu     p    += count;
870f68a32c8SEmil Constantinescu     *p++  = ' ';
871f68a32c8SEmil Constantinescu   }
872f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
873f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
874f68a32c8SEmil Constantinescu }
875f68a32c8SEmil Constantinescu 
876f68a32c8SEmil Constantinescu #undef __FUNCT__
8775f70b478SJed Brown #define __FUNCT__ "TSView_RK"
8785f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
879e4dd521cSBarry Smith {
8805f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8818370ee3bSLisandro Dalcin   PetscBool      iascii;
882dfbe8321SBarry Smith   PetscErrorCode ierr;
8832cdf8120Sasbjorn 
884e4dd521cSBarry Smith   PetscFunctionBegin;
885251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8868370ee3bSLisandro Dalcin   if (iascii) {
8879c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
888f68a32c8SEmil Constantinescu     TSRKType  rktype;
889f68a32c8SEmil Constantinescu     char      buf[512];
890f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
891f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
892f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
893f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
894d760c35bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
8958370ee3bSLisandro Dalcin   }
8969c334d8fSLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
897f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
898f68a32c8SEmil Constantinescu }
899f68a32c8SEmil Constantinescu 
900f68a32c8SEmil Constantinescu #undef __FUNCT__
901f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK"
902f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
903f68a32c8SEmil Constantinescu {
904f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
9059c334d8fSLisandro Dalcin   TSAdapt        adapt;
906f68a32c8SEmil Constantinescu 
907f68a32c8SEmil Constantinescu   PetscFunctionBegin;
9089c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
9099c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
910f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
911f68a32c8SEmil Constantinescu }
912f68a32c8SEmil Constantinescu 
913f68a32c8SEmil Constantinescu #undef __FUNCT__
914f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType"
915f68a32c8SEmil Constantinescu /*@C
916f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
917f68a32c8SEmil Constantinescu 
918f68a32c8SEmil Constantinescu   Logically collective
919f68a32c8SEmil Constantinescu 
920f68a32c8SEmil Constantinescu   Input Parameter:
921f68a32c8SEmil Constantinescu +  ts - timestepping context
922f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
923f68a32c8SEmil Constantinescu 
924f68a32c8SEmil Constantinescu   Level: intermediate
925f68a32c8SEmil Constantinescu 
926f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
927f68a32c8SEmil Constantinescu @*/
928f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
929f68a32c8SEmil Constantinescu {
930f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
931f68a32c8SEmil Constantinescu 
932f68a32c8SEmil Constantinescu   PetscFunctionBegin;
933f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
934f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
935f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
936f68a32c8SEmil Constantinescu }
937f68a32c8SEmil Constantinescu 
938f68a32c8SEmil Constantinescu #undef __FUNCT__
939f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType"
940f68a32c8SEmil Constantinescu /*@C
941f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
942f68a32c8SEmil Constantinescu 
943f68a32c8SEmil Constantinescu   Logically collective
944f68a32c8SEmil Constantinescu 
945f68a32c8SEmil Constantinescu   Input Parameter:
946f68a32c8SEmil Constantinescu .  ts - timestepping context
947f68a32c8SEmil Constantinescu 
948f68a32c8SEmil Constantinescu   Output Parameter:
949f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
950f68a32c8SEmil Constantinescu 
951f68a32c8SEmil Constantinescu   Level: intermediate
952f68a32c8SEmil Constantinescu 
953f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
954f68a32c8SEmil Constantinescu @*/
955f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
956f68a32c8SEmil Constantinescu {
957f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
958f68a32c8SEmil Constantinescu 
959f68a32c8SEmil Constantinescu   PetscFunctionBegin;
960f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
961f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
962f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
963f68a32c8SEmil Constantinescu }
964f68a32c8SEmil Constantinescu 
965f68a32c8SEmil Constantinescu #undef __FUNCT__
966f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK"
967560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
968f68a32c8SEmil Constantinescu {
969f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
970f68a32c8SEmil Constantinescu 
971f68a32c8SEmil Constantinescu   PetscFunctionBegin;
972f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
973f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
974f68a32c8SEmil Constantinescu }
975f68a32c8SEmil Constantinescu #undef __FUNCT__
976f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK"
977560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
978f68a32c8SEmil Constantinescu {
979f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
980f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
981f68a32c8SEmil Constantinescu   PetscBool      match;
982f68a32c8SEmil Constantinescu   RKTableauLink  link;
983f68a32c8SEmil Constantinescu 
984f68a32c8SEmil Constantinescu   PetscFunctionBegin;
985f68a32c8SEmil Constantinescu   if (rk->tableau) {
986f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
987f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
988f68a32c8SEmil Constantinescu   }
989f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
990f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
991f68a32c8SEmil Constantinescu     if (match) {
992be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
993f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
994be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
995f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
996f68a32c8SEmil Constantinescu     }
997f68a32c8SEmil Constantinescu   }
998f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
999e4dd521cSBarry Smith   PetscFunctionReturn(0);
1000e4dd521cSBarry Smith }
1001e4dd521cSBarry Smith 
1002ff22ae23SHong Zhang #undef __FUNCT__
1003ff22ae23SHong Zhang #define __FUNCT__ "TSGetStages_RK"
1004ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
1005ff22ae23SHong Zhang {
1006ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
1007ff22ae23SHong Zhang 
1008ff22ae23SHong Zhang   PetscFunctionBegin;
1009ff22ae23SHong Zhang   *ns = rk->tableau->s;
1010d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
1011ff22ae23SHong Zhang   PetscFunctionReturn(0);
1012ff22ae23SHong Zhang }
1013ff22ae23SHong Zhang 
1014ff22ae23SHong Zhang 
1015e4dd521cSBarry Smith /* ------------------------------------------------------------ */
1016ebe3b25bSBarry Smith /*MC
1017f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
1018ebe3b25bSBarry Smith 
1019f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
1020f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
1021ebe3b25bSBarry Smith 
1022f68a32c8SEmil Constantinescu   Notes:
1023f68a32c8SEmil Constantinescu   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
1024ebe3b25bSBarry Smith 
1025d41222bbSBarry Smith   Level: beginner
1026d41222bbSBarry Smith 
1027f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
1028f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
1029ebe3b25bSBarry Smith 
1030ebe3b25bSBarry Smith M*/
1031e4dd521cSBarry Smith #undef __FUNCT__
10325f70b478SJed Brown #define __FUNCT__ "TSCreate_RK"
10338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1034e4dd521cSBarry Smith {
10351566a47fSLisandro Dalcin   TS_RK          *rk;
1036dfbe8321SBarry Smith   PetscErrorCode ierr;
1037e4dd521cSBarry Smith 
1038e4dd521cSBarry Smith   PetscFunctionBegin;
1039f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1040f68a32c8SEmil Constantinescu 
1041f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10425f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10435f70b478SJed Brown   ts->ops->view           = TSView_RK;
1044f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1045f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
104642f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1047f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1048f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1049f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1050fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1051f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1052ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
105342f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1054e4dd521cSBarry Smith 
10550f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10560f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10570f7a1166SHong Zhang 
10581566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10591566a47fSLisandro Dalcin   ts->data = (void*)rk;
1060e4dd521cSBarry Smith 
1061f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1062f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1063be5899b3SLisandro Dalcin 
1064be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1065e4dd521cSBarry Smith   PetscFunctionReturn(0);
1066e4dd521cSBarry Smith }
1067