12302aa1bSDebojyoti Ghosh /* 22302aa1bSDebojyoti Ghosh Code for time stepping with the General Linear with Error Estimation method 32302aa1bSDebojyoti Ghosh 42302aa1bSDebojyoti Ghosh Notes: 52302aa1bSDebojyoti Ghosh The general system is written as 62302aa1bSDebojyoti Ghosh 72302aa1bSDebojyoti Ghosh Udot = F(t,U) 82302aa1bSDebojyoti Ghosh 92302aa1bSDebojyoti Ghosh */ 102302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 112302aa1bSDebojyoti Ghosh #include <petscdm.h> 122302aa1bSDebojyoti Ghosh 13642f3702SEmil Constantinescu static PetscBool cited = PETSC_FALSE; 14642f3702SEmil Constantinescu static const char citation[] = 15642f3702SEmil Constantinescu "@ARTICLE{Constantinescu_TR2016b,\n" 16642f3702SEmil Constantinescu " author = {Constantinescu, E.M.},\n" 17642f3702SEmil Constantinescu " title = {Estimating Global Errors in Time Stepping},\n" 18642f3702SEmil Constantinescu " journal = {ArXiv e-prints},\n" 19642f3702SEmil Constantinescu " year = 2016,\n" 20642f3702SEmil Constantinescu " adsurl = {http://adsabs.harvard.edu/abs/2015arXiv150305166C}\n}\n"; 21642f3702SEmil Constantinescu 222302aa1bSDebojyoti Ghosh static TSGLEEType TSGLEEDefaultType = TSGLEE35; 232302aa1bSDebojyoti Ghosh static PetscBool TSGLEERegisterAllCalled; 242302aa1bSDebojyoti Ghosh static PetscBool TSGLEEPackageInitialized; 252302aa1bSDebojyoti Ghosh static PetscInt explicit_stage_time_id; 262302aa1bSDebojyoti Ghosh 272302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau; 282302aa1bSDebojyoti Ghosh struct _GLEETableau { 292302aa1bSDebojyoti Ghosh char *name; 302302aa1bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i*/ 312302aa1bSDebojyoti Ghosh PetscInt s; /* Number of stages */ 322302aa1bSDebojyoti Ghosh PetscInt r; /* Number of steps */ 332302aa1bSDebojyoti Ghosh PetscReal gamma; /* LTE ratio */ 342302aa1bSDebojyoti Ghosh PetscReal *A,*B,*U,*V,*S,*F,*c; /* Tableau */ 352302aa1bSDebojyoti Ghosh PetscReal *Fembed; /* Embedded final method coefficients */ 36fd96ebc2SDebojyoti Ghosh PetscReal *Ferror; /* Coefficients for computing error */ 3757df6a1bSDebojyoti Ghosh PetscReal *Serror; /* Coefficients for initializing the error */ 382302aa1bSDebojyoti Ghosh PetscInt pinterp; /* Interpolation order */ 392302aa1bSDebojyoti Ghosh PetscReal *binterp; /* Interpolation coefficients */ 402302aa1bSDebojyoti Ghosh PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 412302aa1bSDebojyoti Ghosh }; 422302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink; 432302aa1bSDebojyoti Ghosh struct _GLEETableauLink { 442302aa1bSDebojyoti Ghosh struct _GLEETableau tab; 452302aa1bSDebojyoti Ghosh GLEETableauLink next; 462302aa1bSDebojyoti Ghosh }; 472302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList; 482302aa1bSDebojyoti Ghosh 492302aa1bSDebojyoti Ghosh typedef struct { 502302aa1bSDebojyoti Ghosh GLEETableau tableau; 512302aa1bSDebojyoti Ghosh Vec *Y; /* Solution vector (along with auxiliary solution y~ or eps) */ 522302aa1bSDebojyoti Ghosh Vec *X; /* Temporary solution vector */ 532302aa1bSDebojyoti Ghosh Vec *YStage; /* Stage values */ 542302aa1bSDebojyoti Ghosh Vec *YdotStage; /* Stage right hand side */ 552302aa1bSDebojyoti Ghosh Vec W; /* Right-hand-side for implicit stage solve */ 562302aa1bSDebojyoti Ghosh Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 570a01e1b2SEmil Constantinescu Vec yGErr; /* Vector holding the global error after a step is completed */ 58ec0e3ee2SEmil Constantinescu PetscScalar *swork; /* Scalar work (size of the number of stages)*/ 59ec0e3ee2SEmil Constantinescu PetscScalar *rwork; /* Scalar work (size of the number of steps)*/ 602302aa1bSDebojyoti Ghosh PetscReal scoeff; /* shift = scoeff/dt */ 612302aa1bSDebojyoti Ghosh PetscReal stage_time; 622302aa1bSDebojyoti Ghosh TSStepStatus status; 632302aa1bSDebojyoti Ghosh } TS_GLEE; 642302aa1bSDebojyoti Ghosh 652302aa1bSDebojyoti Ghosh /*MC 662302aa1bSDebojyoti Ghosh TSGLEE23 - Second order three stage GLEE method 672302aa1bSDebojyoti Ghosh 682302aa1bSDebojyoti Ghosh This method has three stages. 692302aa1bSDebojyoti Ghosh s = 3, r = 2 702302aa1bSDebojyoti Ghosh 712302aa1bSDebojyoti Ghosh Level: advanced 722302aa1bSDebojyoti Ghosh 732302aa1bSDebojyoti Ghosh .seealso: TSGLEE 7436c34d14SEmil Constantinescu M*/ 752302aa1bSDebojyoti Ghosh /*MC 762302aa1bSDebojyoti Ghosh TSGLEE24 - Second order four stage GLEE method 772302aa1bSDebojyoti Ghosh 782302aa1bSDebojyoti Ghosh This method has four stages. 792302aa1bSDebojyoti Ghosh s = 4, r = 2 802302aa1bSDebojyoti Ghosh 812302aa1bSDebojyoti Ghosh Level: advanced 822302aa1bSDebojyoti Ghosh 832302aa1bSDebojyoti Ghosh .seealso: TSGLEE 8436c34d14SEmil Constantinescu M*/ 852302aa1bSDebojyoti Ghosh /*MC 862302aa1bSDebojyoti Ghosh TSGLEE25i - Second order five stage GLEE method 872302aa1bSDebojyoti Ghosh 882302aa1bSDebojyoti Ghosh This method has five stages. 892302aa1bSDebojyoti Ghosh s = 5, r = 2 902302aa1bSDebojyoti Ghosh 912302aa1bSDebojyoti Ghosh Level: advanced 922302aa1bSDebojyoti Ghosh 932302aa1bSDebojyoti Ghosh .seealso: TSGLEE 9436c34d14SEmil Constantinescu M*/ 952302aa1bSDebojyoti Ghosh /*MC 962302aa1bSDebojyoti Ghosh TSGLEE35 - Third order five stage GLEE method 972302aa1bSDebojyoti Ghosh 982302aa1bSDebojyoti Ghosh This method has five stages. 992302aa1bSDebojyoti Ghosh s = 5, r = 2 1002302aa1bSDebojyoti Ghosh 1012302aa1bSDebojyoti Ghosh Level: advanced 1022302aa1bSDebojyoti Ghosh 1032302aa1bSDebojyoti Ghosh .seealso: TSGLEE 10436c34d14SEmil Constantinescu M*/ 1052302aa1bSDebojyoti Ghosh /*MC 1062302aa1bSDebojyoti Ghosh TSGLEEEXRK2A - Second order six stage GLEE method 1072302aa1bSDebojyoti Ghosh 1082302aa1bSDebojyoti Ghosh This method has six stages. 1092302aa1bSDebojyoti Ghosh s = 6, r = 2 1102302aa1bSDebojyoti Ghosh 1112302aa1bSDebojyoti Ghosh Level: advanced 1122302aa1bSDebojyoti Ghosh 1132302aa1bSDebojyoti Ghosh .seealso: TSGLEE 11436c34d14SEmil Constantinescu M*/ 1152302aa1bSDebojyoti Ghosh /*MC 1162302aa1bSDebojyoti Ghosh TSGLEERK32G1 - Third order eight stage GLEE method 1172302aa1bSDebojyoti Ghosh 1182302aa1bSDebojyoti Ghosh This method has eight stages. 1192302aa1bSDebojyoti Ghosh s = 8, r = 2 1202302aa1bSDebojyoti Ghosh 1212302aa1bSDebojyoti Ghosh Level: advanced 1222302aa1bSDebojyoti Ghosh 1232302aa1bSDebojyoti Ghosh .seealso: TSGLEE 12436c34d14SEmil Constantinescu M*/ 1252302aa1bSDebojyoti Ghosh /*MC 1262302aa1bSDebojyoti Ghosh TSGLEERK285EX - Second order nine stage GLEE method 1272302aa1bSDebojyoti Ghosh 1282302aa1bSDebojyoti Ghosh This method has nine stages. 1292302aa1bSDebojyoti Ghosh s = 9, r = 2 1302302aa1bSDebojyoti Ghosh 1312302aa1bSDebojyoti Ghosh Level: advanced 1322302aa1bSDebojyoti Ghosh 1332302aa1bSDebojyoti Ghosh .seealso: TSGLEE 13436c34d14SEmil Constantinescu M*/ 1352302aa1bSDebojyoti Ghosh 1362302aa1bSDebojyoti Ghosh /*@C 1372302aa1bSDebojyoti Ghosh TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE 1382302aa1bSDebojyoti Ghosh 1392302aa1bSDebojyoti Ghosh Not Collective, but should be called by all processes which will need the schemes to be registered 1402302aa1bSDebojyoti Ghosh 1412302aa1bSDebojyoti Ghosh Level: advanced 1422302aa1bSDebojyoti Ghosh 1432302aa1bSDebojyoti Ghosh .seealso: TSGLEERegisterDestroy() 1442302aa1bSDebojyoti Ghosh @*/ 1452302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void) 1462302aa1bSDebojyoti Ghosh { 1472302aa1bSDebojyoti Ghosh PetscFunctionBegin; 1482302aa1bSDebojyoti Ghosh if (TSGLEERegisterAllCalled) PetscFunctionReturn(0); 1492302aa1bSDebojyoti Ghosh TSGLEERegisterAllCalled = PETSC_TRUE; 1502302aa1bSDebojyoti Ghosh 1512302aa1bSDebojyoti Ghosh { 152b1fbfd33SSatish Balay #define GAMMA 0.5 1532302aa1bSDebojyoti Ghosh /* y-eps form */ 1542302aa1bSDebojyoti Ghosh const PetscInt 155ab8c5912SEmil Constantinescu p = 1, 156ab8c5912SEmil Constantinescu s = 3, 157ab8c5912SEmil Constantinescu r = 2; 158ab8c5912SEmil Constantinescu const PetscReal 159ab8c5912SEmil Constantinescu A[3][3] = {{1.0,0,0},{0,0.5,0},{0,0.5,0.5}}, 160ab8c5912SEmil Constantinescu B[2][3] = {{1.0,0,0},{-2.0,1.0,1.0}}, 161ab8c5912SEmil Constantinescu U[3][2] = {{1.0,0},{1.0,0.5},{1.0,0.5}}, 162ab8c5912SEmil Constantinescu V[2][2] = {{1,0},{0,1}}, 163ab8c5912SEmil Constantinescu S[2] = {1,0}, 164ab8c5912SEmil Constantinescu F[2] = {1,0}, 165b1fbfd33SSatish Balay Fembed[2] = {1,1-GAMMA}, 16657df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 16757df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 1689566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEEi1,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 169ab8c5912SEmil Constantinescu } 170ab8c5912SEmil Constantinescu { 171b1fbfd33SSatish Balay #undef GAMMA 172b1fbfd33SSatish Balay #define GAMMA 0.0 173ab8c5912SEmil Constantinescu /* y-eps form */ 174ab8c5912SEmil Constantinescu const PetscInt 1752302aa1bSDebojyoti Ghosh p = 2, 1762302aa1bSDebojyoti Ghosh s = 3, 1772302aa1bSDebojyoti Ghosh r = 2; 1782302aa1bSDebojyoti Ghosh const PetscReal 1792302aa1bSDebojyoti Ghosh A[3][3] = {{0,0,0},{1,0,0},{0.25,0.25,0}}, 1802302aa1bSDebojyoti Ghosh B[2][3] = {{1.0/12.0,1.0/12.0,5.0/6.0},{1.0/12.0,1.0/12.0,-1.0/6.0}}, 1812302aa1bSDebojyoti Ghosh U[3][2] = {{1,0},{1,10},{1,-1}}, 1822302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 1832302aa1bSDebojyoti Ghosh S[2] = {1,0}, 1842302aa1bSDebojyoti Ghosh F[2] = {1,0}, 185b1fbfd33SSatish Balay Fembed[2] = {1,1-GAMMA}, 18657df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 18757df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 1889566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEE23,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 1892302aa1bSDebojyoti Ghosh } 1902302aa1bSDebojyoti Ghosh { 191b1fbfd33SSatish Balay #undef GAMMA 192b1fbfd33SSatish Balay #define GAMMA 0.0 1932302aa1bSDebojyoti Ghosh /* y-y~ form */ 1942302aa1bSDebojyoti Ghosh const PetscInt 1952302aa1bSDebojyoti Ghosh p = 2, 1962302aa1bSDebojyoti Ghosh s = 4, 1972302aa1bSDebojyoti Ghosh r = 2; 1982302aa1bSDebojyoti Ghosh const PetscReal 1992302aa1bSDebojyoti Ghosh A[4][4] = {{0,0,0,0},{0.75,0,0,0},{0.25,29.0/60.0,0,0},{-21.0/44.0,145.0/44.0,-20.0/11.0,0}}, 2002302aa1bSDebojyoti Ghosh B[2][4] = {{109.0/275.0,58.0/75.0,-37.0/110.0,1.0/6.0},{3.0/11.0,0,75.0/88.0,-1.0/8.0}}, 2012302aa1bSDebojyoti Ghosh U[4][2] = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}}, 2022302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2032302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2042302aa1bSDebojyoti Ghosh F[2] = {1,0}, 205fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 206b1fbfd33SSatish Balay Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)}, 207b1fbfd33SSatish Balay Serror[2] = {1.0-GAMMA,1.0}; 2089566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEE24,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 2092302aa1bSDebojyoti Ghosh } 2102302aa1bSDebojyoti Ghosh { 211b1fbfd33SSatish Balay #undef GAMMA 212b1fbfd33SSatish Balay #define GAMMA 0.0 2132302aa1bSDebojyoti Ghosh /* y-y~ form */ 2142302aa1bSDebojyoti Ghosh const PetscInt 2152302aa1bSDebojyoti Ghosh p = 2, 2162302aa1bSDebojyoti Ghosh s = 5, 2172302aa1bSDebojyoti Ghosh r = 2; 2182302aa1bSDebojyoti Ghosh const PetscReal 2192302aa1bSDebojyoti Ghosh A[5][5] = {{0,0,0,0,0}, 2202302aa1bSDebojyoti Ghosh {-0.94079244066783383269,0,0,0,0}, 2212302aa1bSDebojyoti Ghosh {0.64228187778301907108,0.10915356933958500042,0,0,0}, 2222302aa1bSDebojyoti Ghosh {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0}, 2232302aa1bSDebojyoti Ghosh {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}}, 2242302aa1bSDebojyoti Ghosh B[2][5] = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276}, 2252302aa1bSDebojyoti Ghosh {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}}, 2262302aa1bSDebojyoti Ghosh U[5][2] = {{0.16911424754448327735,0.83088575245551672265}, 2272302aa1bSDebojyoti Ghosh {0.53638465733199574340,0.46361534266800425660}, 2282302aa1bSDebojyoti Ghosh {0.39901579167169582526,0.60098420832830417474}, 2292302aa1bSDebojyoti Ghosh {0.87689005530618575480,0.12310994469381424520}, 2302302aa1bSDebojyoti Ghosh {0.99056100455550913009,0.0094389954444908699092}}, 2312302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2322302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2332302aa1bSDebojyoti Ghosh F[2] = {1,0}, 234fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 235b1fbfd33SSatish Balay Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)}, 236b1fbfd33SSatish Balay Serror[2] = {1.0-GAMMA,1.0}; 2379566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEE25I,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 2382302aa1bSDebojyoti Ghosh } 2392302aa1bSDebojyoti Ghosh { 240b1fbfd33SSatish Balay #undef GAMMA 241b1fbfd33SSatish Balay #define GAMMA 0.0 2422302aa1bSDebojyoti Ghosh /* y-y~ form */ 2432302aa1bSDebojyoti Ghosh const PetscInt 2442302aa1bSDebojyoti Ghosh p = 3, 2452302aa1bSDebojyoti Ghosh s = 5, 2462302aa1bSDebojyoti Ghosh r = 2; 2472302aa1bSDebojyoti Ghosh const PetscReal 2482302aa1bSDebojyoti Ghosh A[5][5] = {{0,0,0,0,0}, 2492302aa1bSDebojyoti Ghosh {- 2169604947363702313.0 / 24313474998937147335.0,0,0,0,0}, 2502302aa1bSDebojyoti Ghosh {46526746497697123895.0 / 94116917485856474137.0,-10297879244026594958.0 / 49199457603717988219.0,0,0,0}, 2512302aa1bSDebojyoti Ghosh {23364788935845982499.0 / 87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 / 36487615018004984309.0,0,0}, 25217d8b14cSSatish Balay {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 / 79257927066651693390.0,0}}, 2532302aa1bSDebojyoti Ghosh B[2][5] = {{61546696837458703723.0 / 56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 / 7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0}, 2542302aa1bSDebojyoti Ghosh {- 9738262186984159168.0 / 99299082461487742983.0,-32797097931948613195.0 / 61521565616362163366.0,42895514606418420631.0 / 71714201188501437336.0,22608567633166065068.0 / 55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}}, 2552302aa1bSDebojyoti Ghosh U[5][2] = {{70820309139834661559.0 / 80863923579509469826.0,10043614439674808267.0 / 80863923579509469826.0}, 2562302aa1bSDebojyoti Ghosh {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0}, 2572302aa1bSDebojyoti Ghosh {78486094644566264568.0 / 88171030896733822981.0,9684936252167558413.0 / 88171030896733822981.0}, 2582302aa1bSDebojyoti Ghosh {65394922146334854435.0 / 84570853840405479554.0,19175931694070625119.0 / 84570853840405479554.0}, 2592302aa1bSDebojyoti Ghosh {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}}, 2602302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2612302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2622302aa1bSDebojyoti Ghosh F[2] = {1,0}, 263fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 264b1fbfd33SSatish Balay Ferror[2] = {-1.0/(1.0-GAMMA),1.0/(1.0-GAMMA)}, 265b1fbfd33SSatish Balay Serror[2] = {1.0-GAMMA,1.0}; 2669566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEE35,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 2672302aa1bSDebojyoti Ghosh } 2682302aa1bSDebojyoti Ghosh { 269b1fbfd33SSatish Balay #undef GAMMA 270b1fbfd33SSatish Balay #define GAMMA 0.25 2712302aa1bSDebojyoti Ghosh /* y-eps form */ 2722302aa1bSDebojyoti Ghosh const PetscInt 2732302aa1bSDebojyoti Ghosh p = 2, 2742302aa1bSDebojyoti Ghosh s = 6, 2752302aa1bSDebojyoti Ghosh r = 2; 2762302aa1bSDebojyoti Ghosh const PetscReal 2772302aa1bSDebojyoti Ghosh A[6][6] = {{0,0,0,0,0,0}, 2782302aa1bSDebojyoti Ghosh {1,0,0,0,0,0}, 2792302aa1bSDebojyoti Ghosh {0,0,0,0,0,0}, 2802302aa1bSDebojyoti Ghosh {0,0,0.5,0,0,0}, 2812302aa1bSDebojyoti Ghosh {0,0,0.25,0.25,0,0}, 2822302aa1bSDebojyoti Ghosh {0,0,0.25,0.25,0.5,0}}, 2832302aa1bSDebojyoti Ghosh B[2][6] = {{0.5,0.5,0,0,0,0}, 2842302aa1bSDebojyoti Ghosh {-2.0/3.0,-2.0/3.0,1.0/3.0,1.0/3.0,1.0/3.0,1.0/3.0}}, 2852302aa1bSDebojyoti Ghosh U[6][2] = {{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75}}, 2862302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2872302aa1bSDebojyoti Ghosh S[2] = {1,0}, 2882302aa1bSDebojyoti Ghosh F[2] = {1,0}, 289b1fbfd33SSatish Balay Fembed[2] = {1,1-GAMMA}, 29057df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 29157df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 2929566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEEEXRK2A,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 2932302aa1bSDebojyoti Ghosh } 2942302aa1bSDebojyoti Ghosh { 295b1fbfd33SSatish Balay #undef GAMMA 296b1fbfd33SSatish Balay #define GAMMA 0.0 2972302aa1bSDebojyoti Ghosh /* y-eps form */ 2982302aa1bSDebojyoti Ghosh const PetscInt 2992302aa1bSDebojyoti Ghosh p = 3, 3002302aa1bSDebojyoti Ghosh s = 8, 3012302aa1bSDebojyoti Ghosh r = 2; 3022302aa1bSDebojyoti Ghosh const PetscReal 3032302aa1bSDebojyoti Ghosh A[8][8] = {{0,0,0,0,0,0,0,0}, 3042302aa1bSDebojyoti Ghosh {0.5,0,0,0,0,0,0,0}, 3052302aa1bSDebojyoti Ghosh {-1,2,0,0,0,0,0,0}, 3062302aa1bSDebojyoti Ghosh {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0}, 3072302aa1bSDebojyoti Ghosh {0,0,0,0,0,0,0,0}, 3082302aa1bSDebojyoti Ghosh {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0}, 3092302aa1bSDebojyoti Ghosh {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0}, 3102302aa1bSDebojyoti Ghosh {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}}, 3112302aa1bSDebojyoti Ghosh B[2][8] = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0}, 3122302aa1bSDebojyoti Ghosh {-1.0/6.0,-2.0/3.0,-1.0/6.0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}}, 3132302aa1bSDebojyoti Ghosh U[8][2] = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}}, 3142302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 3152302aa1bSDebojyoti Ghosh S[2] = {1,0}, 3162302aa1bSDebojyoti Ghosh F[2] = {1,0}, 317b1fbfd33SSatish Balay Fembed[2] = {1,1-GAMMA}, 31857df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 31957df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 3209566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEERK32G1,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 3212302aa1bSDebojyoti Ghosh } 3222302aa1bSDebojyoti Ghosh { 323b1fbfd33SSatish Balay #undef GAMMA 324b1fbfd33SSatish Balay #define GAMMA 0.25 3252302aa1bSDebojyoti Ghosh /* y-eps form */ 3262302aa1bSDebojyoti Ghosh const PetscInt 3272302aa1bSDebojyoti Ghosh p = 2, 3282302aa1bSDebojyoti Ghosh s = 9, 3292302aa1bSDebojyoti Ghosh r = 2; 3302302aa1bSDebojyoti Ghosh const PetscReal 3312302aa1bSDebojyoti Ghosh A[9][9] = {{0,0,0,0,0,0,0,0,0}, 3322302aa1bSDebojyoti Ghosh {0.585786437626904966,0,0,0,0,0,0,0,0}, 3332302aa1bSDebojyoti Ghosh {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0}, 3342302aa1bSDebojyoti Ghosh {0,0,0,0,0,0,0,0,0}, 3352302aa1bSDebojyoti Ghosh {0,0,0,0.292893218813452483,0,0,0,0,0}, 3362302aa1bSDebojyoti Ghosh {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0}, 3372302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0}, 3382302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0}, 3392302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}}, 3402302aa1bSDebojyoti Ghosh B[2][9] = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0}, 3412302aa1bSDebojyoti Ghosh {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}}, 3422302aa1bSDebojyoti Ghosh U[9][2] = {{1,0},{1,0},{1,0},{1,0.75},{1,0.75},{1,0.75},{1,0.75},{1,0.75},{1,0.75}}, 3432302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 3442302aa1bSDebojyoti Ghosh S[2] = {1,0}, 3452302aa1bSDebojyoti Ghosh F[2] = {1,0}, 346b1fbfd33SSatish Balay Fembed[2] = {1,1-GAMMA}, 34757df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 34857df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 3499566063dSJacob Faibussowitsch PetscCall(TSGLEERegister(TSGLEERK285EX,p,s,r,GAMMA,&A[0][0],&B[0][0],&U[0][0],&V[0][0],S,F,NULL,Fembed,Ferror,Serror,0,NULL)); 3502302aa1bSDebojyoti Ghosh } 3512302aa1bSDebojyoti Ghosh 3522302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 3532302aa1bSDebojyoti Ghosh } 3542302aa1bSDebojyoti Ghosh 3552302aa1bSDebojyoti Ghosh /*@C 3562302aa1bSDebojyoti Ghosh TSGLEERegisterDestroy - Frees the list of schemes that were registered by TSGLEERegister(). 3572302aa1bSDebojyoti Ghosh 3582302aa1bSDebojyoti Ghosh Not Collective 3592302aa1bSDebojyoti Ghosh 3602302aa1bSDebojyoti Ghosh Level: advanced 3612302aa1bSDebojyoti Ghosh 3622302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll() 3632302aa1bSDebojyoti Ghosh @*/ 3642302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void) 3652302aa1bSDebojyoti Ghosh { 3662302aa1bSDebojyoti Ghosh GLEETableauLink link; 3672302aa1bSDebojyoti Ghosh 3682302aa1bSDebojyoti Ghosh PetscFunctionBegin; 3692302aa1bSDebojyoti Ghosh while ((link = GLEETableauList)) { 3702302aa1bSDebojyoti Ghosh GLEETableau t = &link->tab; 3712302aa1bSDebojyoti Ghosh GLEETableauList = link->next; 3729566063dSJacob Faibussowitsch PetscCall(PetscFree5(t->A,t->B,t->U,t->V,t->c)); 373*d0609cedSBarry Smith PetscCall(PetscFree2(t->S,t->F)); 374*d0609cedSBarry Smith PetscCall(PetscFree (t->Fembed)); 375*d0609cedSBarry Smith PetscCall(PetscFree (t->Ferror)); 376*d0609cedSBarry Smith PetscCall(PetscFree (t->Serror)); 377*d0609cedSBarry Smith PetscCall(PetscFree (t->binterp)); 378*d0609cedSBarry Smith PetscCall(PetscFree (t->name)); 379*d0609cedSBarry Smith PetscCall(PetscFree (link)); 3802302aa1bSDebojyoti Ghosh } 3812302aa1bSDebojyoti Ghosh TSGLEERegisterAllCalled = PETSC_FALSE; 3822302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 3832302aa1bSDebojyoti Ghosh } 3842302aa1bSDebojyoti Ghosh 3852302aa1bSDebojyoti Ghosh /*@C 3862302aa1bSDebojyoti Ghosh TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called 3878a690491SBarry Smith from TSInitializePackage(). 3882302aa1bSDebojyoti Ghosh 3892302aa1bSDebojyoti Ghosh Level: developer 3902302aa1bSDebojyoti Ghosh 3912302aa1bSDebojyoti Ghosh .seealso: PetscInitialize() 3922302aa1bSDebojyoti Ghosh @*/ 3932302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void) 3942302aa1bSDebojyoti Ghosh { 3952302aa1bSDebojyoti Ghosh PetscFunctionBegin; 3962302aa1bSDebojyoti Ghosh if (TSGLEEPackageInitialized) PetscFunctionReturn(0); 3972302aa1bSDebojyoti Ghosh TSGLEEPackageInitialized = PETSC_TRUE; 3989566063dSJacob Faibussowitsch PetscCall(TSGLEERegisterAll()); 3999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataRegister(&explicit_stage_time_id)); 4009566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSGLEEFinalizePackage)); 4012302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 4022302aa1bSDebojyoti Ghosh } 4032302aa1bSDebojyoti Ghosh 4042302aa1bSDebojyoti Ghosh /*@C 4052302aa1bSDebojyoti Ghosh TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is 4062302aa1bSDebojyoti Ghosh called from PetscFinalize(). 4072302aa1bSDebojyoti Ghosh 4082302aa1bSDebojyoti Ghosh Level: developer 4092302aa1bSDebojyoti Ghosh 4102302aa1bSDebojyoti Ghosh .seealso: PetscFinalize() 4112302aa1bSDebojyoti Ghosh @*/ 4122302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void) 4132302aa1bSDebojyoti Ghosh { 4142302aa1bSDebojyoti Ghosh PetscFunctionBegin; 4152302aa1bSDebojyoti Ghosh TSGLEEPackageInitialized = PETSC_FALSE; 4169566063dSJacob Faibussowitsch PetscCall(TSGLEERegisterDestroy()); 4172302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 4182302aa1bSDebojyoti Ghosh } 4192302aa1bSDebojyoti Ghosh 4202302aa1bSDebojyoti Ghosh /*@C 4212302aa1bSDebojyoti Ghosh TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau 4222302aa1bSDebojyoti Ghosh 4232302aa1bSDebojyoti Ghosh Not Collective, but the same schemes should be registered on all processes on which they will be used 4242302aa1bSDebojyoti Ghosh 4252302aa1bSDebojyoti Ghosh Input Parameters: 4262302aa1bSDebojyoti Ghosh + name - identifier for method 4272302aa1bSDebojyoti Ghosh . order - order of method 4282302aa1bSDebojyoti Ghosh . s - number of stages 4292302aa1bSDebojyoti Ghosh . r - number of steps 4302302aa1bSDebojyoti Ghosh . gamma - LTE ratio 4312302aa1bSDebojyoti Ghosh . A - stage coefficients (dimension s*s, row-major) 4322302aa1bSDebojyoti Ghosh . B - step completion coefficients (dimension r*s, row-major) 4332302aa1bSDebojyoti Ghosh . U - method coefficients (dimension s*r, row-major) 4342302aa1bSDebojyoti Ghosh . V - method coefficients (dimension r*r, row-major) 4352302aa1bSDebojyoti Ghosh . S - starting coefficients 4362302aa1bSDebojyoti Ghosh . F - finishing coefficients 4372302aa1bSDebojyoti Ghosh . c - abscissa (dimension s; NULL to use row sums of A) 4382302aa1bSDebojyoti Ghosh . Fembed - step completion coefficients for embedded method 439fd96ebc2SDebojyoti Ghosh . Ferror - error computation coefficients 44057df6a1bSDebojyoti Ghosh . Serror - error initialization coefficients 4412302aa1bSDebojyoti Ghosh . pinterp - order of interpolation (0 if unavailable) 4422302aa1bSDebojyoti Ghosh - binterp - array of interpolation coefficients (NULL if unavailable) 4432302aa1bSDebojyoti Ghosh 4442302aa1bSDebojyoti Ghosh Notes: 4452302aa1bSDebojyoti Ghosh Several GLEE methods are provided, this function is only needed to create new methods. 4462302aa1bSDebojyoti Ghosh 4472302aa1bSDebojyoti Ghosh Level: advanced 4482302aa1bSDebojyoti Ghosh 4492302aa1bSDebojyoti Ghosh .seealso: TSGLEE 4502302aa1bSDebojyoti Ghosh @*/ 4512302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r, 4522302aa1bSDebojyoti Ghosh PetscReal gamma, 4532302aa1bSDebojyoti Ghosh const PetscReal A[],const PetscReal B[], 4542302aa1bSDebojyoti Ghosh const PetscReal U[],const PetscReal V[], 4552302aa1bSDebojyoti Ghosh const PetscReal S[],const PetscReal F[], 456fd96ebc2SDebojyoti Ghosh const PetscReal c[], 457fd96ebc2SDebojyoti Ghosh const PetscReal Fembed[],const PetscReal Ferror[], 45857df6a1bSDebojyoti Ghosh const PetscReal Serror[], 4592302aa1bSDebojyoti Ghosh PetscInt pinterp, const PetscReal binterp[]) 4602302aa1bSDebojyoti Ghosh { 4612302aa1bSDebojyoti Ghosh GLEETableauLink link; 4622302aa1bSDebojyoti Ghosh GLEETableau t; 4632302aa1bSDebojyoti Ghosh PetscInt i,j; 4642302aa1bSDebojyoti Ghosh 4652302aa1bSDebojyoti Ghosh PetscFunctionBegin; 4669566063dSJacob Faibussowitsch PetscCall(TSGLEEInitializePackage()); 4679566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 4682302aa1bSDebojyoti Ghosh t = &link->tab; 4699566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&t->name)); 4702302aa1bSDebojyoti Ghosh t->order = order; 4712302aa1bSDebojyoti Ghosh t->s = s; 4722302aa1bSDebojyoti Ghosh t->r = r; 4732302aa1bSDebojyoti Ghosh t->gamma = gamma; 4749566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U)); 4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(r,&t->S,r,&t->F)); 4769566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->A,A,s*s)); 4779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->B,B,r*s)); 4789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->U,U,s*r)); 4799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->V,V,r*r)); 4809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->S,S,r)); 4819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->F,F,r)); 4822302aa1bSDebojyoti Ghosh if (c) { 4839566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->c,c,s)); 4842302aa1bSDebojyoti Ghosh } else { 4852302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 4862302aa1bSDebojyoti Ghosh } 4879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(r,&t->Fembed)); 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(r,&t->Ferror)); 4899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(r,&t->Serror)); 4909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Fembed,Fembed,r)); 4919566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Ferror,Ferror,r)); 4929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Serror,Serror,r)); 4932302aa1bSDebojyoti Ghosh t->pinterp = pinterp; 4949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s*pinterp,&t->binterp)); 4959566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->binterp,binterp,s*pinterp)); 4962302aa1bSDebojyoti Ghosh 4972302aa1bSDebojyoti Ghosh link->next = GLEETableauList; 4982302aa1bSDebojyoti Ghosh GLEETableauList = link; 4992302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5002302aa1bSDebojyoti Ghosh } 5012302aa1bSDebojyoti Ghosh 5022302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done) 5032302aa1bSDebojyoti Ghosh { 5042302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*) ts->data; 5052302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 5062302aa1bSDebojyoti Ghosh PetscReal h, *B = tab->B, *V = tab->V, 507fd96ebc2SDebojyoti Ghosh *F = tab->F, 508fd96ebc2SDebojyoti Ghosh *Fembed = tab->Fembed; 5092302aa1bSDebojyoti Ghosh PetscInt s = tab->s, r = tab->r, i, j; 5102302aa1bSDebojyoti Ghosh Vec *Y = glee->Y, *YdotStage = glee->YdotStage; 511ec0e3ee2SEmil Constantinescu PetscScalar *ws = glee->swork, *wr = glee->rwork; 5122302aa1bSDebojyoti Ghosh 5132302aa1bSDebojyoti Ghosh PetscFunctionBegin; 5142302aa1bSDebojyoti Ghosh 5152302aa1bSDebojyoti Ghosh switch (glee->status) { 5162302aa1bSDebojyoti Ghosh case TS_STEP_INCOMPLETE: 5172302aa1bSDebojyoti Ghosh case TS_STEP_PENDING: 5182302aa1bSDebojyoti Ghosh h = ts->time_step; break; 5192302aa1bSDebojyoti Ghosh case TS_STEP_COMPLETE: 520ed450d42SEmil Constantinescu h = ts->ptime - ts->ptime_prev; break; 5212302aa1bSDebojyoti Ghosh default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 5222302aa1bSDebojyoti Ghosh } 5232302aa1bSDebojyoti Ghosh 5242302aa1bSDebojyoti Ghosh if (order == tab->order) { 5252302aa1bSDebojyoti Ghosh 5262302aa1bSDebojyoti Ghosh /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE 527fd96ebc2SDebojyoti Ghosh or TS_STEP_COMPLETE, glee->X has the solution at the 5282302aa1bSDebojyoti Ghosh beginning of the time step. So no need to roll-back. 5292302aa1bSDebojyoti Ghosh */ 5302302aa1bSDebojyoti Ghosh if (glee->status == TS_STEP_INCOMPLETE) { 5312302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 5329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y[i])); 533ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = V[i*r+j]; 5349566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],r,wr,glee->X)); 535ec0e3ee2SEmil Constantinescu for (j=0; j<s; j++) ws[j] = h*B[i*s+j]; 5369566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],s,ws,YdotStage)); 5372302aa1bSDebojyoti Ghosh } 5389566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 539ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = F[j]; 5409566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,r,wr,Y)); 5419566063dSJacob Faibussowitsch } else PetscCall(VecCopy(ts->vec_sol,X)); 5422302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5432302aa1bSDebojyoti Ghosh 5442302aa1bSDebojyoti Ghosh } else if (order == tab->order-1) { 5452302aa1bSDebojyoti Ghosh 5462302aa1bSDebojyoti Ghosh /* Complete with the embedded method (Fembed) */ 5472302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 5489566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Y[i])); 549ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = V[i*r+j]; 5509566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],r,wr,glee->X)); 551ec0e3ee2SEmil Constantinescu for (j=0; j<s; j++) ws[j] = h*B[i*s+j]; 5529566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y[i],s,ws,YdotStage)); 5532302aa1bSDebojyoti Ghosh } 5549566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(X)); 555ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = Fembed[j]; 5569566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,r,wr,Y)); 557fd96ebc2SDebojyoti Ghosh 5582302aa1bSDebojyoti Ghosh if (done) *done = PETSC_TRUE; 5592302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5602302aa1bSDebojyoti Ghosh } 5612302aa1bSDebojyoti Ghosh if (done) *done = PETSC_FALSE; 56298921bdaSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"GLEE '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 5632302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5642302aa1bSDebojyoti Ghosh } 5652302aa1bSDebojyoti Ghosh 5662302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts) 5672302aa1bSDebojyoti Ghosh { 5682302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 5692302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 5702302aa1bSDebojyoti Ghosh const PetscInt s = tab->s, r = tab->r; 5712302aa1bSDebojyoti Ghosh PetscReal *A = tab->A, *U = tab->U, 5722302aa1bSDebojyoti Ghosh *F = tab->F, 5732302aa1bSDebojyoti Ghosh *c = tab->c; 5742302aa1bSDebojyoti Ghosh Vec *Y = glee->Y, *X = glee->X, 5752302aa1bSDebojyoti Ghosh *YStage = glee->YStage, 5762302aa1bSDebojyoti Ghosh *YdotStage = glee->YdotStage, 5772302aa1bSDebojyoti Ghosh W = glee->W; 5782302aa1bSDebojyoti Ghosh SNES snes; 579ec0e3ee2SEmil Constantinescu PetscScalar *ws = glee->swork, *wr = glee->rwork; 5802302aa1bSDebojyoti Ghosh TSAdapt adapt; 5812302aa1bSDebojyoti Ghosh PetscInt i,j,reject,next_scheme,its,lits; 5822302aa1bSDebojyoti Ghosh PetscReal next_time_step; 5832302aa1bSDebojyoti Ghosh PetscReal t; 5842302aa1bSDebojyoti Ghosh PetscBool accept; 5852302aa1bSDebojyoti Ghosh 5862302aa1bSDebojyoti Ghosh PetscFunctionBegin; 5879566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation,&cited)); 588642f3702SEmil Constantinescu 5899566063dSJacob Faibussowitsch for (i=0; i<r; i++) PetscCall(VecCopy(Y[i],X[i])); 5902302aa1bSDebojyoti Ghosh 5919566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 5922302aa1bSDebojyoti Ghosh next_time_step = ts->time_step; 5932302aa1bSDebojyoti Ghosh t = ts->ptime; 5942302aa1bSDebojyoti Ghosh accept = PETSC_TRUE; 5952302aa1bSDebojyoti Ghosh glee->status = TS_STEP_INCOMPLETE; 5962302aa1bSDebojyoti Ghosh 5972302aa1bSDebojyoti Ghosh for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 5982302aa1bSDebojyoti Ghosh 5992302aa1bSDebojyoti Ghosh PetscReal h = ts->time_step; 6009566063dSJacob Faibussowitsch PetscCall(TSPreStep(ts)); 6012302aa1bSDebojyoti Ghosh 6022302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) { 6032302aa1bSDebojyoti Ghosh 6042302aa1bSDebojyoti Ghosh glee->stage_time = t + h*c[i]; 6059566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,glee->stage_time)); 6062302aa1bSDebojyoti Ghosh 6072302aa1bSDebojyoti Ghosh if (A[i*s+i] == 0) { /* Explicit stage */ 6089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(YStage[i])); 609ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = U[i*r+j]; 6109566063dSJacob Faibussowitsch PetscCall(VecMAXPY(YStage[i],r,wr,X)); 611ec0e3ee2SEmil Constantinescu for (j=0; j<i; j++) ws[j] = h*A[i*s+j]; 6129566063dSJacob Faibussowitsch PetscCall(VecMAXPY(YStage[i],i,ws,YdotStage)); 6132302aa1bSDebojyoti Ghosh } else { /* Implicit stage */ 6142302aa1bSDebojyoti Ghosh glee->scoeff = 1.0/A[i*s+i]; 6152302aa1bSDebojyoti Ghosh /* compute right-hand-side */ 6169566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(W)); 617ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = U[i*r+j]; 6189566063dSJacob Faibussowitsch PetscCall(VecMAXPY(W,r,wr,X)); 619ec0e3ee2SEmil Constantinescu for (j=0; j<i; j++) ws[j] = h*A[i*s+j]; 6209566063dSJacob Faibussowitsch PetscCall(VecMAXPY(W,i,ws,YdotStage)); 6219566063dSJacob Faibussowitsch PetscCall(VecScale(W,glee->scoeff/h)); 6222302aa1bSDebojyoti Ghosh /* set initial guess */ 6239566063dSJacob Faibussowitsch PetscCall(VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i])); 6242302aa1bSDebojyoti Ghosh /* solve for this stage */ 6259566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,W,YStage[i])); 6269566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 6279566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(snes,&lits)); 6282302aa1bSDebojyoti Ghosh ts->snes_its += its; ts->ksp_its += lits; 6292302aa1bSDebojyoti Ghosh } 6309566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 6319566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(adapt,ts,glee->stage_time,YStage[i],&accept)); 6322302aa1bSDebojyoti Ghosh if (!accept) goto reject_step; 6339566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,glee->stage_time,i,YStage)); 6349566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i])); 6352302aa1bSDebojyoti Ghosh } 6369566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL)); 6372302aa1bSDebojyoti Ghosh glee->status = TS_STEP_PENDING; 6382302aa1bSDebojyoti Ghosh 6392302aa1bSDebojyoti Ghosh /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 6409566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 6419566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(adapt)); 6429566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,(PetscReal)tab->s,PETSC_TRUE)); 6439566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept)); 6442302aa1bSDebojyoti Ghosh if (accept) { 6452302aa1bSDebojyoti Ghosh /* ignore next_scheme for now */ 6462302aa1bSDebojyoti Ghosh ts->ptime += ts->time_step; 6472302aa1bSDebojyoti Ghosh ts->time_step = next_time_step; 6482302aa1bSDebojyoti Ghosh glee->status = TS_STEP_COMPLETE; 6490a01e1b2SEmil Constantinescu /* compute and store the global error */ 6500a01e1b2SEmil Constantinescu /* Note: this is not needed if TSAdaptGLEE is not used */ 6519566063dSJacob Faibussowitsch PetscCall(TSGetTimeError(ts,0,&(glee->yGErr))); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime)); 6532302aa1bSDebojyoti Ghosh break; 6542302aa1bSDebojyoti Ghosh } else { /* Roll back the current step */ 655ec0e3ee2SEmil Constantinescu for (j=0; j<r; j++) wr[j] = F[j]; 6569566063dSJacob Faibussowitsch PetscCall(VecMAXPY(ts->vec_sol,r,wr,X)); 6572302aa1bSDebojyoti Ghosh ts->time_step = next_time_step; 6582302aa1bSDebojyoti Ghosh glee->status = TS_STEP_INCOMPLETE; 6592302aa1bSDebojyoti Ghosh } 6602302aa1bSDebojyoti Ghosh reject_step: continue; 6612302aa1bSDebojyoti Ghosh } 6622302aa1bSDebojyoti Ghosh if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 6632302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 6642302aa1bSDebojyoti Ghosh } 6652302aa1bSDebojyoti Ghosh 6662302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X) 6672302aa1bSDebojyoti Ghosh { 6682302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 6692302aa1bSDebojyoti Ghosh PetscInt s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j; 6702302aa1bSDebojyoti Ghosh PetscReal h,tt,t; 6712302aa1bSDebojyoti Ghosh PetscScalar *b; 6722302aa1bSDebojyoti Ghosh const PetscReal *B = glee->tableau->binterp; 6732302aa1bSDebojyoti Ghosh 6742302aa1bSDebojyoti Ghosh PetscFunctionBegin; 6753c633725SBarry Smith PetscCheck(B,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name); 6762302aa1bSDebojyoti Ghosh switch (glee->status) { 6772302aa1bSDebojyoti Ghosh case TS_STEP_INCOMPLETE: 6782302aa1bSDebojyoti Ghosh case TS_STEP_PENDING: 6792302aa1bSDebojyoti Ghosh h = ts->time_step; 6802302aa1bSDebojyoti Ghosh t = (itime - ts->ptime)/h; 6812302aa1bSDebojyoti Ghosh break; 6822302aa1bSDebojyoti Ghosh case TS_STEP_COMPLETE: 683ed450d42SEmil Constantinescu h = ts->ptime - ts->ptime_prev; 6842302aa1bSDebojyoti Ghosh t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 6852302aa1bSDebojyoti Ghosh break; 6862302aa1bSDebojyoti Ghosh default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 6872302aa1bSDebojyoti Ghosh } 6889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s,&b)); 6892302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) b[i] = 0; 6902302aa1bSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 6912302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) { 6922302aa1bSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 6932302aa1bSDebojyoti Ghosh } 6942302aa1bSDebojyoti Ghosh } 6959566063dSJacob Faibussowitsch PetscCall(VecCopy(glee->YStage[0],X)); 6969566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,b,glee->YdotStage)); 6979566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 6982302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 6992302aa1bSDebojyoti Ghosh } 7002302aa1bSDebojyoti Ghosh 7012302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/ 7022302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts) 7032302aa1bSDebojyoti Ghosh { 7042302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7052302aa1bSDebojyoti Ghosh PetscInt s, r; 7062302aa1bSDebojyoti Ghosh 7072302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7082302aa1bSDebojyoti Ghosh if (!glee->tableau) PetscFunctionReturn(0); 7092302aa1bSDebojyoti Ghosh s = glee->tableau->s; 7102302aa1bSDebojyoti Ghosh r = glee->tableau->r; 7119566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(r,&glee->Y)); 7129566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(r,&glee->X)); 7139566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(s,&glee->YStage)); 7149566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(s,&glee->YdotStage)); 7159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&glee->Ydot)); 7169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&glee->yGErr)); 7179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&glee->W)); 7189566063dSJacob Faibussowitsch PetscCall(PetscFree2(glee->swork,glee->rwork)); 7192302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7202302aa1bSDebojyoti Ghosh } 7212302aa1bSDebojyoti Ghosh 7222302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot) 7232302aa1bSDebojyoti Ghosh { 7242302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7252302aa1bSDebojyoti Ghosh 7262302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7272302aa1bSDebojyoti Ghosh if (Ydot) { 7282302aa1bSDebojyoti Ghosh if (dm && dm != ts->dm) { 7299566063dSJacob Faibussowitsch PetscCall(DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot)); 7302302aa1bSDebojyoti Ghosh } else *Ydot = glee->Ydot; 7312302aa1bSDebojyoti Ghosh } 7322302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7332302aa1bSDebojyoti Ghosh } 7342302aa1bSDebojyoti Ghosh 7352302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot) 7362302aa1bSDebojyoti Ghosh { 7372302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7382302aa1bSDebojyoti Ghosh if (Ydot) { 7392302aa1bSDebojyoti Ghosh if (dm && dm != ts->dm) { 7409566063dSJacob Faibussowitsch PetscCall(DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot)); 7412302aa1bSDebojyoti Ghosh } 7422302aa1bSDebojyoti Ghosh } 7432302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7442302aa1bSDebojyoti Ghosh } 7452302aa1bSDebojyoti Ghosh 7462302aa1bSDebojyoti Ghosh /* 7472302aa1bSDebojyoti Ghosh This defines the nonlinear equation that is to be solved with SNES 7482302aa1bSDebojyoti Ghosh */ 7492302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts) 7502302aa1bSDebojyoti Ghosh { 7512302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7522302aa1bSDebojyoti Ghosh DM dm,dmsave; 7532302aa1bSDebojyoti Ghosh Vec Ydot; 7542302aa1bSDebojyoti Ghosh PetscReal shift = glee->scoeff / ts->time_step; 7552302aa1bSDebojyoti Ghosh 7562302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7579566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 7589566063dSJacob Faibussowitsch PetscCall(TSGLEEGetVecs(ts,dm,&Ydot)); 7592302aa1bSDebojyoti Ghosh /* Set Ydot = shift*X */ 7609566063dSJacob Faibussowitsch PetscCall(VecCopy(X,Ydot)); 7619566063dSJacob Faibussowitsch PetscCall(VecScale(Ydot,shift)); 7622302aa1bSDebojyoti Ghosh dmsave = ts->dm; 7632302aa1bSDebojyoti Ghosh ts->dm = dm; 7642302aa1bSDebojyoti Ghosh 7659566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE)); 7662302aa1bSDebojyoti Ghosh 7672302aa1bSDebojyoti Ghosh ts->dm = dmsave; 7689566063dSJacob Faibussowitsch PetscCall(TSGLEERestoreVecs(ts,dm,&Ydot)); 7692302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7702302aa1bSDebojyoti Ghosh } 7712302aa1bSDebojyoti Ghosh 7722302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts) 7732302aa1bSDebojyoti Ghosh { 7742302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7752302aa1bSDebojyoti Ghosh DM dm,dmsave; 7762302aa1bSDebojyoti Ghosh Vec Ydot; 7772302aa1bSDebojyoti Ghosh PetscReal shift = glee->scoeff / ts->time_step; 7782302aa1bSDebojyoti Ghosh 7792302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7809566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 7819566063dSJacob Faibussowitsch PetscCall(TSGLEEGetVecs(ts,dm,&Ydot)); 7822302aa1bSDebojyoti Ghosh /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */ 7832302aa1bSDebojyoti Ghosh dmsave = ts->dm; 7842302aa1bSDebojyoti Ghosh ts->dm = dm; 7852302aa1bSDebojyoti Ghosh 7869566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE)); 7872302aa1bSDebojyoti Ghosh 7882302aa1bSDebojyoti Ghosh ts->dm = dmsave; 7899566063dSJacob Faibussowitsch PetscCall(TSGLEERestoreVecs(ts,dm,&Ydot)); 7902302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7912302aa1bSDebojyoti Ghosh } 7922302aa1bSDebojyoti Ghosh 7932302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx) 7942302aa1bSDebojyoti Ghosh { 7952302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7962302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7972302aa1bSDebojyoti Ghosh } 7982302aa1bSDebojyoti Ghosh 7992302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 8002302aa1bSDebojyoti Ghosh { 8012302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8022302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8032302aa1bSDebojyoti Ghosh } 8042302aa1bSDebojyoti Ghosh 8052302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx) 8062302aa1bSDebojyoti Ghosh { 8072302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8082302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8092302aa1bSDebojyoti Ghosh } 8102302aa1bSDebojyoti Ghosh 8112302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 8122302aa1bSDebojyoti Ghosh { 8132302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8142302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8152302aa1bSDebojyoti Ghosh } 8162302aa1bSDebojyoti Ghosh 8172302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts) 8182302aa1bSDebojyoti Ghosh { 8192302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 8202302aa1bSDebojyoti Ghosh GLEETableau tab; 8212302aa1bSDebojyoti Ghosh PetscInt s,r; 8222302aa1bSDebojyoti Ghosh DM dm; 8232302aa1bSDebojyoti Ghosh 8242302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8252302aa1bSDebojyoti Ghosh if (!glee->tableau) { 8269566063dSJacob Faibussowitsch PetscCall(TSGLEESetType(ts,TSGLEEDefaultType)); 8272302aa1bSDebojyoti Ghosh } 8282302aa1bSDebojyoti Ghosh tab = glee->tableau; 8292302aa1bSDebojyoti Ghosh s = tab->s; 8302302aa1bSDebojyoti Ghosh r = tab->r; 8319566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,r,&glee->Y)); 8329566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,r,&glee->X)); 8339566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,s,&glee->YStage)); 8349566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage)); 8359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&glee->Ydot)); 8369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&glee->yGErr)); 8379566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(glee->yGErr)); 8389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&glee->W)); 8399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(s,&glee->swork,r,&glee->rwork)); 8409566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 8419566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts)); 8429566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts)); 8432302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8442302aa1bSDebojyoti Ghosh } 8452302aa1bSDebojyoti Ghosh 8462302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts) 8472302aa1bSDebojyoti Ghosh { 8482302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 8496e2e232bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 8506e2e232bSDebojyoti Ghosh PetscInt r=tab->r,i; 8516e2e232bSDebojyoti Ghosh PetscReal *S=tab->S; 8522302aa1bSDebojyoti Ghosh 8532302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8542302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 8559566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(glee->Y[i])); 8569566063dSJacob Faibussowitsch PetscCall(VecAXPY(glee->Y[i],S[i],ts->vec_sol)); 8572302aa1bSDebojyoti Ghosh } 8582302aa1bSDebojyoti Ghosh 8592302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8602302aa1bSDebojyoti Ghosh } 8612302aa1bSDebojyoti Ghosh 8622302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/ 8632302aa1bSDebojyoti Ghosh 8646b234c65SEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TS ts) 8652302aa1bSDebojyoti Ghosh { 8662302aa1bSDebojyoti Ghosh char gleetype[256]; 8672302aa1bSDebojyoti Ghosh 8682302aa1bSDebojyoti Ghosh PetscFunctionBegin; 869*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"GLEE ODE solver options"); 8702302aa1bSDebojyoti Ghosh { 8712302aa1bSDebojyoti Ghosh GLEETableauLink link; 8722302aa1bSDebojyoti Ghosh PetscInt count,choice; 8732302aa1bSDebojyoti Ghosh PetscBool flg; 8742302aa1bSDebojyoti Ghosh const char **namelist; 8752302aa1bSDebojyoti Ghosh 8769566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype))); 8772302aa1bSDebojyoti Ghosh for (link=GLEETableauList,count=0; link; link=link->next,count++) ; 8789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count,(char***)&namelist)); 8792302aa1bSDebojyoti Ghosh for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 8809566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg)); 8819566063dSJacob Faibussowitsch PetscCall(TSGLEESetType(ts,flg ? namelist[choice] : gleetype)); 8829566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 8832302aa1bSDebojyoti Ghosh } 884*d0609cedSBarry Smith PetscOptionsHeadEnd(); 8852302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8862302aa1bSDebojyoti Ghosh } 8872302aa1bSDebojyoti Ghosh 8882302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer) 8892302aa1bSDebojyoti Ghosh { 8902302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 8912302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 8922302aa1bSDebojyoti Ghosh PetscBool iascii; 8932302aa1bSDebojyoti Ghosh 8942302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 8962302aa1bSDebojyoti Ghosh if (iascii) { 8972302aa1bSDebojyoti Ghosh TSGLEEType gleetype; 8982302aa1bSDebojyoti Ghosh char buf[512]; 8999566063dSJacob Faibussowitsch PetscCall(TSGLEEGetType(ts,&gleetype)); 9009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," GLEE type %s\n",gleetype)); 9019566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c)); 9029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf)); 9032302aa1bSDebojyoti Ghosh /* Note: print out r as well */ 9042302aa1bSDebojyoti Ghosh } 9052302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9062302aa1bSDebojyoti Ghosh } 9072302aa1bSDebojyoti Ghosh 9082302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer) 9092302aa1bSDebojyoti Ghosh { 9102302aa1bSDebojyoti Ghosh SNES snes; 9112302aa1bSDebojyoti Ghosh TSAdapt tsadapt; 9122302aa1bSDebojyoti Ghosh 9132302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9149566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&tsadapt)); 9159566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(tsadapt,viewer)); 9169566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 9179566063dSJacob Faibussowitsch PetscCall(SNESLoad(snes,viewer)); 9182302aa1bSDebojyoti Ghosh /* function and Jacobian context for SNES when used with TS is always ts object */ 9199566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,NULL,NULL,ts)); 9209566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,NULL,NULL,NULL,ts)); 9212302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9222302aa1bSDebojyoti Ghosh } 9232302aa1bSDebojyoti Ghosh 9242302aa1bSDebojyoti Ghosh /*@C 9252302aa1bSDebojyoti Ghosh TSGLEESetType - Set the type of GLEE scheme 9262302aa1bSDebojyoti Ghosh 9272302aa1bSDebojyoti Ghosh Logically collective 9282302aa1bSDebojyoti Ghosh 929d8d19677SJose E. Roman Input Parameters: 9302302aa1bSDebojyoti Ghosh + ts - timestepping context 9312302aa1bSDebojyoti Ghosh - gleetype - type of GLEE-scheme 9322302aa1bSDebojyoti Ghosh 9332302aa1bSDebojyoti Ghosh Level: intermediate 9342302aa1bSDebojyoti Ghosh 9352302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE 9362302aa1bSDebojyoti Ghosh @*/ 9372302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype) 9382302aa1bSDebojyoti Ghosh { 9392302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9402302aa1bSDebojyoti Ghosh PetscValidHeaderSpecific(ts,TS_CLASSID,1); 941b92453a8SLisandro Dalcin PetscValidCharPointer(gleetype,2); 942cac4c232SBarry Smith PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype)); 9432302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9442302aa1bSDebojyoti Ghosh } 9452302aa1bSDebojyoti Ghosh 9462302aa1bSDebojyoti Ghosh /*@C 9472302aa1bSDebojyoti Ghosh TSGLEEGetType - Get the type of GLEE scheme 9482302aa1bSDebojyoti Ghosh 9492302aa1bSDebojyoti Ghosh Logically collective 9502302aa1bSDebojyoti Ghosh 9512302aa1bSDebojyoti Ghosh Input Parameter: 9522302aa1bSDebojyoti Ghosh . ts - timestepping context 9532302aa1bSDebojyoti Ghosh 9542302aa1bSDebojyoti Ghosh Output Parameter: 9552302aa1bSDebojyoti Ghosh . gleetype - type of GLEE-scheme 9562302aa1bSDebojyoti Ghosh 9572302aa1bSDebojyoti Ghosh Level: intermediate 9582302aa1bSDebojyoti Ghosh 9592302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType() 9602302aa1bSDebojyoti Ghosh @*/ 9612302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype) 9622302aa1bSDebojyoti Ghosh { 9632302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9642302aa1bSDebojyoti Ghosh PetscValidHeaderSpecific(ts,TS_CLASSID,1); 965cac4c232SBarry Smith PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype)); 9662302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9672302aa1bSDebojyoti Ghosh } 9682302aa1bSDebojyoti Ghosh 9692302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype) 9702302aa1bSDebojyoti Ghosh { 9712302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 9722302aa1bSDebojyoti Ghosh 9732302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9742302aa1bSDebojyoti Ghosh if (!glee->tableau) { 9759566063dSJacob Faibussowitsch PetscCall(TSGLEESetType(ts,TSGLEEDefaultType)); 9762302aa1bSDebojyoti Ghosh } 9772302aa1bSDebojyoti Ghosh *gleetype = glee->tableau->name; 9782302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9792302aa1bSDebojyoti Ghosh } 9802302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype) 9812302aa1bSDebojyoti Ghosh { 9822302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 9832302aa1bSDebojyoti Ghosh PetscBool match; 9842302aa1bSDebojyoti Ghosh GLEETableauLink link; 9852302aa1bSDebojyoti Ghosh 9862302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9872302aa1bSDebojyoti Ghosh if (glee->tableau) { 9889566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(glee->tableau->name,gleetype,&match)); 9892302aa1bSDebojyoti Ghosh if (match) PetscFunctionReturn(0); 9902302aa1bSDebojyoti Ghosh } 9912302aa1bSDebojyoti Ghosh for (link = GLEETableauList; link; link=link->next) { 9929566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name,gleetype,&match)); 9932302aa1bSDebojyoti Ghosh if (match) { 9949566063dSJacob Faibussowitsch PetscCall(TSReset_GLEE(ts)); 9952302aa1bSDebojyoti Ghosh glee->tableau = &link->tab; 9962302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9972302aa1bSDebojyoti Ghosh } 9982302aa1bSDebojyoti Ghosh } 99998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype); 10002302aa1bSDebojyoti Ghosh } 10012302aa1bSDebojyoti Ghosh 10022302aa1bSDebojyoti Ghosh static PetscErrorCode TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y) 10032302aa1bSDebojyoti Ghosh { 10042302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10052302aa1bSDebojyoti Ghosh 10062302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10070429704eSStefano Zampini if (ns) *ns = glee->tableau->s; 1008d5509560SEmil Constantinescu if (Y) *Y = glee->YStage; 10092302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10102302aa1bSDebojyoti Ghosh } 10112302aa1bSDebojyoti Ghosh 1012b2bf4f3aSDebojyoti Ghosh PetscErrorCode TSGetSolutionComponents_GLEE(TS ts,PetscInt *n,Vec *Y) 10132302aa1bSDebojyoti Ghosh { 10142302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10152302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 10162302aa1bSDebojyoti Ghosh 10172302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10184cdd57e5SDebojyoti Ghosh if (!Y) *n = tab->r; 10192302aa1bSDebojyoti Ghosh else { 10204cdd57e5SDebojyoti Ghosh if ((*n >= 0) && (*n < tab->r)) { 10219566063dSJacob Faibussowitsch PetscCall(VecCopy(glee->Y[*n],*Y)); 102298921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Second argument (%d) out of range[%d,%d].",*n,0,tab->r-1); 10232302aa1bSDebojyoti Ghosh } 10242302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10252302aa1bSDebojyoti Ghosh } 10262302aa1bSDebojyoti Ghosh 10274cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,Vec *X) 10284cdd57e5SDebojyoti Ghosh { 10294cdd57e5SDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10304cdd57e5SDebojyoti Ghosh GLEETableau tab = glee->tableau; 10314cdd57e5SDebojyoti Ghosh PetscReal *F = tab->Fembed; 10324cdd57e5SDebojyoti Ghosh PetscInt r = tab->r; 10334cdd57e5SDebojyoti Ghosh Vec *Y = glee->Y; 1034ec0e3ee2SEmil Constantinescu PetscScalar *wr = glee->rwork; 1035ec0e3ee2SEmil Constantinescu PetscInt i; 10364cdd57e5SDebojyoti Ghosh 10374cdd57e5SDebojyoti Ghosh PetscFunctionBegin; 10389566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(*X)); 1039ec0e3ee2SEmil Constantinescu for (i=0; i<r; i++) wr[i] = F[i]; 10409566063dSJacob Faibussowitsch PetscCall(VecMAXPY((*X),r,wr,Y)); 10414cdd57e5SDebojyoti Ghosh PetscFunctionReturn(0); 10424cdd57e5SDebojyoti Ghosh } 10434cdd57e5SDebojyoti Ghosh 10440a01e1b2SEmil Constantinescu PetscErrorCode TSGetTimeError_GLEE(TS ts,PetscInt n,Vec *X) 10454cdd57e5SDebojyoti Ghosh { 10464cdd57e5SDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10474cdd57e5SDebojyoti Ghosh GLEETableau tab = glee->tableau; 10484cdd57e5SDebojyoti Ghosh PetscReal *F = tab->Ferror; 10494cdd57e5SDebojyoti Ghosh PetscInt r = tab->r; 10504cdd57e5SDebojyoti Ghosh Vec *Y = glee->Y; 1051ec0e3ee2SEmil Constantinescu PetscScalar *wr = glee->rwork; 1052ec0e3ee2SEmil Constantinescu PetscInt i; 10534cdd57e5SDebojyoti Ghosh 10544cdd57e5SDebojyoti Ghosh PetscFunctionBegin; 10559566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(*X)); 1056657c1e31SEmil Constantinescu if (n==0) { 1057ec0e3ee2SEmil Constantinescu for (i=0; i<r; i++) wr[i] = F[i]; 10589566063dSJacob Faibussowitsch PetscCall(VecMAXPY((*X),r,wr,Y)); 10590a01e1b2SEmil Constantinescu } else if (n==-1) { 1060657c1e31SEmil Constantinescu *X=glee->yGErr; 10610a01e1b2SEmil Constantinescu } 10624cdd57e5SDebojyoti Ghosh PetscFunctionReturn(0); 10634cdd57e5SDebojyoti Ghosh } 10644cdd57e5SDebojyoti Ghosh 106557df6a1bSDebojyoti Ghosh PetscErrorCode TSSetTimeError_GLEE(TS ts,Vec X) 106657df6a1bSDebojyoti Ghosh { 106757df6a1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 106857df6a1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 106957df6a1bSDebojyoti Ghosh PetscReal *S = tab->Serror; 107057df6a1bSDebojyoti Ghosh PetscInt r = tab->r,i; 107157df6a1bSDebojyoti Ghosh Vec *Y = glee->Y; 107257df6a1bSDebojyoti Ghosh 107357df6a1bSDebojyoti Ghosh PetscFunctionBegin; 10743c633725SBarry Smith PetscCheck(r == 2,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSSetTimeError_GLEE not supported for '%s' with r=%D.",tab->name,tab->r); 107557df6a1bSDebojyoti Ghosh for (i=1; i<r; i++) { 10769566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,Y[i])); 10779566063dSJacob Faibussowitsch PetscCall(VecAXPBY(Y[i],S[0],S[1],X)); 10789566063dSJacob Faibussowitsch PetscCall(VecCopy(X,glee->yGErr)); 107957df6a1bSDebojyoti Ghosh } 108057df6a1bSDebojyoti Ghosh PetscFunctionReturn(0); 108157df6a1bSDebojyoti Ghosh } 108257df6a1bSDebojyoti Ghosh 1083b3a6b972SJed Brown static PetscErrorCode TSDestroy_GLEE(TS ts) 1084b3a6b972SJed Brown { 1085b3a6b972SJed Brown PetscFunctionBegin; 10869566063dSJacob Faibussowitsch PetscCall(TSReset_GLEE(ts)); 1087b3a6b972SJed Brown if (ts->dm) { 10889566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts)); 10899566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts)); 1090b3a6b972SJed Brown } 10919566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 10929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL)); 10939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL)); 1094b3a6b972SJed Brown PetscFunctionReturn(0); 1095b3a6b972SJed Brown } 1096b3a6b972SJed Brown 10972302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */ 10982302aa1bSDebojyoti Ghosh /*MC 10992302aa1bSDebojyoti Ghosh TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes 11002302aa1bSDebojyoti Ghosh 11012302aa1bSDebojyoti Ghosh The user should provide the right hand side of the equation 11022302aa1bSDebojyoti Ghosh using TSSetRHSFunction(). 11032302aa1bSDebojyoti Ghosh 11042302aa1bSDebojyoti Ghosh Notes: 11052302aa1bSDebojyoti Ghosh The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type 11062302aa1bSDebojyoti Ghosh 11072302aa1bSDebojyoti Ghosh Level: beginner 11082302aa1bSDebojyoti Ghosh 11092302aa1bSDebojyoti Ghosh .seealso: TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(), 11102302aa1bSDebojyoti Ghosh TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A, 11112302aa1bSDebojyoti Ghosh TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister() 11122302aa1bSDebojyoti Ghosh 11132302aa1bSDebojyoti Ghosh M*/ 11142302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts) 11152302aa1bSDebojyoti Ghosh { 11162302aa1bSDebojyoti Ghosh TS_GLEE *th; 11172302aa1bSDebojyoti Ghosh 11182302aa1bSDebojyoti Ghosh PetscFunctionBegin; 11199566063dSJacob Faibussowitsch PetscCall(TSGLEEInitializePackage()); 11202302aa1bSDebojyoti Ghosh 11212302aa1bSDebojyoti Ghosh ts->ops->reset = TSReset_GLEE; 11222302aa1bSDebojyoti Ghosh ts->ops->destroy = TSDestroy_GLEE; 11232302aa1bSDebojyoti Ghosh ts->ops->view = TSView_GLEE; 11242302aa1bSDebojyoti Ghosh ts->ops->load = TSLoad_GLEE; 11252302aa1bSDebojyoti Ghosh ts->ops->setup = TSSetUp_GLEE; 11262302aa1bSDebojyoti Ghosh ts->ops->step = TSStep_GLEE; 11272302aa1bSDebojyoti Ghosh ts->ops->interpolate = TSInterpolate_GLEE; 11282302aa1bSDebojyoti Ghosh ts->ops->evaluatestep = TSEvaluateStep_GLEE; 11292302aa1bSDebojyoti Ghosh ts->ops->setfromoptions = TSSetFromOptions_GLEE; 11302302aa1bSDebojyoti Ghosh ts->ops->getstages = TSGetStages_GLEE; 11312302aa1bSDebojyoti Ghosh ts->ops->snesfunction = SNESTSFormFunction_GLEE; 11322302aa1bSDebojyoti Ghosh ts->ops->snesjacobian = SNESTSFormJacobian_GLEE; 1133b2bf4f3aSDebojyoti Ghosh ts->ops->getsolutioncomponents = TSGetSolutionComponents_GLEE; 11344cdd57e5SDebojyoti Ghosh ts->ops->getauxsolution = TSGetAuxSolution_GLEE; 11354cdd57e5SDebojyoti Ghosh ts->ops->gettimeerror = TSGetTimeError_GLEE; 113657df6a1bSDebojyoti Ghosh ts->ops->settimeerror = TSSetTimeError_GLEE; 11372302aa1bSDebojyoti Ghosh ts->ops->startingmethod = TSStartingMethod_GLEE; 1138b80d5313SLisandro Dalcin ts->default_adapt_type = TSADAPTGLEE; 11392302aa1bSDebojyoti Ghosh 1140efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 1141efd4aadfSBarry Smith 11429566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&th)); 11432302aa1bSDebojyoti Ghosh ts->data = (void*)th; 11442302aa1bSDebojyoti Ghosh 11459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE)); 11469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE)); 11472302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 11482302aa1bSDebojyoti Ghosh } 1149