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 52*287c2655SBarry Smith TSRK1FE - First order forward Euler scheme. 53262ff7bbSBarry Smith 54f68a32c8SEmil Constantinescu This method has one stage. 55f68a32c8SEmil Constantinescu 56*287c2655SBarry Smith Options database: 57*287c2655SBarry Smith . -ts_tk_type 1fe 58*287c2655SBarry Smith 59f68a32c8SEmil Constantinescu Level: advanced 60f68a32c8SEmil Constantinescu 61*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 62f68a32c8SEmil Constantinescu M*/ 63f68a32c8SEmil Constantinescu /*MC 642109b73fSDebojyoti Ghosh TSRK2A - Second order RK scheme. 65f68a32c8SEmil Constantinescu 66f68a32c8SEmil Constantinescu This method has two stages. 67f68a32c8SEmil Constantinescu 68*287c2655SBarry Smith Options database: 69*287c2655SBarry Smith . -ts_tk_type 2a 70*287c2655SBarry Smith 71f68a32c8SEmil Constantinescu Level: advanced 72f68a32c8SEmil Constantinescu 73*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 74f68a32c8SEmil Constantinescu M*/ 75f68a32c8SEmil Constantinescu /*MC 76f68a32c8SEmil Constantinescu TSRK3 - Third order RK scheme. 77f68a32c8SEmil Constantinescu 78f68a32c8SEmil Constantinescu This method has three stages. 79f68a32c8SEmil Constantinescu 80*287c2655SBarry Smith Options database: 81*287c2655SBarry Smith . -ts_tk_type 3 82*287c2655SBarry Smith 83f68a32c8SEmil Constantinescu Level: advanced 84f68a32c8SEmil Constantinescu 85*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 86f68a32c8SEmil Constantinescu M*/ 87f68a32c8SEmil Constantinescu /*MC 882109b73fSDebojyoti Ghosh TSRK3BS - Third order RK scheme of Bogacki-Shampine with 2nd order embedded method. 892109b73fSDebojyoti Ghosh 90d1905564SLisandro Dalcin This method has four stages with the First Same As Last (FSAL) property. 912109b73fSDebojyoti Ghosh 92*287c2655SBarry Smith Options database: 93*287c2655SBarry Smith . -ts_tk_type 3bs 94*287c2655SBarry Smith 952109b73fSDebojyoti Ghosh Level: advanced 962109b73fSDebojyoti Ghosh 9798164e13SLisandro Dalcin References: https://doi.org/10.1016/0893-9659(89)90079-7 98d1905564SLisandro Dalcin 99*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1002109b73fSDebojyoti Ghosh M*/ 1012109b73fSDebojyoti Ghosh /*MC 102f68a32c8SEmil Constantinescu TSRK4 - Fourth order RK scheme. 103f68a32c8SEmil Constantinescu 1042e698f8bSDebojyoti Ghosh This is the classical Runge-Kutta method with four stages. 105f68a32c8SEmil Constantinescu 106*287c2655SBarry Smith Options database: 107*287c2655SBarry Smith . -ts_tk_type 4 108*287c2655SBarry Smith 109f68a32c8SEmil Constantinescu Level: advanced 110f68a32c8SEmil Constantinescu 111*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 112f68a32c8SEmil Constantinescu M*/ 113f68a32c8SEmil Constantinescu /*MC 1142e698f8bSDebojyoti Ghosh TSRK5F - Fifth order Fehlberg RK scheme with a 4th order embedded method. 115f68a32c8SEmil Constantinescu 116f68a32c8SEmil Constantinescu This method has six stages. 117f68a32c8SEmil Constantinescu 118*287c2655SBarry Smith Options database: 119*287c2655SBarry Smith . -ts_tk_type 5f 120*287c2655SBarry Smith 121f68a32c8SEmil Constantinescu Level: advanced 122f68a32c8SEmil Constantinescu 123*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 124f68a32c8SEmil Constantinescu M*/ 1252109b73fSDebojyoti Ghosh /*MC 1262e698f8bSDebojyoti Ghosh TSRK5DP - Fifth order Dormand-Prince RK scheme with the 4th order embedded method. 1272109b73fSDebojyoti Ghosh 128d1905564SLisandro Dalcin This method has seven stages with the First Same As Last (FSAL) property. 1292109b73fSDebojyoti Ghosh 130*287c2655SBarry Smith Options database: 131*287c2655SBarry Smith . -ts_tk_type 5dp 132*287c2655SBarry Smith 1332109b73fSDebojyoti Ghosh Level: advanced 1342109b73fSDebojyoti Ghosh 13598164e13SLisandro Dalcin References: https://doi.org/10.1016/0771-050X(80)90013-3 136d1905564SLisandro Dalcin 137*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 1382109b73fSDebojyoti Ghosh M*/ 13905e23783SLisandro Dalcin /*MC 14005e23783SLisandro Dalcin TSRK5BS - Fifth order Bogacki-Shampine RK scheme with 4th order embedded method. 14105e23783SLisandro Dalcin 14205e23783SLisandro Dalcin This method has eight stages with the First Same As Last (FSAL) property. 14305e23783SLisandro Dalcin 144*287c2655SBarry Smith Options database: 145*287c2655SBarry Smith . -ts_tk_type 5bs 146*287c2655SBarry Smith 14705e23783SLisandro Dalcin Level: advanced 14805e23783SLisandro Dalcin 14998164e13SLisandro Dalcin References: https://doi.org/10.1016/0898-1221(96)00141-1 15005e23783SLisandro Dalcin 151*287c2655SBarry Smith .seealso: TSRK, TSRKType, TSRKSetType() 15205e23783SLisandro Dalcin M*/ 153262ff7bbSBarry Smith 154f68a32c8SEmil Constantinescu /*@C 155f68a32c8SEmil Constantinescu TSRKRegisterAll - Registers all of the Runge-Kutta explicit methods in TSRK 156262ff7bbSBarry Smith 157f68a32c8SEmil Constantinescu Not Collective, but should be called by all processes which will need the schemes to be registered 158262ff7bbSBarry Smith 159f68a32c8SEmil Constantinescu Level: advanced 160262ff7bbSBarry Smith 161f68a32c8SEmil Constantinescu .keywords: TS, TSRK, register, all 162262ff7bbSBarry Smith 163f68a32c8SEmil Constantinescu .seealso: TSRKRegisterDestroy() 164262ff7bbSBarry Smith @*/ 165f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterAll(void) 166262ff7bbSBarry Smith { 1674ac538c5SBarry Smith PetscErrorCode ierr; 168262ff7bbSBarry Smith 169262ff7bbSBarry Smith PetscFunctionBegin; 170f68a32c8SEmil Constantinescu if (TSRKRegisterAllCalled) PetscFunctionReturn(0); 171f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_TRUE; 172e4dd521cSBarry Smith 173e4dd521cSBarry Smith { 174f68a32c8SEmil Constantinescu const PetscReal 175f68a32c8SEmil Constantinescu A[1][1] = {{0.0}}, 176f68a32c8SEmil Constantinescu b[1] = {1.0}; 1773a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK1FE,1,1,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 178e4dd521cSBarry Smith } 179db046809SBarry Smith { 180f68a32c8SEmil Constantinescu const PetscReal 181f68a32c8SEmil Constantinescu A[2][2] = {{0.0,0.0}, 182f68a32c8SEmil Constantinescu {1.0,0.0}}, 183f68a32c8SEmil Constantinescu b[2] = {0.5,0.5}, 184f68a32c8SEmil Constantinescu bembed[2] = {1.0,0}; 1853a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK2A,2,2,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 186db046809SBarry Smith } 187f68a32c8SEmil Constantinescu { 188f68a32c8SEmil Constantinescu const PetscReal 18917f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 19017f6384fSBarry Smith A[3][3] = {{0,0,0}, 19117f6384fSBarry Smith {2.0q/3.0q,0,0}, 19217f6384fSBarry Smith {-1.0q/3.0q,1.0,0}}, 19317f6384fSBarry Smith #else 194f68a32c8SEmil Constantinescu A[3][3] = {{0,0,0}, 195f68a32c8SEmil Constantinescu {2.0/3.0,0,0}, 196f68a32c8SEmil Constantinescu {-1.0/3.0,1.0,0}}, 19717f6384fSBarry Smith #endif 198f68a32c8SEmil Constantinescu b[3] = {0.25,0.5,0.25}; 1993a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3,3,3,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 200fdefd5e5SDebojyoti Ghosh } 201fdefd5e5SDebojyoti Ghosh { 202fdefd5e5SDebojyoti Ghosh const PetscReal 20317f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 20417f6384fSBarry Smith A[4][4] = {{0,0,0,0}, 20517f6384fSBarry Smith {1.0q/2.0q,0,0,0}, 20617f6384fSBarry Smith {0,3.0q/4.0q,0,0}, 20717f6384fSBarry Smith {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}}, 20817f6384fSBarry Smith b[4] = {2.0q/9.0q,1.0q/3.0q,4.0q/9.0q,0}, 20917f6384fSBarry Smith bembed[4] = {7.0q/24.0q,1.0q/4.0q,1.0q/3.0q,1.0q/8.0q}; 21017f6384fSBarry Smith #else 211fdefd5e5SDebojyoti Ghosh A[4][4] = {{0,0,0,0}, 212fdefd5e5SDebojyoti Ghosh {1.0/2.0,0,0,0}, 213fdefd5e5SDebojyoti Ghosh {0,3.0/4.0,0,0}, 214fdefd5e5SDebojyoti Ghosh {2.0/9.0,1.0/3.0,4.0/9.0,0}}, 215fdefd5e5SDebojyoti Ghosh b[4] = {2.0/9.0,1.0/3.0,4.0/9.0,0}, 216fdefd5e5SDebojyoti Ghosh bembed[4] = {7.0/24.0,1.0/4.0,1.0/3.0,1.0/8.0}; 21717f6384fSBarry Smith #endif 2183a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK3BS,3,4,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 219db046809SBarry Smith } 220f68a32c8SEmil Constantinescu { 221f68a32c8SEmil Constantinescu const PetscReal 222f68a32c8SEmil Constantinescu A[4][4] = {{0,0,0,0}, 223f68a32c8SEmil Constantinescu {0.5,0,0,0}, 224f68a32c8SEmil Constantinescu {0,0.5,0,0}, 225f68a32c8SEmil Constantinescu {0,0,1.0,0}}, 22617f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 22717f6384fSBarry Smith b[4] = {1.0q/6.0q,1.0q/3.0q,1.0q/3.0q,1.0q/6.0q}; 22817f6384fSBarry Smith #else 229f68a32c8SEmil Constantinescu b[4] = {1.0/6.0,1.0/3.0,1.0/3.0,1.0/6.0}; 23017f6384fSBarry Smith #endif 2313a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK4,4,4,&A[0][0],b,NULL,NULL,0,NULL);CHKERRQ(ierr); 232db046809SBarry Smith } 233f68a32c8SEmil Constantinescu { 234f68a32c8SEmil Constantinescu const PetscReal 23517f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 23617f6384fSBarry Smith A[6][6] = {{0,0,0,0,0,0}, 23717f6384fSBarry Smith {0.25q,0,0,0,0,0}, 23817f6384fSBarry Smith {3.0q/32.0q,9.0q/32.0q,0,0,0,0}, 23917f6384fSBarry Smith {1932.0q/2197.0q,-7200.0q/2197.0q,7296.0q/2197.0q,0,0,0}, 24017f6384fSBarry Smith {439.0q/216.0q,-8.0q,3680.0q/513.0q,-845.0q/4104.0q,0,0}, 24117f6384fSBarry Smith {-8.0q/27.0q,2.0q,-3544.0q/2565.0q,1859.0q/4104.0q,-11.0q/40.0q,0}}, 24217f6384fSBarry Smith b[6] = {16.0q/135.0q,0,6656.0q/12825.0q,28561.0q/56430.0q,-9.0q/50.0q,2.0q/55.0q}, 24317f6384fSBarry Smith bembed[6] = {25.0q/216.0q,0,1408.0q/2565.0q,2197.0q/4104.0q,-1.0q/5.0q,0}; 24417f6384fSBarry Smith #else 245f68a32c8SEmil Constantinescu A[6][6] = {{0,0,0,0,0,0}, 246f68a32c8SEmil Constantinescu {0.25,0,0,0,0,0}, 247f68a32c8SEmil Constantinescu {3.0/32.0,9.0/32.0,0,0,0,0}, 248f68a32c8SEmil Constantinescu {1932.0/2197.0,-7200.0/2197.0,7296.0/2197.0,0,0,0}, 249f68a32c8SEmil Constantinescu {439.0/216.0,-8.0,3680.0/513.0,-845.0/4104.0,0,0}, 250f68a32c8SEmil Constantinescu {-8.0/27.0,2.0,-3544.0/2565.0,1859.0/4104.0,-11.0/40.0,0}}, 251f68a32c8SEmil Constantinescu b[6] = {16.0/135.0,0,6656.0/12825.0,28561.0/56430.0,-9.0/50.0,2.0/55.0}, 252f68a32c8SEmil Constantinescu bembed[6] = {25.0/216.0,0,1408.0/2565.0,2197.0/4104.0,-1.0/5.0,0}; 25317f6384fSBarry Smith #endif 2543a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5F,5,6,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 255fdefd5e5SDebojyoti Ghosh } 256fdefd5e5SDebojyoti Ghosh { 257fdefd5e5SDebojyoti Ghosh const PetscReal 25817f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 25917f6384fSBarry Smith A[7][7] = {{0,0,0,0,0,0,0}, 26017f6384fSBarry Smith {1.0q/5.0q,0,0,0,0,0,0}, 26117f6384fSBarry Smith {3.0q/40.0q,9.0q/40.0q,0,0,0,0,0}, 26217f6384fSBarry Smith {44.0q/45.0q,-56.0q/15.0q,32.0q/9.0q,0,0,0,0}, 26317f6384fSBarry Smith {19372.0q/6561.0q,-25360.0q/2187.0q,64448.0q/6561.0q,-212.0q/729.0q,0,0,0}, 26417f6384fSBarry Smith {9017.0q/3168.0q,-355.0q/33.0q,46732.0q/5247.0q,49.0q/176.0q,-5103.0q/18656.0q,0,0}, 26517f6384fSBarry Smith {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}}, 26617f6384fSBarry Smith b[7] = {35.0q/384.0q,0,500.0q/1113.0q,125.0q/192.0q,-2187.0q/6784.0q,11.0q/84.0q,0}, 26717f6384fSBarry Smith bembed[7] = {5179.0q/57600.0q,0,7571.0q/16695.0q,393.0q/640.0q,-92097.0q/339200.0q,187.0q/2100.0q,1.0q/40.0q}; 26817f6384fSBarry Smith #else 269fdefd5e5SDebojyoti Ghosh A[7][7] = {{0,0,0,0,0,0,0}, 270fdefd5e5SDebojyoti Ghosh {1.0/5.0,0,0,0,0,0,0}, 271fdefd5e5SDebojyoti Ghosh {3.0/40.0,9.0/40.0,0,0,0,0,0}, 272fdefd5e5SDebojyoti Ghosh {44.0/45.0,-56.0/15.0,32.0/9.0,0,0,0,0}, 273fdefd5e5SDebojyoti Ghosh {19372.0/6561.0,-25360.0/2187.0,64448.0/6561.0,-212.0/729.0,0,0,0}, 274fdefd5e5SDebojyoti Ghosh {9017.0/3168.0,-355.0/33.0,46732.0/5247.0,49.0/176.0,-5103.0/18656.0,0,0}, 275fdefd5e5SDebojyoti Ghosh {35.0/384.0,0,500.0/1113.0,125.0/192.0,-2187.0/6784.0,11.0/84.0,0}}, 276fdefd5e5SDebojyoti 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}, 277fdefd5e5SDebojyoti 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}; 27817f6384fSBarry Smith #endif 2793a8a9803SLisandro Dalcin ierr = TSRKRegister(TSRK5DP,5,7,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 280f68a32c8SEmil Constantinescu } 28105e23783SLisandro Dalcin { 28205e23783SLisandro Dalcin const PetscReal 28317f6384fSBarry Smith #if defined(PETSC_USE_REAL___FLOAT128) 28417f6384fSBarry Smith A[8][8] = {{0,0,0,0,0,0,0,0}, 28517f6384fSBarry Smith {1.0q/6.0q,0,0,0,0,0,0,0}, 28617f6384fSBarry Smith {2.0q/27.0q,4.0q/27.0q,0,0,0,0,0,0}, 28717f6384fSBarry Smith {183.0q/1372.0q,-162.0q/343.0q,1053.0q/1372.0q,0,0,0,0,0}, 28817f6384fSBarry Smith {68.0q/297.0q,-4.0q/11.0q,42.0q/143.0q,1960.0q/3861.0q,0,0,0,0}, 28917f6384fSBarry Smith {597.0q/22528.0q,81.0q/352.0q,63099.0q/585728.0q,58653.0q/366080.0q,4617.0q/20480.0q,0,0,0}, 29017f6384fSBarry Smith {174197.0q/959244.0q,-30942.0q/79937.0q,8152137.0q/19744439.0q,666106.0q/1039181.0q,-29421.0q/29068.0q,482048.0q/414219.0q,0,0}, 29117f6384fSBarry Smith {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}}, 29217f6384fSBarry Smith b[8] = {587.0q/8064.0q,0,4440339.0q/15491840.0q,24353.0q/124800.0q,387.0q/44800.0q,2152.0q/5985.0q,7267.0q/94080.0q,0}, 29317f6384fSBarry Smith bembed[8] = {2479.0q/34992.0q,0,123.0q/416.0q,612941.0q/3411720.0q,43.0q/1440.0q,2272.0q/6561.0q,79937.0q/1113912.0q,3293.0q/556956.0q}; 29417f6384fSBarry Smith #else 29505e23783SLisandro Dalcin A[8][8] = {{0,0,0,0,0,0,0,0}, 29605e23783SLisandro Dalcin {1.0/6.0,0,0,0,0,0,0,0}, 29705e23783SLisandro Dalcin {2.0/27.0,4.0/27.0,0,0,0,0,0,0}, 29805e23783SLisandro Dalcin {183.0/1372.0,-162.0/343.0,1053.0/1372.0,0,0,0,0,0}, 29905e23783SLisandro Dalcin {68.0/297.0,-4.0/11.0,42.0/143.0,1960.0/3861.0,0,0,0,0}, 30005e23783SLisandro Dalcin {597.0/22528.0,81.0/352.0,63099.0/585728.0,58653.0/366080.0,4617.0/20480.0,0,0,0}, 30105e23783SLisandro Dalcin {174197.0/959244.0,-30942.0/79937.0,8152137.0/19744439.0,666106.0/1039181.0,-29421.0/29068.0,482048.0/414219.0,0,0}, 30205e23783SLisandro Dalcin {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}}, 30305e23783SLisandro Dalcin b[8] = {587.0/8064.0,0,4440339.0/15491840.0,24353.0/124800.0,387.0/44800.0,2152.0/5985.0,7267.0/94080.0,0}, 30405e23783SLisandro Dalcin bembed[8] = {2479.0/34992.0,0,123.0/416.0,612941.0/3411720.0,43.0/1440.0,2272.0/6561.0,79937.0/1113912.0,3293.0/556956.0}; 30517f6384fSBarry Smith #endif 30605e23783SLisandro Dalcin ierr = TSRKRegister(TSRK5BS,5,8,&A[0][0],b,NULL,bembed,0,NULL);CHKERRQ(ierr); 30705e23783SLisandro Dalcin } 308db046809SBarry Smith PetscFunctionReturn(0); 309db046809SBarry Smith } 310db046809SBarry Smith 311f68a32c8SEmil Constantinescu /*@C 312f68a32c8SEmil Constantinescu TSRKRegisterDestroy - Frees the list of schemes that were registered by TSRKRegister(). 313f68a32c8SEmil Constantinescu 314f68a32c8SEmil Constantinescu Not Collective 315f68a32c8SEmil Constantinescu 316f68a32c8SEmil Constantinescu Level: advanced 317f68a32c8SEmil Constantinescu 318f68a32c8SEmil Constantinescu .keywords: TSRK, register, destroy 319f68a32c8SEmil Constantinescu .seealso: TSRKRegister(), TSRKRegisterAll() 320f68a32c8SEmil Constantinescu @*/ 321f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegisterDestroy(void) 322e4dd521cSBarry Smith { 323f68a32c8SEmil Constantinescu PetscErrorCode ierr; 324f68a32c8SEmil Constantinescu RKTableauLink link; 325f68a32c8SEmil Constantinescu 326f68a32c8SEmil Constantinescu PetscFunctionBegin; 327f68a32c8SEmil Constantinescu while ((link = RKTableauList)) { 328f68a32c8SEmil Constantinescu RKTableau t = &link->tab; 329f68a32c8SEmil Constantinescu RKTableauList = link->next; 330f68a32c8SEmil Constantinescu ierr = PetscFree3(t->A,t->b,t->c); CHKERRQ(ierr); 331f68a32c8SEmil Constantinescu ierr = PetscFree (t->bembed); CHKERRQ(ierr); 332f68a32c8SEmil Constantinescu ierr = PetscFree (t->binterp); CHKERRQ(ierr); 333f68a32c8SEmil Constantinescu ierr = PetscFree (t->name); CHKERRQ(ierr); 334f68a32c8SEmil Constantinescu ierr = PetscFree (link); CHKERRQ(ierr); 335f68a32c8SEmil Constantinescu } 336f68a32c8SEmil Constantinescu TSRKRegisterAllCalled = PETSC_FALSE; 337f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 338f68a32c8SEmil Constantinescu } 339f68a32c8SEmil Constantinescu 340f68a32c8SEmil Constantinescu /*@C 341f68a32c8SEmil Constantinescu TSRKInitializePackage - This function initializes everything in the TSRK package. It is called 342f68a32c8SEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RK() 343f68a32c8SEmil Constantinescu when using static libraries. 344f68a32c8SEmil Constantinescu 345f68a32c8SEmil Constantinescu Level: developer 346f68a32c8SEmil Constantinescu 347f68a32c8SEmil Constantinescu .keywords: TS, TSRK, initialize, package 348f68a32c8SEmil Constantinescu .seealso: PetscInitialize() 349f68a32c8SEmil Constantinescu @*/ 350f68a32c8SEmil Constantinescu PetscErrorCode TSRKInitializePackage(void) 351f68a32c8SEmil Constantinescu { 352186e87acSLisandro Dalcin PetscErrorCode ierr; 353e4dd521cSBarry Smith 354e4dd521cSBarry Smith PetscFunctionBegin; 355f68a32c8SEmil Constantinescu if (TSRKPackageInitialized) PetscFunctionReturn(0); 356f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_TRUE; 357f68a32c8SEmil Constantinescu ierr = TSRKRegisterAll();CHKERRQ(ierr); 358f68a32c8SEmil Constantinescu ierr = PetscRegisterFinalize(TSRKFinalizePackage);CHKERRQ(ierr); 359f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 360f68a32c8SEmil Constantinescu } 361186e87acSLisandro Dalcin 362f68a32c8SEmil Constantinescu /*@C 363f68a32c8SEmil Constantinescu TSRKFinalizePackage - This function destroys everything in the TSRK package. It is 364f68a32c8SEmil Constantinescu called from PetscFinalize(). 365186e87acSLisandro Dalcin 366f68a32c8SEmil Constantinescu Level: developer 367f68a32c8SEmil Constantinescu 368f68a32c8SEmil Constantinescu .keywords: Petsc, destroy, package 369f68a32c8SEmil Constantinescu .seealso: PetscFinalize() 370f68a32c8SEmil Constantinescu @*/ 371f68a32c8SEmil Constantinescu PetscErrorCode TSRKFinalizePackage(void) 372f68a32c8SEmil Constantinescu { 373f68a32c8SEmil Constantinescu PetscErrorCode ierr; 374f68a32c8SEmil Constantinescu 375f68a32c8SEmil Constantinescu PetscFunctionBegin; 376f68a32c8SEmil Constantinescu TSRKPackageInitialized = PETSC_FALSE; 377f68a32c8SEmil Constantinescu ierr = TSRKRegisterDestroy();CHKERRQ(ierr); 378f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 379f68a32c8SEmil Constantinescu } 380f68a32c8SEmil Constantinescu 381f68a32c8SEmil Constantinescu /*@C 382f68a32c8SEmil Constantinescu TSRKRegister - register an RK scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 383f68a32c8SEmil Constantinescu 384f68a32c8SEmil Constantinescu Not Collective, but the same schemes should be registered on all processes on which they will be used 385f68a32c8SEmil Constantinescu 386f68a32c8SEmil Constantinescu Input Parameters: 387f68a32c8SEmil Constantinescu + name - identifier for method 388f68a32c8SEmil Constantinescu . order - approximation order of method 389f68a32c8SEmil Constantinescu . s - number of stages, this is the dimension of the matrices below 390f68a32c8SEmil Constantinescu . A - stage coefficients (dimension s*s, row-major) 391f68a32c8SEmil Constantinescu . b - step completion table (dimension s; NULL to use last row of A) 392f68a32c8SEmil Constantinescu . c - abscissa (dimension s; NULL to use row sums of A) 393f68a32c8SEmil Constantinescu . bembed - completion table for embedded method (dimension s; NULL if not available) 3943a8a9803SLisandro Dalcin . p - Order of the interpolation scheme, equal to the number of columns of binterp 3953a8a9803SLisandro Dalcin - binterp - Coefficients of the interpolation formula (dimension s*p; NULL to reuse b with p=1) 396f68a32c8SEmil Constantinescu 397f68a32c8SEmil Constantinescu Notes: 398f68a32c8SEmil Constantinescu Several RK methods are provided, this function is only needed to create new methods. 399f68a32c8SEmil Constantinescu 400f68a32c8SEmil Constantinescu Level: advanced 401f68a32c8SEmil Constantinescu 402f68a32c8SEmil Constantinescu .keywords: TS, register 403f68a32c8SEmil Constantinescu 404f68a32c8SEmil Constantinescu .seealso: TSRK 405f68a32c8SEmil Constantinescu @*/ 406f68a32c8SEmil Constantinescu PetscErrorCode TSRKRegister(TSRKType name,PetscInt order,PetscInt s, 407f68a32c8SEmil Constantinescu const PetscReal A[],const PetscReal b[],const PetscReal c[], 4083a8a9803SLisandro Dalcin const PetscReal bembed[],PetscInt p,const PetscReal binterp[]) 409f68a32c8SEmil Constantinescu { 410f68a32c8SEmil Constantinescu PetscErrorCode ierr; 411f68a32c8SEmil Constantinescu RKTableauLink link; 412f68a32c8SEmil Constantinescu RKTableau t; 413f68a32c8SEmil Constantinescu PetscInt i,j; 414f68a32c8SEmil Constantinescu 415f68a32c8SEmil Constantinescu PetscFunctionBegin; 4163a8a9803SLisandro Dalcin PetscValidCharPointer(name,1); 4173a8a9803SLisandro Dalcin PetscValidRealPointer(A,4); 4183a8a9803SLisandro Dalcin if (b) PetscValidRealPointer(b,5); 4193a8a9803SLisandro Dalcin if (c) PetscValidRealPointer(c,6); 4203a8a9803SLisandro Dalcin if (bembed) PetscValidRealPointer(bembed,7); 4213a8a9803SLisandro Dalcin if (binterp || p > 1) PetscValidRealPointer(binterp,9); 4223a8a9803SLisandro Dalcin 42395dccacaSBarry Smith ierr = PetscNew(&link);CHKERRQ(ierr); 424f68a32c8SEmil Constantinescu t = &link->tab; 4253a8a9803SLisandro Dalcin 426f68a32c8SEmil Constantinescu ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 427f68a32c8SEmil Constantinescu t->order = order; 428f68a32c8SEmil Constantinescu t->s = s; 429dcca6d9dSJed Brown ierr = PetscMalloc3(s*s,&t->A,s,&t->b,s,&t->c);CHKERRQ(ierr); 430f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 431f68a32c8SEmil Constantinescu if (b) { ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr); } 432f68a32c8SEmil Constantinescu else for (i=0; i<s; i++) t->b[i] = A[(s-1)*s+i]; 433f68a32c8SEmil Constantinescu if (c) { ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr); } 434f68a32c8SEmil 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]; 435d760c35bSDebojyoti Ghosh t->FSAL = PETSC_TRUE; 436d760c35bSDebojyoti Ghosh for (i=0; i<s; i++) if (t->A[(s-1)*s+i] != t->b[i]) t->FSAL = PETSC_FALSE; 4373a8a9803SLisandro Dalcin 438f68a32c8SEmil Constantinescu if (bembed) { 439785e854fSJed Brown ierr = PetscMalloc1(s,&t->bembed);CHKERRQ(ierr); 440f68a32c8SEmil Constantinescu ierr = PetscMemcpy(t->bembed,bembed,s*sizeof(bembed[0]));CHKERRQ(ierr); 441f68a32c8SEmil Constantinescu } 442f68a32c8SEmil Constantinescu 4433a8a9803SLisandro Dalcin if (!binterp) { p = 1; binterp = t->b; } 4443a8a9803SLisandro Dalcin t->p = p; 4453a8a9803SLisandro Dalcin ierr = PetscMalloc1(s*p,&t->binterp);CHKERRQ(ierr); 4463a8a9803SLisandro Dalcin ierr = PetscMemcpy(t->binterp,binterp,s*p*sizeof(binterp[0]));CHKERRQ(ierr); 4473a8a9803SLisandro Dalcin 448f68a32c8SEmil Constantinescu link->next = RKTableauList; 449f68a32c8SEmil Constantinescu RKTableauList = link; 450f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 451f68a32c8SEmil Constantinescu } 452f68a32c8SEmil Constantinescu 453e4dd521cSBarry Smith /* 454f68a32c8SEmil Constantinescu The step completion formula is 455e4dd521cSBarry Smith 456f68a32c8SEmil Constantinescu x1 = x0 + h b^T YdotRHS 457f68a32c8SEmil Constantinescu 458f68a32c8SEmil Constantinescu This function can be called before or after ts->vec_sol has been updated. 459f68a32c8SEmil Constantinescu Suppose we have a completion formula (b) and an embedded formula (be) of different order. 460f68a32c8SEmil Constantinescu We can write 461f68a32c8SEmil Constantinescu 462f68a32c8SEmil Constantinescu x1e = x0 + h be^T YdotRHS 463f68a32c8SEmil Constantinescu = x1 - h b^T YdotRHS + h be^T YdotRHS 464f68a32c8SEmil Constantinescu = x1 + h (be - b)^T YdotRHS 465f68a32c8SEmil Constantinescu 466f68a32c8SEmil Constantinescu so we can evaluate the method with different order even after the step has been optimistically completed. 467e4dd521cSBarry Smith */ 468f68a32c8SEmil Constantinescu static PetscErrorCode TSEvaluateStep_RK(TS ts,PetscInt order,Vec X,PetscBool *done) 469f68a32c8SEmil Constantinescu { 470f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 471f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 472f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 473f68a32c8SEmil Constantinescu PetscReal h; 474f68a32c8SEmil Constantinescu PetscInt s = tab->s,j; 475f68a32c8SEmil Constantinescu PetscErrorCode ierr; 476f68a32c8SEmil Constantinescu 477f68a32c8SEmil Constantinescu PetscFunctionBegin; 478f68a32c8SEmil Constantinescu switch (rk->status) { 479f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 480f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 481f68a32c8SEmil Constantinescu h = ts->time_step; break; 482f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 483be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 484f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 485f68a32c8SEmil Constantinescu } 486f68a32c8SEmil Constantinescu if (order == tab->order) { 487f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { 488f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 489f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->b[j]; 490f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 491f68a32c8SEmil Constantinescu } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 492f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 493f68a32c8SEmil Constantinescu } else if (order == tab->order-1) { 494f68a32c8SEmil Constantinescu if (!tab->bembed) goto unavailable; 495f68a32c8SEmil Constantinescu if (rk->status == TS_STEP_INCOMPLETE) { /* Complete with the embedded method (be) */ 496f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 497f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*tab->bembed[j]; 498f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 499f68a32c8SEmil Constantinescu } else { /* Rollback and re-complete using (be-b) */ 500f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 501f68a32c8SEmil Constantinescu for (j=0; j<s; j++) w[j] = h*(tab->bembed[j] - tab->b[j]); 502f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,w,rk->YdotRHS);CHKERRQ(ierr); 5030f7a1166SHong Zhang if (ts->vec_costintegral && ts->costintegralfwd) { 5040f7a1166SHong Zhang ierr = VecCopy(rk->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 5050f7a1166SHong Zhang } 506f68a32c8SEmil Constantinescu } 507f68a32c8SEmil Constantinescu if (done) *done = PETSC_TRUE; 508f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 509f68a32c8SEmil Constantinescu } 510f68a32c8SEmil Constantinescu unavailable: 511f68a32c8SEmil Constantinescu if (done) *done = PETSC_FALSE; 512a7fac7c2SEmil 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); 513f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 514f68a32c8SEmil Constantinescu } 515f68a32c8SEmil Constantinescu 5160f7a1166SHong Zhang static PetscErrorCode TSForwardCostIntegral_RK(TS ts) 5170f7a1166SHong Zhang { 5180f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5190f7a1166SHong Zhang RKTableau tab = rk->tableau; 5200f7a1166SHong Zhang const PetscInt s = tab->s; 5210f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5220f7a1166SHong Zhang Vec *Y = rk->Y; 5230f7a1166SHong Zhang PetscInt i; 5240f7a1166SHong Zhang PetscErrorCode ierr; 5250f7a1166SHong Zhang 5260f7a1166SHong Zhang PetscFunctionBegin; 5270f7a1166SHong Zhang /* backup cost integral */ 5280f7a1166SHong Zhang ierr = VecCopy(ts->vec_costintegral,rk->VecCostIntegral0);CHKERRQ(ierr); 5290f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5300f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 5310f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,rk->ptime+rk->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5320f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,rk->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5330f7a1166SHong Zhang } 5340f7a1166SHong Zhang PetscFunctionReturn(0); 5350f7a1166SHong Zhang } 5360f7a1166SHong Zhang 5370f7a1166SHong Zhang static PetscErrorCode TSAdjointCostIntegral_RK(TS ts) 5380f7a1166SHong Zhang { 5390f7a1166SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 5400f7a1166SHong Zhang RKTableau tab = rk->tableau; 5410f7a1166SHong Zhang const PetscInt s = tab->s; 5420f7a1166SHong Zhang const PetscReal *b = tab->b,*c = tab->c; 5430f7a1166SHong Zhang Vec *Y = rk->Y; 5440f7a1166SHong Zhang PetscInt i; 5450f7a1166SHong Zhang PetscErrorCode ierr; 5460f7a1166SHong Zhang 5470f7a1166SHong Zhang PetscFunctionBegin; 5480f7a1166SHong Zhang for (i=s-1; i>=0; i--) { 5490f7a1166SHong Zhang /* Evolve ts->vec_costintegral to compute integrals */ 5500f7a1166SHong Zhang ierr = TSAdjointComputeCostIntegrand(ts,ts->ptime-ts->time_step*(1.0-c[i]),Y[i],ts->vec_costintegrand);CHKERRQ(ierr); 5510f7a1166SHong Zhang ierr = VecAXPY(ts->vec_costintegral,-ts->time_step*b[i],ts->vec_costintegrand);CHKERRQ(ierr); 5520f7a1166SHong Zhang } 5530f7a1166SHong Zhang PetscFunctionReturn(0); 5540f7a1166SHong Zhang } 5550f7a1166SHong Zhang 556fecfb714SLisandro Dalcin static PetscErrorCode TSRollBack_RK(TS ts) 557fecfb714SLisandro Dalcin { 558fecfb714SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 559fecfb714SLisandro Dalcin RKTableau tab = rk->tableau; 560fecfb714SLisandro Dalcin const PetscInt s = tab->s; 561fecfb714SLisandro Dalcin const PetscReal *b = tab->b; 562fecfb714SLisandro Dalcin PetscScalar *w = rk->work; 563fecfb714SLisandro Dalcin Vec *YdotRHS = rk->YdotRHS; 564fecfb714SLisandro Dalcin PetscInt j; 565fecfb714SLisandro Dalcin PetscReal h; 566fecfb714SLisandro Dalcin PetscErrorCode ierr; 567fecfb714SLisandro Dalcin 568fecfb714SLisandro Dalcin PetscFunctionBegin; 569fecfb714SLisandro Dalcin switch (rk->status) { 570fecfb714SLisandro Dalcin case TS_STEP_INCOMPLETE: 571fecfb714SLisandro Dalcin case TS_STEP_PENDING: 572fecfb714SLisandro Dalcin h = ts->time_step; break; 573fecfb714SLisandro Dalcin case TS_STEP_COMPLETE: 574fecfb714SLisandro Dalcin h = ts->ptime - ts->ptime_prev; break; 575fecfb714SLisandro Dalcin default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 576fecfb714SLisandro Dalcin } 577fecfb714SLisandro Dalcin for (j=0; j<s; j++) w[j] = -h*b[j]; 578fecfb714SLisandro Dalcin ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 579fecfb714SLisandro Dalcin PetscFunctionReturn(0); 580fecfb714SLisandro Dalcin } 581fecfb714SLisandro Dalcin 582f68a32c8SEmil Constantinescu static PetscErrorCode TSStep_RK(TS ts) 583f68a32c8SEmil Constantinescu { 584f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 585f68a32c8SEmil Constantinescu RKTableau tab = rk->tableau; 586f68a32c8SEmil Constantinescu const PetscInt s = tab->s; 587fecfb714SLisandro Dalcin const PetscReal *A = tab->A,*c = tab->c; 588f68a32c8SEmil Constantinescu PetscScalar *w = rk->work; 589f68a32c8SEmil Constantinescu Vec *Y = rk->Y,*YdotRHS = rk->YdotRHS; 590d1905564SLisandro Dalcin PetscBool FSAL = tab->FSAL; 591f68a32c8SEmil Constantinescu TSAdapt adapt; 592fecfb714SLisandro Dalcin PetscInt i,j; 593be5899b3SLisandro Dalcin PetscInt rejections = 0; 594be5899b3SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 595be5899b3SLisandro Dalcin PetscReal next_time_step = ts->time_step; 596f68a32c8SEmil Constantinescu PetscErrorCode ierr; 597f68a32c8SEmil Constantinescu 598f68a32c8SEmil Constantinescu PetscFunctionBegin; 599d1905564SLisandro Dalcin if (ts->steprollback || ts->steprestart) FSAL = PETSC_FALSE; 600d1905564SLisandro Dalcin if (FSAL) { ierr = VecCopy(YdotRHS[s-1],YdotRHS[0]);CHKERRQ(ierr); } 601d1905564SLisandro Dalcin 602f68a32c8SEmil Constantinescu rk->status = TS_STEP_INCOMPLETE; 603be5899b3SLisandro Dalcin while (!ts->reason && rk->status != TS_STEP_COMPLETE) { 604be5899b3SLisandro Dalcin PetscReal t = ts->ptime; 605f68a32c8SEmil Constantinescu PetscReal h = ts->time_step; 606f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 6079be3e283SDebojyoti Ghosh rk->stage_time = t + h*c[i]; 6089be3e283SDebojyoti Ghosh ierr = TSPreStage(ts,rk->stage_time); CHKERRQ(ierr); 609f68a32c8SEmil Constantinescu ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 610f68a32c8SEmil Constantinescu for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 611f68a32c8SEmil Constantinescu ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 6129be3e283SDebojyoti Ghosh ierr = TSPostStage(ts,rk->stage_time,i,Y); CHKERRQ(ierr); 613f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 614be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(adapt,ts,rk->stage_time,Y[i],&stageok);CHKERRQ(ierr); 615fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 6168f738a7cSLisandro Dalcin if (FSAL && !i) continue; 617f68a32c8SEmil Constantinescu ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 618f68a32c8SEmil Constantinescu } 619be5899b3SLisandro Dalcin 620be5899b3SLisandro Dalcin rk->status = TS_STEP_INCOMPLETE; 621f68a32c8SEmil Constantinescu ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 622f68a32c8SEmil Constantinescu rk->status = TS_STEP_PENDING; 623f68a32c8SEmil Constantinescu ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 624f68a32c8SEmil Constantinescu ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 6251917a363SLisandro Dalcin ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE);CHKERRQ(ierr); 626fecfb714SLisandro Dalcin ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 627be5899b3SLisandro Dalcin rk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 628be5899b3SLisandro Dalcin if (!accept) { /* Roll back the current step */ 629fecfb714SLisandro Dalcin ierr = TSRollBack_RK(ts);CHKERRQ(ierr); 630be5899b3SLisandro Dalcin ts->time_step = next_time_step; 631be5899b3SLisandro Dalcin goto reject_step; 632be5899b3SLisandro Dalcin } 633be5899b3SLisandro Dalcin 6340f7a1166SHong Zhang if (ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation*/ 6350f7a1166SHong Zhang rk->ptime = ts->ptime; 6360f7a1166SHong Zhang rk->time_step = ts->time_step; 637493ed95fSHong Zhang } 638be5899b3SLisandro Dalcin 639f68a32c8SEmil Constantinescu ts->ptime += ts->time_step; 640f68a32c8SEmil Constantinescu ts->time_step = next_time_step; 641f68a32c8SEmil Constantinescu break; 642be5899b3SLisandro Dalcin 643be5899b3SLisandro Dalcin reject_step: 644fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 645be5899b3SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 646be5899b3SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 647be5899b3SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 648f68a32c8SEmil Constantinescu } 649f68a32c8SEmil Constantinescu } 650f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 651e4dd521cSBarry Smith } 652e4dd521cSBarry Smith 653f6a906c0SBarry Smith static PetscErrorCode TSAdjointSetUp_RK(TS ts) 654f6a906c0SBarry Smith { 655f6a906c0SBarry Smith TS_RK *rk = (TS_RK*)ts->data; 656be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 657be5899b3SLisandro Dalcin PetscInt s = tab->s; 658f6a906c0SBarry Smith PetscErrorCode ierr; 659f6a906c0SBarry Smith 660f6a906c0SBarry Smith PetscFunctionBegin; 661f6a906c0SBarry Smith if (ts->adjointsetupcalled++) PetscFunctionReturn(0); 662abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 663f6a906c0SBarry Smith if(ts->vecs_sensip) { 664abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensip[0],s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 665f6a906c0SBarry Smith } 666abc2977eSBarry Smith ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 667f6a906c0SBarry Smith PetscFunctionReturn(0); 668f6a906c0SBarry Smith } 669f6a906c0SBarry Smith 67042f2b339SBarry Smith static PetscErrorCode TSAdjointStep_RK(TS ts) 671d2daff3dSHong Zhang { 672c235aa8dSHong Zhang TS_RK *rk = (TS_RK*)ts->data; 673c235aa8dSHong Zhang RKTableau tab = rk->tableau; 674c235aa8dSHong Zhang const PetscInt s = tab->s; 675c235aa8dSHong Zhang const PetscReal *A = tab->A,*b = tab->b,*c = tab->c; 676c235aa8dSHong Zhang PetscScalar *w = rk->work; 677ad8e2604SHong Zhang Vec *Y = rk->Y,*VecDeltaLam = rk->VecDeltaLam,*VecDeltaMu = rk->VecDeltaMu,*VecSensiTemp = rk->VecSensiTemp; 678b18ea86cSHong Zhang PetscInt i,j,nadj; 6793d3eaba7SBarry Smith PetscReal t = ts->ptime; 680d2daff3dSHong Zhang PetscErrorCode ierr; 681cef76868SBarry Smith PetscReal h = ts->time_step; 682cef76868SBarry Smith Mat J,Jp; 683c235aa8dSHong Zhang 684d2daff3dSHong Zhang PetscFunctionBegin; 685c235aa8dSHong Zhang rk->status = TS_STEP_INCOMPLETE; 686c235aa8dSHong Zhang for (i=s-1; i>=0; i--) { 687c235aa8dSHong Zhang rk->stage_time = t + h*(1.0-c[i]); 688abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 689b18ea86cSHong Zhang ierr = VecCopy(ts->vecs_sensi[nadj],VecSensiTemp[nadj]);CHKERRQ(ierr); 690302440fdSBarry Smith ierr = VecScale(VecSensiTemp[nadj],-h*b[i]);CHKERRQ(ierr); 691c235aa8dSHong Zhang for (j=i+1; j<s; j++) { 692302440fdSBarry Smith ierr = VecAXPY(VecSensiTemp[nadj],-h*A[j*s+i],VecDeltaLam[nadj*s+j]);CHKERRQ(ierr); 693b18ea86cSHong Zhang } 694c235aa8dSHong Zhang } 695ad8e2604SHong Zhang /* Stage values of lambda */ 696c235aa8dSHong Zhang ierr = TSGetRHSJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 697c235aa8dSHong Zhang ierr = TSComputeRHSJacobian(ts,rk->stage_time,Y[i],J,Jp);CHKERRQ(ierr); 698493ed95fSHong Zhang if (ts->vec_costintegral) { 699493ed95fSHong Zhang ierr = TSAdjointComputeDRDYFunction(ts,rk->stage_time,Y[i],ts->vecs_drdy);CHKERRQ(ierr); 700493ed95fSHong Zhang } 701abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 702b18ea86cSHong Zhang ierr = MatMultTranspose(J,VecSensiTemp[nadj],VecDeltaLam[nadj*s+i]);CHKERRQ(ierr); 703493ed95fSHong Zhang if (ts->vec_costintegral) { 704493ed95fSHong Zhang ierr = VecAXPY(VecDeltaLam[nadj*s+i],-h*b[i],ts->vecs_drdy[nadj]);CHKERRQ(ierr); 705493ed95fSHong Zhang } 706b18ea86cSHong Zhang } 7076497ce14SHong Zhang 708ad8e2604SHong Zhang /* Stage values of mu */ 7096497ce14SHong Zhang if(ts->vecs_sensip) { 7105bf8c567SBarry Smith ierr = TSAdjointComputeRHSJacobian(ts,rk->stage_time,Y[i],ts->Jacp);CHKERRQ(ierr); 711493ed95fSHong Zhang if (ts->vec_costintegral) { 712493ed95fSHong Zhang ierr = TSAdjointComputeDRDPFunction(ts,rk->stage_time,Y[i],ts->vecs_drdp);CHKERRQ(ierr); 713493ed95fSHong Zhang } 714493ed95fSHong Zhang 715abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 716ad8e2604SHong Zhang ierr = MatMultTranspose(ts->Jacp,VecSensiTemp[nadj],VecDeltaMu[nadj*s+i]);CHKERRQ(ierr); 717493ed95fSHong Zhang if (ts->vec_costintegral) { 718493ed95fSHong Zhang ierr = VecAXPY(VecDeltaMu[nadj*s+i],-h*b[i],ts->vecs_drdp[nadj]);CHKERRQ(ierr); 719493ed95fSHong Zhang } 720ad8e2604SHong Zhang } 721c235aa8dSHong Zhang } 7226497ce14SHong Zhang } 723c235aa8dSHong Zhang 724c235aa8dSHong Zhang for (j=0; j<s; j++) w[j] = 1.0; 725abc2977eSBarry Smith for (nadj=0; nadj<ts->numcost; nadj++) { 726b18ea86cSHong Zhang ierr = VecMAXPY(ts->vecs_sensi[nadj],s,w,&VecDeltaLam[nadj*s]);CHKERRQ(ierr); 7276497ce14SHong Zhang if(ts->vecs_sensip) { 728ad8e2604SHong Zhang ierr = VecMAXPY(ts->vecs_sensip[nadj],s,w,&VecDeltaMu[nadj*s]);CHKERRQ(ierr); 729b18ea86cSHong Zhang } 7306497ce14SHong Zhang } 731c235aa8dSHong Zhang rk->status = TS_STEP_COMPLETE; 732d2daff3dSHong Zhang PetscFunctionReturn(0); 733d2daff3dSHong Zhang } 734d2daff3dSHong Zhang 735f68a32c8SEmil Constantinescu static PetscErrorCode TSInterpolate_RK(TS ts,PetscReal itime,Vec X) 736f68a32c8SEmil Constantinescu { 737f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 7383a8a9803SLisandro Dalcin PetscInt s = rk->tableau->s,p = rk->tableau->p,i,j; 739f68a32c8SEmil Constantinescu PetscReal h; 740f68a32c8SEmil Constantinescu PetscReal tt,t; 741f68a32c8SEmil Constantinescu PetscScalar *b; 742f68a32c8SEmil Constantinescu const PetscReal *B = rk->tableau->binterp; 743f68a32c8SEmil Constantinescu PetscErrorCode ierr; 744e4dd521cSBarry Smith 745f68a32c8SEmil Constantinescu PetscFunctionBegin; 746f68a32c8SEmil Constantinescu if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSRK %s does not have an interpolation formula",rk->tableau->name); 747be5899b3SLisandro Dalcin 748f68a32c8SEmil Constantinescu switch (rk->status) { 749f68a32c8SEmil Constantinescu case TS_STEP_INCOMPLETE: 750f68a32c8SEmil Constantinescu case TS_STEP_PENDING: 751f68a32c8SEmil Constantinescu h = ts->time_step; 752f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h; 753f68a32c8SEmil Constantinescu break; 754f68a32c8SEmil Constantinescu case TS_STEP_COMPLETE: 755be5899b3SLisandro Dalcin h = ts->ptime - ts->ptime_prev; 756f68a32c8SEmil Constantinescu t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 757f68a32c8SEmil Constantinescu break; 758f68a32c8SEmil Constantinescu default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 759e4dd521cSBarry Smith } 760785e854fSJed Brown ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 761f68a32c8SEmil Constantinescu for (i=0; i<s; i++) b[i] = 0; 7623a8a9803SLisandro Dalcin for (j=0,tt=t; j<p; j++,tt*=t) { 763f68a32c8SEmil Constantinescu for (i=0; i<s; i++) { 7643a8a9803SLisandro Dalcin b[i] += h * B[i*p+j] * tt; 765e4dd521cSBarry Smith } 766f68a32c8SEmil Constantinescu } 767be5899b3SLisandro Dalcin 768f68a32c8SEmil Constantinescu ierr = VecCopy(rk->Y[0],X);CHKERRQ(ierr); 769f68a32c8SEmil Constantinescu ierr = VecMAXPY(X,s,b,rk->YdotRHS);CHKERRQ(ierr); 770be5899b3SLisandro Dalcin 771f68a32c8SEmil Constantinescu ierr = PetscFree(b);CHKERRQ(ierr); 772e4dd521cSBarry Smith PetscFunctionReturn(0); 773e4dd521cSBarry Smith } 774e4dd521cSBarry Smith 775e4dd521cSBarry Smith /*------------------------------------------------------------*/ 776be5899b3SLisandro Dalcin 777be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauReset(TS ts) 778be5899b3SLisandro Dalcin { 779be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 780be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 781be5899b3SLisandro Dalcin PetscErrorCode ierr; 782be5899b3SLisandro Dalcin 783be5899b3SLisandro Dalcin PetscFunctionBegin; 784be5899b3SLisandro Dalcin if (!tab) PetscFunctionReturn(0); 785be5899b3SLisandro Dalcin ierr = PetscFree(rk->work);CHKERRQ(ierr); 786be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->Y);CHKERRQ(ierr); 787be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s,&rk->YdotRHS);CHKERRQ(ierr); 788be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaLam);CHKERRQ(ierr); 789be5899b3SLisandro Dalcin ierr = VecDestroyVecs(tab->s*ts->numcost,&rk->VecDeltaMu);CHKERRQ(ierr); 790be5899b3SLisandro Dalcin PetscFunctionReturn(0); 791be5899b3SLisandro Dalcin } 792be5899b3SLisandro Dalcin 793277b19d0SLisandro Dalcin static PetscErrorCode TSReset_RK(TS ts) 794e4dd521cSBarry Smith { 7955f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 7966849ba73SBarry Smith PetscErrorCode ierr; 797e4dd521cSBarry Smith 798e4dd521cSBarry Smith PetscFunctionBegin; 799be5899b3SLisandro Dalcin ierr = TSRKTableauReset(ts);CHKERRQ(ierr); 8000f7a1166SHong Zhang ierr = VecDestroy(&rk->VecCostIntegral0);CHKERRQ(ierr); 801abc2977eSBarry Smith ierr = VecDestroyVecs(ts->numcost,&rk->VecSensiTemp);CHKERRQ(ierr); 802277b19d0SLisandro Dalcin PetscFunctionReturn(0); 803e4dd521cSBarry Smith } 804277b19d0SLisandro Dalcin 805277b19d0SLisandro Dalcin static PetscErrorCode TSDestroy_RK(TS ts) 806277b19d0SLisandro Dalcin { 807277b19d0SLisandro Dalcin PetscErrorCode ierr; 808277b19d0SLisandro Dalcin 809277b19d0SLisandro Dalcin PetscFunctionBegin; 810277b19d0SLisandro Dalcin ierr = TSReset_RK(ts);CHKERRQ(ierr); 811277b19d0SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 812f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",NULL);CHKERRQ(ierr); 813f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",NULL);CHKERRQ(ierr); 814f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 815f68a32c8SEmil Constantinescu } 816f68a32c8SEmil Constantinescu 817f68a32c8SEmil Constantinescu 818f68a32c8SEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSRK(DM fine,DM coarse,void *ctx) 819f68a32c8SEmil Constantinescu { 820f68a32c8SEmil Constantinescu PetscFunctionBegin; 821f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 822f68a32c8SEmil Constantinescu } 823f68a32c8SEmil Constantinescu 824f68a32c8SEmil Constantinescu static PetscErrorCode DMRestrictHook_TSRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 825f68a32c8SEmil Constantinescu { 826f68a32c8SEmil Constantinescu PetscFunctionBegin; 827f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 828f68a32c8SEmil Constantinescu } 829f68a32c8SEmil Constantinescu 830f68a32c8SEmil Constantinescu 831f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSRK(DM dm,DM subdm,void *ctx) 832f68a32c8SEmil Constantinescu { 833f68a32c8SEmil Constantinescu PetscFunctionBegin; 834f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 835f68a32c8SEmil Constantinescu } 836f68a32c8SEmil Constantinescu 837f68a32c8SEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 838f68a32c8SEmil Constantinescu { 839f68a32c8SEmil Constantinescu 840f68a32c8SEmil Constantinescu PetscFunctionBegin; 841f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 842f68a32c8SEmil Constantinescu } 843c235aa8dSHong Zhang /* 844d2daff3dSHong Zhang static PetscErrorCode RKSetAdjCoe(RKTableau tab) 845d2daff3dSHong Zhang { 846d2daff3dSHong Zhang PetscReal *A,*b; 847d2daff3dSHong Zhang PetscInt s,i,j; 848d2daff3dSHong Zhang PetscErrorCode ierr; 849d2daff3dSHong Zhang 850d2daff3dSHong Zhang PetscFunctionBegin; 851d2daff3dSHong Zhang s = tab->s; 852d2daff3dSHong Zhang ierr = PetscMalloc2(s*s,&A,s,&b);CHKERRQ(ierr); 853d2daff3dSHong Zhang 854d2daff3dSHong Zhang for (i=0; i<s; i++) 855d2daff3dSHong Zhang for (j=0; j<s; j++) { 856d2daff3dSHong 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]; 857d2daff3dSHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"Coefficients: A[%D][%D]=%.6f\n",i,j,A[i*s+j]);CHKERRQ(ierr); 858d2daff3dSHong Zhang } 859d2daff3dSHong Zhang 860d2daff3dSHong Zhang for (i=0; i<s; i++) b[i] = (tab->b[s-1-i]==0)? 0: -tab->b[s-1-i]; 861d2daff3dSHong Zhang 862d2daff3dSHong Zhang ierr = PetscMemcpy(tab->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 863d2daff3dSHong Zhang ierr = PetscMemcpy(tab->b,b,s*sizeof(b[0]));CHKERRQ(ierr); 864d2daff3dSHong Zhang ierr = PetscFree2(A,b);CHKERRQ(ierr); 865d2daff3dSHong Zhang PetscFunctionReturn(0); 866d2daff3dSHong Zhang } 867c235aa8dSHong Zhang */ 868be5899b3SLisandro Dalcin 869be5899b3SLisandro Dalcin static PetscErrorCode TSRKTableauSetUp(TS ts) 870be5899b3SLisandro Dalcin { 871be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 872be5899b3SLisandro Dalcin RKTableau tab = rk->tableau; 873be5899b3SLisandro Dalcin PetscErrorCode ierr; 874be5899b3SLisandro Dalcin 875be5899b3SLisandro Dalcin PetscFunctionBegin; 876be5899b3SLisandro Dalcin ierr = PetscMalloc1(tab->s,&rk->work);CHKERRQ(ierr); 877be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->Y);CHKERRQ(ierr); 878be5899b3SLisandro Dalcin ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&rk->YdotRHS);CHKERRQ(ierr); 879be5899b3SLisandro Dalcin PetscFunctionReturn(0); 880be5899b3SLisandro Dalcin } 881be5899b3SLisandro Dalcin 882be5899b3SLisandro Dalcin 883f68a32c8SEmil Constantinescu static PetscErrorCode TSSetUp_RK(TS ts) 884f68a32c8SEmil Constantinescu { 885f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 886f68a32c8SEmil Constantinescu PetscErrorCode ierr; 887f68a32c8SEmil Constantinescu DM dm; 888f68a32c8SEmil Constantinescu 889f68a32c8SEmil Constantinescu PetscFunctionBegin; 8903dd83b38SBarry Smith ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 891be5899b3SLisandro Dalcin ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr); 8920f7a1166SHong Zhang if (!rk->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 8930f7a1166SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&rk->VecCostIntegral0);CHKERRQ(ierr); 8940f7a1166SHong Zhang } 895f68a32c8SEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 896f68a32c8SEmil Constantinescu if (dm) { 897f68a32c8SEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSRK,DMRestrictHook_TSRK,ts);CHKERRQ(ierr); 898f68a32c8SEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSRK,DMSubDomainRestrictHook_TSRK,ts);CHKERRQ(ierr); 899f68a32c8SEmil Constantinescu } 900e4dd521cSBarry Smith PetscFunctionReturn(0); 901e4dd521cSBarry Smith } 902d2daff3dSHong Zhang 903f6a906c0SBarry Smith 904e4dd521cSBarry Smith /*------------------------------------------------------------*/ 905e4dd521cSBarry Smith 9064416b707SBarry Smith static PetscErrorCode TSSetFromOptions_RK(PetscOptionItems *PetscOptionsObject,TS ts) 907e4dd521cSBarry Smith { 908be5899b3SLisandro Dalcin TS_RK *rk = (TS_RK*)ts->data; 909dfbe8321SBarry Smith PetscErrorCode ierr; 910262ff7bbSBarry Smith 911e4dd521cSBarry Smith PetscFunctionBegin; 912e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"RK ODE solver options");CHKERRQ(ierr); 913f68a32c8SEmil Constantinescu { 914f68a32c8SEmil Constantinescu RKTableauLink link; 915f68a32c8SEmil Constantinescu PetscInt count,choice; 916f68a32c8SEmil Constantinescu PetscBool flg; 917f68a32c8SEmil Constantinescu const char **namelist; 918f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) ; 919785e854fSJed Brown ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 920f68a32c8SEmil Constantinescu for (link=RKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 921be5899b3SLisandro Dalcin ierr = PetscOptionsEList("-ts_rk_type","Family of RK method","TSRKSetType",(const char*const*)namelist,count,rk->tableau->name,&choice,&flg);CHKERRQ(ierr); 922be5899b3SLisandro Dalcin if (flg) {ierr = TSRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 923f68a32c8SEmil Constantinescu ierr = PetscFree(namelist);CHKERRQ(ierr); 924f68a32c8SEmil Constantinescu } 925262ff7bbSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 926e4dd521cSBarry Smith PetscFunctionReturn(0); 927e4dd521cSBarry Smith } 928e4dd521cSBarry Smith 929f68a32c8SEmil Constantinescu static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 930f68a32c8SEmil Constantinescu { 931f68a32c8SEmil Constantinescu PetscErrorCode ierr; 932f68a32c8SEmil Constantinescu PetscInt i; 933f68a32c8SEmil Constantinescu size_t left,count; 934f68a32c8SEmil Constantinescu char *p; 935f68a32c8SEmil Constantinescu 936f68a32c8SEmil Constantinescu PetscFunctionBegin; 937f68a32c8SEmil Constantinescu for (i=0,p=buf,left=len; i<n; i++) { 938f68a32c8SEmil Constantinescu ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 939f68a32c8SEmil Constantinescu if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 940f68a32c8SEmil Constantinescu left -= count; 941f68a32c8SEmil Constantinescu p += count; 942f68a32c8SEmil Constantinescu *p++ = ' '; 943f68a32c8SEmil Constantinescu } 944f68a32c8SEmil Constantinescu p[i ? 0 : -1] = 0; 945f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 946f68a32c8SEmil Constantinescu } 947f68a32c8SEmil Constantinescu 9485f70b478SJed Brown static PetscErrorCode TSView_RK(TS ts,PetscViewer viewer) 949e4dd521cSBarry Smith { 9505f70b478SJed Brown TS_RK *rk = (TS_RK*)ts->data; 9518370ee3bSLisandro Dalcin PetscBool iascii; 952dfbe8321SBarry Smith PetscErrorCode ierr; 9532cdf8120Sasbjorn 954e4dd521cSBarry Smith PetscFunctionBegin; 955251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9568370ee3bSLisandro Dalcin if (iascii) { 9579c334d8fSLisandro Dalcin RKTableau tab = rk->tableau; 958f68a32c8SEmil Constantinescu TSRKType rktype; 959f68a32c8SEmil Constantinescu char buf[512]; 960f68a32c8SEmil Constantinescu ierr = TSRKGetType(ts,&rktype);CHKERRQ(ierr); 9610c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," RK: %s\n",rktype);CHKERRQ(ierr); 9620c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 9630c37eb8eSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," FSAL property: %s\n",tab->FSAL ? "yes" : "no");CHKERRQ(ierr); 964f68a32c8SEmil Constantinescu ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 965f68a32c8SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 9668370ee3bSLisandro Dalcin } 9679c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 968f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 969f68a32c8SEmil Constantinescu } 970f68a32c8SEmil Constantinescu 971f68a32c8SEmil Constantinescu static PetscErrorCode TSLoad_RK(TS ts,PetscViewer viewer) 972f68a32c8SEmil Constantinescu { 973f68a32c8SEmil Constantinescu PetscErrorCode ierr; 9749c334d8fSLisandro Dalcin TSAdapt adapt; 975f68a32c8SEmil Constantinescu 976f68a32c8SEmil Constantinescu PetscFunctionBegin; 9779c334d8fSLisandro Dalcin ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 9789c334d8fSLisandro Dalcin ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 979f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 980f68a32c8SEmil Constantinescu } 981f68a32c8SEmil Constantinescu 982f68a32c8SEmil Constantinescu /*@C 983f68a32c8SEmil Constantinescu TSRKSetType - Set the type of RK scheme 984f68a32c8SEmil Constantinescu 985f68a32c8SEmil Constantinescu Logically collective 986f68a32c8SEmil Constantinescu 987f68a32c8SEmil Constantinescu Input Parameter: 988f68a32c8SEmil Constantinescu + ts - timestepping context 989f68a32c8SEmil Constantinescu - rktype - type of RK-scheme 990f68a32c8SEmil Constantinescu 991*287c2655SBarry Smith Options Database: 992*287c2655SBarry Smith . -ts_tk_type - <1fe,2a,3,3bs,4,5f,5dp,5bs> 993*287c2655SBarry Smith 994f68a32c8SEmil Constantinescu Level: intermediate 995f68a32c8SEmil Constantinescu 996*287c2655SBarry Smith .seealso: TSRKGetType(), TSRK, TSRKType, TSRK1FE, TSRK2A, TSRK3, TSRK3BS, TSRK4, TSRK5F, TSRK5DP, TSRK5BS 997f68a32c8SEmil Constantinescu @*/ 998f68a32c8SEmil Constantinescu PetscErrorCode TSRKSetType(TS ts,TSRKType rktype) 999f68a32c8SEmil Constantinescu { 1000f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1001f68a32c8SEmil Constantinescu 1002f68a32c8SEmil Constantinescu PetscFunctionBegin; 1003f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1004b92453a8SLisandro Dalcin PetscValidCharPointer(rktype,2); 1005f68a32c8SEmil Constantinescu ierr = PetscTryMethod(ts,"TSRKSetType_C",(TS,TSRKType),(ts,rktype));CHKERRQ(ierr); 1006f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1007f68a32c8SEmil Constantinescu } 1008f68a32c8SEmil Constantinescu 1009f68a32c8SEmil Constantinescu /*@C 1010f68a32c8SEmil Constantinescu TSRKGetType - Get the type of RK scheme 1011f68a32c8SEmil Constantinescu 1012f68a32c8SEmil Constantinescu Logically collective 1013f68a32c8SEmil Constantinescu 1014f68a32c8SEmil Constantinescu Input Parameter: 1015f68a32c8SEmil Constantinescu . ts - timestepping context 1016f68a32c8SEmil Constantinescu 1017f68a32c8SEmil Constantinescu Output Parameter: 1018f68a32c8SEmil Constantinescu . rktype - type of RK-scheme 1019f68a32c8SEmil Constantinescu 1020f68a32c8SEmil Constantinescu Level: intermediate 1021f68a32c8SEmil Constantinescu 1022f68a32c8SEmil Constantinescu .seealso: TSRKGetType() 1023f68a32c8SEmil Constantinescu @*/ 1024f68a32c8SEmil Constantinescu PetscErrorCode TSRKGetType(TS ts,TSRKType *rktype) 1025f68a32c8SEmil Constantinescu { 1026f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1027f68a32c8SEmil Constantinescu 1028f68a32c8SEmil Constantinescu PetscFunctionBegin; 1029f68a32c8SEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1030f68a32c8SEmil Constantinescu ierr = PetscUseMethod(ts,"TSRKGetType_C",(TS,TSRKType*),(ts,rktype));CHKERRQ(ierr); 1031f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1032f68a32c8SEmil Constantinescu } 1033f68a32c8SEmil Constantinescu 1034560360afSLisandro Dalcin static PetscErrorCode TSRKGetType_RK(TS ts,TSRKType *rktype) 1035f68a32c8SEmil Constantinescu { 1036f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1037f68a32c8SEmil Constantinescu 1038f68a32c8SEmil Constantinescu PetscFunctionBegin; 1039f68a32c8SEmil Constantinescu *rktype = rk->tableau->name; 1040f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1041f68a32c8SEmil Constantinescu } 1042560360afSLisandro Dalcin static PetscErrorCode TSRKSetType_RK(TS ts,TSRKType rktype) 1043f68a32c8SEmil Constantinescu { 1044f68a32c8SEmil Constantinescu TS_RK *rk = (TS_RK*)ts->data; 1045f68a32c8SEmil Constantinescu PetscErrorCode ierr; 1046f68a32c8SEmil Constantinescu PetscBool match; 1047f68a32c8SEmil Constantinescu RKTableauLink link; 1048f68a32c8SEmil Constantinescu 1049f68a32c8SEmil Constantinescu PetscFunctionBegin; 1050f68a32c8SEmil Constantinescu if (rk->tableau) { 1051f68a32c8SEmil Constantinescu ierr = PetscStrcmp(rk->tableau->name,rktype,&match);CHKERRQ(ierr); 1052f68a32c8SEmil Constantinescu if (match) PetscFunctionReturn(0); 1053f68a32c8SEmil Constantinescu } 1054f68a32c8SEmil Constantinescu for (link = RKTableauList; link; link=link->next) { 1055f68a32c8SEmil Constantinescu ierr = PetscStrcmp(link->tab.name,rktype,&match);CHKERRQ(ierr); 1056f68a32c8SEmil Constantinescu if (match) { 1057be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauReset(ts);CHKERRQ(ierr);} 1058f68a32c8SEmil Constantinescu rk->tableau = &link->tab; 1059be5899b3SLisandro Dalcin if (ts->setupcalled) {ierr = TSRKTableauSetUp(ts);CHKERRQ(ierr);} 10602ffb9264SLisandro Dalcin ts->default_adapt_type = rk->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE; 1061f68a32c8SEmil Constantinescu PetscFunctionReturn(0); 1062f68a32c8SEmil Constantinescu } 1063f68a32c8SEmil Constantinescu } 1064f68a32c8SEmil Constantinescu SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",rktype); 1065e4dd521cSBarry Smith PetscFunctionReturn(0); 1066e4dd521cSBarry Smith } 1067e4dd521cSBarry Smith 1068ff22ae23SHong Zhang static PetscErrorCode TSGetStages_RK(TS ts,PetscInt *ns,Vec **Y) 1069ff22ae23SHong Zhang { 1070ff22ae23SHong Zhang TS_RK *rk = (TS_RK*)ts->data; 1071ff22ae23SHong Zhang 1072ff22ae23SHong Zhang PetscFunctionBegin; 1073ff22ae23SHong Zhang *ns = rk->tableau->s; 1074d2daff3dSHong Zhang if(Y) *Y = rk->Y; 1075ff22ae23SHong Zhang PetscFunctionReturn(0); 1076ff22ae23SHong Zhang } 1077ff22ae23SHong Zhang 1078ff22ae23SHong Zhang 1079e4dd521cSBarry Smith /* ------------------------------------------------------------ */ 1080ebe3b25bSBarry Smith /*MC 1081f68a32c8SEmil Constantinescu TSRK - ODE and DAE solver using Runge-Kutta schemes 1082ebe3b25bSBarry Smith 1083f68a32c8SEmil Constantinescu The user should provide the right hand side of the equation 1084f68a32c8SEmil Constantinescu using TSSetRHSFunction(). 1085ebe3b25bSBarry Smith 1086f68a32c8SEmil Constantinescu Notes: 108798164e13SLisandro Dalcin The default is TSRK3BS, it can be changed with TSRKSetType() or -ts_rk_type 1088ebe3b25bSBarry Smith 1089d41222bbSBarry Smith Level: beginner 1090d41222bbSBarry Smith 1091f68a32c8SEmil Constantinescu .seealso: TSCreate(), TS, TSSetType(), TSRKSetType(), TSRKGetType(), TSRKSetFullyImplicit(), TSRK2D, TTSRK2E, TSRK3, 1092f68a32c8SEmil Constantinescu TSRK4, TSRK5, TSRKPRSSP2, TSRKBPR3, TSRKType, TSRKRegister() 1093ebe3b25bSBarry Smith 1094ebe3b25bSBarry Smith M*/ 10958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_RK(TS ts) 1096e4dd521cSBarry Smith { 10971566a47fSLisandro Dalcin TS_RK *rk; 1098dfbe8321SBarry Smith PetscErrorCode ierr; 1099e4dd521cSBarry Smith 1100e4dd521cSBarry Smith PetscFunctionBegin; 1101f68a32c8SEmil Constantinescu ierr = TSRKInitializePackage();CHKERRQ(ierr); 1102f68a32c8SEmil Constantinescu 1103f68a32c8SEmil Constantinescu ts->ops->reset = TSReset_RK; 11045f70b478SJed Brown ts->ops->destroy = TSDestroy_RK; 11055f70b478SJed Brown ts->ops->view = TSView_RK; 1106f68a32c8SEmil Constantinescu ts->ops->load = TSLoad_RK; 1107f68a32c8SEmil Constantinescu ts->ops->setup = TSSetUp_RK; 110842f2b339SBarry Smith ts->ops->adjointsetup = TSAdjointSetUp_RK; 1109f68a32c8SEmil Constantinescu ts->ops->step = TSStep_RK; 1110f68a32c8SEmil Constantinescu ts->ops->interpolate = TSInterpolate_RK; 1111f68a32c8SEmil Constantinescu ts->ops->evaluatestep = TSEvaluateStep_RK; 1112fecfb714SLisandro Dalcin ts->ops->rollback = TSRollBack_RK; 1113f68a32c8SEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_RK; 1114ff22ae23SHong Zhang ts->ops->getstages = TSGetStages_RK; 111542f2b339SBarry Smith ts->ops->adjointstep = TSAdjointStep_RK; 1116e4dd521cSBarry Smith 11170f7a1166SHong Zhang ts->ops->adjointintegral = TSAdjointCostIntegral_RK; 11180f7a1166SHong Zhang ts->ops->forwardintegral = TSForwardCostIntegral_RK; 11190f7a1166SHong Zhang 11201566a47fSLisandro Dalcin ierr = PetscNewLog(ts,&rk);CHKERRQ(ierr); 11211566a47fSLisandro Dalcin ts->data = (void*)rk; 1122e4dd521cSBarry Smith 1123f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKGetType_C",TSRKGetType_RK);CHKERRQ(ierr); 1124f68a32c8SEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSRKSetType_C",TSRKSetType_RK);CHKERRQ(ierr); 1125be5899b3SLisandro Dalcin 1126be5899b3SLisandro Dalcin ierr = TSRKSetType(ts,TSRKDefault);CHKERRQ(ierr); 1127e4dd521cSBarry Smith PetscFunctionReturn(0); 1128e4dd521cSBarry Smith } 1129