xref: /petsc/src/ts/impls/explicit/rk/rk.c (revision d190556429ad48cdb35b176ba6712769ba23bffc)
1e4dd521cSBarry Smith /*
22e698f8bSDebojyoti Ghosh   Code for time stepping with the Runge-Kutta method
3f68a32c8SEmil Constantinescu 
4f68a32c8SEmil Constantinescu   Notes:
5f68a32c8SEmil Constantinescu   The general system is written as
6f68a32c8SEmil Constantinescu 
72e698f8bSDebojyoti Ghosh   Udot = F(t,U)
8f68a32c8SEmil Constantinescu 
9e4dd521cSBarry Smith */
10af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
11f68a32c8SEmil Constantinescu #include <petscdm.h>
12f68a32c8SEmil Constantinescu 
13484bcdc7SDebojyoti Ghosh static TSRKType  TSRKDefault = TSRK3BS;
14f68a32c8SEmil Constantinescu static PetscBool TSRKRegisterAllCalled;
15f68a32c8SEmil Constantinescu static PetscBool TSRKPackageInitialized;
16f68a32c8SEmil Constantinescu 
17f68a32c8SEmil Constantinescu typedef struct _RKTableau *RKTableau;
18f68a32c8SEmil Constantinescu struct _RKTableau {
19f68a32c8SEmil Constantinescu   char      *name;
20d760c35bSDebojyoti Ghosh   PetscInt   order;               /* Classical approximation order of the method i              */
21f68a32c8SEmil Constantinescu   PetscInt   s;                   /* Number of stages                                           */
223a8a9803SLisandro Dalcin   PetscInt   p;                   /* Interpolation order                                        */
23d760c35bSDebojyoti Ghosh   PetscBool  FSAL;                /* flag to indicate if tableau is FSAL                        */
24f68a32c8SEmil Constantinescu   PetscReal *A,*b,*c;             /* Tableau                                                    */
25f68a32c8SEmil Constantinescu   PetscReal *bembed;              /* Embedded formula of order one less (order-1)               */
26f68a32c8SEmil Constantinescu   PetscReal *binterp;             /* Dense output formula                                       */
27f68a32c8SEmil Constantinescu   PetscReal  ccfl;                /* Placeholder for CFL coefficient relative to forward Euler  */
28f68a32c8SEmil Constantinescu };
29f68a32c8SEmil Constantinescu typedef struct _RKTableauLink *RKTableauLink;
30f68a32c8SEmil Constantinescu struct _RKTableauLink {
31f68a32c8SEmil Constantinescu   struct _RKTableau tab;
32f68a32c8SEmil Constantinescu   RKTableauLink     next;
33f68a32c8SEmil Constantinescu };
34f68a32c8SEmil Constantinescu static RKTableauLink RKTableauList;
35e4dd521cSBarry Smith 
36e4dd521cSBarry Smith typedef struct {
37f68a32c8SEmil Constantinescu   RKTableau    tableau;
38f68a32c8SEmil Constantinescu   Vec          *Y;               /* States computed during the step */
39f68a32c8SEmil Constantinescu   Vec          *YdotRHS;         /* Function evaluations for the non-stiff part */
40ad8e2604SHong Zhang   Vec          *VecDeltaLam;     /* Increment of the adjoint sensitivity w.r.t IC at stage */
41ad8e2604SHong Zhang   Vec          *VecDeltaMu;      /* Increment of the adjoint sensitivity w.r.t P at stage */
42b18ea86cSHong Zhang   Vec          *VecSensiTemp;    /* Vector to be timed with Jacobian transpose */
430f7a1166SHong Zhang   Vec          VecCostIntegral0; /* backup for roll-backs due to events */
44f68a32c8SEmil Constantinescu   PetscScalar  *work;            /* Scalar work */
45f68a32c8SEmil Constantinescu   PetscReal    stage_time;
46f68a32c8SEmil Constantinescu   TSStepStatus status;
470f7a1166SHong Zhang   PetscReal    ptime;
480f7a1166SHong Zhang   PetscReal    time_step;
495f70b478SJed Brown } TS_RK;
50e4dd521cSBarry Smith 
51f68a32c8SEmil Constantinescu /*MC
52f68a32c8SEmil Constantinescu      TSRK1 - First order forward Euler scheme.
53262ff7bbSBarry Smith 
54f68a32c8SEmil Constantinescu      This method has one stage.
55f68a32c8SEmil Constantinescu 
56f68a32c8SEmil Constantinescu      Level: advanced
57f68a32c8SEmil Constantinescu 
58f68a32c8SEmil Constantinescu .seealso: TSRK
59f68a32c8SEmil Constantinescu M*/
60f68a32c8SEmil Constantinescu /*MC
612109b73fSDebojyoti Ghosh      TSRK2A - Second order RK scheme.
62f68a32c8SEmil Constantinescu 
63f68a32c8SEmil Constantinescu      This method has two stages.
64f68a32c8SEmil Constantinescu 
65f68a32c8SEmil Constantinescu      Level: advanced
66f68a32c8SEmil Constantinescu 
67f68a32c8SEmil Constantinescu .seealso: TSRK
68f68a32c8SEmil Constantinescu M*/
69f68a32c8SEmil Constantinescu /*MC
70f68a32c8SEmil Constantinescu      TSRK3 - Third order RK scheme.
71f68a32c8SEmil Constantinescu 
72f68a32c8SEmil Constantinescu      This method has three stages.
73f68a32c8SEmil Constantinescu 
74f68a32c8SEmil Constantinescu      Level: advanced
75f68a32c8SEmil Constantinescu 
76f68a32c8SEmil Constantinescu .seealso: TSRK
77f68a32c8SEmil Constantinescu M*/
78f68a32c8SEmil Constantinescu /*MC
792109b73fSDebojyoti Ghosh      TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method.
802109b73fSDebojyoti Ghosh 
81*d1905564SLisandro Dalcin      This method has four stages with the First Same As Last (FSAL) property.
822109b73fSDebojyoti Ghosh 
832109b73fSDebojyoti Ghosh      Level: advanced
842109b73fSDebojyoti Ghosh 
85*d1905564SLisandro Dalcin      Notes: https://doi.org/10.1016/0893-9659(89)90079-7
86*d1905564SLisandro Dalcin 
872109b73fSDebojyoti Ghosh .seealso: TSRK
882109b73fSDebojyoti Ghosh M*/
892109b73fSDebojyoti Ghosh /*MC
90f68a32c8SEmil Constantinescu      TSRK4 - Fourth order RK scheme.
91f68a32c8SEmil Constantinescu 
922e698f8bSDebojyoti Ghosh      This is the classical Runge-Kutta method with four stages.
93f68a32c8SEmil Constantinescu 
94f68a32c8SEmil Constantinescu      Level: advanced
95f68a32c8SEmil Constantinescu 
96f68a32c8SEmil Constantinescu .seealso: TSRK
97f68a32c8SEmil Constantinescu M*/
98f68a32c8SEmil Constantinescu /*MC
992e698f8bSDebojyoti Ghosh      TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method.
100f68a32c8SEmil Constantinescu 
101f68a32c8SEmil Constantinescu      This method has six stages.
102f68a32c8SEmil Constantinescu 
103f68a32c8SEmil Constantinescu      Level: advanced
104f68a32c8SEmil Constantinescu 
105f68a32c8SEmil Constantinescu .seealso: TSRK
106f68a32c8SEmil Constantinescu M*/
1072109b73fSDebojyoti Ghosh /*MC
1082e698f8bSDebojyoti Ghosh      TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method.
1092109b73fSDebojyoti Ghosh 
110*d1905564SLisandro Dalcin      This method has seven stages with the First Same As Last (FSAL) property.
1112109b73fSDebojyoti Ghosh 
1122109b73fSDebojyoti Ghosh      Level: advanced
1132109b73fSDebojyoti Ghosh 
114*d1905564SLisandro Dalcin      Notes: https://doi.org/10.1016/0771-050X(80)90013-3
115*d1905564SLisandro Dalcin 
1162109b73fSDebojyoti Ghosh .seealso: TSRK
1172109b73fSDebojyoti Ghosh M*/
118262ff7bbSBarry Smith 
119f68a32c8SEmil Constantinescu /*@C
120f68a32c8SEmil Constantinescu   TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK
121262ff7bbSBarry Smith 
122f68a32c8SEmil Constantinescu   Not Collective, but should be called by all processes which will need the schemes to be registered
123262ff7bbSBarry Smith 
124f68a32c8SEmil Constantinescu   Level: advanced
125262ff7bbSBarry Smith 
126f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all
127262ff7bbSBarry Smith 
128f68a32c8SEmil Constantinescu .seealso:  TSRKRegisterDestroy()
129262ff7bbSBarry Smith @*/
130f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void)
131262ff7bbSBarry Smith {
1324ac538c5SBarry Smith   PetscErrorCode ierr;
133262ff7bbSBarry Smith 
134262ff7bbSBarry Smith   PetscFunctionBegin;
135f68a32c8SEmil Constantinescu   if (TSRKRegisterAllCalled) PetscFunctionReturn(0);
136f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_TRUE;
137e4dd521cSBarry Smith 
138e4dd521cSBarry Smith   {
139f68a32c8SEmil Constantinescu     const PetscReal
140f68a32c8SEmil Constantinescu       A[1][1] = {{0.0}},
141f68a32c8SEmil Constantinescu       b[1]    = {1.0};
1423a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
143e4dd521cSBarry Smith   }
144db046809SBarry Smith   {
145f68a32c8SEmil Constantinescu     const PetscReal
146f68a32c8SEmil Constantinescu       A[2][2]     = {{0.0,0.0},
147f68a32c8SEmil Constantinescu                     {1.0,0.0}},
148f68a32c8SEmil Constantinescu       b[2]        = {0.5,0.5},
149f68a32c8SEmil Constantinescu       bembed[2]   = {1.0,0};
1503a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
151db046809SBarry Smith   }
152f68a32c8SEmil Constantinescu   {
153f68a32c8SEmil Constantinescu     const PetscReal
154f68a32c8SEmil Constantinescu       A[3][3] = {{0,0,0},
155f68a32c8SEmil Constantinescu                  {2.0/3.0,0,0},
156f68a32c8SEmil Constantinescu                  {-1.0/3.0,1.0,0}},
157f68a32c8SEmil Constantinescu       b[3]    = {0.25,0.5,0.25};
1583a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
159fdefd5e5SDebojyoti Ghosh   }
160fdefd5e5SDebojyoti Ghosh   {
161fdefd5e5SDebojyoti Ghosh     const PetscReal
162fdefd5e5SDebojyoti Ghosh       A[4][4] = {{0,0,0,0},
163fdefd5e5SDebojyoti Ghosh                  {1.0/2.0,0,0,0},
164fdefd5e5SDebojyoti Ghosh                  {0,3.0/4.0,0,0},
165fdefd5e5SDebojyoti Ghosh                  {2.0/9.0,1.0/3.0,4.0/9.0,0}},
166fdefd5e5SDebojyoti Ghosh       b[4]    = {2.0/9.0,1.0/3.0,4.0/9.0,0},
167fdefd5e5SDebojyoti Ghosh       bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0};
1683a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
169db046809SBarry Smith   }
170f68a32c8SEmil Constantinescu   {
171f68a32c8SEmil Constantinescu     const PetscReal
172f68a32c8SEmil Constantinescu       A[4][4] = {{0,0,0,0},
173f68a32c8SEmil Constantinescu                  {0.5,0,0,0},
174f68a32c8SEmil Constantinescu                  {0,0.5,0,0},
175f68a32c8SEmil Constantinescu                  {0,0,1.0,0}},
176f68a32c8SEmil Constantinescu       b[4]    = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0};
1773a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr);
178db046809SBarry Smith   }
179f68a32c8SEmil Constantinescu   {
180f68a32c8SEmil Constantinescu     const PetscReal
181f68a32c8SEmil Constantinescu       A[6][6]   = {{0,0,0,0,0,0},
182f68a32c8SEmil Constantinescu                    {0.25,0,0,0,0,0},
183f68a32c8SEmil Constantinescu                    {3.0/32.0,9.0/32.0,0,0,0,0},
184f68a32c8SEmil Constantinescu                    {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0},
185f68a32c8SEmil Constantinescu                    {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0},
186f68a32c8SEmil Constantinescu                    {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}},
187f68a32c8SEmil Constantinescu       b[6]      = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0},
188f68a32c8SEmil Constantinescu       bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0};
1893a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
190fdefd5e5SDebojyoti Ghosh   }
191fdefd5e5SDebojyoti Ghosh   {
192fdefd5e5SDebojyoti Ghosh     const PetscReal
193fdefd5e5SDebojyoti Ghosh       A[7][7]   = {{0,0,0,0,0,0,0},
194fdefd5e5SDebojyoti Ghosh                    {1.0/5.0,0,0,0,0,0,0},
195fdefd5e5SDebojyoti Ghosh                    {3.0/40.0,9.0/40.0,0,0,0,0,0},
196fdefd5e5SDebojyoti Ghosh                    {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0},
197fdefd5e5SDebojyoti Ghosh                    {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0},
198fdefd5e5SDebojyoti Ghosh                    {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0},
199fdefd5e5SDebojyoti Ghosh                    {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}},
200fdefd5e5SDebojyoti 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},
201fdefd5e5SDebojyoti 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};
2023a8a9803SLisandro Dalcin     ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr);
203f68a32c8SEmil Constantinescu   }
204db046809SBarry Smith   PetscFunctionReturn(0);
205db046809SBarry Smith }
206db046809SBarry Smith 
207f68a32c8SEmil Constantinescu /*@C
208f68a32c8SEmil Constantinescu    TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister().
209f68a32c8SEmil Constantinescu 
210f68a32c8SEmil Constantinescu    Not Collective
211f68a32c8SEmil Constantinescu 
212f68a32c8SEmil Constantinescu    Level: advanced
213f68a32c8SEmil Constantinescu 
214f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy
215f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll()
216f68a32c8SEmil Constantinescu @*/
217f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void)
218e4dd521cSBarry Smith {
219f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
220f68a32c8SEmil Constantinescu   RKTableauLink link;
221f68a32c8SEmil Constantinescu 
222f68a32c8SEmil Constantinescu   PetscFunctionBegin;
223f68a32c8SEmil Constantinescu   while ((link = RKTableauList)) {
224f68a32c8SEmil Constantinescu     RKTableau t = &link->tab;
225f68a32c8SEmil Constantinescu     RKTableauList = link->next;
226f68a32c8SEmil Constantinescu     ierr = PetscFree3(t->A,t->b,t->c);  CHKERRQ(ierr);
227f68a32c8SEmil Constantinescu     ierr = PetscFree (t->bembed);       CHKERRQ(ierr);
228f68a32c8SEmil Constantinescu     ierr = PetscFree (t->binterp);      CHKERRQ(ierr);
229f68a32c8SEmil Constantinescu     ierr = PetscFree (t->name);         CHKERRQ(ierr);
230f68a32c8SEmil Constantinescu     ierr = PetscFree (link);            CHKERRQ(ierr);
231f68a32c8SEmil Constantinescu   }
232f68a32c8SEmil Constantinescu   TSRKRegisterAllCalled = PETSC_FALSE;
233f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
234f68a32c8SEmil Constantinescu }
235f68a32c8SEmil Constantinescu 
236f68a32c8SEmil Constantinescu /*@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 = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr);
255f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
256f68a32c8SEmil Constantinescu }
257186e87acSLisandro Dalcin 
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 /*@C
278f68a32c8SEmil Constantinescu    TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
279f68a32c8SEmil Constantinescu 
280f68a32c8SEmil Constantinescu    Not Collective, but the same schemes should be registered on all processes on which they will be used
281f68a32c8SEmil Constantinescu 
282f68a32c8SEmil Constantinescu    Input Parameters:
283f68a32c8SEmil Constantinescu +  name - identifier for method
284f68a32c8SEmil Constantinescu .  order - approximation order of method
285f68a32c8SEmil Constantinescu .  s - number of stages, this is the dimension of the matrices below
286f68a32c8SEmil Constantinescu .  A - stage coefficients (dimension s*s, row-major)
287f68a32c8SEmil Constantinescu .  b - step completion table (dimension s; NULL to use last row of A)
288f68a32c8SEmil Constantinescu .  c - abscissa (dimension s; NULL to use row sums of A)
289f68a32c8SEmil Constantinescu .  bembed - completion table for embedded method (dimension s; NULL if not available)
2903a8a9803SLisandro Dalcin .  p - Order of the interpolation scheme, equal to the number of columns of binterp
2913a8a9803SLisandro Dalcin -  binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1)
292f68a32c8SEmil Constantinescu 
293f68a32c8SEmil Constantinescu    Notes:
294f68a32c8SEmil Constantinescu    Several RK methods are provided, this function is only needed to create new methods.
295f68a32c8SEmil Constantinescu 
296f68a32c8SEmil Constantinescu    Level: advanced
297f68a32c8SEmil Constantinescu 
298f68a32c8SEmil Constantinescu .keywords: TS, register
299f68a32c8SEmil Constantinescu 
300f68a32c8SEmil Constantinescu .seealso: TSRK
301f68a32c8SEmil Constantinescu @*/
302f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s,
303f68a32c8SEmil Constantinescu                             const PetscReal A[],const PetscReal b[],const PetscReal c[],
3043a8a9803SLisandro Dalcin                             const PetscReal bembed[],PetscInt p,const PetscReal binterp[])
305f68a32c8SEmil Constantinescu {
306f68a32c8SEmil Constantinescu   PetscErrorCode  ierr;
307f68a32c8SEmil Constantinescu   RKTableauLink   link;
308f68a32c8SEmil Constantinescu   RKTableau       t;
309f68a32c8SEmil Constantinescu   PetscInt        i,j;
310f68a32c8SEmil Constantinescu 
311f68a32c8SEmil Constantinescu   PetscFunctionBegin;
3123a8a9803SLisandro Dalcin   PetscValidCharPointer(name,1);
3133a8a9803SLisandro Dalcin   PetscValidRealPointer(A,4);
3143a8a9803SLisandro Dalcin   if (b) PetscValidRealPointer(b,5);
3153a8a9803SLisandro Dalcin   if (c) PetscValidRealPointer(c,6);
3163a8a9803SLisandro Dalcin   if (bembed) PetscValidRealPointer(bembed,7);
3173a8a9803SLisandro Dalcin   if (binterp || p > 1) PetscValidRealPointer(binterp,9);
3183a8a9803SLisandro Dalcin 
31995dccacaSBarry Smith   ierr = PetscNew(&link);CHKERRQ(ierr);
320f68a32c8SEmil Constantinescu   t = &link->tab;
3213a8a9803SLisandro Dalcin 
322f68a32c8SEmil Constantinescu   ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr);
323f68a32c8SEmil Constantinescu   t->order = order;
324f68a32c8SEmil Constantinescu   t->s = s;
325dcca6d9dSJed Brown   ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr);
326f68a32c8SEmil Constantinescu   ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
327f68a32c8SEmil Constantinescu   if (b)  { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); }
328f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i];
329f68a32c8SEmil Constantinescu   if (c)  { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); }
330f68a32c8SEmil Constantinescu   else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j];
331d760c35bSDebojyoti Ghosh   t->FSAL = PETSC_TRUE;
332d760c35bSDebojyoti Ghosh   for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE;
3333a8a9803SLisandro Dalcin 
334f68a32c8SEmil Constantinescu   if (bembed) {
335785e854fSJed Brown     ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr);
336f68a32c8SEmil Constantinescu     ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr);
337f68a32c8SEmil Constantinescu   }
338f68a32c8SEmil Constantinescu 
3393a8a9803SLisandro Dalcin   if (!binterp) { p = 1; binterp = t->b; }
3403a8a9803SLisandro Dalcin   t->p = p;
3413a8a9803SLisandro Dalcin   ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr);
3423a8a9803SLisandro Dalcin   ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr);
3433a8a9803SLisandro Dalcin 
344f68a32c8SEmil Constantinescu   link->next = RKTableauList;
345f68a32c8SEmil Constantinescu   RKTableauList = link;
346f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
347f68a32c8SEmil Constantinescu }
348f68a32c8SEmil Constantinescu 
349e4dd521cSBarry Smith /*
350f68a32c8SEmil Constantinescu  The step completion formula is
351e4dd521cSBarry Smith 
352f68a32c8SEmil Constantinescu  x1 = x0 + h b^T YdotRHS
353f68a32c8SEmil Constantinescu 
354f68a32c8SEmil Constantinescu  This function can be called before or after ts->vec_sol has been updated.
355f68a32c8SEmil Constantinescu  Suppose we have a completion formula (b) and an embedded formula (be) of different order.
356f68a32c8SEmil Constantinescu  We can write
357f68a32c8SEmil Constantinescu 
358f68a32c8SEmil Constantinescu  x1e = x0 + h be^T YdotRHS
359f68a32c8SEmil Constantinescu      = x1 - h b^T YdotRHS + h be^T YdotRHS
360f68a32c8SEmil Constantinescu      = x1 + h (be - b)^T YdotRHS
361f68a32c8SEmil Constantinescu 
362f68a32c8SEmil Constantinescu  so we can evaluate the method with different order even after the step has been optimistically completed.
363e4dd521cSBarry Smith */
364f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done)
365f68a32c8SEmil Constantinescu {
366f68a32c8SEmil Constantinescu   TS_RK         *rk   = (TS_RK*)ts->data;
367f68a32c8SEmil Constantinescu   RKTableau      tab  = rk->tableau;
368f68a32c8SEmil Constantinescu   PetscScalar   *w    = rk->work;
369f68a32c8SEmil Constantinescu   PetscReal      h;
370f68a32c8SEmil Constantinescu   PetscInt       s    = tab->s,j;
371f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
372f68a32c8SEmil Constantinescu 
373f68a32c8SEmil Constantinescu   PetscFunctionBegin;
374f68a32c8SEmil Constantinescu   switch (rk->status) {
375f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
376f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
377f68a32c8SEmil Constantinescu     h = ts->time_step; break;
378f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
379be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
380f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
381f68a32c8SEmil Constantinescu   }
382f68a32c8SEmil Constantinescu   if (order == tab->order) {
383f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) {
384f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
385f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->b[j];
386f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
387f68a32c8SEmil Constantinescu     } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);}
388f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
389f68a32c8SEmil Constantinescu   } else if (order == tab->order-1) {
390f68a32c8SEmil Constantinescu     if (!tab->bembed) goto unavailable;
391f68a32c8SEmil Constantinescu     if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */
392f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
393f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*tab->bembed[j];
394f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
395f68a32c8SEmil Constantinescu     } else { /* Rollback and re-complete using (be-b) */
396f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
397f68a32c8SEmil Constantinescu       for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]);
398f68a32c8SEmil Constantinescu       ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr);
3990f7a1166SHong Zhang       if (ts->vec_costintegral && ts->costintegralfwd) {
4000f7a1166SHong Zhang         ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr);
4010f7a1166SHong Zhang       }
402f68a32c8SEmil Constantinescu     }
403f68a32c8SEmil Constantinescu     if (done) *done = PETSC_TRUE;
404f68a32c8SEmil Constantinescu     PetscFunctionReturn(0);
405f68a32c8SEmil Constantinescu   }
406f68a32c8SEmil Constantinescu unavailable:
407f68a32c8SEmil Constantinescu   if (done) *done = PETSC_FALSE;
408a7fac7c2SEmil Constantinescu   else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"RK '%s' of order %D cannot evaluate step at order %D. Consider using -ts_adapt_type none or a different method that has an embedded estimate.",tab->name,tab->order,order);
409f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
410f68a32c8SEmil Constantinescu }
411f68a32c8SEmil Constantinescu 
4120f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts)
4130f7a1166SHong Zhang {
4140f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4150f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4160f7a1166SHong Zhang   const PetscInt  s = tab->s;
4170f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4180f7a1166SHong Zhang   Vec             *Y = rk->Y;
4190f7a1166SHong Zhang   PetscInt        i;
4200f7a1166SHong Zhang   PetscErrorCode  ierr;
4210f7a1166SHong Zhang 
4220f7a1166SHong Zhang   PetscFunctionBegin;
4230f7a1166SHong Zhang   /* backup cost integral */
4240f7a1166SHong Zhang   ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr);
4250f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4260f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4270f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4280f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4290f7a1166SHong Zhang   }
4300f7a1166SHong Zhang   PetscFunctionReturn(0);
4310f7a1166SHong Zhang }
4320f7a1166SHong Zhang 
4330f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts)
4340f7a1166SHong Zhang {
4350f7a1166SHong Zhang   TS_RK           *rk = (TS_RK*)ts->data;
4360f7a1166SHong Zhang   RKTableau       tab = rk->tableau;
4370f7a1166SHong Zhang   const PetscInt  s = tab->s;
4380f7a1166SHong Zhang   const PetscReal *b = tab->b,*c = tab->c;
4390f7a1166SHong Zhang   Vec             *Y = rk->Y;
4400f7a1166SHong Zhang   PetscInt        i;
4410f7a1166SHong Zhang   PetscErrorCode  ierr;
4420f7a1166SHong Zhang 
4430f7a1166SHong Zhang   PetscFunctionBegin;
4440f7a1166SHong Zhang   for (i=s-1; i>=0; i--) {
4450f7a1166SHong Zhang     /* Evolve ts->vec_costintegral to compute integrals */
4460f7a1166SHong Zhang     ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr);
4470f7a1166SHong Zhang     ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr);
4480f7a1166SHong Zhang   }
4490f7a1166SHong Zhang   PetscFunctionReturn(0);
4500f7a1166SHong Zhang }
4510f7a1166SHong Zhang 
452fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts)
453fecfb714SLisandro Dalcin {
454fecfb714SLisandro Dalcin   TS_RK           *rk = (TS_RK*)ts->data;
455fecfb714SLisandro Dalcin   RKTableau       tab = rk->tableau;
456fecfb714SLisandro Dalcin   const PetscInt  s  = tab->s;
457fecfb714SLisandro Dalcin   const PetscReal *b = tab->b;
458fecfb714SLisandro Dalcin   PetscScalar     *w = rk->work;
459fecfb714SLisandro Dalcin   Vec             *YdotRHS = rk->YdotRHS;
460fecfb714SLisandro Dalcin   PetscInt        j;
461fecfb714SLisandro Dalcin   PetscReal       h;
462fecfb714SLisandro Dalcin   PetscErrorCode  ierr;
463fecfb714SLisandro Dalcin 
464fecfb714SLisandro Dalcin   PetscFunctionBegin;
465fecfb714SLisandro Dalcin   switch (rk->status) {
466fecfb714SLisandro Dalcin   case TS_STEP_INCOMPLETE:
467fecfb714SLisandro Dalcin   case TS_STEP_PENDING:
468fecfb714SLisandro Dalcin     h = ts->time_step; break;
469fecfb714SLisandro Dalcin   case TS_STEP_COMPLETE:
470fecfb714SLisandro Dalcin     h = ts->ptime - ts->ptime_prev; break;
471fecfb714SLisandro Dalcin   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
472fecfb714SLisandro Dalcin   }
473fecfb714SLisandro Dalcin   for (j=0; j<s; j++) w[j] = -h*b[j];
474fecfb714SLisandro Dalcin   ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr);
475fecfb714SLisandro Dalcin   PetscFunctionReturn(0);
476fecfb714SLisandro Dalcin }
477fecfb714SLisandro Dalcin 
478f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts)
479f68a32c8SEmil Constantinescu {
480f68a32c8SEmil Constantinescu   TS_RK           *rk   = (TS_RK*)ts->data;
481f68a32c8SEmil Constantinescu   RKTableau        tab  = rk->tableau;
482f68a32c8SEmil Constantinescu   const PetscInt   s = tab->s;
483fecfb714SLisandro Dalcin   const PetscReal *A = tab->A,*c = tab->c;
484f68a32c8SEmil Constantinescu   PetscScalar     *w = rk->work;
485f68a32c8SEmil Constantinescu   Vec             *Y = rk->Y,*YdotRHS = rk->YdotRHS;
486*d1905564SLisandro Dalcin   PetscBool        FSAL = tab->FSAL;
487f68a32c8SEmil Constantinescu   TSAdapt          adapt;
488fecfb714SLisandro Dalcin   PetscInt         i,j;
489be5899b3SLisandro Dalcin   PetscInt         rejections = 0;
490be5899b3SLisandro Dalcin   PetscBool        stageok,accept = PETSC_TRUE;
491be5899b3SLisandro Dalcin   PetscReal        next_time_step = ts->time_step;
492f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
493f68a32c8SEmil Constantinescu 
494f68a32c8SEmil Constantinescu   PetscFunctionBegin;
495*d1905564SLisandro Dalcin   if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE;
496*d1905564SLisandro Dalcin   if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); }
497*d1905564SLisandro Dalcin 
498f68a32c8SEmil Constantinescu   rk->status = TS_STEP_INCOMPLETE;
499be5899b3SLisandro Dalcin   while (!ts->reason && rk->status != TS_STEP_COMPLETE) {
500be5899b3SLisandro Dalcin     PetscReal t = ts->ptime;
501f68a32c8SEmil Constantinescu     PetscReal h = ts->time_step;
502f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
503*d1905564SLisandro Dalcin       if (FSAL && !i) continue;
5049be3e283SDebojyoti Ghosh       rk->stage_time = t + h*c[i];
5059be3e283SDebojyoti Ghosh       ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr);
506f68a32c8SEmil Constantinescu       ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr);
507f68a32c8SEmil Constantinescu       for (j=0; j<i; j++) w[j] = h*A[i*s+j];
508f68a32c8SEmil Constantinescu       ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr);
5099be3e283SDebojyoti Ghosh       ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr);
510f68a32c8SEmil Constantinescu       ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
511be5899b3SLisandro Dalcin       ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr);
512fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
513f68a32c8SEmil Constantinescu       ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr);
514f68a32c8SEmil Constantinescu     }
515be5899b3SLisandro Dalcin 
516be5899b3SLisandro Dalcin     rk->status = TS_STEP_INCOMPLETE;
517f68a32c8SEmil Constantinescu     ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr);
518f68a32c8SEmil Constantinescu     rk->status = TS_STEP_PENDING;
519f68a32c8SEmil Constantinescu     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
520f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr);
521f68a32c8SEmil Constantinescu     ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr);
522fecfb714SLisandro Dalcin     ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
523be5899b3SLisandro Dalcin     rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
524be5899b3SLisandro Dalcin     if (!accept) { /* Roll back the current step */
525fecfb714SLisandro Dalcin       ierr = TSRollBack_RK(ts);CHKERRQ(ierr);
526be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
527be5899b3SLisandro Dalcin       goto reject_step;
528be5899b3SLisandro Dalcin     }
529be5899b3SLisandro Dalcin 
5300f7a1166SHong Zhang     if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/
5310f7a1166SHong Zhang       rk->ptime     = ts->ptime;
5320f7a1166SHong Zhang       rk->time_step = ts->time_step;
533493ed95fSHong Zhang     }
534be5899b3SLisandro Dalcin 
535f68a32c8SEmil Constantinescu     ts->ptime += ts->time_step;
536f68a32c8SEmil Constantinescu     ts->time_step = next_time_step;
537f68a32c8SEmil Constantinescu     break;
538be5899b3SLisandro Dalcin 
539be5899b3SLisandro Dalcin   reject_step:
540fecfb714SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
541be5899b3SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
542be5899b3SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
543be5899b3SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
544f68a32c8SEmil Constantinescu     }
545f68a32c8SEmil Constantinescu   }
546f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
547e4dd521cSBarry Smith }
548e4dd521cSBarry Smith 
549f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts)
550f6a906c0SBarry Smith {
551f6a906c0SBarry Smith   TS_RK         *rk  = (TS_RK*)ts->data;
552be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
553be5899b3SLisandro Dalcin   PetscInt       s   = tab->s;
554f6a906c0SBarry Smith   PetscErrorCode ierr;
555f6a906c0SBarry Smith 
556f6a906c0SBarry Smith   PetscFunctionBegin;
557f6a906c0SBarry Smith   if (ts->adjointsetupcalled++) PetscFunctionReturn(0);
558abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
559f6a906c0SBarry Smith   if(ts->vecs_sensip) {
560abc2977eSBarry Smith     ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
561f6a906c0SBarry Smith   }
562abc2977eSBarry Smith   ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
563f6a906c0SBarry Smith   PetscFunctionReturn(0);
564f6a906c0SBarry Smith }
565f6a906c0SBarry Smith 
56642f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts)
567d2daff3dSHong Zhang {
568c235aa8dSHong Zhang   TS_RK           *rk   = (TS_RK*)ts->data;
569c235aa8dSHong Zhang   RKTableau        tab  = rk->tableau;
570c235aa8dSHong Zhang   const PetscInt   s    = tab->s;
571c235aa8dSHong Zhang   const PetscReal *A = tab->A,*b = tab->b,*c = tab->c;
572c235aa8dSHong Zhang   PetscScalar     *w    = rk->work;
573ad8e2604SHong Zhang   Vec             *Y    = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp;
574b18ea86cSHong Zhang   PetscInt         i,j,nadj;
5753d3eaba7SBarry Smith   PetscReal        t = ts->ptime;
576d2daff3dSHong Zhang   PetscErrorCode   ierr;
577cef76868SBarry Smith   PetscReal        h = ts->time_step;
578cef76868SBarry Smith   Mat              J,Jp;
579c235aa8dSHong Zhang 
580d2daff3dSHong Zhang   PetscFunctionBegin;
581c235aa8dSHong Zhang   rk->status = TS_STEP_INCOMPLETE;
582c235aa8dSHong Zhang   for (i=s-1; i>=0; i--) {
583c235aa8dSHong Zhang     rk->stage_time = t + h*(1.0-c[i]);
584abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
585b18ea86cSHong Zhang       ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr);
586302440fdSBarry Smith       ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr);
587c235aa8dSHong Zhang       for (j=i+1; j<s; j++) {
588302440fdSBarry Smith         ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr);
589b18ea86cSHong Zhang       }
590c235aa8dSHong Zhang     }
591ad8e2604SHong Zhang     /* Stage values of lambda */
592c235aa8dSHong Zhang     ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr);
593c235aa8dSHong Zhang     ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr);
594493ed95fSHong Zhang     if (ts->vec_costintegral) {
595493ed95fSHong Zhang       ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr);
596493ed95fSHong Zhang     }
597abc2977eSBarry Smith     for (nadj=0; nadj<ts->numcost; nadj++) {
598b18ea86cSHong Zhang       ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr);
599493ed95fSHong Zhang       if (ts->vec_costintegral) {
600493ed95fSHong Zhang         ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr);
601493ed95fSHong Zhang       }
602b18ea86cSHong Zhang     }
6036497ce14SHong Zhang 
604ad8e2604SHong Zhang     /* Stage values of mu */
6056497ce14SHong Zhang     if(ts->vecs_sensip) {
6065bf8c567SBarry Smith       ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr);
607493ed95fSHong Zhang       if (ts->vec_costintegral) {
608493ed95fSHong Zhang         ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr);
609493ed95fSHong Zhang       }
610493ed95fSHong Zhang 
611abc2977eSBarry Smith       for (nadj=0; nadj<ts->numcost; nadj++) {
612ad8e2604SHong Zhang         ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr);
613493ed95fSHong Zhang         if (ts->vec_costintegral) {
614493ed95fSHong Zhang           ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr);
615493ed95fSHong Zhang         }
616ad8e2604SHong Zhang       }
617c235aa8dSHong Zhang     }
6186497ce14SHong Zhang   }
619c235aa8dSHong Zhang 
620c235aa8dSHong Zhang   for (j=0; j<s; j++) w[j] = 1.0;
621abc2977eSBarry Smith   for (nadj=0; nadj<ts->numcost; nadj++) {
622b18ea86cSHong Zhang     ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr);
6236497ce14SHong Zhang     if(ts->vecs_sensip) {
624ad8e2604SHong Zhang       ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr);
625b18ea86cSHong Zhang     }
6266497ce14SHong Zhang   }
627c235aa8dSHong Zhang   rk->status = TS_STEP_COMPLETE;
628d2daff3dSHong Zhang   PetscFunctionReturn(0);
629d2daff3dSHong Zhang }
630d2daff3dSHong Zhang 
631f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X)
632f68a32c8SEmil Constantinescu {
633f68a32c8SEmil Constantinescu   TS_RK           *rk = (TS_RK*)ts->data;
6343a8a9803SLisandro Dalcin   PetscInt         s  = rk->tableau->s,p = rk->tableau->p,i,j;
635f68a32c8SEmil Constantinescu   PetscReal        h;
636f68a32c8SEmil Constantinescu   PetscReal        tt,t;
637f68a32c8SEmil Constantinescu   PetscScalar     *b;
638f68a32c8SEmil Constantinescu   const PetscReal *B = rk->tableau->binterp;
639f68a32c8SEmil Constantinescu   PetscErrorCode   ierr;
640e4dd521cSBarry Smith 
641f68a32c8SEmil Constantinescu   PetscFunctionBegin;
642f68a32c8SEmil Constantinescu   if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name);
643be5899b3SLisandro Dalcin 
644f68a32c8SEmil Constantinescu   switch (rk->status) {
645f68a32c8SEmil Constantinescu   case TS_STEP_INCOMPLETE:
646f68a32c8SEmil Constantinescu   case TS_STEP_PENDING:
647f68a32c8SEmil Constantinescu     h = ts->time_step;
648f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h;
649f68a32c8SEmil Constantinescu     break;
650f68a32c8SEmil Constantinescu   case TS_STEP_COMPLETE:
651be5899b3SLisandro Dalcin     h = ts->ptime - ts->ptime_prev;
652f68a32c8SEmil Constantinescu     t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */
653f68a32c8SEmil Constantinescu     break;
654f68a32c8SEmil Constantinescu   default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus");
655e4dd521cSBarry Smith   }
656785e854fSJed Brown   ierr = PetscMalloc1(s,&b);CHKERRQ(ierr);
657f68a32c8SEmil Constantinescu   for (i=0; i<s; i++) b[i] = 0;
6583a8a9803SLisandro Dalcin   for (j=0,tt=t; j<p; j++,tt*=t) {
659f68a32c8SEmil Constantinescu     for (i=0; i<s; i++) {
6603a8a9803SLisandro Dalcin       b[i]  += h * B[i*p+j] * tt;
661e4dd521cSBarry Smith     }
662f68a32c8SEmil Constantinescu   }
663be5899b3SLisandro Dalcin 
664f68a32c8SEmil Constantinescu   ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr);
665f68a32c8SEmil Constantinescu   ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr);
666be5899b3SLisandro Dalcin 
667f68a32c8SEmil Constantinescu   ierr = PetscFree(b);CHKERRQ(ierr);
668e4dd521cSBarry Smith   PetscFunctionReturn(0);
669e4dd521cSBarry Smith }
670e4dd521cSBarry Smith 
671e4dd521cSBarry Smith /*------------------------------------------------------------*/
672be5899b3SLisandro Dalcin 
673be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts)
674be5899b3SLisandro Dalcin {
675be5899b3SLisandro Dalcin   TS_RK          *rk = (TS_RK*)ts->data;
676be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
677be5899b3SLisandro Dalcin   PetscErrorCode ierr;
678be5899b3SLisandro Dalcin 
679be5899b3SLisandro Dalcin   PetscFunctionBegin;
680be5899b3SLisandro Dalcin   if (!tab) PetscFunctionReturn(0);
681be5899b3SLisandro Dalcin   ierr = PetscFree(rk->work);CHKERRQ(ierr);
682be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr);
683be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr);
684be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr);
685be5899b3SLisandro Dalcin   ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr);
686be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
687be5899b3SLisandro Dalcin }
688be5899b3SLisandro Dalcin 
689277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts)
690e4dd521cSBarry Smith {
6915f70b478SJed Brown   TS_RK         *rk = (TS_RK*)ts->data;
6926849ba73SBarry Smith   PetscErrorCode ierr;
693e4dd521cSBarry Smith 
694e4dd521cSBarry Smith   PetscFunctionBegin;
695be5899b3SLisandro Dalcin   ierr = TSRKTableauReset(ts);CHKERRQ(ierr);
6960f7a1166SHong Zhang   ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr);
697abc2977eSBarry Smith   ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr);
698277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
699e4dd521cSBarry Smith }
700277b19d0SLisandro Dalcin 
701277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts)
702277b19d0SLisandro Dalcin {
703277b19d0SLisandro Dalcin   PetscErrorCode ierr;
704277b19d0SLisandro Dalcin 
705277b19d0SLisandro Dalcin   PetscFunctionBegin;
706277b19d0SLisandro Dalcin   ierr = TSReset_RK(ts);CHKERRQ(ierr);
707277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
708f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr);
709f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr);
710f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
711f68a32c8SEmil Constantinescu }
712f68a32c8SEmil Constantinescu 
713f68a32c8SEmil Constantinescu 
714f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx)
715f68a32c8SEmil Constantinescu {
716f68a32c8SEmil Constantinescu   PetscFunctionBegin;
717f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
718f68a32c8SEmil Constantinescu }
719f68a32c8SEmil Constantinescu 
720f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
721f68a32c8SEmil Constantinescu {
722f68a32c8SEmil Constantinescu   PetscFunctionBegin;
723f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
724f68a32c8SEmil Constantinescu }
725f68a32c8SEmil Constantinescu 
726f68a32c8SEmil Constantinescu 
727f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx)
728f68a32c8SEmil Constantinescu {
729f68a32c8SEmil Constantinescu   PetscFunctionBegin;
730f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
731f68a32c8SEmil Constantinescu }
732f68a32c8SEmil Constantinescu 
733f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx)
734f68a32c8SEmil Constantinescu {
735f68a32c8SEmil Constantinescu 
736f68a32c8SEmil Constantinescu   PetscFunctionBegin;
737f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
738f68a32c8SEmil Constantinescu }
739c235aa8dSHong Zhang /*
740d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab)
741d2daff3dSHong Zhang {
742d2daff3dSHong Zhang   PetscReal *A,*b;
743d2daff3dSHong Zhang   PetscInt        s,i,j;
744d2daff3dSHong Zhang   PetscErrorCode  ierr;
745d2daff3dSHong Zhang 
746d2daff3dSHong Zhang   PetscFunctionBegin;
747d2daff3dSHong Zhang   s    = tab->s;
748d2daff3dSHong Zhang   ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr);
749d2daff3dSHong Zhang 
750d2daff3dSHong Zhang   for (i=0; i<s; i++)
751d2daff3dSHong Zhang     for (j=0; j<s; j++) {
752d2daff3dSHong 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];
753d2daff3dSHong Zhang       ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr);
754d2daff3dSHong Zhang     }
755d2daff3dSHong Zhang 
756d2daff3dSHong Zhang   for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i];
757d2daff3dSHong Zhang 
758d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr);
759d2daff3dSHong Zhang   ierr  = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr);
760d2daff3dSHong Zhang   ierr  = PetscFree2(A,b);CHKERRQ(ierr);
761d2daff3dSHong Zhang   PetscFunctionReturn(0);
762d2daff3dSHong Zhang }
763c235aa8dSHong Zhang */
764be5899b3SLisandro Dalcin 
765be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts)
766be5899b3SLisandro Dalcin {
767be5899b3SLisandro Dalcin   TS_RK         *rk  = (TS_RK*)ts->data;
768be5899b3SLisandro Dalcin   RKTableau      tab = rk->tableau;
769be5899b3SLisandro Dalcin   PetscErrorCode ierr;
770be5899b3SLisandro Dalcin 
771be5899b3SLisandro Dalcin   PetscFunctionBegin;
772be5899b3SLisandro Dalcin   ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr);
773be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr);
774be5899b3SLisandro Dalcin   ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr);
775be5899b3SLisandro Dalcin   PetscFunctionReturn(0);
776be5899b3SLisandro Dalcin }
777be5899b3SLisandro Dalcin 
778be5899b3SLisandro Dalcin 
779f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts)
780f68a32c8SEmil Constantinescu {
781f68a32c8SEmil Constantinescu   TS_RK         *rk = (TS_RK*)ts->data;
782f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
783f68a32c8SEmil Constantinescu   DM             dm;
784f68a32c8SEmil Constantinescu 
785f68a32c8SEmil Constantinescu   PetscFunctionBegin;
7863dd83b38SBarry Smith   ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr);
787be5899b3SLisandro Dalcin   ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);
7880f7a1166SHong Zhang   if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */
7890f7a1166SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr);
7900f7a1166SHong Zhang   }
791f68a32c8SEmil Constantinescu   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
792f68a32c8SEmil Constantinescu   if (dm) {
793f68a32c8SEmil Constantinescu     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr);
794f68a32c8SEmil Constantinescu     ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr);
795f68a32c8SEmil Constantinescu   }
796e4dd521cSBarry Smith   PetscFunctionReturn(0);
797e4dd521cSBarry Smith }
798d2daff3dSHong Zhang 
799f6a906c0SBarry Smith 
800e4dd521cSBarry Smith /*------------------------------------------------------------*/
801e4dd521cSBarry Smith 
8024416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts)
803e4dd521cSBarry Smith {
804be5899b3SLisandro Dalcin   TS_RK         *rk = (TS_RK*)ts->data;
805dfbe8321SBarry Smith   PetscErrorCode ierr;
806262ff7bbSBarry Smith 
807e4dd521cSBarry Smith   PetscFunctionBegin;
808e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr);
809f68a32c8SEmil Constantinescu   {
810f68a32c8SEmil Constantinescu     RKTableauLink  link;
811f68a32c8SEmil Constantinescu     PetscInt       count,choice;
812f68a32c8SEmil Constantinescu     PetscBool      flg;
813f68a32c8SEmil Constantinescu     const char   **namelist;
814f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) ;
815785e854fSJed Brown     ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr);
816f68a32c8SEmil Constantinescu     for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name;
817be5899b3SLisandro Dalcin     ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr);
818be5899b3SLisandro Dalcin     if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);}
819f68a32c8SEmil Constantinescu     ierr = PetscFree(namelist);CHKERRQ(ierr);
820f68a32c8SEmil Constantinescu   }
821262ff7bbSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
822e4dd521cSBarry Smith   PetscFunctionReturn(0);
823e4dd521cSBarry Smith }
824e4dd521cSBarry Smith 
825f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[])
826f68a32c8SEmil Constantinescu {
827f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
828f68a32c8SEmil Constantinescu   PetscInt       i;
829f68a32c8SEmil Constantinescu   size_t         left,count;
830f68a32c8SEmil Constantinescu   char           *p;
831f68a32c8SEmil Constantinescu 
832f68a32c8SEmil Constantinescu   PetscFunctionBegin;
833f68a32c8SEmil Constantinescu   for (i=0,p=buf,left=len; i<n; i++) {
834f68a32c8SEmil Constantinescu     ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr);
835f68a32c8SEmil Constantinescu     if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer");
836f68a32c8SEmil Constantinescu     left -= count;
837f68a32c8SEmil Constantinescu     p    += count;
838f68a32c8SEmil Constantinescu     *p++  = ' ';
839f68a32c8SEmil Constantinescu   }
840f68a32c8SEmil Constantinescu   p[i ? 0 : -1] = 0;
841f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
842f68a32c8SEmil Constantinescu }
843f68a32c8SEmil Constantinescu 
8445f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer)
845e4dd521cSBarry Smith {
8465f70b478SJed Brown   TS_RK          *rk = (TS_RK*)ts->data;
8478370ee3bSLisandro Dalcin   PetscBool      iascii;
848dfbe8321SBarry Smith   PetscErrorCode ierr;
8492cdf8120Sasbjorn 
850e4dd521cSBarry Smith   PetscFunctionBegin;
851251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
8528370ee3bSLisandro Dalcin   if (iascii) {
8539c334d8fSLisandro Dalcin     RKTableau tab  = rk->tableau;
854f68a32c8SEmil Constantinescu     TSRKType  rktype;
855f68a32c8SEmil Constantinescu     char      buf[512];
856f68a32c8SEmil Constantinescu     ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr);
8570c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  RK: %s\n",rktype);CHKERRQ(ierr);
8580c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order: %D\n",tab->order);CHKERRQ(ierr);
8590c37eb8eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr);
860f68a32c8SEmil Constantinescu     ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr);
861f68a32c8SEmil Constantinescu     ierr = PetscViewerASCIIPrintf(viewer,"  Abscissa c = %s\n",buf);CHKERRQ(ierr);
8628370ee3bSLisandro Dalcin   }
8639c334d8fSLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
864f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
865f68a32c8SEmil Constantinescu }
866f68a32c8SEmil Constantinescu 
867f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer)
868f68a32c8SEmil Constantinescu {
869f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
8709c334d8fSLisandro Dalcin   TSAdapt        adapt;
871f68a32c8SEmil Constantinescu 
872f68a32c8SEmil Constantinescu   PetscFunctionBegin;
8739c334d8fSLisandro Dalcin   ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
8749c334d8fSLisandro Dalcin   ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr);
875f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
876f68a32c8SEmil Constantinescu }
877f68a32c8SEmil Constantinescu 
878f68a32c8SEmil Constantinescu /*@C
879f68a32c8SEmil Constantinescu   TSRKSetType - Set the type of RK scheme
880f68a32c8SEmil Constantinescu 
881f68a32c8SEmil Constantinescu   Logically collective
882f68a32c8SEmil Constantinescu 
883f68a32c8SEmil Constantinescu   Input Parameter:
884f68a32c8SEmil Constantinescu +  ts - timestepping context
885f68a32c8SEmil Constantinescu -  rktype - type of RK-scheme
886f68a32c8SEmil Constantinescu 
887f68a32c8SEmil Constantinescu   Level: intermediate
888f68a32c8SEmil Constantinescu 
889f68a32c8SEmil Constantinescu .seealso: TSRKGetType(), TSRK, TSRK2, TSRK3, TSRKPRSSP2, TSRK4, TSRK5
890f68a32c8SEmil Constantinescu @*/
891f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype)
892f68a32c8SEmil Constantinescu {
893f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
894f68a32c8SEmil Constantinescu 
895f68a32c8SEmil Constantinescu   PetscFunctionBegin;
896f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
897f68a32c8SEmil Constantinescu   ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr);
898f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
899f68a32c8SEmil Constantinescu }
900f68a32c8SEmil Constantinescu 
901f68a32c8SEmil Constantinescu /*@C
902f68a32c8SEmil Constantinescu   TSRKGetType - Get the type of RK scheme
903f68a32c8SEmil Constantinescu 
904f68a32c8SEmil Constantinescu   Logically collective
905f68a32c8SEmil Constantinescu 
906f68a32c8SEmil Constantinescu   Input Parameter:
907f68a32c8SEmil Constantinescu .  ts - timestepping context
908f68a32c8SEmil Constantinescu 
909f68a32c8SEmil Constantinescu   Output Parameter:
910f68a32c8SEmil Constantinescu .  rktype - type of RK-scheme
911f68a32c8SEmil Constantinescu 
912f68a32c8SEmil Constantinescu   Level: intermediate
913f68a32c8SEmil Constantinescu 
914f68a32c8SEmil Constantinescu .seealso: TSRKGetType()
915f68a32c8SEmil Constantinescu @*/
916f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype)
917f68a32c8SEmil Constantinescu {
918f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
919f68a32c8SEmil Constantinescu 
920f68a32c8SEmil Constantinescu   PetscFunctionBegin;
921f68a32c8SEmil Constantinescu   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
922f68a32c8SEmil Constantinescu   ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr);
923f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
924f68a32c8SEmil Constantinescu }
925f68a32c8SEmil Constantinescu 
926560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype)
927f68a32c8SEmil Constantinescu {
928f68a32c8SEmil Constantinescu   TS_RK *rk = (TS_RK*)ts->data;
929f68a32c8SEmil Constantinescu 
930f68a32c8SEmil Constantinescu   PetscFunctionBegin;
931f68a32c8SEmil Constantinescu   *rktype = rk->tableau->name;
932f68a32c8SEmil Constantinescu   PetscFunctionReturn(0);
933f68a32c8SEmil Constantinescu }
934560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype)
935f68a32c8SEmil Constantinescu {
936f68a32c8SEmil Constantinescu   TS_RK          *rk = (TS_RK*)ts->data;
937f68a32c8SEmil Constantinescu   PetscErrorCode ierr;
938f68a32c8SEmil Constantinescu   PetscBool      match;
939f68a32c8SEmil Constantinescu   RKTableauLink  link;
940f68a32c8SEmil Constantinescu 
941f68a32c8SEmil Constantinescu   PetscFunctionBegin;
942f68a32c8SEmil Constantinescu   if (rk->tableau) {
943f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr);
944f68a32c8SEmil Constantinescu     if (match) PetscFunctionReturn(0);
945f68a32c8SEmil Constantinescu   }
946f68a32c8SEmil Constantinescu   for (link = RKTableauList; link; link=link->next) {
947f68a32c8SEmil Constantinescu     ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr);
948f68a32c8SEmil Constantinescu     if (match) {
949be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);}
950f68a32c8SEmil Constantinescu       rk->tableau = &link->tab;
951be5899b3SLisandro Dalcin       if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);}
952f68a32c8SEmil Constantinescu       PetscFunctionReturn(0);
953f68a32c8SEmil Constantinescu     }
954f68a32c8SEmil Constantinescu   }
955f68a32c8SEmil Constantinescu   SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype);
956e4dd521cSBarry Smith   PetscFunctionReturn(0);
957e4dd521cSBarry Smith }
958e4dd521cSBarry Smith 
959ff22ae23SHong Zhang static PetscErrorCode  TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y)
960ff22ae23SHong Zhang {
961ff22ae23SHong Zhang   TS_RK          *rk = (TS_RK*)ts->data;
962ff22ae23SHong Zhang 
963ff22ae23SHong Zhang   PetscFunctionBegin;
964ff22ae23SHong Zhang   *ns = rk->tableau->s;
965d2daff3dSHong Zhang   if(Y) *Y  = rk->Y;
966ff22ae23SHong Zhang   PetscFunctionReturn(0);
967ff22ae23SHong Zhang }
968ff22ae23SHong Zhang 
969ff22ae23SHong Zhang 
970e4dd521cSBarry Smith /* ------------------------------------------------------------ */
971ebe3b25bSBarry Smith /*MC
972f68a32c8SEmil Constantinescu       TSRK - ODE and DAE solver using Runge-Kutta schemes
973ebe3b25bSBarry Smith 
974f68a32c8SEmil Constantinescu   The user should provide the right hand side of the equation
975f68a32c8SEmil Constantinescu   using TSSetRHSFunction().
976ebe3b25bSBarry Smith 
977f68a32c8SEmil Constantinescu   Notes:
978f68a32c8SEmil Constantinescu   The default is TSRK3, it can be changed with TSRKSetType() or -ts_rk_type
979ebe3b25bSBarry Smith 
980d41222bbSBarry Smith   Level: beginner
981d41222bbSBarry Smith 
982f68a32c8SEmil Constantinescu .seealso:  TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3,
983f68a32c8SEmil Constantinescu            TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister()
984ebe3b25bSBarry Smith 
985ebe3b25bSBarry Smith M*/
9868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts)
987e4dd521cSBarry Smith {
9881566a47fSLisandro Dalcin   TS_RK          *rk;
989dfbe8321SBarry Smith   PetscErrorCode ierr;
990e4dd521cSBarry Smith 
991e4dd521cSBarry Smith   PetscFunctionBegin;
992f68a32c8SEmil Constantinescu   ierr = TSRKInitializePackage();CHKERRQ(ierr);
993f68a32c8SEmil Constantinescu 
994f68a32c8SEmil Constantinescu   ts->ops->reset          = TSReset_RK;
9955f70b478SJed Brown   ts->ops->destroy        = TSDestroy_RK;
9965f70b478SJed Brown   ts->ops->view           = TSView_RK;
997f68a32c8SEmil Constantinescu   ts->ops->load           = TSLoad_RK;
998f68a32c8SEmil Constantinescu   ts->ops->setup          = TSSetUp_RK;
99942f2b339SBarry Smith   ts->ops->adjointsetup   = TSAdjointSetUp_RK;
1000f68a32c8SEmil Constantinescu   ts->ops->step           = TSStep_RK;
1001f68a32c8SEmil Constantinescu   ts->ops->interpolate    = TSInterpolate_RK;
1002f68a32c8SEmil Constantinescu   ts->ops->evaluatestep   = TSEvaluateStep_RK;
1003fecfb714SLisandro Dalcin   ts->ops->rollback       = TSRollBack_RK;
1004f68a32c8SEmil Constantinescu   ts->ops->setfromoptions = TSSetFromOptions_RK;
1005ff22ae23SHong Zhang   ts->ops->getstages      = TSGetStages_RK;
100642f2b339SBarry Smith   ts->ops->adjointstep    = TSAdjointStep_RK;
1007e4dd521cSBarry Smith 
10080f7a1166SHong Zhang   ts->ops->adjointintegral = TSAdjointCostIntegral_RK;
10090f7a1166SHong Zhang   ts->ops->forwardintegral = TSForwardCostIntegral_RK;
10100f7a1166SHong Zhang 
10111566a47fSLisandro Dalcin   ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr);
10121566a47fSLisandro Dalcin   ts->data = (void*)rk;
1013e4dd521cSBarry Smith 
1014f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr);
1015f68a32c8SEmil Constantinescu   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr);
1016be5899b3SLisandro Dalcin 
1017be5899b3SLisandro Dalcin   ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr);
1018e4dd521cSBarry Smith   PetscFunctionReturn(0);
1019e4dd521cSBarry Smith }
1020