xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 0f7a116688c2ba2d67457c3de78440cacd479e7e)
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 static PetscInt  explicit_stage_time_id;
17f68a32c8SEmil Constantinescu 
18f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau;
19f68a32c8SEmil Constantinescu struct _RKTableau {
20f68a32c8SEmil Constantinescu   char      *name;
21d760c35bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i              */
22f68a32c8SEmil Constantinescu   PetscInt   s;                   /* Number of stages                                           */
23d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
24f68a32c8SEmil Constantinescu   PetscInt   pinterp;             /* Interpolation order                                        */
25f68a32c8SEmil Constantinescu   PetscReal *A,*b,*c;             /* Tableau                                                    */
26f68a32c8SEmil Constantinescu   PetscReal *bembed;              /* Embedded formula of order one less (order-1)               */
27f68a32c8SEmil Constantinescu   PetscReal *binterp;             /* Dense output formula                                       */
28f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
29f68a32c8SEmil Constantinescu };
30f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
31f68a32c8SEmil Constantinescu struct _RKTableauLink {
32f68a32c8SEmil Constantinescu   struct _RKTableau tab;
33f68a32c8SEmil Constantinescu   RKTableauLink     next;
34f68a32c8SEmil Constantinescu };
35f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
36e4dd521cSBarry Smith 
37e4dd521cSBarry Smith typedef struct {
38f68a32c8SEmil Constantinescu   RKTableau   tableau;
39f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step */
40f68a32c8SEmil Constantinescu   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
41ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage*/
42ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage*/
43b18ea86cSHong Zhang   Vec          *VecSensiTemp;    /* Vector to be timed with Jacobian transpose*/
44*0f7a1166SHong Zhang   Vec          VecCostIntegral0;          /* backup for roll-backs due to events */
45f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
46f68a32c8SEmil Constantinescu   PetscReal    stage_time;
47f68a32c8SEmil Constantinescu   TSStepStatus status;
48*0f7a1166SHong Zhang   PetscReal    ptime;
49*0f7a1166SHong Zhang   PetscReal    time_step;
505f70b478SJed Brown } TS_RK;
51e4dd521cSBarry Smith 
52f68a32c8SEmil Constantinescu /*MC
53f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
54262ff7bbSBarry Smith 
55f68a32c8SEmil Constantinescu      This method has one stage.
56f68a32c8SEmil Constantinescu 
57f68a32c8SEmil Constantinescu      Level: advanced
58f68a32c8SEmil Constantinescu 
59f68a32c8SEmil Constantinescu .seealso: TSRK
60f68a32c8SEmil Constantinescu M*/
61f68a32c8SEmil Constantinescu /*MC
622109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
63f68a32c8SEmil Constantinescu 
64f68a32c8SEmil Constantinescu      This method has two stages.
65f68a32c8SEmil Constantinescu 
66f68a32c8SEmil Constantinescu      Level: advanced
67f68a32c8SEmil Constantinescu 
68f68a32c8SEmil Constantinescu .seealso: TSRK
69f68a32c8SEmil Constantinescu M*/
70f68a32c8SEmil Constantinescu /*MC
71f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
72f68a32c8SEmil Constantinescu 
73f68a32c8SEmil Constantinescu      This method has three stages.
74f68a32c8SEmil Constantinescu 
75f68a32c8SEmil Constantinescu      Level: advanced
76f68a32c8SEmil Constantinescu 
77f68a32c8SEmil Constantinescu .seealso: TSRK
78f68a32c8SEmil Constantinescu M*/
79f68a32c8SEmil Constantinescu /*MC
802109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
812109b73fSDebojyoti Ghosh 
822109b73fSDebojyoti Ghosh      This method has four stages.
832109b73fSDebojyoti Ghosh 
842109b73fSDebojyoti Ghosh      Level: advanced
852109b73fSDebojyoti Ghosh 
862109b73fSDebojyoti Ghosh .seealso: TSRK
872109b73fSDebojyoti Ghosh M*/
882109b73fSDebojyoti Ghosh /*MC
89f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
90f68a32c8SEmil Constantinescu 
912e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
92f68a32c8SEmil Constantinescu 
93f68a32c8SEmil Constantinescu      Level: advanced
94f68a32c8SEmil Constantinescu 
95f68a32c8SEmil Constantinescu .seealso: TSRK
96f68a32c8SEmil Constantinescu M*/
97f68a32c8SEmil Constantinescu /*MC
982e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
99f68a32c8SEmil Constantinescu 
100f68a32c8SEmil Constantinescu      This method has six stages.
101f68a32c8SEmil Constantinescu 
102f68a32c8SEmil Constantinescu      Level: advanced
103f68a32c8SEmil Constantinescu 
104f68a32c8SEmil Constantinescu .seealso: TSRK
105f68a32c8SEmil Constantinescu M*/
1062109b73fSDebojyoti Ghosh /*MC
1072e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1082109b73fSDebojyoti Ghosh 
1092109b73fSDebojyoti Ghosh      This method has seven stages.
1102109b73fSDebojyoti Ghosh 
1112109b73fSDebojyoti Ghosh      Level: advanced
1122109b73fSDebojyoti Ghosh 
1132109b73fSDebojyoti Ghosh .seealso: TSRK
1142109b73fSDebojyoti Ghosh M*/
115262ff7bbSBarry Smith 
116262ff7bbSBarry Smith #undef __FUNCT__
117f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll"
118f68a32c8SEmil Constantinescu /*@C
119f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
120262ff7bbSBarry Smith 
121f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
122262ff7bbSBarry Smith 
123f68a32c8SEmil Constantinescu   Level: advanced
124262ff7bbSBarry Smith 
125f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
126262ff7bbSBarry Smith 
127f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
128262ff7bbSBarry Smith @*/
129f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
130262ff7bbSBarry Smith {
1314ac538c5SBarry Smith   PetscErrorCode ierr;
132262ff7bbSBarry Smith 
133262ff7bbSBarry Smith   PetscFunctionBegin;
134f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
135f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
136e4dd521cSBarry Smith 
137e4dd521cSBarry Smith   {
138f68a32c8SEmil Constantinescu     const PetscReal
139f68a32c8SEmil Constantinescu       A[1][1] = {{0.0}},
140f68a32c8SEmil Constantinescu       b[1]    = {1.0};
141f68a32c8SEmil Constantinescu     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
142e4dd521cSBarry Smith   }
143db046809SBarry Smith   {
144f68a32c8SEmil Constantinescu     const PetscReal
145f68a32c8SEmil Constantinescu       A[2][2]     = {{0.0,0.0},
146f68a32c8SEmil Constantinescu                     {1.0,0.0}},
147f68a32c8SEmil Constantinescu       b[2]        = {0.5,0.5},
148f68a32c8SEmil Constantinescu       bembed[2]   = {1.0,0};
149fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr);
150db046809SBarry Smith   }
151f68a32c8SEmil Constantinescu   {
152f68a32c8SEmil Constantinescu     const PetscReal
153f68a32c8SEmil Constantinescu       A[3][3] = {{0,0,0},
154f68a32c8SEmil Constantinescu                  {2.0/3.0,0,0},
155f68a32c8SEmil Constantinescu                  {-1.0/3.0,1.0,0}},
156f68a32c8SEmil Constantinescu       b[3]    = {0.25,0.5,0.25};
157fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr);
158fdefd5e5SDebojyoti Ghosh   }
159fdefd5e5SDebojyoti Ghosh   {
160fdefd5e5SDebojyoti Ghosh     const PetscReal
161fdefd5e5SDebojyoti Ghosh       A[4][4] = {{0,0,0,0},
162fdefd5e5SDebojyoti Ghosh                  {1.0/2.0,0,0,0},
163fdefd5e5SDebojyoti Ghosh                  {0,3.0/4.0,0,0},
164fdefd5e5SDebojyoti Ghosh                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
165fdefd5e5SDebojyoti Ghosh       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
166fdefd5e5SDebojyoti Ghosh       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
167fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr);
168db046809SBarry Smith   }
169f68a32c8SEmil Constantinescu   {
170f68a32c8SEmil Constantinescu     const PetscReal
171f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
172f68a32c8SEmil Constantinescu                  {0.5,0,0,0},
173f68a32c8SEmil Constantinescu                  {0,0.5,0,0},
174f68a32c8SEmil Constantinescu                  {0,0,1.0,0}},
175f68a32c8SEmil Constantinescu       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
176fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr);
177db046809SBarry Smith   }
178f68a32c8SEmil Constantinescu   {
179f68a32c8SEmil Constantinescu     const PetscReal
180f68a32c8SEmil Constantinescu       A[6][6]   = {{0,0,0,0,0,0},
181f68a32c8SEmil Constantinescu                    {0.25,0,0,0,0,0},
182f68a32c8SEmil Constantinescu                    {3.0/32.0,9.0/32.0,0,0,0,0},
183f68a32c8SEmil Constantinescu                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
184f68a32c8SEmil Constantinescu                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
185f68a32c8SEmil Constantinescu                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
186f68a32c8SEmil Constantinescu       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
187f68a32c8SEmil Constantinescu       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
188fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
189fdefd5e5SDebojyoti Ghosh   }
190fdefd5e5SDebojyoti Ghosh   {
191fdefd5e5SDebojyoti Ghosh     const PetscReal
192fdefd5e5SDebojyoti Ghosh       A[7][7]   = {{0,0,0,0,0,0,0},
193fdefd5e5SDebojyoti Ghosh                    {1.0/5.0,0,0,0,0,0,0},
194fdefd5e5SDebojyoti Ghosh                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
195fdefd5e5SDebojyoti Ghosh                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
196fdefd5e5SDebojyoti Ghosh                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
197fdefd5e5SDebojyoti Ghosh                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
198fdefd5e5SDebojyoti Ghosh                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
199fdefd5e5SDebojyoti 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},
200fdefd5e5SDebojyoti 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};
201fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
202f68a32c8SEmil Constantinescu   }
203db046809SBarry Smith   PetscFunctionReturn(0);
204db046809SBarry Smith }
205db046809SBarry Smith 
206e4dd521cSBarry Smith #undef __FUNCT__
207f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy"
208f68a32c8SEmil Constantinescu /*@C
209f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
210f68a32c8SEmil Constantinescu 
211f68a32c8SEmil Constantinescu    Not Collective
212f68a32c8SEmil Constantinescu 
213f68a32c8SEmil Constantinescu    Level: advanced
214f68a32c8SEmil Constantinescu 
215f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
216f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
217f68a32c8SEmil Constantinescu @*/
218f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
219e4dd521cSBarry Smith {
220f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
221f68a32c8SEmil Constantinescu   RKTableauLink link;
222f68a32c8SEmil Constantinescu 
223f68a32c8SEmil Constantinescu   PetscFunctionBegin;
224f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
225f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
226f68a32c8SEmil Constantinescu     RKTableauList = link->next;
227f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
228f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
229f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
230f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
231f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
232f68a32c8SEmil Constantinescu   }
233f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
234f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
235f68a32c8SEmil Constantinescu }
236f68a32c8SEmil Constantinescu 
237f68a32c8SEmil Constantinescu #undef __FUNCT__
238f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage"
239f68a32c8SEmil Constantinescu /*@C
240f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
241f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
242f68a32c8SEmil Constantinescu   when using static libraries.
243f68a32c8SEmil Constantinescu 
244f68a32c8SEmil Constantinescu   Level: developer
245f68a32c8SEmil Constantinescu 
246f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
247f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
248f68a32c8SEmil Constantinescu @*/
249f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
250f68a32c8SEmil Constantinescu {
251186e87acSLisandro Dalcin   PetscErrorCode ierr;
252e4dd521cSBarry Smith 
253e4dd521cSBarry Smith   PetscFunctionBegin;
254f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
255f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
256f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
257f68a32c8SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
258f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
259f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
260f68a32c8SEmil Constantinescu }
261186e87acSLisandro Dalcin 
262f68a32c8SEmil Constantinescu #undef __FUNCT__
263f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage"
264f68a32c8SEmil Constantinescu /*@C
265f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
266f68a32c8SEmil Constantinescu   called from PetscFinalize().
267186e87acSLisandro Dalcin 
268f68a32c8SEmil Constantinescu   Level: developer
269f68a32c8SEmil Constantinescu 
270f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
271f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
272f68a32c8SEmil Constantinescu @*/
273f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
274f68a32c8SEmil Constantinescu {
275f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
276f68a32c8SEmil Constantinescu 
277f68a32c8SEmil Constantinescu   PetscFunctionBegin;
278f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
279f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
280f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
281f68a32c8SEmil Constantinescu }
282f68a32c8SEmil Constantinescu 
283f68a32c8SEmil Constantinescu #undef __FUNCT__
284f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister"
285f68a32c8SEmil Constantinescu /*@C
286f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
287f68a32c8SEmil Constantinescu 
288f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
289f68a32c8SEmil Constantinescu 
290f68a32c8SEmil Constantinescu    Input Parameters:
291f68a32c8SEmil Constantinescu +  name - identifier for method
292f68a32c8SEmil Constantinescu .  order - approximation order of method
293f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
294f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
295f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
296f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
297f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
298f68a32c8SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterp
299f68a32c8SEmil Constantinescu -  binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt)
300f68a32c8SEmil Constantinescu 
301f68a32c8SEmil Constantinescu    Notes:
302f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
303f68a32c8SEmil Constantinescu 
304f68a32c8SEmil Constantinescu    Level: advanced
305f68a32c8SEmil Constantinescu 
306f68a32c8SEmil Constantinescu .keywords: TS, register
307f68a32c8SEmil Constantinescu 
308f68a32c8SEmil Constantinescu .seealso: TSRK
309f68a32c8SEmil Constantinescu @*/
310f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
311f68a32c8SEmil Constantinescu                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
312f68a32c8SEmil Constantinescu                                  const PetscReal bembed[],
313f68a32c8SEmil Constantinescu                                  PetscInt pinterp,const PetscReal binterp[])
314f68a32c8SEmil Constantinescu {
315f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
316f68a32c8SEmil Constantinescu   RKTableauLink   link;
317f68a32c8SEmil Constantinescu   RKTableau       t;
318f68a32c8SEmil Constantinescu   PetscInt        i,j;
319f68a32c8SEmil Constantinescu 
320f68a32c8SEmil Constantinescu   PetscFunctionBegin;
321f68a32c8SEmil Constantinescu   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
322f68a32c8SEmil Constantinescu   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
323f68a32c8SEmil Constantinescu   t        = &link->tab;
324f68a32c8SEmil Constantinescu   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
325f68a32c8SEmil Constantinescu   t->order = order;
326f68a32c8SEmil Constantinescu   t->s     = s;
327dcca6d9dSJed Brown   ierr     = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
328f68a32c8SEmil Constantinescu   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
329f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
330f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
331f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
332f68a32c8SEmil 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];
333d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
334d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
335f68a32c8SEmil Constantinescu   if (bembed) {
336785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
337f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
338f68a32c8SEmil Constantinescu   }
339f68a32c8SEmil Constantinescu 
340f68a32c8SEmil Constantinescu   t->pinterp     = pinterp;
341785e854fSJed Brown   ierr           = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
342f68a32c8SEmil Constantinescu   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
343f68a32c8SEmil Constantinescu   link->next     = RKTableauList;
344f68a32c8SEmil Constantinescu   RKTableauList = link;
345f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
346f68a32c8SEmil Constantinescu }
347f68a32c8SEmil Constantinescu 
348f68a32c8SEmil Constantinescu #undef __FUNCT__
349f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK"
350e4dd521cSBarry Smith /*
351f68a32c8SEmil Constantinescu  The step completion formula is
352e4dd521cSBarry Smith 
353f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
354f68a32c8SEmil Constantinescu 
355f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
356f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
357f68a32c8SEmil Constantinescu  We can write
358f68a32c8SEmil Constantinescu 
359f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
360f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
361f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
362f68a32c8SEmil Constantinescu 
363f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
364e4dd521cSBarry Smith */
365f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
366f68a32c8SEmil Constantinescu {
367f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
368f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
369f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
370f68a32c8SEmil Constantinescu   PetscReal      h;
371f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
372f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
373f68a32c8SEmil Constantinescu 
374f68a32c8SEmil Constantinescu   PetscFunctionBegin;
375f68a32c8SEmil Constantinescu   switch (rk->status) {
376f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
377f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
378f68a32c8SEmil Constantinescu     h = ts->time_step; break;
379f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
380f68a32c8SEmil Constantinescu     h = ts->time_step_prev; break;
381f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
382f68a32c8SEmil Constantinescu   }
383f68a32c8SEmil Constantinescu   if (order == tab->order) {
384f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
385f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
386f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
387f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
388f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
389f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
390f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
391f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
392f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
393f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
394f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
395f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
396f68a32c8SEmil Constantinescu     } else {                    /* Rollback and re-complete using (be-b) */
397f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
398f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
399f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
400*0f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
401*0f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
402*0f7a1166SHong Zhang       }
403f68a32c8SEmil Constantinescu     }
404f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
405f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
406f68a32c8SEmil Constantinescu   }
407f68a32c8SEmil Constantinescu unavailable:
408f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
409f68a32c8SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order);
410f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
411f68a32c8SEmil Constantinescu }
412f68a32c8SEmil Constantinescu 
413f68a32c8SEmil Constantinescu #undef __FUNCT__
414*0f7a1166SHong Zhang #define __FUNCT__ "TSForwardCostIntegral_RK"
415*0f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
416*0f7a1166SHong Zhang {
417*0f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
418*0f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
419*0f7a1166SHong Zhang   const PetscInt  s = tab->s;
420*0f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
421*0f7a1166SHong Zhang   Vec             *Y = rk->Y;
422*0f7a1166SHong Zhang   PetscInt        i;
423*0f7a1166SHong Zhang   PetscErrorCode  ierr;
424*0f7a1166SHong Zhang 
425*0f7a1166SHong Zhang   PetscFunctionBegin;
426*0f7a1166SHong Zhang   /* backup cost integral */
427*0f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
428*0f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
429*0f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
430*0f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
431*0f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
432*0f7a1166SHong Zhang   }
433*0f7a1166SHong Zhang   PetscFunctionReturn(0);
434*0f7a1166SHong Zhang }
435*0f7a1166SHong Zhang 
436*0f7a1166SHong Zhang #undef __FUNCT__
437*0f7a1166SHong Zhang #define __FUNCT__ "TSAdjointCostIntegral_RK"
438*0f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
439*0f7a1166SHong Zhang {
440*0f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
441*0f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
442*0f7a1166SHong Zhang   const PetscInt  s = tab->s;
443*0f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
444*0f7a1166SHong Zhang   Vec             *Y = rk->Y;
445*0f7a1166SHong Zhang   PetscInt        i;
446*0f7a1166SHong Zhang   PetscErrorCode  ierr;
447*0f7a1166SHong Zhang 
448*0f7a1166SHong Zhang   PetscFunctionBegin;
449*0f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
450*0f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
451*0f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
452*0f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
453*0f7a1166SHong Zhang   }
454*0f7a1166SHong Zhang   PetscFunctionReturn(0);
455*0f7a1166SHong Zhang }
456*0f7a1166SHong Zhang 
457*0f7a1166SHong Zhang #undef __FUNCT__
458f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK"
459f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
460f68a32c8SEmil Constantinescu {
461f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
462f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
463f68a32c8SEmil Constantinescu   const PetscInt   s    = tab->s;
464f68a32c8SEmil Constantinescu   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
465f68a32c8SEmil Constantinescu   PetscScalar     *w    = rk->work;
466f68a32c8SEmil Constantinescu   Vec             *Y    = rk->Y,*YdotRHS = rk->YdotRHS;
467f68a32c8SEmil Constantinescu   TSAdapt          adapt;
468f68a32c8SEmil Constantinescu   PetscInt         i,j,reject,next_scheme;
469f68a32c8SEmil Constantinescu   PetscReal        next_time_step;
470f68a32c8SEmil Constantinescu   PetscReal        t;
471f68a32c8SEmil Constantinescu   PetscBool        accept;
472f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
473f68a32c8SEmil Constantinescu 
474f68a32c8SEmil Constantinescu   PetscFunctionBegin;
475f68a32c8SEmil Constantinescu 
476f68a32c8SEmil Constantinescu   next_time_step = ts->time_step;
477f68a32c8SEmil Constantinescu   t              = ts->ptime;
478f68a32c8SEmil Constantinescu   accept         = PETSC_TRUE;
479f68a32c8SEmil Constantinescu   rk->status     = TS_STEP_INCOMPLETE;
480f68a32c8SEmil Constantinescu 
481f68a32c8SEmil Constantinescu 
482f68a32c8SEmil Constantinescu   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
483f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
484f68a32c8SEmil Constantinescu     ierr = TSPreStep(ts);CHKERRQ(ierr);
485f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
4869be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
4879be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
488f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
489f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
490f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
4919be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
492f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
493b295832fSPierre Barbier de Reuille       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&accept);CHKERRQ(ierr);
494cb9d8021SPierre Barbier de Reuille       if (!accept) break;
495f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
496f68a32c8SEmil Constantinescu     }
497cb9d8021SPierre Barbier de Reuille     if(!accept) continue;
498f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
499f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
500f68a32c8SEmil Constantinescu 
501f68a32c8SEmil Constantinescu     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
502f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
503f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
504f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
505f68a32c8SEmil Constantinescu     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
506f68a32c8SEmil Constantinescu     if (accept) {
507*0f7a1166SHong Zhang       if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
508*0f7a1166SHong Zhang         rk->ptime     = ts->ptime;
509*0f7a1166SHong Zhang         rk->time_step = ts->time_step;
510493ed95fSHong Zhang       }
511f68a32c8SEmil Constantinescu       /* ignore next_scheme for now */
512f68a32c8SEmil Constantinescu       ts->ptime    += ts->time_step;
513f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
514e3caeda6SBarry Smith       ts->steps++;
515f68a32c8SEmil Constantinescu       rk->status = TS_STEP_COMPLETE;
516f68a32c8SEmil Constantinescu       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
517f68a32c8SEmil Constantinescu       break;
518f68a32c8SEmil Constantinescu     } else {                    /* Roll back the current step */
519f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
520f68a32c8SEmil Constantinescu       ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr);
521f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
522f68a32c8SEmil Constantinescu       rk->status   = TS_STEP_INCOMPLETE;
523f68a32c8SEmil Constantinescu     }
524f68a32c8SEmil Constantinescu   }
525f68a32c8SEmil Constantinescu   if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
526f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
527e4dd521cSBarry Smith }
528e4dd521cSBarry Smith 
529f68a32c8SEmil Constantinescu #undef __FUNCT__
530f6a906c0SBarry Smith #define __FUNCT__ "TSAdjointSetUp_RK"
531f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
532f6a906c0SBarry Smith {
533f6a906c0SBarry Smith   TS_RK         *rk = (TS_RK*)ts->data;
534f6a906c0SBarry Smith   RKTableau      tab;
535f6a906c0SBarry Smith   PetscInt       s;
536f6a906c0SBarry Smith   PetscErrorCode ierr;
537f6a906c0SBarry Smith 
538f6a906c0SBarry Smith   PetscFunctionBegin;
539f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
540f6a906c0SBarry Smith   tab  = rk->tableau;
541f6a906c0SBarry Smith   s    = tab->s;
542abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
543f6a906c0SBarry Smith   if(ts->vecs_sensip) {
544abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
545f6a906c0SBarry Smith   }
546abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
547f6a906c0SBarry Smith   PetscFunctionReturn(0);
548f6a906c0SBarry Smith }
549f6a906c0SBarry Smith 
550f6a906c0SBarry Smith #undef __FUNCT__
55142f2b339SBarry Smith #define __FUNCT__ "TSAdjointStep_RK"
55242f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
553d2daff3dSHong Zhang {
554c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
555c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
556c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
557c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
558c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
559ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
560b18ea86cSHong Zhang   PetscInt         i,j,nadj;
561c235aa8dSHong Zhang   PetscReal        t;
562d2daff3dSHong Zhang   PetscErrorCode   ierr;
563cef76868SBarry Smith   PetscReal        h = ts->time_step;
564cef76868SBarry Smith   Mat              J,Jp;
565c235aa8dSHong Zhang 
566d2daff3dSHong Zhang   PetscFunctionBegin;
567c235aa8dSHong Zhang   t          = ts->ptime;
568c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
569cef76868SBarry Smith   h = ts->time_step;
570c235aa8dSHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
571c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
572c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
573abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
574b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
575302440fdSBarry Smith       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
576c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
577302440fdSBarry Smith         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
578b18ea86cSHong Zhang       }
579c235aa8dSHong Zhang     }
580ad8e2604SHong Zhang     /* Stage values of lambda */
581c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
582c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
583493ed95fSHong Zhang     if (ts->vec_costintegral) {
584493ed95fSHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
585493ed95fSHong Zhang     }
586abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
587b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
588493ed95fSHong Zhang       if (ts->vec_costintegral) {
589493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
590493ed95fSHong Zhang       }
591b18ea86cSHong Zhang     }
5926497ce14SHong Zhang 
593ad8e2604SHong Zhang     /* Stage values of mu */
5946497ce14SHong Zhang     if(ts->vecs_sensip) {
5955bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
596493ed95fSHong Zhang       if (ts->vec_costintegral) {
597493ed95fSHong Zhang         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
598493ed95fSHong Zhang       }
599493ed95fSHong Zhang 
600abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
601ad8e2604SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
602493ed95fSHong Zhang         if (ts->vec_costintegral) {
603493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
604493ed95fSHong Zhang         }
605ad8e2604SHong Zhang       }
606c235aa8dSHong Zhang     }
6076497ce14SHong Zhang   }
608c235aa8dSHong Zhang 
609c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
610abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
611b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6126497ce14SHong Zhang     if(ts->vecs_sensip) {
613ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
614b18ea86cSHong Zhang     }
6156497ce14SHong Zhang   }
616c235aa8dSHong Zhang   ts->steps++;
617c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
618d2daff3dSHong Zhang   PetscFunctionReturn(0);
619d2daff3dSHong Zhang }
620d2daff3dSHong Zhang 
621d2daff3dSHong Zhang #undef __FUNCT__
622f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK"
623f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
624f68a32c8SEmil Constantinescu {
625f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
626f68a32c8SEmil Constantinescu   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
627f68a32c8SEmil Constantinescu   PetscReal        h;
628f68a32c8SEmil Constantinescu   PetscReal        tt,t;
629f68a32c8SEmil Constantinescu   PetscScalar     *b;
630f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
631f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
632e4dd521cSBarry Smith 
633f68a32c8SEmil Constantinescu   PetscFunctionBegin;
634f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
635f68a32c8SEmil Constantinescu   switch (rk->status) {
636f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
637f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
638f68a32c8SEmil Constantinescu     h = ts->time_step;
639f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
640f68a32c8SEmil Constantinescu     break;
641f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
642f68a32c8SEmil Constantinescu     h = ts->time_step_prev;
643f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
644f68a32c8SEmil Constantinescu     break;
645f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
646e4dd521cSBarry Smith   }
647785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
648f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
649f68a32c8SEmil Constantinescu   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
650f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
651f68a32c8SEmil Constantinescu       b[i]  += h * B[i*pinterp+j] * tt;
652e4dd521cSBarry Smith     }
653f68a32c8SEmil Constantinescu   }
654f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
655f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
656f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
657e4dd521cSBarry Smith   PetscFunctionReturn(0);
658e4dd521cSBarry Smith }
659e4dd521cSBarry Smith 
660e4dd521cSBarry Smith /*------------------------------------------------------------*/
661e4dd521cSBarry Smith #undef __FUNCT__
662277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK"
663277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
664e4dd521cSBarry Smith {
6655f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
666f68a32c8SEmil Constantinescu   PetscInt       s;
6676849ba73SBarry Smith   PetscErrorCode ierr;
668e4dd521cSBarry Smith 
669e4dd521cSBarry Smith   PetscFunctionBegin;
670f68a32c8SEmil Constantinescu   if (!rk->tableau) PetscFunctionReturn(0);
671f68a32c8SEmil Constantinescu   s    = rk->tableau->s;
672f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr);
673f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr);
674*0f7a1166SHong Zhang   if (rk->VecCostIntegral0) {
675*0f7a1166SHong Zhang     ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
676*0f7a1166SHong Zhang   }
677abc2977eSBarry Smith   ierr = VecDestroyVecs(s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
678abc2977eSBarry Smith   ierr = VecDestroyVecs(s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
679abc2977eSBarry Smith   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
680f68a32c8SEmil Constantinescu   ierr = PetscFree(rk->work);CHKERRQ(ierr);
681277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
682e4dd521cSBarry Smith }
683277b19d0SLisandro Dalcin 
684277b19d0SLisandro Dalcin #undef __FUNCT__
685277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK"
686277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
687277b19d0SLisandro Dalcin {
688277b19d0SLisandro Dalcin   PetscErrorCode ierr;
689277b19d0SLisandro Dalcin 
690277b19d0SLisandro Dalcin   PetscFunctionBegin;
691277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
692277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
693f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
694f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
695f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
696f68a32c8SEmil Constantinescu }
697f68a32c8SEmil Constantinescu 
698f68a32c8SEmil Constantinescu 
699f68a32c8SEmil Constantinescu #undef __FUNCT__
700f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK"
701f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
702f68a32c8SEmil Constantinescu {
703f68a32c8SEmil Constantinescu   PetscFunctionBegin;
704f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
705f68a32c8SEmil Constantinescu }
706f68a32c8SEmil Constantinescu 
707f68a32c8SEmil Constantinescu #undef __FUNCT__
708f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK"
709f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
710f68a32c8SEmil Constantinescu {
711f68a32c8SEmil Constantinescu   PetscFunctionBegin;
712f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
713f68a32c8SEmil Constantinescu }
714f68a32c8SEmil Constantinescu 
715f68a32c8SEmil Constantinescu 
716f68a32c8SEmil Constantinescu #undef __FUNCT__
717f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK"
718f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
719f68a32c8SEmil Constantinescu {
720f68a32c8SEmil Constantinescu   PetscFunctionBegin;
721f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
722f68a32c8SEmil Constantinescu }
723f68a32c8SEmil Constantinescu 
724f68a32c8SEmil Constantinescu #undef __FUNCT__
725f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
726f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
727f68a32c8SEmil Constantinescu {
728f68a32c8SEmil Constantinescu 
729f68a32c8SEmil Constantinescu   PetscFunctionBegin;
730f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
731f68a32c8SEmil Constantinescu }
732c235aa8dSHong Zhang /*
733f68a32c8SEmil Constantinescu #undef __FUNCT__
734d2daff3dSHong Zhang #define __FUNCT__ "RKSetAdjCoe"
735d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
736d2daff3dSHong Zhang {
737d2daff3dSHong Zhang   PetscReal *A,*b;
738d2daff3dSHong Zhang   PetscInt        s,i,j;
739d2daff3dSHong Zhang   PetscErrorCode  ierr;
740d2daff3dSHong Zhang 
741d2daff3dSHong Zhang   PetscFunctionBegin;
742d2daff3dSHong Zhang   s    = tab->s;
743d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
744d2daff3dSHong Zhang 
745d2daff3dSHong Zhang   for (i=0; i<s; i++)
746d2daff3dSHong Zhang     for (j=0; j<s; j++) {
747d2daff3dSHong 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];
748d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
749d2daff3dSHong Zhang     }
750d2daff3dSHong Zhang 
751d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
752d2daff3dSHong Zhang 
753d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
754d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
755d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
756d2daff3dSHong Zhang   PetscFunctionReturn(0);
757d2daff3dSHong Zhang }
758c235aa8dSHong Zhang */
759d2daff3dSHong Zhang #undef __FUNCT__
760f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK"
761f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
762f68a32c8SEmil Constantinescu {
763f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
764f68a32c8SEmil Constantinescu   RKTableau      tab;
765f68a32c8SEmil Constantinescu   PetscInt       s;
766f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
767f68a32c8SEmil Constantinescu   DM             dm;
768f68a32c8SEmil Constantinescu 
769f68a32c8SEmil Constantinescu   PetscFunctionBegin;
770*0f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
771*0f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
772*0f7a1166SHong Zhang   }
773f68a32c8SEmil Constantinescu   if (!rk->tableau) {
774f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
775f68a32c8SEmil Constantinescu   }
776f68a32c8SEmil Constantinescu   tab  = rk->tableau;
777f68a32c8SEmil Constantinescu   s    = tab->s;
778f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr);
779f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr);
780785e854fSJed Brown   ierr = PetscMalloc1(s,&rk->work);CHKERRQ(ierr);
781f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
782f68a32c8SEmil Constantinescu   if (dm) {
783f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
784f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
785f68a32c8SEmil Constantinescu   }
786e4dd521cSBarry Smith   PetscFunctionReturn(0);
787e4dd521cSBarry Smith }
788d2daff3dSHong Zhang 
789f6a906c0SBarry Smith 
790e4dd521cSBarry Smith /*------------------------------------------------------------*/
791e4dd521cSBarry Smith 
792e4dd521cSBarry Smith #undef __FUNCT__
7935f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK"
7944416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
795e4dd521cSBarry Smith {
796dfbe8321SBarry Smith   PetscErrorCode ierr;
797f68a32c8SEmil Constantinescu   char           rktype[256];
798262ff7bbSBarry Smith 
799e4dd521cSBarry Smith   PetscFunctionBegin;
800e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
801f68a32c8SEmil Constantinescu   {
802f68a32c8SEmil Constantinescu     RKTableauLink  link;
803f68a32c8SEmil Constantinescu     PetscInt       count,choice;
804f68a32c8SEmil Constantinescu     PetscBool      flg;
805f68a32c8SEmil Constantinescu     const char   **namelist;
806f68a32c8SEmil Constantinescu     ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr);
807f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
808785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
809f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
810f68a32c8SEmil Constantinescu     ierr      = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr);
811f68a32c8SEmil Constantinescu     ierr      = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr);
812f68a32c8SEmil Constantinescu     ierr      = PetscFree(namelist);CHKERRQ(ierr);
813f68a32c8SEmil Constantinescu   }
814262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
815e4dd521cSBarry Smith   PetscFunctionReturn(0);
816e4dd521cSBarry Smith }
817e4dd521cSBarry Smith 
818e4dd521cSBarry Smith #undef __FUNCT__
819f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray"
820f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
821f68a32c8SEmil Constantinescu {
822f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
823f68a32c8SEmil Constantinescu   PetscInt       i;
824f68a32c8SEmil Constantinescu   size_t         left,count;
825f68a32c8SEmil Constantinescu   char           *p;
826f68a32c8SEmil Constantinescu 
827f68a32c8SEmil Constantinescu   PetscFunctionBegin;
828f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
829f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
830f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
831f68a32c8SEmil Constantinescu     left -= count;
832f68a32c8SEmil Constantinescu     p    += count;
833f68a32c8SEmil Constantinescu     *p++  = ' ';
834f68a32c8SEmil Constantinescu   }
835f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
836f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
837f68a32c8SEmil Constantinescu }
838f68a32c8SEmil Constantinescu 
839f68a32c8SEmil Constantinescu #undef __FUNCT__
8405f70b478SJed Brown #define __FUNCT__ "TSView_RK"
8415f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
842e4dd521cSBarry Smith {
8435f70b478SJed Brown   TS_RK         *rk   = (TS_RK*)ts->data;
844f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
8458370ee3bSLisandro Dalcin   PetscBool      iascii;
846dfbe8321SBarry Smith   PetscErrorCode ierr;
847f68a32c8SEmil Constantinescu   TSAdapt        adapt;
8482cdf8120Sasbjorn 
849e4dd521cSBarry Smith   PetscFunctionBegin;
850251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8518370ee3bSLisandro Dalcin   if (iascii) {
852f68a32c8SEmil Constantinescu     TSRKType rktype;
853f68a32c8SEmil Constantinescu     char     buf[512];
854f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
855f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
856f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
857f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
858d760c35bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
8598370ee3bSLisandro Dalcin   }
860f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
861f68a32c8SEmil Constantinescu   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
862f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
863f68a32c8SEmil Constantinescu }
864f68a32c8SEmil Constantinescu 
865f68a32c8SEmil Constantinescu #undef __FUNCT__
866f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK"
867f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
868f68a32c8SEmil Constantinescu {
869f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
870f68a32c8SEmil Constantinescu   TSAdapt        tsadapt;
871f68a32c8SEmil Constantinescu 
872f68a32c8SEmil Constantinescu   PetscFunctionBegin;
873f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
874f68a32c8SEmil Constantinescu   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
875f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
876f68a32c8SEmil Constantinescu }
877f68a32c8SEmil Constantinescu 
878f68a32c8SEmil Constantinescu #undef __FUNCT__
879f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType"
880f68a32c8SEmil Constantinescu /*@C
881f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
882f68a32c8SEmil Constantinescu 
883f68a32c8SEmil Constantinescu   Logically collective
884f68a32c8SEmil Constantinescu 
885f68a32c8SEmil Constantinescu   Input Parameter:
886f68a32c8SEmil Constantinescu +  ts - timestepping context
887f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
888f68a32c8SEmil Constantinescu 
889f68a32c8SEmil Constantinescu   Level: intermediate
890f68a32c8SEmil Constantinescu 
891f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
892f68a32c8SEmil Constantinescu @*/
893f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
894f68a32c8SEmil Constantinescu {
895f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
896f68a32c8SEmil Constantinescu 
897f68a32c8SEmil Constantinescu   PetscFunctionBegin;
898f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
899f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
900f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
901f68a32c8SEmil Constantinescu }
902f68a32c8SEmil Constantinescu 
903f68a32c8SEmil Constantinescu #undef __FUNCT__
904f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType"
905f68a32c8SEmil Constantinescu /*@C
906f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
907f68a32c8SEmil Constantinescu 
908f68a32c8SEmil Constantinescu   Logically collective
909f68a32c8SEmil Constantinescu 
910f68a32c8SEmil Constantinescu   Input Parameter:
911f68a32c8SEmil Constantinescu .  ts - timestepping context
912f68a32c8SEmil Constantinescu 
913f68a32c8SEmil Constantinescu   Output Parameter:
914f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
915f68a32c8SEmil Constantinescu 
916f68a32c8SEmil Constantinescu   Level: intermediate
917f68a32c8SEmil Constantinescu 
918f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
919f68a32c8SEmil Constantinescu @*/
920f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
921f68a32c8SEmil Constantinescu {
922f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
923f68a32c8SEmil Constantinescu 
924f68a32c8SEmil Constantinescu   PetscFunctionBegin;
925f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
926f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
927f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
928f68a32c8SEmil Constantinescu }
929f68a32c8SEmil Constantinescu 
930f68a32c8SEmil Constantinescu #undef __FUNCT__
931f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK"
932f68a32c8SEmil Constantinescu PetscErrorCode  TSRKGetType_RK(TS ts,TSRKType *rktype)
933f68a32c8SEmil Constantinescu {
934f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
935f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
936f68a32c8SEmil Constantinescu 
937f68a32c8SEmil Constantinescu   PetscFunctionBegin;
938f68a32c8SEmil Constantinescu   if (!rk->tableau) {
939f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
940f68a32c8SEmil Constantinescu   }
941f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
942f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
943f68a32c8SEmil Constantinescu }
944f68a32c8SEmil Constantinescu #undef __FUNCT__
945f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK"
946f68a32c8SEmil Constantinescu PetscErrorCode  TSRKSetType_RK(TS ts,TSRKType rktype)
947f68a32c8SEmil Constantinescu {
948f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
949f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
950f68a32c8SEmil Constantinescu   PetscBool      match;
951f68a32c8SEmil Constantinescu   RKTableauLink link;
952f68a32c8SEmil Constantinescu 
953f68a32c8SEmil Constantinescu   PetscFunctionBegin;
954f68a32c8SEmil Constantinescu   if (rk->tableau) {
955f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
956f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
957f68a32c8SEmil Constantinescu   }
958f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
959f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
960f68a32c8SEmil Constantinescu     if (match) {
961f68a32c8SEmil Constantinescu       ierr = TSReset_RK(ts);CHKERRQ(ierr);
962f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
963f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
964f68a32c8SEmil Constantinescu     }
965f68a32c8SEmil Constantinescu   }
966f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
967e4dd521cSBarry Smith   PetscFunctionReturn(0);
968e4dd521cSBarry Smith }
969e4dd521cSBarry Smith 
970ff22ae23SHong Zhang #undef __FUNCT__
971ff22ae23SHong Zhang #define __FUNCT__ "TSGetStages_RK"
972ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
973ff22ae23SHong Zhang {
974ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
975ff22ae23SHong Zhang 
976ff22ae23SHong Zhang   PetscFunctionBegin;
977ff22ae23SHong Zhang   *ns = rk->tableau->s;
978d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
979ff22ae23SHong Zhang   PetscFunctionReturn(0);
980ff22ae23SHong Zhang }
981ff22ae23SHong Zhang 
982ff22ae23SHong Zhang 
983e4dd521cSBarry Smith /* ------------------------------------------------------------ */
984ebe3b25bSBarry Smith /*MC
985f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
986ebe3b25bSBarry Smith 
987f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
988f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
989ebe3b25bSBarry Smith 
990f68a32c8SEmil Constantinescu   Notes:
991f68a32c8SEmil Constantinescu   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
992ebe3b25bSBarry Smith 
993d41222bbSBarry Smith   Level: beginner
994d41222bbSBarry Smith 
995f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
996f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
997ebe3b25bSBarry Smith 
998ebe3b25bSBarry Smith M*/
999e4dd521cSBarry Smith #undef __FUNCT__
10005f70b478SJed Brown #define __FUNCT__ "TSCreate_RK"
10018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
1002e4dd521cSBarry Smith {
1003f68a32c8SEmil Constantinescu   TS_RK     *th;
1004dfbe8321SBarry Smith   PetscErrorCode ierr;
1005e4dd521cSBarry Smith 
1006e4dd521cSBarry Smith   PetscFunctionBegin;
1007f68a32c8SEmil Constantinescu #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1008f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
1009f68a32c8SEmil Constantinescu #endif
1010f68a32c8SEmil Constantinescu 
1011f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
10125f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
10135f70b478SJed Brown   ts->ops->view           = TSView_RK;
1014f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
1015f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
101642f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1017f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1018f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1019f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1020f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1021ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
102242f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1023e4dd521cSBarry Smith 
1024*0f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
1025*0f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
1026*0f7a1166SHong Zhang 
1027b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
1028f68a32c8SEmil Constantinescu   ts->data = (void*)th;
1029e4dd521cSBarry Smith 
1030f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1031f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1032e4dd521cSBarry Smith   PetscFunctionReturn(0);
1033e4dd521cSBarry Smith }
1034