xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision ad8e26048d5e2f552b5f15f71fe8cd76ee6c3899)
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 */
10b45d2f2cSJed Brown #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 */
41*ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage*/
42*ad8e2604SHong 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*/
44f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
45f68a32c8SEmil Constantinescu   PetscReal    stage_time;
46f68a32c8SEmil Constantinescu   TSStepStatus status;
475f70b478SJed Brown } TS_RK;
48e4dd521cSBarry Smith 
49f68a32c8SEmil Constantinescu /*MC
50f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
51262ff7bbSBarry Smith 
52f68a32c8SEmil Constantinescu      This method has one stage.
53f68a32c8SEmil Constantinescu 
54f68a32c8SEmil Constantinescu      Level: advanced
55f68a32c8SEmil Constantinescu 
56f68a32c8SEmil Constantinescu .seealso: TSRK
57f68a32c8SEmil Constantinescu M*/
58f68a32c8SEmil Constantinescu /*MC
592109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
60f68a32c8SEmil Constantinescu 
61f68a32c8SEmil Constantinescu      This method has two stages.
62f68a32c8SEmil Constantinescu 
63f68a32c8SEmil Constantinescu      Level: advanced
64f68a32c8SEmil Constantinescu 
65f68a32c8SEmil Constantinescu .seealso: TSRK
66f68a32c8SEmil Constantinescu M*/
67f68a32c8SEmil Constantinescu /*MC
68f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
69f68a32c8SEmil Constantinescu 
70f68a32c8SEmil Constantinescu      This method has three stages.
71f68a32c8SEmil Constantinescu 
72f68a32c8SEmil Constantinescu      Level: advanced
73f68a32c8SEmil Constantinescu 
74f68a32c8SEmil Constantinescu .seealso: TSRK
75f68a32c8SEmil Constantinescu M*/
76f68a32c8SEmil Constantinescu /*MC
772109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
782109b73fSDebojyoti Ghosh 
792109b73fSDebojyoti Ghosh      This method has four stages.
802109b73fSDebojyoti Ghosh 
812109b73fSDebojyoti Ghosh      Level: advanced
822109b73fSDebojyoti Ghosh 
832109b73fSDebojyoti Ghosh .seealso: TSRK
842109b73fSDebojyoti Ghosh M*/
852109b73fSDebojyoti Ghosh /*MC
86f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
87f68a32c8SEmil Constantinescu 
882e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
89f68a32c8SEmil Constantinescu 
90f68a32c8SEmil Constantinescu      Level: advanced
91f68a32c8SEmil Constantinescu 
92f68a32c8SEmil Constantinescu .seealso: TSRK
93f68a32c8SEmil Constantinescu M*/
94f68a32c8SEmil Constantinescu /*MC
952e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
96f68a32c8SEmil Constantinescu 
97f68a32c8SEmil Constantinescu      This method has six stages.
98f68a32c8SEmil Constantinescu 
99f68a32c8SEmil Constantinescu      Level: advanced
100f68a32c8SEmil Constantinescu 
101f68a32c8SEmil Constantinescu .seealso: TSRK
102f68a32c8SEmil Constantinescu M*/
1032109b73fSDebojyoti Ghosh /*MC
1042e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1052109b73fSDebojyoti Ghosh 
1062109b73fSDebojyoti Ghosh      This method has seven stages.
1072109b73fSDebojyoti Ghosh 
1082109b73fSDebojyoti Ghosh      Level: advanced
1092109b73fSDebojyoti Ghosh 
1102109b73fSDebojyoti Ghosh .seealso: TSRK
1112109b73fSDebojyoti Ghosh M*/
112262ff7bbSBarry Smith 
113262ff7bbSBarry Smith #undef __FUNCT__
114f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll"
115f68a32c8SEmil Constantinescu /*@C
116f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
117262ff7bbSBarry Smith 
118f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
119262ff7bbSBarry Smith 
120f68a32c8SEmil Constantinescu   Level: advanced
121262ff7bbSBarry Smith 
122f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
123262ff7bbSBarry Smith 
124f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
125262ff7bbSBarry Smith @*/
126f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
127262ff7bbSBarry Smith {
1284ac538c5SBarry Smith   PetscErrorCode ierr;
129262ff7bbSBarry Smith 
130262ff7bbSBarry Smith   PetscFunctionBegin;
131f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
132f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
133e4dd521cSBarry Smith 
134e4dd521cSBarry Smith   {
135f68a32c8SEmil Constantinescu     const PetscReal
136f68a32c8SEmil Constantinescu       A[1][1] = {{0.0}},
137f68a32c8SEmil Constantinescu       b[1]    = {1.0};
138f68a32c8SEmil Constantinescu     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
139e4dd521cSBarry Smith   }
140db046809SBarry Smith   {
141f68a32c8SEmil Constantinescu     const PetscReal
142f68a32c8SEmil Constantinescu       A[2][2]     = {{0.0,0.0},
143f68a32c8SEmil Constantinescu                     {1.0,0.0}},
144f68a32c8SEmil Constantinescu       b[2]        = {0.5,0.5},
145f68a32c8SEmil Constantinescu       bembed[2]   = {1.0,0};
146fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr);
147db046809SBarry Smith   }
148f68a32c8SEmil Constantinescu   {
149f68a32c8SEmil Constantinescu     const PetscReal
150f68a32c8SEmil Constantinescu       A[3][3] = {{0,0,0},
151f68a32c8SEmil Constantinescu                  {2.0/3.0,0,0},
152f68a32c8SEmil Constantinescu                  {-1.0/3.0,1.0,0}},
153f68a32c8SEmil Constantinescu       b[3]    = {0.25,0.5,0.25};
154fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr);
155fdefd5e5SDebojyoti Ghosh   }
156fdefd5e5SDebojyoti Ghosh   {
157fdefd5e5SDebojyoti Ghosh     const PetscReal
158fdefd5e5SDebojyoti Ghosh       A[4][4] = {{0,0,0,0},
159fdefd5e5SDebojyoti Ghosh                  {1.0/2.0,0,0,0},
160fdefd5e5SDebojyoti Ghosh                  {0,3.0/4.0,0,0},
161fdefd5e5SDebojyoti Ghosh                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
162fdefd5e5SDebojyoti Ghosh       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
163fdefd5e5SDebojyoti Ghosh       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
164fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr);
165db046809SBarry Smith   }
166f68a32c8SEmil Constantinescu   {
167f68a32c8SEmil Constantinescu     const PetscReal
168f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
169f68a32c8SEmil Constantinescu                  {0.5,0,0,0},
170f68a32c8SEmil Constantinescu                  {0,0.5,0,0},
171f68a32c8SEmil Constantinescu                  {0,0,1.0,0}},
172f68a32c8SEmil Constantinescu       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
173fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr);
174db046809SBarry Smith   }
175f68a32c8SEmil Constantinescu   {
176f68a32c8SEmil Constantinescu     const PetscReal
177f68a32c8SEmil Constantinescu       A[6][6]   = {{0,0,0,0,0,0},
178f68a32c8SEmil Constantinescu                    {0.25,0,0,0,0,0},
179f68a32c8SEmil Constantinescu                    {3.0/32.0,9.0/32.0,0,0,0,0},
180f68a32c8SEmil Constantinescu                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
181f68a32c8SEmil Constantinescu                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
182f68a32c8SEmil Constantinescu                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
183f68a32c8SEmil Constantinescu       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
184f68a32c8SEmil Constantinescu       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
185fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
186fdefd5e5SDebojyoti Ghosh   }
187fdefd5e5SDebojyoti Ghosh   {
188fdefd5e5SDebojyoti Ghosh     const PetscReal
189fdefd5e5SDebojyoti Ghosh       A[7][7]   = {{0,0,0,0,0,0,0},
190fdefd5e5SDebojyoti Ghosh                    {1.0/5.0,0,0,0,0,0,0},
191fdefd5e5SDebojyoti Ghosh                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
192fdefd5e5SDebojyoti Ghosh                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
193fdefd5e5SDebojyoti Ghosh                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
194fdefd5e5SDebojyoti Ghosh                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
195fdefd5e5SDebojyoti Ghosh                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
196fdefd5e5SDebojyoti 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},
197fdefd5e5SDebojyoti 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};
198fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
199f68a32c8SEmil Constantinescu   }
200db046809SBarry Smith   PetscFunctionReturn(0);
201db046809SBarry Smith }
202db046809SBarry Smith 
203e4dd521cSBarry Smith #undef __FUNCT__
204f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy"
205f68a32c8SEmil Constantinescu /*@C
206f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
207f68a32c8SEmil Constantinescu 
208f68a32c8SEmil Constantinescu    Not Collective
209f68a32c8SEmil Constantinescu 
210f68a32c8SEmil Constantinescu    Level: advanced
211f68a32c8SEmil Constantinescu 
212f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
213f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
214f68a32c8SEmil Constantinescu @*/
215f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
216e4dd521cSBarry Smith {
217f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
218f68a32c8SEmil Constantinescu   RKTableauLink link;
219f68a32c8SEmil Constantinescu 
220f68a32c8SEmil Constantinescu   PetscFunctionBegin;
221f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
222f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
223f68a32c8SEmil Constantinescu     RKTableauList = link->next;
224f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
225f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
226f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
227f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
228f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
229f68a32c8SEmil Constantinescu   }
230f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
231f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
232f68a32c8SEmil Constantinescu }
233f68a32c8SEmil Constantinescu 
234f68a32c8SEmil Constantinescu #undef __FUNCT__
235f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage"
236f68a32c8SEmil Constantinescu /*@C
237f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
238f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
239f68a32c8SEmil Constantinescu   when using static libraries.
240f68a32c8SEmil Constantinescu 
241f68a32c8SEmil Constantinescu   Level: developer
242f68a32c8SEmil Constantinescu 
243f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
244f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
245f68a32c8SEmil Constantinescu @*/
246f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
247f68a32c8SEmil Constantinescu {
248186e87acSLisandro Dalcin   PetscErrorCode ierr;
249e4dd521cSBarry Smith 
250e4dd521cSBarry Smith   PetscFunctionBegin;
251f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
252f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
253f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
254f68a32c8SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
255f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
256f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
257f68a32c8SEmil Constantinescu }
258186e87acSLisandro Dalcin 
259f68a32c8SEmil Constantinescu #undef __FUNCT__
260f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage"
261f68a32c8SEmil Constantinescu /*@C
262f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
263f68a32c8SEmil Constantinescu   called from PetscFinalize().
264186e87acSLisandro Dalcin 
265f68a32c8SEmil Constantinescu   Level: developer
266f68a32c8SEmil Constantinescu 
267f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
268f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
269f68a32c8SEmil Constantinescu @*/
270f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
271f68a32c8SEmil Constantinescu {
272f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
273f68a32c8SEmil Constantinescu 
274f68a32c8SEmil Constantinescu   PetscFunctionBegin;
275f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
276f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
277f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
278f68a32c8SEmil Constantinescu }
279f68a32c8SEmil Constantinescu 
280f68a32c8SEmil Constantinescu #undef __FUNCT__
281f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister"
282f68a32c8SEmil Constantinescu /*@C
283f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
284f68a32c8SEmil Constantinescu 
285f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
286f68a32c8SEmil Constantinescu 
287f68a32c8SEmil Constantinescu    Input Parameters:
288f68a32c8SEmil Constantinescu +  name - identifier for method
289f68a32c8SEmil Constantinescu .  order - approximation order of method
290f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
291f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
292f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
293f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
294f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
295f68a32c8SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterp
296f68a32c8SEmil Constantinescu -  binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt)
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu    Notes:
299f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
300f68a32c8SEmil Constantinescu 
301f68a32c8SEmil Constantinescu    Level: advanced
302f68a32c8SEmil Constantinescu 
303f68a32c8SEmil Constantinescu .keywords: TS, register
304f68a32c8SEmil Constantinescu 
305f68a32c8SEmil Constantinescu .seealso: TSRK
306f68a32c8SEmil Constantinescu @*/
307f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
308f68a32c8SEmil Constantinescu                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
309f68a32c8SEmil Constantinescu                                  const PetscReal bembed[],
310f68a32c8SEmil Constantinescu                                  PetscInt pinterp,const PetscReal binterp[])
311f68a32c8SEmil Constantinescu {
312f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
313f68a32c8SEmil Constantinescu   RKTableauLink   link;
314f68a32c8SEmil Constantinescu   RKTableau       t;
315f68a32c8SEmil Constantinescu   PetscInt        i,j;
316f68a32c8SEmil Constantinescu 
317f68a32c8SEmil Constantinescu   PetscFunctionBegin;
318f68a32c8SEmil Constantinescu   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
319f68a32c8SEmil Constantinescu   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
320f68a32c8SEmil Constantinescu   t        = &link->tab;
321f68a32c8SEmil Constantinescu   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
322f68a32c8SEmil Constantinescu   t->order = order;
323f68a32c8SEmil Constantinescu   t->s     = s;
324dcca6d9dSJed Brown   ierr     = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
325f68a32c8SEmil Constantinescu   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
326f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
327f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
328f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
329f68a32c8SEmil 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];
330d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
331d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
332f68a32c8SEmil Constantinescu   if (bembed) {
333785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
334f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
335f68a32c8SEmil Constantinescu   }
336f68a32c8SEmil Constantinescu 
337f68a32c8SEmil Constantinescu   t->pinterp     = pinterp;
338785e854fSJed Brown   ierr           = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr);
339f68a32c8SEmil Constantinescu   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
340f68a32c8SEmil Constantinescu   link->next     = RKTableauList;
341f68a32c8SEmil Constantinescu   RKTableauList = link;
342f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
343f68a32c8SEmil Constantinescu }
344f68a32c8SEmil Constantinescu 
345f68a32c8SEmil Constantinescu #undef __FUNCT__
346f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK"
347e4dd521cSBarry Smith /*
348f68a32c8SEmil Constantinescu  The step completion formula is
349e4dd521cSBarry Smith 
350f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
351f68a32c8SEmil Constantinescu 
352f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
353f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
354f68a32c8SEmil Constantinescu  We can write
355f68a32c8SEmil Constantinescu 
356f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
357f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
358f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
359f68a32c8SEmil Constantinescu 
360f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
361e4dd521cSBarry Smith */
362f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
363f68a32c8SEmil Constantinescu {
364f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
365f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
366f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
367f68a32c8SEmil Constantinescu   PetscReal      h;
368f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
369f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
370f68a32c8SEmil Constantinescu 
371f68a32c8SEmil Constantinescu   PetscFunctionBegin;
372f68a32c8SEmil Constantinescu   switch (rk->status) {
373f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
374f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
375f68a32c8SEmil Constantinescu     h = ts->time_step; break;
376f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
377f68a32c8SEmil Constantinescu     h = ts->time_step_prev; break;
378f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
379f68a32c8SEmil Constantinescu   }
380f68a32c8SEmil Constantinescu   if (order == tab->order) {
381f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
382f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
383f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
384f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
385f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
386f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
387f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
388f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
389f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
390f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
391f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
392f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
393f68a32c8SEmil Constantinescu     } else {                    /* Rollback and re-complete using (be-b) */
394f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
395f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
396f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
397f68a32c8SEmil Constantinescu     }
398f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
399f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
400f68a32c8SEmil Constantinescu   }
401f68a32c8SEmil Constantinescu unavailable:
402f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
403f68a32c8SEmil 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);
404f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
405f68a32c8SEmil Constantinescu }
406f68a32c8SEmil Constantinescu 
407f68a32c8SEmil Constantinescu #undef __FUNCT__
408f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK"
409f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
410f68a32c8SEmil Constantinescu {
411f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
412f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
413f68a32c8SEmil Constantinescu   const PetscInt   s    = tab->s;
414f68a32c8SEmil Constantinescu   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
415f68a32c8SEmil Constantinescu   PetscScalar     *w    = rk->work;
416f68a32c8SEmil Constantinescu   Vec             *Y    = rk->Y,*YdotRHS = rk->YdotRHS;
417f68a32c8SEmil Constantinescu   TSAdapt          adapt;
418f68a32c8SEmil Constantinescu   PetscInt         i,j,reject,next_scheme;
419f68a32c8SEmil Constantinescu   PetscReal        next_time_step;
420f68a32c8SEmil Constantinescu   PetscReal        t;
421f68a32c8SEmil Constantinescu   PetscBool        accept;
422f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu   PetscFunctionBegin;
425f68a32c8SEmil Constantinescu 
426f68a32c8SEmil Constantinescu   next_time_step = ts->time_step;
427f68a32c8SEmil Constantinescu   t              = ts->ptime;
428f68a32c8SEmil Constantinescu   accept         = PETSC_TRUE;
429f68a32c8SEmil Constantinescu   rk->status     = TS_STEP_INCOMPLETE;
430f68a32c8SEmil Constantinescu 
431f68a32c8SEmil Constantinescu 
432f68a32c8SEmil Constantinescu   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
433f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
434f68a32c8SEmil Constantinescu     ierr = TSPreStep(ts);CHKERRQ(ierr);
435f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
4369be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
4379be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
438f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
439f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
440f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
4419be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
442f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
443f68a32c8SEmil Constantinescu       ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
444f68a32c8SEmil Constantinescu       if (!accept) goto reject_step;
445f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
446f68a32c8SEmil Constantinescu     }
447f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
448f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
449f68a32c8SEmil Constantinescu 
450f68a32c8SEmil Constantinescu     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
451f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
452f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
453f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
454f68a32c8SEmil Constantinescu     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
455f68a32c8SEmil Constantinescu     if (accept) {
456f68a32c8SEmil Constantinescu       /* ignore next_scheme for now */
457f68a32c8SEmil Constantinescu       ts->ptime    += ts->time_step;
458f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
459e3caeda6SBarry Smith       ts->steps++;
460f68a32c8SEmil Constantinescu       rk->status = TS_STEP_COMPLETE;
461f68a32c8SEmil Constantinescu       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
462f68a32c8SEmil Constantinescu       break;
463f68a32c8SEmil Constantinescu     } else {                    /* Roll back the current step */
464f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
465f68a32c8SEmil Constantinescu       ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr);
466f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
467f68a32c8SEmil Constantinescu       rk->status   = TS_STEP_INCOMPLETE;
468f68a32c8SEmil Constantinescu     }
469f68a32c8SEmil Constantinescu reject_step: continue;
470f68a32c8SEmil Constantinescu   }
471f68a32c8SEmil Constantinescu   if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
472f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
473e4dd521cSBarry Smith }
474e4dd521cSBarry Smith 
475f68a32c8SEmil Constantinescu #undef __FUNCT__
476d2daff3dSHong Zhang #define __FUNCT__ "TSStepAdj_RK"
477c235aa8dSHong Zhang static PetscErrorCode TSStepAdj_RK(TS ts)
478d2daff3dSHong Zhang {
479c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
480c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
481c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
482c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
483c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
484*ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
485b18ea86cSHong Zhang   PetscInt         i,j,nadj;
486c235aa8dSHong Zhang   PetscReal        t;
487d2daff3dSHong Zhang   PetscErrorCode   ierr;
488c235aa8dSHong Zhang 
489d2daff3dSHong Zhang   PetscFunctionBegin;
490d2daff3dSHong Zhang 
491b18ea86cSHong Zhang   /*ierr = TSStep_RK(ts);CHKERRQ(ierr); // reuse TSStep */
492d2daff3dSHong Zhang 
493c235aa8dSHong Zhang   t          = ts->ptime;
494c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
495b18ea86cSHong Zhang   PetscReal h = ts->time_step;
496c235aa8dSHong Zhang   ierr = TSPreStep(ts);CHKERRQ(ierr);
497c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
498c235aa8dSHong Zhang     Mat J,Jp;
499c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
500b18ea86cSHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
501b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
502b18ea86cSHong Zhang       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);
503c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
504b18ea86cSHong Zhang         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);
505b18ea86cSHong Zhang       }
506c235aa8dSHong Zhang     }
507*ad8e2604SHong Zhang     /* Stage values of lambda */
508c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
509c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
510b18ea86cSHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
511b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
512b18ea86cSHong Zhang     }
513*ad8e2604SHong Zhang     /* Stage values of mu */
514*ad8e2604SHong Zhang     ierr = TSRHSJacobianP(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
515*ad8e2604SHong Zhang     for (nadj=0; nadj<ts->numberadjs; nadj++) {
516*ad8e2604SHong Zhang       ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
517*ad8e2604SHong Zhang     }
518c235aa8dSHong Zhang   }
519c235aa8dSHong Zhang 
520c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
521b18ea86cSHong Zhang   for (nadj=0; nadj<ts->numberadjs; nadj++) {
522b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
523*ad8e2604SHong Zhang     ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
524b18ea86cSHong Zhang   }
525c235aa8dSHong Zhang   ts->ptime += ts->time_step;
526c235aa8dSHong Zhang   ts->steps++;
527c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
528b18ea86cSHong Zhang   ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vecs_sensi,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
529*ad8e2604SHong Zhang   ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vecs_sensip,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
530d2daff3dSHong Zhang   PetscFunctionReturn(0);
531d2daff3dSHong Zhang }
532d2daff3dSHong Zhang 
533d2daff3dSHong Zhang #undef __FUNCT__
534f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK"
535f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
536f68a32c8SEmil Constantinescu {
537f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
538f68a32c8SEmil Constantinescu   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
539f68a32c8SEmil Constantinescu   PetscReal        h;
540f68a32c8SEmil Constantinescu   PetscReal        tt,t;
541f68a32c8SEmil Constantinescu   PetscScalar     *b;
542f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
543f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
544e4dd521cSBarry Smith 
545f68a32c8SEmil Constantinescu   PetscFunctionBegin;
546f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
547f68a32c8SEmil Constantinescu   switch (rk->status) {
548f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
549f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
550f68a32c8SEmil Constantinescu     h = ts->time_step;
551f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
552f68a32c8SEmil Constantinescu     break;
553f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
554f68a32c8SEmil Constantinescu     h = ts->time_step_prev;
555f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
556f68a32c8SEmil Constantinescu     break;
557f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
558e4dd521cSBarry Smith   }
559785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
560f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
561f68a32c8SEmil Constantinescu   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
562f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
563f68a32c8SEmil Constantinescu       b[i]  += h * B[i*pinterp+j] * tt;
564e4dd521cSBarry Smith     }
565f68a32c8SEmil Constantinescu   }
566f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
567f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
568f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
569e4dd521cSBarry Smith   PetscFunctionReturn(0);
570e4dd521cSBarry Smith }
571e4dd521cSBarry Smith 
572e4dd521cSBarry Smith /*------------------------------------------------------------*/
573e4dd521cSBarry Smith #undef __FUNCT__
574277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK"
575277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
576e4dd521cSBarry Smith {
5775f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
578f68a32c8SEmil Constantinescu   PetscInt       s;
5796849ba73SBarry Smith   PetscErrorCode ierr;
580e4dd521cSBarry Smith 
581e4dd521cSBarry Smith   PetscFunctionBegin;
582f68a32c8SEmil Constantinescu   if (!rk->tableau) PetscFunctionReturn(0);
583f68a32c8SEmil Constantinescu   s    = rk->tableau->s;
584f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr);
585f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr);
586c235aa8dSHong Zhang   if(ts->reverse_mode) {
587b18ea86cSHong Zhang     ierr = VecDestroyVecs(s*ts->numberadjs,&rk->VecDeltaLam);CHKERRQ(ierr);
588*ad8e2604SHong Zhang     if(rk->VecDeltaMu) {
589*ad8e2604SHong Zhang       ierr = VecDestroyVecs(s*ts->numberadjs,&rk->VecDeltaMu);CHKERRQ(ierr);
590*ad8e2604SHong Zhang     }
591b18ea86cSHong Zhang     ierr = VecDestroyVecs(ts->numberadjs,&rk->VecSensiTemp);CHKERRQ(ierr);
592c235aa8dSHong Zhang   }
593f68a32c8SEmil Constantinescu   ierr = PetscFree(rk->work);CHKERRQ(ierr);
594277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
595e4dd521cSBarry Smith }
596277b19d0SLisandro Dalcin 
597277b19d0SLisandro Dalcin #undef __FUNCT__
598277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK"
599277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
600277b19d0SLisandro Dalcin {
601277b19d0SLisandro Dalcin   PetscErrorCode ierr;
602277b19d0SLisandro Dalcin 
603277b19d0SLisandro Dalcin   PetscFunctionBegin;
604277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
605277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
606f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
607f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
608f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
609f68a32c8SEmil Constantinescu }
610f68a32c8SEmil Constantinescu 
611f68a32c8SEmil Constantinescu 
612f68a32c8SEmil Constantinescu #undef __FUNCT__
613f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK"
614f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
615f68a32c8SEmil Constantinescu {
616f68a32c8SEmil Constantinescu   PetscFunctionBegin;
617f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
618f68a32c8SEmil Constantinescu }
619f68a32c8SEmil Constantinescu 
620f68a32c8SEmil Constantinescu #undef __FUNCT__
621f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK"
622f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
623f68a32c8SEmil Constantinescu {
624f68a32c8SEmil Constantinescu   PetscFunctionBegin;
625f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
626f68a32c8SEmil Constantinescu }
627f68a32c8SEmil Constantinescu 
628f68a32c8SEmil Constantinescu 
629f68a32c8SEmil Constantinescu #undef __FUNCT__
630f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK"
631f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
632f68a32c8SEmil Constantinescu {
633f68a32c8SEmil Constantinescu   PetscFunctionBegin;
634f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
635f68a32c8SEmil Constantinescu }
636f68a32c8SEmil Constantinescu 
637f68a32c8SEmil Constantinescu #undef __FUNCT__
638f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
639f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
640f68a32c8SEmil Constantinescu {
641f68a32c8SEmil Constantinescu 
642f68a32c8SEmil Constantinescu   PetscFunctionBegin;
643f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
644f68a32c8SEmil Constantinescu }
645c235aa8dSHong Zhang /*
646f68a32c8SEmil Constantinescu #undef __FUNCT__
647d2daff3dSHong Zhang #define __FUNCT__ "RKSetAdjCoe"
648d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
649d2daff3dSHong Zhang {
650d2daff3dSHong Zhang   PetscReal *A,*b;
651d2daff3dSHong Zhang   PetscInt        s,i,j;
652d2daff3dSHong Zhang   PetscErrorCode  ierr;
653d2daff3dSHong Zhang 
654d2daff3dSHong Zhang   PetscFunctionBegin;
655d2daff3dSHong Zhang   s    = tab->s;
656d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
657d2daff3dSHong Zhang 
658d2daff3dSHong Zhang   for (i=0; i<s; i++)
659d2daff3dSHong Zhang     for (j=0; j<s; j++) {
660d2daff3dSHong 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];
661d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
662d2daff3dSHong Zhang     }
663d2daff3dSHong Zhang 
664d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
665d2daff3dSHong Zhang 
666d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
667d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
668d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
669d2daff3dSHong Zhang   PetscFunctionReturn(0);
670d2daff3dSHong Zhang }
671c235aa8dSHong Zhang */
672d2daff3dSHong Zhang #undef __FUNCT__
673f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK"
674f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
675f68a32c8SEmil Constantinescu {
676f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
677f68a32c8SEmil Constantinescu   RKTableau      tab;
678f68a32c8SEmil Constantinescu   PetscInt       s;
679f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
680f68a32c8SEmil Constantinescu   DM             dm;
681f68a32c8SEmil Constantinescu 
682f68a32c8SEmil Constantinescu   PetscFunctionBegin;
683f68a32c8SEmil Constantinescu   if (!rk->tableau) {
684f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
685f68a32c8SEmil Constantinescu   }
686f68a32c8SEmil Constantinescu   tab  = rk->tableau;
687f68a32c8SEmil Constantinescu   s    = tab->s;
688*ad8e2604SHong Zhang   /* old */
689*ad8e2604SHong Zhang   /*if (ts->reverse_mode) {
690*ad8e2604SHong Zhang     ierr = RKSetAdjCoe(tab);CHKERRQ(ierr);
691*ad8e2604SHong Zhang   }*/
692f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr);
693f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr);
694c235aa8dSHong Zhang   if (ts->reverse_mode) {
695b18ea86cSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numberadjs,&rk->VecDeltaLam);CHKERRQ(ierr);
696*ad8e2604SHong Zhang     if(ts->vecs_sensip) {
697*ad8e2604SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numberadjs,&rk->VecDeltaMu);CHKERRQ(ierr);
698*ad8e2604SHong Zhang     }
699b18ea86cSHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numberadjs,&rk->VecSensiTemp);CHKERRQ(ierr);
700c235aa8dSHong Zhang   }
701785e854fSJed Brown   ierr = PetscMalloc1(s,&rk->work);CHKERRQ(ierr);
702f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
703f68a32c8SEmil Constantinescu   if (dm) {
704f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
705f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
706f68a32c8SEmil Constantinescu   }
707e4dd521cSBarry Smith   PetscFunctionReturn(0);
708e4dd521cSBarry Smith }
709d2daff3dSHong Zhang 
710e4dd521cSBarry Smith /*------------------------------------------------------------*/
711e4dd521cSBarry Smith 
712e4dd521cSBarry Smith #undef __FUNCT__
7135f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK"
7148c34d3f5SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptions *PetscOptionsObject,TS ts)
715e4dd521cSBarry Smith {
716dfbe8321SBarry Smith   PetscErrorCode ierr;
717f68a32c8SEmil Constantinescu   char           rktype[256];
718262ff7bbSBarry Smith 
719e4dd521cSBarry Smith   PetscFunctionBegin;
720e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
721f68a32c8SEmil Constantinescu   {
722f68a32c8SEmil Constantinescu     RKTableauLink  link;
723f68a32c8SEmil Constantinescu     PetscInt       count,choice;
724f68a32c8SEmil Constantinescu     PetscBool      flg;
725f68a32c8SEmil Constantinescu     const char   **namelist;
726f68a32c8SEmil Constantinescu     ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr);
727f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
728785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
729f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
730f68a32c8SEmil Constantinescu     ierr      = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr);
731f68a32c8SEmil Constantinescu     ierr      = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr);
732f68a32c8SEmil Constantinescu     ierr      = PetscFree(namelist);CHKERRQ(ierr);
733f68a32c8SEmil Constantinescu   }
734262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
735e4dd521cSBarry Smith   PetscFunctionReturn(0);
736e4dd521cSBarry Smith }
737e4dd521cSBarry Smith 
738e4dd521cSBarry Smith #undef __FUNCT__
739f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray"
740f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
741f68a32c8SEmil Constantinescu {
742f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
743f68a32c8SEmil Constantinescu   PetscInt       i;
744f68a32c8SEmil Constantinescu   size_t         left,count;
745f68a32c8SEmil Constantinescu   char           *p;
746f68a32c8SEmil Constantinescu 
747f68a32c8SEmil Constantinescu   PetscFunctionBegin;
748f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
749f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
750f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
751f68a32c8SEmil Constantinescu     left -= count;
752f68a32c8SEmil Constantinescu     p    += count;
753f68a32c8SEmil Constantinescu     *p++  = ' ';
754f68a32c8SEmil Constantinescu   }
755f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
756f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
757f68a32c8SEmil Constantinescu }
758f68a32c8SEmil Constantinescu 
759f68a32c8SEmil Constantinescu #undef __FUNCT__
7605f70b478SJed Brown #define __FUNCT__ "TSView_RK"
7615f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
762e4dd521cSBarry Smith {
7635f70b478SJed Brown   TS_RK         *rk   = (TS_RK*)ts->data;
764f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
7658370ee3bSLisandro Dalcin   PetscBool      iascii;
766dfbe8321SBarry Smith   PetscErrorCode ierr;
767f68a32c8SEmil Constantinescu   TSAdapt        adapt;
7682cdf8120Sasbjorn 
769e4dd521cSBarry Smith   PetscFunctionBegin;
770251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
7718370ee3bSLisandro Dalcin   if (iascii) {
772f68a32c8SEmil Constantinescu     TSRKType rktype;
773f68a32c8SEmil Constantinescu     char     buf[512];
774f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
775f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
776f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
777f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
778d760c35bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
7798370ee3bSLisandro Dalcin   }
780f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
781f68a32c8SEmil Constantinescu   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
782f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
783f68a32c8SEmil Constantinescu }
784f68a32c8SEmil Constantinescu 
785f68a32c8SEmil Constantinescu #undef __FUNCT__
786f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK"
787f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
788f68a32c8SEmil Constantinescu {
789f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
790f68a32c8SEmil Constantinescu   TSAdapt        tsadapt;
791f68a32c8SEmil Constantinescu 
792f68a32c8SEmil Constantinescu   PetscFunctionBegin;
793f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
794f68a32c8SEmil Constantinescu   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
795f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
796f68a32c8SEmil Constantinescu }
797f68a32c8SEmil Constantinescu 
798f68a32c8SEmil Constantinescu #undef __FUNCT__
799f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType"
800f68a32c8SEmil Constantinescu /*@C
801f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
802f68a32c8SEmil Constantinescu 
803f68a32c8SEmil Constantinescu   Logically collective
804f68a32c8SEmil Constantinescu 
805f68a32c8SEmil Constantinescu   Input Parameter:
806f68a32c8SEmil Constantinescu +  ts - timestepping context
807f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
808f68a32c8SEmil Constantinescu 
809f68a32c8SEmil Constantinescu   Level: intermediate
810f68a32c8SEmil Constantinescu 
811f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
812f68a32c8SEmil Constantinescu @*/
813f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
814f68a32c8SEmil Constantinescu {
815f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
816f68a32c8SEmil Constantinescu 
817f68a32c8SEmil Constantinescu   PetscFunctionBegin;
818f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
819f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
820f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
821f68a32c8SEmil Constantinescu }
822f68a32c8SEmil Constantinescu 
823f68a32c8SEmil Constantinescu #undef __FUNCT__
824f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType"
825f68a32c8SEmil Constantinescu /*@C
826f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
827f68a32c8SEmil Constantinescu 
828f68a32c8SEmil Constantinescu   Logically collective
829f68a32c8SEmil Constantinescu 
830f68a32c8SEmil Constantinescu   Input Parameter:
831f68a32c8SEmil Constantinescu .  ts - timestepping context
832f68a32c8SEmil Constantinescu 
833f68a32c8SEmil Constantinescu   Output Parameter:
834f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
835f68a32c8SEmil Constantinescu 
836f68a32c8SEmil Constantinescu   Level: intermediate
837f68a32c8SEmil Constantinescu 
838f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
839f68a32c8SEmil Constantinescu @*/
840f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
841f68a32c8SEmil Constantinescu {
842f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
843f68a32c8SEmil Constantinescu 
844f68a32c8SEmil Constantinescu   PetscFunctionBegin;
845f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
846f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
847f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
848f68a32c8SEmil Constantinescu }
849f68a32c8SEmil Constantinescu 
850f68a32c8SEmil Constantinescu #undef __FUNCT__
851f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK"
852f68a32c8SEmil Constantinescu PetscErrorCode  TSRKGetType_RK(TS ts,TSRKType *rktype)
853f68a32c8SEmil Constantinescu {
854f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
855f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
856f68a32c8SEmil Constantinescu 
857f68a32c8SEmil Constantinescu   PetscFunctionBegin;
858f68a32c8SEmil Constantinescu   if (!rk->tableau) {
859f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
860f68a32c8SEmil Constantinescu   }
861f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
862f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
863f68a32c8SEmil Constantinescu }
864f68a32c8SEmil Constantinescu #undef __FUNCT__
865f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK"
866f68a32c8SEmil Constantinescu PetscErrorCode  TSRKSetType_RK(TS ts,TSRKType rktype)
867f68a32c8SEmil Constantinescu {
868f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
869f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
870f68a32c8SEmil Constantinescu   PetscBool      match;
871f68a32c8SEmil Constantinescu   RKTableauLink link;
872f68a32c8SEmil Constantinescu 
873f68a32c8SEmil Constantinescu   PetscFunctionBegin;
874f68a32c8SEmil Constantinescu   if (rk->tableau) {
875f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
876f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
877f68a32c8SEmil Constantinescu   }
878f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
879f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
880f68a32c8SEmil Constantinescu     if (match) {
881f68a32c8SEmil Constantinescu       ierr = TSReset_RK(ts);CHKERRQ(ierr);
882f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
883f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
884f68a32c8SEmil Constantinescu     }
885f68a32c8SEmil Constantinescu   }
886f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
887e4dd521cSBarry Smith   PetscFunctionReturn(0);
888e4dd521cSBarry Smith }
889e4dd521cSBarry Smith 
890ff22ae23SHong Zhang #undef __FUNCT__
891ff22ae23SHong Zhang #define __FUNCT__ "TSGetStages_RK"
892ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
893ff22ae23SHong Zhang {
894ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
895ff22ae23SHong Zhang 
896ff22ae23SHong Zhang   PetscFunctionBegin;
897ff22ae23SHong Zhang   *ns = rk->tableau->s;
898d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
899ff22ae23SHong Zhang   PetscFunctionReturn(0);
900ff22ae23SHong Zhang }
901ff22ae23SHong Zhang 
902ff22ae23SHong Zhang 
903e4dd521cSBarry Smith /* ------------------------------------------------------------ */
904ebe3b25bSBarry Smith /*MC
905f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
906ebe3b25bSBarry Smith 
907f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
908f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
909ebe3b25bSBarry Smith 
910f68a32c8SEmil Constantinescu   Notes:
911f68a32c8SEmil Constantinescu   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
912ebe3b25bSBarry Smith 
913d41222bbSBarry Smith   Level: beginner
914d41222bbSBarry Smith 
915f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
916f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
917ebe3b25bSBarry Smith 
918ebe3b25bSBarry Smith M*/
919e4dd521cSBarry Smith #undef __FUNCT__
9205f70b478SJed Brown #define __FUNCT__ "TSCreate_RK"
9218cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
922e4dd521cSBarry Smith {
923f68a32c8SEmil Constantinescu   TS_RK     *th;
924dfbe8321SBarry Smith   PetscErrorCode ierr;
925e4dd521cSBarry Smith 
926e4dd521cSBarry Smith   PetscFunctionBegin;
927f68a32c8SEmil Constantinescu #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
928f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
929f68a32c8SEmil Constantinescu #endif
930f68a32c8SEmil Constantinescu 
931f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
9325f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
9335f70b478SJed Brown   ts->ops->view           = TSView_RK;
934f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
935f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
936f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
937f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
938f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
939f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
940ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
941d2daff3dSHong Zhang   ts->ops->stepadj        = TSStepAdj_RK;
942e4dd521cSBarry Smith 
943b00a9115SJed Brown   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
944f68a32c8SEmil Constantinescu   ts->data = (void*)th;
945e4dd521cSBarry Smith 
946f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
947f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
948e4dd521cSBarry Smith   PetscFunctionReturn(0);
949e4dd521cSBarry Smith }
950