xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision 484bcdc76a0107e2a88fdaf8f414b7c9da036582)
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 
13*484bcdc7SDebojyoti 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 */
41f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
42f68a32c8SEmil Constantinescu   PetscReal    stage_time;
43f68a32c8SEmil Constantinescu   TSStepStatus status;
445f70b478SJed Brown } TS_RK;
45e4dd521cSBarry Smith 
46f68a32c8SEmil Constantinescu /*MC
47f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
48262ff7bbSBarry Smith 
49f68a32c8SEmil Constantinescu      This method has one stage.
50f68a32c8SEmil Constantinescu 
51f68a32c8SEmil Constantinescu      Level: advanced
52f68a32c8SEmil Constantinescu 
53f68a32c8SEmil Constantinescu .seealso: TSRK
54f68a32c8SEmil Constantinescu M*/
55f68a32c8SEmil Constantinescu /*MC
562109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
57f68a32c8SEmil Constantinescu 
58f68a32c8SEmil Constantinescu      This method has two stages.
59f68a32c8SEmil Constantinescu 
60f68a32c8SEmil Constantinescu      Level: advanced
61f68a32c8SEmil Constantinescu 
62f68a32c8SEmil Constantinescu .seealso: TSRK
63f68a32c8SEmil Constantinescu M*/
64f68a32c8SEmil Constantinescu /*MC
65f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
66f68a32c8SEmil Constantinescu 
67f68a32c8SEmil Constantinescu      This method has three stages.
68f68a32c8SEmil Constantinescu 
69f68a32c8SEmil Constantinescu      Level: advanced
70f68a32c8SEmil Constantinescu 
71f68a32c8SEmil Constantinescu .seealso: TSRK
72f68a32c8SEmil Constantinescu M*/
73f68a32c8SEmil Constantinescu /*MC
742109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
752109b73fSDebojyoti Ghosh 
762109b73fSDebojyoti Ghosh      This method has four stages.
772109b73fSDebojyoti Ghosh 
782109b73fSDebojyoti Ghosh      Level: advanced
792109b73fSDebojyoti Ghosh 
802109b73fSDebojyoti Ghosh .seealso: TSRK
812109b73fSDebojyoti Ghosh M*/
822109b73fSDebojyoti Ghosh /*MC
83f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
84f68a32c8SEmil Constantinescu 
852e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
86f68a32c8SEmil Constantinescu 
87f68a32c8SEmil Constantinescu      Level: advanced
88f68a32c8SEmil Constantinescu 
89f68a32c8SEmil Constantinescu .seealso: TSRK
90f68a32c8SEmil Constantinescu M*/
91f68a32c8SEmil Constantinescu /*MC
922e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
93f68a32c8SEmil Constantinescu 
94f68a32c8SEmil Constantinescu      This method has six stages.
95f68a32c8SEmil Constantinescu 
96f68a32c8SEmil Constantinescu      Level: advanced
97f68a32c8SEmil Constantinescu 
98f68a32c8SEmil Constantinescu .seealso: TSRK
99f68a32c8SEmil Constantinescu M*/
1002109b73fSDebojyoti Ghosh /*MC
1012e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1022109b73fSDebojyoti Ghosh 
1032109b73fSDebojyoti Ghosh      This method has seven stages.
1042109b73fSDebojyoti Ghosh 
1052109b73fSDebojyoti Ghosh      Level: advanced
1062109b73fSDebojyoti Ghosh 
1072109b73fSDebojyoti Ghosh .seealso: TSRK
1082109b73fSDebojyoti Ghosh M*/
109262ff7bbSBarry Smith 
110262ff7bbSBarry Smith #undef __FUNCT__
111f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterAll"
112f68a32c8SEmil Constantinescu /*@C
113f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
114262ff7bbSBarry Smith 
115f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
116262ff7bbSBarry Smith 
117f68a32c8SEmil Constantinescu   Level: advanced
118262ff7bbSBarry Smith 
119f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
120262ff7bbSBarry Smith 
121f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
122262ff7bbSBarry Smith @*/
123f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
124262ff7bbSBarry Smith {
1254ac538c5SBarry Smith   PetscErrorCode ierr;
126262ff7bbSBarry Smith 
127262ff7bbSBarry Smith   PetscFunctionBegin;
128f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
129f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
130e4dd521cSBarry Smith 
131e4dd521cSBarry Smith   {
132f68a32c8SEmil Constantinescu     const PetscReal
133f68a32c8SEmil Constantinescu       A[1][1] = {{0.0}},
134f68a32c8SEmil Constantinescu       b[1]    = {1.0};
135f68a32c8SEmil Constantinescu     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,1,b);CHKERRQ(ierr);
136e4dd521cSBarry Smith   }
137db046809SBarry Smith   {
138f68a32c8SEmil Constantinescu     const PetscReal
139f68a32c8SEmil Constantinescu       A[2][2]     = {{0.0,0.0},
140f68a32c8SEmil Constantinescu                     {1.0,0.0}},
141f68a32c8SEmil Constantinescu       b[2]        = {0.5,0.5},
142f68a32c8SEmil Constantinescu       bembed[2]   = {1.0,0};
143fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,2,b);CHKERRQ(ierr);
144db046809SBarry Smith   }
145f68a32c8SEmil Constantinescu   {
146f68a32c8SEmil Constantinescu     const PetscReal
147f68a32c8SEmil Constantinescu       A[3][3] = {{0,0,0},
148f68a32c8SEmil Constantinescu                  {2.0/3.0,0,0},
149f68a32c8SEmil Constantinescu                  {-1.0/3.0,1.0,0}},
150f68a32c8SEmil Constantinescu       b[3]    = {0.25,0.5,0.25};
151fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,3,b);CHKERRQ(ierr);
152fdefd5e5SDebojyoti Ghosh   }
153fdefd5e5SDebojyoti Ghosh   {
154fdefd5e5SDebojyoti Ghosh     const PetscReal
155fdefd5e5SDebojyoti Ghosh       A[4][4] = {{0,0,0,0},
156fdefd5e5SDebojyoti Ghosh                  {1.0/2.0,0,0,0},
157fdefd5e5SDebojyoti Ghosh                  {0,3.0/4.0,0,0},
158fdefd5e5SDebojyoti Ghosh                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
159fdefd5e5SDebojyoti Ghosh       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
160fdefd5e5SDebojyoti Ghosh       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
161fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,3,b);CHKERRQ(ierr);
162db046809SBarry Smith   }
163f68a32c8SEmil Constantinescu   {
164f68a32c8SEmil Constantinescu     const PetscReal
165f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
166f68a32c8SEmil Constantinescu                  {0.5,0,0,0},
167f68a32c8SEmil Constantinescu                  {0,0.5,0,0},
168f68a32c8SEmil Constantinescu                  {0,0,1.0,0}},
169f68a32c8SEmil Constantinescu       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
170fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,4,b);CHKERRQ(ierr);
171db046809SBarry Smith   }
172f68a32c8SEmil Constantinescu   {
173f68a32c8SEmil Constantinescu     const PetscReal
174f68a32c8SEmil Constantinescu       A[6][6]   = {{0,0,0,0,0,0},
175f68a32c8SEmil Constantinescu                    {0.25,0,0,0,0,0},
176f68a32c8SEmil Constantinescu                    {3.0/32.0,9.0/32.0,0,0,0,0},
177f68a32c8SEmil Constantinescu                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
178f68a32c8SEmil Constantinescu                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
179f68a32c8SEmil Constantinescu                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
180f68a32c8SEmil Constantinescu       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
181f68a32c8SEmil Constantinescu       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
182fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
183fdefd5e5SDebojyoti Ghosh   }
184fdefd5e5SDebojyoti Ghosh   {
185fdefd5e5SDebojyoti Ghosh     const PetscReal
186fdefd5e5SDebojyoti Ghosh       A[7][7]   = {{0,0,0,0,0,0,0},
187fdefd5e5SDebojyoti Ghosh                    {1.0/5.0,0,0,0,0,0,0},
188fdefd5e5SDebojyoti Ghosh                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
189fdefd5e5SDebojyoti Ghosh                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
190fdefd5e5SDebojyoti Ghosh                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
191fdefd5e5SDebojyoti Ghosh                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
192fdefd5e5SDebojyoti Ghosh                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
193fdefd5e5SDebojyoti 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},
194fdefd5e5SDebojyoti 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};
195fdefd5e5SDebojyoti Ghosh     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,5,b);CHKERRQ(ierr);
196f68a32c8SEmil Constantinescu   }
197db046809SBarry Smith   PetscFunctionReturn(0);
198db046809SBarry Smith }
199db046809SBarry Smith 
200e4dd521cSBarry Smith #undef __FUNCT__
201f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegisterDestroy"
202f68a32c8SEmil Constantinescu /*@C
203f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
204f68a32c8SEmil Constantinescu 
205f68a32c8SEmil Constantinescu    Not Collective
206f68a32c8SEmil Constantinescu 
207f68a32c8SEmil Constantinescu    Level: advanced
208f68a32c8SEmil Constantinescu 
209f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
210f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
211f68a32c8SEmil Constantinescu @*/
212f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
213e4dd521cSBarry Smith {
214f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
215f68a32c8SEmil Constantinescu   RKTableauLink link;
216f68a32c8SEmil Constantinescu 
217f68a32c8SEmil Constantinescu   PetscFunctionBegin;
218f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
219f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
220f68a32c8SEmil Constantinescu     RKTableauList = link->next;
221f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
222f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
223f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
224f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
225f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
226f68a32c8SEmil Constantinescu   }
227f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
228f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
229f68a32c8SEmil Constantinescu }
230f68a32c8SEmil Constantinescu 
231f68a32c8SEmil Constantinescu #undef __FUNCT__
232f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKInitializePackage"
233f68a32c8SEmil Constantinescu /*@C
234f68a32c8SEmil Constantinescu   TSRKInitializePackage - This function initializes everything in the TSRK package. It is called
235f68a32c8SEmil Constantinescu   from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK()
236f68a32c8SEmil Constantinescu   when using static libraries.
237f68a32c8SEmil Constantinescu 
238f68a32c8SEmil Constantinescu   Level: developer
239f68a32c8SEmil Constantinescu 
240f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package
241f68a32c8SEmil Constantinescu .seealso: PetscInitialize()
242f68a32c8SEmil Constantinescu @*/
243f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void)
244f68a32c8SEmil Constantinescu {
245186e87acSLisandro Dalcin   PetscErrorCode ierr;
246e4dd521cSBarry Smith 
247e4dd521cSBarry Smith   PetscFunctionBegin;
248f68a32c8SEmil Constantinescu   if (TSRKPackageInitialized) PetscFunctionReturn(0);
249f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_TRUE;
250f68a32c8SEmil Constantinescu   ierr = TSRKRegisterAll();CHKERRQ(ierr);
251f68a32c8SEmil Constantinescu   ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr);
252f68a32c8SEmil Constantinescu   ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
253f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
254f68a32c8SEmil Constantinescu }
255186e87acSLisandro Dalcin 
256f68a32c8SEmil Constantinescu #undef __FUNCT__
257f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKFinalizePackage"
258f68a32c8SEmil Constantinescu /*@C
259f68a32c8SEmil Constantinescu   TSRKFinalizePackage - This function destroys everything in the TSRK package. It is
260f68a32c8SEmil Constantinescu   called from PetscFinalize().
261186e87acSLisandro Dalcin 
262f68a32c8SEmil Constantinescu   Level: developer
263f68a32c8SEmil Constantinescu 
264f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package
265f68a32c8SEmil Constantinescu .seealso: PetscFinalize()
266f68a32c8SEmil Constantinescu @*/
267f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void)
268f68a32c8SEmil Constantinescu {
269f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
270f68a32c8SEmil Constantinescu 
271f68a32c8SEmil Constantinescu   PetscFunctionBegin;
272f68a32c8SEmil Constantinescu   TSRKPackageInitialized = PETSC_FALSE;
273f68a32c8SEmil Constantinescu   ierr = TSRKRegisterDestroy();CHKERRQ(ierr);
274f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
275f68a32c8SEmil Constantinescu }
276f68a32c8SEmil Constantinescu 
277f68a32c8SEmil Constantinescu #undef __FUNCT__
278f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKRegister"
279f68a32c8SEmil Constantinescu /*@C
280f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
281f68a32c8SEmil Constantinescu 
282f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
283f68a32c8SEmil Constantinescu 
284f68a32c8SEmil Constantinescu    Input Parameters:
285f68a32c8SEmil Constantinescu +  name - identifier for method
286f68a32c8SEmil Constantinescu .  order - approximation order of method
287f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
288f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
289f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
290f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
291f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
292f68a32c8SEmil Constantinescu .  pinterp - Order of the interpolation scheme, equal to the number of columns of binterp
293f68a32c8SEmil Constantinescu -  binterp - Coefficients of the interpolation formula (dimension s*pinterp; NULL to reuse binterpt)
294f68a32c8SEmil Constantinescu 
295f68a32c8SEmil Constantinescu    Notes:
296f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu    Level: advanced
299f68a32c8SEmil Constantinescu 
300f68a32c8SEmil Constantinescu .keywords: TS, register
301f68a32c8SEmil Constantinescu 
302f68a32c8SEmil Constantinescu .seealso: TSRK
303f68a32c8SEmil Constantinescu @*/
304f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
305f68a32c8SEmil Constantinescu                                  const PetscReal A[],const PetscReal b[],const PetscReal c[],
306f68a32c8SEmil Constantinescu                                  const PetscReal bembed[],
307f68a32c8SEmil Constantinescu                                  PetscInt pinterp,const PetscReal binterp[])
308f68a32c8SEmil Constantinescu {
309f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
310f68a32c8SEmil Constantinescu   RKTableauLink   link;
311f68a32c8SEmil Constantinescu   RKTableau       t;
312f68a32c8SEmil Constantinescu   PetscInt        i,j;
313f68a32c8SEmil Constantinescu 
314f68a32c8SEmil Constantinescu   PetscFunctionBegin;
315f68a32c8SEmil Constantinescu   ierr     = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr);
316f68a32c8SEmil Constantinescu   ierr     = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr);
317f68a32c8SEmil Constantinescu   t        = &link->tab;
318f68a32c8SEmil Constantinescu   ierr     = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
319f68a32c8SEmil Constantinescu   t->order = order;
320f68a32c8SEmil Constantinescu   t->s     = s;
321f68a32c8SEmil Constantinescu   ierr     = PetscMalloc3(s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr);
322f68a32c8SEmil Constantinescu   ierr     = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
323f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
324f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
325f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
326f68a32c8SEmil 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];
327d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
328d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
329f68a32c8SEmil Constantinescu   if (bembed) {
330f68a32c8SEmil Constantinescu     ierr = PetscMalloc(s*sizeof(PetscReal),&t->bembed);CHKERRQ(ierr);
331f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
332f68a32c8SEmil Constantinescu   }
333f68a32c8SEmil Constantinescu 
334f68a32c8SEmil Constantinescu   t->pinterp     = pinterp;
335f68a32c8SEmil Constantinescu   ierr           = PetscMalloc(s*pinterp*sizeof(PetscReal),&t->binterp);CHKERRQ(ierr);
336f68a32c8SEmil Constantinescu   ierr           = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr);
337f68a32c8SEmil Constantinescu   link->next     = RKTableauList;
338f68a32c8SEmil Constantinescu   RKTableauList = link;
339f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
340f68a32c8SEmil Constantinescu }
341f68a32c8SEmil Constantinescu 
342f68a32c8SEmil Constantinescu #undef __FUNCT__
343f68a32c8SEmil Constantinescu #define __FUNCT__ "TSEvaluateStep_RK"
344e4dd521cSBarry Smith /*
345f68a32c8SEmil Constantinescu  The step completion formula is
346e4dd521cSBarry Smith 
347f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
348f68a32c8SEmil Constantinescu 
349f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
350f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
351f68a32c8SEmil Constantinescu  We can write
352f68a32c8SEmil Constantinescu 
353f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
354f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
355f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
356f68a32c8SEmil Constantinescu 
357f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
358e4dd521cSBarry Smith */
359f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
360f68a32c8SEmil Constantinescu {
361f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
362f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
363f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
364f68a32c8SEmil Constantinescu   PetscReal      h;
365f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
366f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
367f68a32c8SEmil Constantinescu 
368f68a32c8SEmil Constantinescu   PetscFunctionBegin;
369f68a32c8SEmil Constantinescu   switch (rk->status) {
370f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
371f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
372f68a32c8SEmil Constantinescu     h = ts->time_step; break;
373f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
374f68a32c8SEmil Constantinescu     h = ts->time_step_prev; break;
375f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
376f68a32c8SEmil Constantinescu   }
377f68a32c8SEmil Constantinescu   if (order == tab->order) {
378f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
379f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
380f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
381f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
382f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
383f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
384f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
385f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
386f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
387f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
388f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
389f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
390f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
391f68a32c8SEmil Constantinescu     } else {                    /* Rollback and re-complete using (be-b) */
392f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
393f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
394f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
395f68a32c8SEmil Constantinescu     }
396f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
397f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
398f68a32c8SEmil Constantinescu   }
399f68a32c8SEmil Constantinescu unavailable:
400f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
401f68a32c8SEmil 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);
402f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
403f68a32c8SEmil Constantinescu }
404f68a32c8SEmil Constantinescu 
405f68a32c8SEmil Constantinescu #undef __FUNCT__
406f68a32c8SEmil Constantinescu #define __FUNCT__ "TSStep_RK"
407f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
408f68a32c8SEmil Constantinescu {
409f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
410f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
411f68a32c8SEmil Constantinescu   const PetscInt   s    = tab->s;
412f68a32c8SEmil Constantinescu   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
413f68a32c8SEmil Constantinescu   PetscScalar     *w    = rk->work;
414f68a32c8SEmil Constantinescu   Vec             *Y    = rk->Y,*YdotRHS = rk->YdotRHS;
415f68a32c8SEmil Constantinescu   TSAdapt          adapt;
416f68a32c8SEmil Constantinescu   PetscInt         i,j,reject,next_scheme;
417f68a32c8SEmil Constantinescu   PetscReal        next_time_step;
418f68a32c8SEmil Constantinescu   PetscReal        t;
419f68a32c8SEmil Constantinescu   PetscBool        accept;
420f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
421f68a32c8SEmil Constantinescu 
422f68a32c8SEmil Constantinescu   PetscFunctionBegin;
423f68a32c8SEmil Constantinescu 
424f68a32c8SEmil Constantinescu   next_time_step = ts->time_step;
425f68a32c8SEmil Constantinescu   t              = ts->ptime;
426f68a32c8SEmil Constantinescu   accept         = PETSC_TRUE;
427f68a32c8SEmil Constantinescu   rk->status     = TS_STEP_INCOMPLETE;
428f68a32c8SEmil Constantinescu 
429f68a32c8SEmil Constantinescu 
430f68a32c8SEmil Constantinescu   for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) {
431f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
432f68a32c8SEmil Constantinescu     ierr = TSPreStep(ts);CHKERRQ(ierr);
433f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
434f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
435f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
436f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
437f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
438f68a32c8SEmil Constantinescu       ierr = TSAdaptCheckStage(adapt,ts,&accept);CHKERRQ(ierr);
439f68a32c8SEmil Constantinescu       if (!accept) goto reject_step;
440f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
441f68a32c8SEmil Constantinescu     }
442f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
443f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
444f68a32c8SEmil Constantinescu 
445f68a32c8SEmil Constantinescu     /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */
446f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
447f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
448f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
449f68a32c8SEmil Constantinescu     ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr);
450f68a32c8SEmil Constantinescu     if (accept) {
451f68a32c8SEmil Constantinescu       /* ignore next_scheme for now */
452f68a32c8SEmil Constantinescu       ts->ptime    += ts->time_step;
453f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
454e3caeda6SBarry Smith       ts->steps++;
455f68a32c8SEmil Constantinescu       rk->status = TS_STEP_COMPLETE;
456f68a32c8SEmil Constantinescu       ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr);
457f68a32c8SEmil Constantinescu       break;
458f68a32c8SEmil Constantinescu     } else {                    /* Roll back the current step */
459f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = -h*b[j];
460f68a32c8SEmil Constantinescu       ierr = VecMAXPY(ts->vec_sol,s,w,rk->YdotRHS);CHKERRQ(ierr);
461f68a32c8SEmil Constantinescu       ts->time_step = next_time_step;
462f68a32c8SEmil Constantinescu       rk->status   = TS_STEP_INCOMPLETE;
463f68a32c8SEmil Constantinescu     }
464f68a32c8SEmil Constantinescu reject_step: continue;
465f68a32c8SEmil Constantinescu   }
466f68a32c8SEmil Constantinescu   if (rk->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
467f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
468e4dd521cSBarry Smith }
469e4dd521cSBarry Smith 
470f68a32c8SEmil Constantinescu #undef __FUNCT__
471f68a32c8SEmil Constantinescu #define __FUNCT__ "TSInterpolate_RK"
472f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
473f68a32c8SEmil Constantinescu {
474f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
475f68a32c8SEmil Constantinescu   PetscInt         s  = rk->tableau->s,pinterp = rk->tableau->pinterp,i,j;
476f68a32c8SEmil Constantinescu   PetscReal        h;
477f68a32c8SEmil Constantinescu   PetscReal        tt,t;
478f68a32c8SEmil Constantinescu   PetscScalar     *b;
479f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
480f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
481e4dd521cSBarry Smith 
482f68a32c8SEmil Constantinescu   PetscFunctionBegin;
483f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
484f68a32c8SEmil Constantinescu   switch (rk->status) {
485f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
486f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
487f68a32c8SEmil Constantinescu     h = ts->time_step;
488f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
489f68a32c8SEmil Constantinescu     break;
490f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
491f68a32c8SEmil Constantinescu     h = ts->time_step_prev;
492f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
493f68a32c8SEmil Constantinescu     break;
494f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
495e4dd521cSBarry Smith   }
496f68a32c8SEmil Constantinescu   ierr = PetscMalloc(s*sizeof(PetscScalar),&b);CHKERRQ(ierr);
497f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
498f68a32c8SEmil Constantinescu   for (j=0,tt=t; j<pinterp; j++,tt*=t) {
499f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
500f68a32c8SEmil Constantinescu       b[i]  += h * B[i*pinterp+j] * tt;
501e4dd521cSBarry Smith     }
502f68a32c8SEmil Constantinescu   }
503f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
504f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
505f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
506e4dd521cSBarry Smith   PetscFunctionReturn(0);
507e4dd521cSBarry Smith }
508e4dd521cSBarry Smith 
509e4dd521cSBarry Smith /*------------------------------------------------------------*/
510e4dd521cSBarry Smith #undef __FUNCT__
511277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_RK"
512277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
513e4dd521cSBarry Smith {
5145f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
515f68a32c8SEmil Constantinescu   PetscInt       s;
5166849ba73SBarry Smith   PetscErrorCode ierr;
517e4dd521cSBarry Smith 
518e4dd521cSBarry Smith   PetscFunctionBegin;
519f68a32c8SEmil Constantinescu   if (!rk->tableau) PetscFunctionReturn(0);
520f68a32c8SEmil Constantinescu   s    = rk->tableau->s;
521f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->Y);CHKERRQ(ierr);
522f68a32c8SEmil Constantinescu   ierr = VecDestroyVecs(s,&rk->YdotRHS);CHKERRQ(ierr);
523f68a32c8SEmil Constantinescu   ierr = PetscFree(rk->work);CHKERRQ(ierr);
524277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
525e4dd521cSBarry Smith }
526277b19d0SLisandro Dalcin 
527277b19d0SLisandro Dalcin #undef __FUNCT__
528277b19d0SLisandro Dalcin #define __FUNCT__ "TSDestroy_RK"
529277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
530277b19d0SLisandro Dalcin {
531277b19d0SLisandro Dalcin   PetscErrorCode ierr;
532277b19d0SLisandro Dalcin 
533277b19d0SLisandro Dalcin   PetscFunctionBegin;
534277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
535277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
536f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
537f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
538f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
539f68a32c8SEmil Constantinescu }
540f68a32c8SEmil Constantinescu 
541f68a32c8SEmil Constantinescu 
542f68a32c8SEmil Constantinescu #undef __FUNCT__
543f68a32c8SEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSRK"
544f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
545f68a32c8SEmil Constantinescu {
546f68a32c8SEmil Constantinescu   PetscFunctionBegin;
547f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
548f68a32c8SEmil Constantinescu }
549f68a32c8SEmil Constantinescu 
550f68a32c8SEmil Constantinescu #undef __FUNCT__
551f68a32c8SEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSRK"
552f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
553f68a32c8SEmil Constantinescu {
554f68a32c8SEmil Constantinescu   PetscFunctionBegin;
555f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
556f68a32c8SEmil Constantinescu }
557f68a32c8SEmil Constantinescu 
558f68a32c8SEmil Constantinescu 
559f68a32c8SEmil Constantinescu #undef __FUNCT__
560f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSRK"
561f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
562f68a32c8SEmil Constantinescu {
563f68a32c8SEmil Constantinescu   PetscFunctionBegin;
564f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
565f68a32c8SEmil Constantinescu }
566f68a32c8SEmil Constantinescu 
567f68a32c8SEmil Constantinescu #undef __FUNCT__
568f68a32c8SEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSRK"
569f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
570f68a32c8SEmil Constantinescu {
571f68a32c8SEmil Constantinescu 
572f68a32c8SEmil Constantinescu   PetscFunctionBegin;
573f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
574f68a32c8SEmil Constantinescu }
575f68a32c8SEmil Constantinescu 
576f68a32c8SEmil Constantinescu #undef __FUNCT__
577f68a32c8SEmil Constantinescu #define __FUNCT__ "TSSetUp_RK"
578f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
579f68a32c8SEmil Constantinescu {
580f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
581f68a32c8SEmil Constantinescu   RKTableau      tab;
582f68a32c8SEmil Constantinescu   PetscInt       s;
583f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
584f68a32c8SEmil Constantinescu   DM             dm;
585f68a32c8SEmil Constantinescu 
586f68a32c8SEmil Constantinescu   PetscFunctionBegin;
587f68a32c8SEmil Constantinescu   if (!rk->tableau) {
588f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
589f68a32c8SEmil Constantinescu   }
590f68a32c8SEmil Constantinescu   tab  = rk->tableau;
591f68a32c8SEmil Constantinescu   s    = tab->s;
592f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->Y);CHKERRQ(ierr);
593f68a32c8SEmil Constantinescu   ierr = VecDuplicateVecs(ts->vec_sol,s,&rk->YdotRHS);CHKERRQ(ierr);
594f68a32c8SEmil Constantinescu   ierr = PetscMalloc(s*sizeof(rk->work[0]),&rk->work);CHKERRQ(ierr);
595f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
596f68a32c8SEmil Constantinescu   if (dm) {
597f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
598f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
599f68a32c8SEmil Constantinescu   }
600e4dd521cSBarry Smith   PetscFunctionReturn(0);
601e4dd521cSBarry Smith }
602e4dd521cSBarry Smith /*------------------------------------------------------------*/
603e4dd521cSBarry Smith 
604e4dd521cSBarry Smith #undef __FUNCT__
6055f70b478SJed Brown #define __FUNCT__ "TSSetFromOptions_RK"
6065f70b478SJed Brown static PetscErrorCode TSSetFromOptions_RK(TS ts)
607e4dd521cSBarry Smith {
608dfbe8321SBarry Smith   PetscErrorCode ierr;
609f68a32c8SEmil Constantinescu   char           rktype[256];
610262ff7bbSBarry Smith 
611e4dd521cSBarry Smith   PetscFunctionBegin;
612262ff7bbSBarry Smith   ierr = PetscOptionsHead("RK ODE solver options");CHKERRQ(ierr);
613f68a32c8SEmil Constantinescu   {
614f68a32c8SEmil Constantinescu     RKTableauLink  link;
615f68a32c8SEmil Constantinescu     PetscInt       count,choice;
616f68a32c8SEmil Constantinescu     PetscBool      flg;
617f68a32c8SEmil Constantinescu     const char   **namelist;
618f68a32c8SEmil Constantinescu     ierr = PetscStrncpy(rktype,TSRKDefault,sizeof(rktype));CHKERRQ(ierr);
619f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
620f68a32c8SEmil Constantinescu     ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr);
621f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
622f68a32c8SEmil Constantinescu     ierr      = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rktype,&choice,&flg);CHKERRQ(ierr);
623f68a32c8SEmil Constantinescu     ierr      = TSRKSetType(ts,flg ? namelist[choice] : rktype);CHKERRQ(ierr);
624f68a32c8SEmil Constantinescu     ierr      = PetscFree(namelist);CHKERRQ(ierr);
625f68a32c8SEmil Constantinescu   }
626262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
627e4dd521cSBarry Smith   PetscFunctionReturn(0);
628e4dd521cSBarry Smith }
629e4dd521cSBarry Smith 
630e4dd521cSBarry Smith #undef __FUNCT__
631f68a32c8SEmil Constantinescu #define __FUNCT__ "PetscFormatRealArray"
632f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
633f68a32c8SEmil Constantinescu {
634f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
635f68a32c8SEmil Constantinescu   PetscInt       i;
636f68a32c8SEmil Constantinescu   size_t         left,count;
637f68a32c8SEmil Constantinescu   char           *p;
638f68a32c8SEmil Constantinescu 
639f68a32c8SEmil Constantinescu   PetscFunctionBegin;
640f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
641f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
642f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
643f68a32c8SEmil Constantinescu     left -= count;
644f68a32c8SEmil Constantinescu     p    += count;
645f68a32c8SEmil Constantinescu     *p++  = ' ';
646f68a32c8SEmil Constantinescu   }
647f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
648f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
649f68a32c8SEmil Constantinescu }
650f68a32c8SEmil Constantinescu 
651f68a32c8SEmil Constantinescu #undef __FUNCT__
6525f70b478SJed Brown #define __FUNCT__ "TSView_RK"
6535f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
654e4dd521cSBarry Smith {
6555f70b478SJed Brown   TS_RK         *rk   = (TS_RK*)ts->data;
656f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
6578370ee3bSLisandro Dalcin   PetscBool      iascii;
658dfbe8321SBarry Smith   PetscErrorCode ierr;
659f68a32c8SEmil Constantinescu   TSAdapt        adapt;
6602cdf8120Sasbjorn 
661e4dd521cSBarry Smith   PetscFunctionBegin;
662251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6638370ee3bSLisandro Dalcin   if (iascii) {
664f68a32c8SEmil Constantinescu     TSRKType rktype;
665f68a32c8SEmil Constantinescu     char     buf[512];
666f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
667f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  RK %s\n",rktype);CHKERRQ(ierr);
668f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
669f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa     c = %s\n",buf);CHKERRQ(ierr);
670d760c35bSDebojyoti Ghosh     ierr = PetscViewerASCIIPrintf(viewer,"FSAL: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
6718370ee3bSLisandro Dalcin   }
672f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
673f68a32c8SEmil Constantinescu   ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr);
674f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
675f68a32c8SEmil Constantinescu }
676f68a32c8SEmil Constantinescu 
677f68a32c8SEmil Constantinescu #undef __FUNCT__
678f68a32c8SEmil Constantinescu #define __FUNCT__ "TSLoad_RK"
679f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
680f68a32c8SEmil Constantinescu {
681f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
682f68a32c8SEmil Constantinescu   TSAdapt        tsadapt;
683f68a32c8SEmil Constantinescu 
684f68a32c8SEmil Constantinescu   PetscFunctionBegin;
685f68a32c8SEmil Constantinescu   ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr);
686f68a32c8SEmil Constantinescu   ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr);
687f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
688f68a32c8SEmil Constantinescu }
689f68a32c8SEmil Constantinescu 
690f68a32c8SEmil Constantinescu #undef __FUNCT__
691f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType"
692f68a32c8SEmil Constantinescu /*@C
693f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
694f68a32c8SEmil Constantinescu 
695f68a32c8SEmil Constantinescu   Logically collective
696f68a32c8SEmil Constantinescu 
697f68a32c8SEmil Constantinescu   Input Parameter:
698f68a32c8SEmil Constantinescu +  ts - timestepping context
699f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
700f68a32c8SEmil Constantinescu 
701f68a32c8SEmil Constantinescu   Level: intermediate
702f68a32c8SEmil Constantinescu 
703f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
704f68a32c8SEmil Constantinescu @*/
705f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
706f68a32c8SEmil Constantinescu {
707f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
708f68a32c8SEmil Constantinescu 
709f68a32c8SEmil Constantinescu   PetscFunctionBegin;
710f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
711f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
712f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
713f68a32c8SEmil Constantinescu }
714f68a32c8SEmil Constantinescu 
715f68a32c8SEmil Constantinescu #undef __FUNCT__
716f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType"
717f68a32c8SEmil Constantinescu /*@C
718f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
719f68a32c8SEmil Constantinescu 
720f68a32c8SEmil Constantinescu   Logically collective
721f68a32c8SEmil Constantinescu 
722f68a32c8SEmil Constantinescu   Input Parameter:
723f68a32c8SEmil Constantinescu .  ts - timestepping context
724f68a32c8SEmil Constantinescu 
725f68a32c8SEmil Constantinescu   Output Parameter:
726f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
727f68a32c8SEmil Constantinescu 
728f68a32c8SEmil Constantinescu   Level: intermediate
729f68a32c8SEmil Constantinescu 
730f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
731f68a32c8SEmil Constantinescu @*/
732f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
733f68a32c8SEmil Constantinescu {
734f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
735f68a32c8SEmil Constantinescu 
736f68a32c8SEmil Constantinescu   PetscFunctionBegin;
737f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
738f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
739f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
740f68a32c8SEmil Constantinescu }
741f68a32c8SEmil Constantinescu 
742f68a32c8SEmil Constantinescu #undef __FUNCT__
743f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKGetType_RK"
744f68a32c8SEmil Constantinescu PetscErrorCode  TSRKGetType_RK(TS ts,TSRKType *rktype)
745f68a32c8SEmil Constantinescu {
746f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
747f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
748f68a32c8SEmil Constantinescu 
749f68a32c8SEmil Constantinescu   PetscFunctionBegin;
750f68a32c8SEmil Constantinescu   if (!rk->tableau) {
751f68a32c8SEmil Constantinescu     ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
752f68a32c8SEmil Constantinescu   }
753f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
754f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
755f68a32c8SEmil Constantinescu }
756f68a32c8SEmil Constantinescu #undef __FUNCT__
757f68a32c8SEmil Constantinescu #define __FUNCT__ "TSRKSetType_RK"
758f68a32c8SEmil Constantinescu PetscErrorCode  TSRKSetType_RK(TS ts,TSRKType rktype)
759f68a32c8SEmil Constantinescu {
760f68a32c8SEmil Constantinescu   TS_RK     *rk = (TS_RK*)ts->data;
761f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
762f68a32c8SEmil Constantinescu   PetscBool      match;
763f68a32c8SEmil Constantinescu   RKTableauLink link;
764f68a32c8SEmil Constantinescu 
765f68a32c8SEmil Constantinescu   PetscFunctionBegin;
766f68a32c8SEmil Constantinescu   if (rk->tableau) {
767f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
768f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
769f68a32c8SEmil Constantinescu   }
770f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
771f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
772f68a32c8SEmil Constantinescu     if (match) {
773f68a32c8SEmil Constantinescu       ierr = TSReset_RK(ts);CHKERRQ(ierr);
774f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
775f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
776f68a32c8SEmil Constantinescu     }
777f68a32c8SEmil Constantinescu   }
778f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
779e4dd521cSBarry Smith   PetscFunctionReturn(0);
780e4dd521cSBarry Smith }
781e4dd521cSBarry Smith 
782e4dd521cSBarry Smith /* ------------------------------------------------------------ */
783ebe3b25bSBarry Smith /*MC
784f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
785ebe3b25bSBarry Smith 
786f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
787f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
788ebe3b25bSBarry Smith 
789f68a32c8SEmil Constantinescu   Notes:
790f68a32c8SEmil Constantinescu   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
791ebe3b25bSBarry Smith 
792d41222bbSBarry Smith   Level: beginner
793d41222bbSBarry Smith 
794f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
795f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
796ebe3b25bSBarry Smith 
797ebe3b25bSBarry Smith M*/
798e4dd521cSBarry Smith #undef __FUNCT__
7995f70b478SJed Brown #define __FUNCT__ "TSCreate_RK"
8008cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
801e4dd521cSBarry Smith {
802f68a32c8SEmil Constantinescu   TS_RK     *th;
803dfbe8321SBarry Smith   PetscErrorCode ierr;
804e4dd521cSBarry Smith 
805e4dd521cSBarry Smith   PetscFunctionBegin;
806f68a32c8SEmil Constantinescu #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
807f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
808f68a32c8SEmil Constantinescu #endif
809f68a32c8SEmil Constantinescu 
810f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
8115f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
8125f70b478SJed Brown   ts->ops->view           = TSView_RK;
813f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
814f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
815f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
816f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
817f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
818f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
819e4dd521cSBarry Smith 
820f68a32c8SEmil Constantinescu   ierr = PetscNewLog(ts,TS_RK,&th);CHKERRQ(ierr);
821f68a32c8SEmil Constantinescu   ts->data = (void*)th;
822e4dd521cSBarry Smith 
823f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
824f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
825e4dd521cSBarry Smith   PetscFunctionReturn(0);
826e4dd521cSBarry Smith }
827