12302aa1bSDebojyoti Ghosh /* 22302aa1bSDebojyoti Ghosh Code for time stepping with the General Linear with Error Estimation method 32302aa1bSDebojyoti Ghosh 4*642f3702SEmil Constantinescu 52302aa1bSDebojyoti Ghosh Notes: 62302aa1bSDebojyoti Ghosh The general system is written as 72302aa1bSDebojyoti Ghosh 82302aa1bSDebojyoti Ghosh Udot = F(t,U) 92302aa1bSDebojyoti Ghosh 102302aa1bSDebojyoti Ghosh */ 112302aa1bSDebojyoti Ghosh #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 122302aa1bSDebojyoti Ghosh #include <petscdm.h> 132302aa1bSDebojyoti Ghosh 14*642f3702SEmil Constantinescu static PetscBool cited = PETSC_FALSE; 15*642f3702SEmil Constantinescu static const char citation[] = 16*642f3702SEmil Constantinescu "@ARTICLE{Constantinescu_TR2016b,\n" 17*642f3702SEmil Constantinescu " author = {Constantinescu, E.M.},\n" 18*642f3702SEmil Constantinescu " title = {Estimating Global Errors in Time Stepping},\n" 19*642f3702SEmil Constantinescu " journal = {ArXiv e-prints},\n" 20*642f3702SEmil Constantinescu " year = 2016,\n" 21*642f3702SEmil Constantinescu " adsurl = {http://adsabs.harvard.edu/abs/2015arXiv150305166C}\n}\n"; 22*642f3702SEmil Constantinescu 232302aa1bSDebojyoti Ghosh static TSGLEEType TSGLEEDefaultType = TSGLEE35; 242302aa1bSDebojyoti Ghosh static PetscBool TSGLEERegisterAllCalled; 252302aa1bSDebojyoti Ghosh static PetscBool TSGLEEPackageInitialized; 262302aa1bSDebojyoti Ghosh static PetscInt explicit_stage_time_id; 272302aa1bSDebojyoti Ghosh 282302aa1bSDebojyoti Ghosh typedef struct _GLEETableau *GLEETableau; 292302aa1bSDebojyoti Ghosh struct _GLEETableau { 302302aa1bSDebojyoti Ghosh char *name; 312302aa1bSDebojyoti Ghosh PetscInt order; /* Classical approximation order of the method i*/ 322302aa1bSDebojyoti Ghosh PetscInt s; /* Number of stages */ 332302aa1bSDebojyoti Ghosh PetscInt r; /* Number of steps */ 342302aa1bSDebojyoti Ghosh PetscReal gamma; /* LTE ratio */ 352302aa1bSDebojyoti Ghosh PetscReal *A,*B,*U,*V,*S,*F,*c; /* Tableau */ 362302aa1bSDebojyoti Ghosh PetscReal *Fembed; /* Embedded final method coefficients */ 37fd96ebc2SDebojyoti Ghosh PetscReal *Ferror; /* Coefficients for computing error */ 3857df6a1bSDebojyoti Ghosh PetscReal *Serror; /* Coefficients for initializing the error */ 392302aa1bSDebojyoti Ghosh PetscInt pinterp; /* Interpolation order */ 402302aa1bSDebojyoti Ghosh PetscReal *binterp; /* Interpolation coefficients */ 412302aa1bSDebojyoti Ghosh PetscReal ccfl; /* Placeholder for CFL coefficient relative to forward Euler */ 422302aa1bSDebojyoti Ghosh }; 432302aa1bSDebojyoti Ghosh typedef struct _GLEETableauLink *GLEETableauLink; 442302aa1bSDebojyoti Ghosh struct _GLEETableauLink { 452302aa1bSDebojyoti Ghosh struct _GLEETableau tab; 462302aa1bSDebojyoti Ghosh GLEETableauLink next; 472302aa1bSDebojyoti Ghosh }; 482302aa1bSDebojyoti Ghosh static GLEETableauLink GLEETableauList; 492302aa1bSDebojyoti Ghosh 502302aa1bSDebojyoti Ghosh typedef struct { 512302aa1bSDebojyoti Ghosh GLEETableau tableau; 522302aa1bSDebojyoti Ghosh Vec *Y; /* Solution vector (along with auxiliary solution y~ or eps) */ 532302aa1bSDebojyoti Ghosh Vec *X; /* Temporary solution vector */ 542302aa1bSDebojyoti Ghosh Vec *YStage; /* Stage values */ 552302aa1bSDebojyoti Ghosh Vec *YdotStage; /* Stage right hand side */ 562302aa1bSDebojyoti Ghosh Vec W; /* Right-hand-side for implicit stage solve */ 572302aa1bSDebojyoti Ghosh Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 582302aa1bSDebojyoti Ghosh PetscScalar *work; /* Scalar work */ 592302aa1bSDebojyoti Ghosh PetscReal scoeff; /* shift = scoeff/dt */ 602302aa1bSDebojyoti Ghosh PetscReal stage_time; 612302aa1bSDebojyoti Ghosh TSStepStatus status; 622302aa1bSDebojyoti Ghosh } TS_GLEE; 632302aa1bSDebojyoti Ghosh 642302aa1bSDebojyoti Ghosh /*MC 652302aa1bSDebojyoti Ghosh TSGLEE23 - Second order three stage GLEE method 662302aa1bSDebojyoti Ghosh 672302aa1bSDebojyoti Ghosh This method has three stages. 682302aa1bSDebojyoti Ghosh s = 3, r = 2 692302aa1bSDebojyoti Ghosh 702302aa1bSDebojyoti Ghosh Level: advanced 712302aa1bSDebojyoti Ghosh 722302aa1bSDebojyoti Ghosh .seealso: TSGLEE 732302aa1bSDebojyoti Ghosh */ 742302aa1bSDebojyoti Ghosh /*MC 752302aa1bSDebojyoti Ghosh TSGLEE24 - Second order four stage GLEE method 762302aa1bSDebojyoti Ghosh 772302aa1bSDebojyoti Ghosh This method has four stages. 782302aa1bSDebojyoti Ghosh s = 4, r = 2 792302aa1bSDebojyoti Ghosh 802302aa1bSDebojyoti Ghosh Level: advanced 812302aa1bSDebojyoti Ghosh 822302aa1bSDebojyoti Ghosh .seealso: TSGLEE 832302aa1bSDebojyoti Ghosh */ 842302aa1bSDebojyoti Ghosh /*MC 852302aa1bSDebojyoti Ghosh TSGLEE25i - Second order five stage GLEE method 862302aa1bSDebojyoti Ghosh 872302aa1bSDebojyoti Ghosh This method has five stages. 882302aa1bSDebojyoti Ghosh s = 5, r = 2 892302aa1bSDebojyoti Ghosh 902302aa1bSDebojyoti Ghosh Level: advanced 912302aa1bSDebojyoti Ghosh 922302aa1bSDebojyoti Ghosh .seealso: TSGLEE 932302aa1bSDebojyoti Ghosh */ 942302aa1bSDebojyoti Ghosh /*MC 952302aa1bSDebojyoti Ghosh TSGLEE35 - Third order five stage GLEE method 962302aa1bSDebojyoti Ghosh 972302aa1bSDebojyoti Ghosh This method has five stages. 982302aa1bSDebojyoti Ghosh s = 5, r = 2 992302aa1bSDebojyoti Ghosh 1002302aa1bSDebojyoti Ghosh Level: advanced 1012302aa1bSDebojyoti Ghosh 1022302aa1bSDebojyoti Ghosh .seealso: TSGLEE 1032302aa1bSDebojyoti Ghosh */ 1042302aa1bSDebojyoti Ghosh /*MC 1052302aa1bSDebojyoti Ghosh TSGLEEEXRK2A - Second order six stage GLEE method 1062302aa1bSDebojyoti Ghosh 1072302aa1bSDebojyoti Ghosh This method has six stages. 1082302aa1bSDebojyoti Ghosh s = 6, r = 2 1092302aa1bSDebojyoti Ghosh 1102302aa1bSDebojyoti Ghosh Level: advanced 1112302aa1bSDebojyoti Ghosh 1122302aa1bSDebojyoti Ghosh .seealso: TSGLEE 1132302aa1bSDebojyoti Ghosh */ 1142302aa1bSDebojyoti Ghosh /*MC 1152302aa1bSDebojyoti Ghosh TSGLEERK32G1 - Third order eight stage GLEE method 1162302aa1bSDebojyoti Ghosh 1172302aa1bSDebojyoti Ghosh This method has eight stages. 1182302aa1bSDebojyoti Ghosh s = 8, r = 2 1192302aa1bSDebojyoti Ghosh 1202302aa1bSDebojyoti Ghosh Level: advanced 1212302aa1bSDebojyoti Ghosh 1222302aa1bSDebojyoti Ghosh .seealso: TSGLEE 1232302aa1bSDebojyoti Ghosh */ 1242302aa1bSDebojyoti Ghosh /*MC 1252302aa1bSDebojyoti Ghosh TSGLEERK285EX - Second order nine stage GLEE method 1262302aa1bSDebojyoti Ghosh 1272302aa1bSDebojyoti Ghosh This method has nine stages. 1282302aa1bSDebojyoti Ghosh s = 9, r = 2 1292302aa1bSDebojyoti Ghosh 1302302aa1bSDebojyoti Ghosh Level: advanced 1312302aa1bSDebojyoti Ghosh 1322302aa1bSDebojyoti Ghosh .seealso: TSGLEE 1332302aa1bSDebojyoti Ghosh */ 1342302aa1bSDebojyoti Ghosh 1352302aa1bSDebojyoti Ghosh #undef __FUNCT__ 1362302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterAll" 1372302aa1bSDebojyoti Ghosh /*@C 1382302aa1bSDebojyoti Ghosh TSGLEERegisterAll - Registers all of the General Linear with Error Estimation methods in TSGLEE 1392302aa1bSDebojyoti Ghosh 1402302aa1bSDebojyoti Ghosh Not Collective, but should be called by all processes which will need the schemes to be registered 1412302aa1bSDebojyoti Ghosh 1422302aa1bSDebojyoti Ghosh Level: advanced 1432302aa1bSDebojyoti Ghosh 1442302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, register, all 1452302aa1bSDebojyoti Ghosh 1462302aa1bSDebojyoti Ghosh .seealso: TSGLEERegisterDestroy() 1472302aa1bSDebojyoti Ghosh @*/ 1482302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterAll(void) 1492302aa1bSDebojyoti Ghosh { 1502302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 1512302aa1bSDebojyoti Ghosh 1522302aa1bSDebojyoti Ghosh PetscFunctionBegin; 1532302aa1bSDebojyoti Ghosh if (TSGLEERegisterAllCalled) PetscFunctionReturn(0); 1542302aa1bSDebojyoti Ghosh TSGLEERegisterAllCalled = PETSC_TRUE; 1552302aa1bSDebojyoti Ghosh 1562302aa1bSDebojyoti Ghosh { 1572302aa1bSDebojyoti Ghosh /* y-eps form */ 1582302aa1bSDebojyoti Ghosh const PetscInt 159ab8c5912SEmil Constantinescu p = 1, 160ab8c5912SEmil Constantinescu s = 3, 161ab8c5912SEmil Constantinescu r = 2; 162ab8c5912SEmil Constantinescu const PetscReal 163ab8c5912SEmil Constantinescu gamma = 0.5, 164ab8c5912SEmil Constantinescu A[3][3] = {{1.0,0,0},{0,0.5,0},{0,0.5,0.5}}, 165ab8c5912SEmil Constantinescu B[2][3] = {{1.0,0,0},{-2.0,1.0,1.0}}, 166ab8c5912SEmil Constantinescu U[3][2] = {{1.0,0},{1.0,0.5},{1.0,0.5}}, 167ab8c5912SEmil Constantinescu V[2][2] = {{1,0},{0,1}}, 168ab8c5912SEmil Constantinescu S[2] = {1,0}, 169ab8c5912SEmil Constantinescu F[2] = {1,0}, 170fd96ebc2SDebojyoti Ghosh Fembed[2] = {1,1-gamma}, 17157df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 17257df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 17357df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 174ab8c5912SEmil Constantinescu } 175ab8c5912SEmil Constantinescu { 176ab8c5912SEmil Constantinescu /* y-eps form */ 177ab8c5912SEmil Constantinescu const PetscInt 1782302aa1bSDebojyoti Ghosh p = 2, 1792302aa1bSDebojyoti Ghosh s = 3, 1802302aa1bSDebojyoti Ghosh r = 2; 1812302aa1bSDebojyoti Ghosh const PetscReal 1822302aa1bSDebojyoti Ghosh gamma = 0, 1832302aa1bSDebojyoti Ghosh A[3][3] = {{0,0,0},{1,0,0},{0.25,0.25,0}}, 1842302aa1bSDebojyoti 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}}, 1852302aa1bSDebojyoti Ghosh U[3][2] = {{1,0},{1,10},{1,-1}}, 1862302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 1872302aa1bSDebojyoti Ghosh S[2] = {1,0}, 1882302aa1bSDebojyoti Ghosh F[2] = {1,0}, 189fd96ebc2SDebojyoti Ghosh Fembed[2] = {1,1-gamma}, 19057df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 19157df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 19257df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 1932302aa1bSDebojyoti Ghosh } 1942302aa1bSDebojyoti Ghosh { 1952302aa1bSDebojyoti Ghosh /* y-y~ form */ 1962302aa1bSDebojyoti Ghosh const PetscInt 1972302aa1bSDebojyoti Ghosh p = 2, 1982302aa1bSDebojyoti Ghosh s = 4, 1992302aa1bSDebojyoti Ghosh r = 2; 2002302aa1bSDebojyoti Ghosh const PetscReal 2012302aa1bSDebojyoti Ghosh gamma = 0, 2022302aa1bSDebojyoti 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}}, 2032302aa1bSDebojyoti 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}}, 2042302aa1bSDebojyoti Ghosh U[4][2] = {{0,1},{75.0/58.0,-17.0/58.0},{0,1},{0,1}}, 2052302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2062302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2072302aa1bSDebojyoti Ghosh F[2] = {1,0}, 208fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 20957df6a1bSDebojyoti Ghosh Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)}, 21057df6a1bSDebojyoti Ghosh Serror[2] = {1.0-gamma,1.0}; 21157df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 2122302aa1bSDebojyoti Ghosh } 2132302aa1bSDebojyoti Ghosh { 2142302aa1bSDebojyoti Ghosh /* y-y~ form */ 2152302aa1bSDebojyoti Ghosh const PetscInt 2162302aa1bSDebojyoti Ghosh p = 2, 2172302aa1bSDebojyoti Ghosh s = 5, 2182302aa1bSDebojyoti Ghosh r = 2; 2192302aa1bSDebojyoti Ghosh const PetscReal 2202302aa1bSDebojyoti Ghosh gamma = 0, 2212302aa1bSDebojyoti Ghosh A[5][5] = {{0,0,0,0,0}, 2222302aa1bSDebojyoti Ghosh {-0.94079244066783383269,0,0,0,0}, 2232302aa1bSDebojyoti Ghosh {0.64228187778301907108,0.10915356933958500042,0,0,0}, 2242302aa1bSDebojyoti Ghosh {-0.51764297742287450812,0.74414270351096040738,-0.71404164927824538121,0,0}, 2252302aa1bSDebojyoti Ghosh {-0.44696561556825969206,-0.76768425657590196518,0.20111608138142987881,0.93828186737840469796,0}}, 2262302aa1bSDebojyoti Ghosh B[2][5] = {{-0.029309178948150356153,-0.49671981884013874923,0.34275801517650053274,0.32941112623949194988,0.85385985637229662276}, 2272302aa1bSDebojyoti Ghosh {0.78133219686062535272,0.074238691892675897635,0.57957363498384957966,-0.24638502829674959968,-0.18875949544040123033}}, 2282302aa1bSDebojyoti Ghosh U[5][2] = {{0.16911424754448327735,0.83088575245551672265}, 2292302aa1bSDebojyoti Ghosh {0.53638465733199574340,0.46361534266800425660}, 2302302aa1bSDebojyoti Ghosh {0.39901579167169582526,0.60098420832830417474}, 2312302aa1bSDebojyoti Ghosh {0.87689005530618575480,0.12310994469381424520}, 2322302aa1bSDebojyoti Ghosh {0.99056100455550913009,0.0094389954444908699092}}, 2332302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2342302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2352302aa1bSDebojyoti Ghosh F[2] = {1,0}, 236fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 23757df6a1bSDebojyoti Ghosh Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)}, 23857df6a1bSDebojyoti Ghosh Serror[2] = {1.0-gamma,1.0}; 23957df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 2402302aa1bSDebojyoti Ghosh } 2412302aa1bSDebojyoti Ghosh { 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 gamma = 0, 2492302aa1bSDebojyoti Ghosh A[5][5] = {{0,0,0,0,0}, 2502302aa1bSDebojyoti Ghosh {- 2169604947363702313.0 / 24313474998937147335.0,0,0,0,0}, 2512302aa1bSDebojyoti Ghosh {46526746497697123895.0 / 94116917485856474137.0,-10297879244026594958.0 / 49199457603717988219.0,0,0,0}, 2522302aa1bSDebojyoti Ghosh {23364788935845982499.0 / 87425311444725389446.0,-79205144337496116638.0 / 148994349441340815519.0,40051189859317443782.0 / 36487615018004984309.0,0,0}, 2532302aa1bSDebojyoti Ghosh {42089522664062539205.0 / 124911313006412840286.0,-15074384760342762939.0 / 137927286865289746282.0,-62274678522253371016.0 / 125918573676298591413.0,13755475729852471739.0 / 79257927066651693390.0}}, 2542302aa1bSDebojyoti Ghosh B[2][5] = {{61546696837458703723.0 / 56982519523786160813.0,-55810892792806293355.0 / 206957624151308356511.0,24061048952676379087.0 / 158739347956038723465.0,3577972206874351339.0 / 7599733370677197135.0,-59449832954780563947.0 / 137360038685338563670.0}, 2552302aa1bSDebojyoti Ghosh {- 9738262186984159168.0 / 99299082461487742983.0,-32797097931948613195.0 / 61521565616362163366.0,42895514606418420631.0 / 71714201188501437336.0,22608567633166065068.0 / 55371917805607957003.0,94655809487476459565.0 / 151517167160302729021.0}}, 2562302aa1bSDebojyoti Ghosh U[5][2] = {{70820309139834661559.0 / 80863923579509469826.0,10043614439674808267.0 / 80863923579509469826.0}, 2572302aa1bSDebojyoti Ghosh {161694774978034105510.0 / 106187653640211060371.0,-55507121337823045139.0 / 106187653640211060371.0}, 2582302aa1bSDebojyoti Ghosh {78486094644566264568.0 / 88171030896733822981.0,9684936252167558413.0 / 88171030896733822981.0}, 2592302aa1bSDebojyoti Ghosh {65394922146334854435.0 / 84570853840405479554.0,19175931694070625119.0 / 84570853840405479554.0}, 2602302aa1bSDebojyoti Ghosh {8607282770183754108.0 / 108658046436496925911.0,100050763666313171803.0 / 108658046436496925911.0}}, 2612302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 2622302aa1bSDebojyoti Ghosh S[2] = {1,1}, 2632302aa1bSDebojyoti Ghosh F[2] = {1,0}, 264fd96ebc2SDebojyoti Ghosh Fembed[2] = {0,1}, 26557df6a1bSDebojyoti Ghosh Ferror[2] = {-1.0/(1.0-gamma),1.0/(1.0-gamma)}, 26657df6a1bSDebojyoti Ghosh Serror[2] = {1.0-gamma,1.0}; 26757df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 2682302aa1bSDebojyoti Ghosh } 2692302aa1bSDebojyoti Ghosh { 2702302aa1bSDebojyoti Ghosh /* y-eps form */ 2712302aa1bSDebojyoti Ghosh const PetscInt 2722302aa1bSDebojyoti Ghosh p = 2, 2732302aa1bSDebojyoti Ghosh s = 6, 2742302aa1bSDebojyoti Ghosh r = 2; 2752302aa1bSDebojyoti Ghosh const PetscReal 2762302aa1bSDebojyoti Ghosh gamma = 0.25, 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}, 289fd96ebc2SDebojyoti Ghosh Fembed[2] = {1,1-gamma}, 29057df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 29157df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 29257df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 2932302aa1bSDebojyoti Ghosh } 2942302aa1bSDebojyoti Ghosh { 2952302aa1bSDebojyoti Ghosh /* y-eps form */ 2962302aa1bSDebojyoti Ghosh const PetscInt 2972302aa1bSDebojyoti Ghosh p = 3, 2982302aa1bSDebojyoti Ghosh s = 8, 2992302aa1bSDebojyoti Ghosh r = 2; 3002302aa1bSDebojyoti Ghosh const PetscReal 3012302aa1bSDebojyoti Ghosh gamma = 0, 3022302aa1bSDebojyoti Ghosh A[8][8] = {{0,0,0,0,0,0,0,0}, 3032302aa1bSDebojyoti Ghosh {0.5,0,0,0,0,0,0,0}, 3042302aa1bSDebojyoti Ghosh {-1,2,0,0,0,0,0,0}, 3052302aa1bSDebojyoti Ghosh {1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0}, 3062302aa1bSDebojyoti Ghosh {0,0,0,0,0,0,0,0}, 3072302aa1bSDebojyoti Ghosh {-7.0/24.0,1.0/3.0,1.0/12.0,-1.0/8.0,0.5,0,0,0}, 3082302aa1bSDebojyoti Ghosh {7.0/6.0,-4.0/3.0,-1.0/3.0,0.5,-1.0,2.0,0,0}, 3092302aa1bSDebojyoti Ghosh {0,0,0,0,1.0/6.0,2.0/3.0,1.0/6.0,0}}, 3102302aa1bSDebojyoti Ghosh B[2][8] = {{1.0/6.0,2.0/3.0,1.0/6.0,0,0,0,0,0}, 3112302aa1bSDebojyoti 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}}, 3122302aa1bSDebojyoti Ghosh U[8][2] = {{1,0},{1,0},{1,0},{1,0},{1,1},{1,1},{1,1},{1,1}}, 3132302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 3142302aa1bSDebojyoti Ghosh S[2] = {1,0}, 3152302aa1bSDebojyoti Ghosh F[2] = {1,0}, 316fd96ebc2SDebojyoti Ghosh Fembed[2] = {1,1-gamma}, 31757df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 31857df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 31957df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 3202302aa1bSDebojyoti Ghosh } 3212302aa1bSDebojyoti Ghosh { 3222302aa1bSDebojyoti Ghosh /* y-eps form */ 3232302aa1bSDebojyoti Ghosh const PetscInt 3242302aa1bSDebojyoti Ghosh p = 2, 3252302aa1bSDebojyoti Ghosh s = 9, 3262302aa1bSDebojyoti Ghosh r = 2; 3272302aa1bSDebojyoti Ghosh const PetscReal 3282302aa1bSDebojyoti Ghosh gamma = 0.25, 3292302aa1bSDebojyoti Ghosh A[9][9] = {{0,0,0,0,0,0,0,0,0}, 3302302aa1bSDebojyoti Ghosh {0.585786437626904966,0,0,0,0,0,0,0,0}, 3312302aa1bSDebojyoti Ghosh {0.149999999999999994,0.849999999999999978,0,0,0,0,0,0,0}, 3322302aa1bSDebojyoti Ghosh {0,0,0,0,0,0,0,0,0}, 3332302aa1bSDebojyoti Ghosh {0,0,0,0.292893218813452483,0,0,0,0,0}, 3342302aa1bSDebojyoti Ghosh {0,0,0,0.0749999999999999972,0.424999999999999989,0,0,0,0}, 3352302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0,0,0}, 3362302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.292893218813452483,0,0}, 3372302aa1bSDebojyoti Ghosh {0,0,0,0.176776695296636893,0.176776695296636893,0.146446609406726241,0.0749999999999999972,0.424999999999999989,0}}, 3382302aa1bSDebojyoti Ghosh B[2][9] = {{0.353553390593273786,0.353553390593273786,0.292893218813452483,0,0,0,0,0,0}, 3392302aa1bSDebojyoti Ghosh {-0.471404520791031678,-0.471404520791031678,-0.390524291751269959,0.235702260395515839,0.235702260395515839,0.195262145875634979,0.235702260395515839,0.235702260395515839,0.195262145875634979}}, 3402302aa1bSDebojyoti 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}}, 3412302aa1bSDebojyoti Ghosh V[2][2] = {{1,0},{0,1}}, 3422302aa1bSDebojyoti Ghosh S[2] = {1,0}, 3432302aa1bSDebojyoti Ghosh F[2] = {1,0}, 344fd96ebc2SDebojyoti Ghosh Fembed[2] = {1,1-gamma}, 34557df6a1bSDebojyoti Ghosh Ferror[2] = {0,1}, 34657df6a1bSDebojyoti Ghosh Serror[2] = {1,0}; 34757df6a1bSDebojyoti Ghosh ierr = 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);CHKERRQ(ierr); 3482302aa1bSDebojyoti Ghosh } 3492302aa1bSDebojyoti Ghosh 3502302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 3512302aa1bSDebojyoti Ghosh } 3522302aa1bSDebojyoti Ghosh 3532302aa1bSDebojyoti Ghosh #undef __FUNCT__ 3542302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegisterDestroy" 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 .keywords: TSGLEE, register, destroy 3632302aa1bSDebojyoti Ghosh .seealso: TSGLEERegister(), TSGLEERegisterAll() 3642302aa1bSDebojyoti Ghosh @*/ 3652302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegisterDestroy(void) 3662302aa1bSDebojyoti Ghosh { 3672302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 3682302aa1bSDebojyoti Ghosh GLEETableauLink link; 3692302aa1bSDebojyoti Ghosh 3702302aa1bSDebojyoti Ghosh PetscFunctionBegin; 3712302aa1bSDebojyoti Ghosh while ((link = GLEETableauList)) { 3722302aa1bSDebojyoti Ghosh GLEETableau t = &link->tab; 3732302aa1bSDebojyoti Ghosh GLEETableauList = link->next; 3742302aa1bSDebojyoti Ghosh ierr = PetscFree5(t->A,t->B,t->U,t->V,t->c); CHKERRQ(ierr); 3752302aa1bSDebojyoti Ghosh ierr = PetscFree2(t->S,t->F); CHKERRQ(ierr); 3762302aa1bSDebojyoti Ghosh ierr = PetscFree (t->Fembed); CHKERRQ(ierr); 377fd96ebc2SDebojyoti Ghosh ierr = PetscFree (t->Ferror); CHKERRQ(ierr); 37857df6a1bSDebojyoti Ghosh ierr = PetscFree (t->Serror); CHKERRQ(ierr); 3792302aa1bSDebojyoti Ghosh ierr = PetscFree (t->binterp); CHKERRQ(ierr); 3802302aa1bSDebojyoti Ghosh ierr = PetscFree (t->name); CHKERRQ(ierr); 3812302aa1bSDebojyoti Ghosh ierr = PetscFree (link); CHKERRQ(ierr); 3822302aa1bSDebojyoti Ghosh } 3832302aa1bSDebojyoti Ghosh TSGLEERegisterAllCalled = PETSC_FALSE; 3842302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 3852302aa1bSDebojyoti Ghosh } 3862302aa1bSDebojyoti Ghosh 3872302aa1bSDebojyoti Ghosh #undef __FUNCT__ 3882302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEInitializePackage" 3892302aa1bSDebojyoti Ghosh /*@C 3902302aa1bSDebojyoti Ghosh TSGLEEInitializePackage - This function initializes everything in the TSGLEE package. It is called 3912302aa1bSDebojyoti Ghosh from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLEE() 3922302aa1bSDebojyoti Ghosh when using static libraries. 3932302aa1bSDebojyoti Ghosh 3942302aa1bSDebojyoti Ghosh Level: developer 3952302aa1bSDebojyoti Ghosh 3962302aa1bSDebojyoti Ghosh .keywords: TS, TSGLEE, initialize, package 3972302aa1bSDebojyoti Ghosh .seealso: PetscInitialize() 3982302aa1bSDebojyoti Ghosh @*/ 3992302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEInitializePackage(void) 4002302aa1bSDebojyoti Ghosh { 4012302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 4022302aa1bSDebojyoti Ghosh 4032302aa1bSDebojyoti Ghosh PetscFunctionBegin; 4042302aa1bSDebojyoti Ghosh if (TSGLEEPackageInitialized) PetscFunctionReturn(0); 4052302aa1bSDebojyoti Ghosh TSGLEEPackageInitialized = PETSC_TRUE; 4062302aa1bSDebojyoti Ghosh ierr = TSGLEERegisterAll();CHKERRQ(ierr); 4072302aa1bSDebojyoti Ghosh ierr = PetscObjectComposedDataRegister(&explicit_stage_time_id);CHKERRQ(ierr); 4082302aa1bSDebojyoti Ghosh ierr = PetscRegisterFinalize(TSGLEEFinalizePackage);CHKERRQ(ierr); 4092302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 4102302aa1bSDebojyoti Ghosh } 4112302aa1bSDebojyoti Ghosh 4122302aa1bSDebojyoti Ghosh #undef __FUNCT__ 4132302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEFinalizePackage" 4142302aa1bSDebojyoti Ghosh /*@C 4152302aa1bSDebojyoti Ghosh TSGLEEFinalizePackage - This function destroys everything in the TSGLEE package. It is 4162302aa1bSDebojyoti Ghosh called from PetscFinalize(). 4172302aa1bSDebojyoti Ghosh 4182302aa1bSDebojyoti Ghosh Level: developer 4192302aa1bSDebojyoti Ghosh 4202302aa1bSDebojyoti Ghosh .keywords: Petsc, destroy, package 4212302aa1bSDebojyoti Ghosh .seealso: PetscFinalize() 4222302aa1bSDebojyoti Ghosh @*/ 4232302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEFinalizePackage(void) 4242302aa1bSDebojyoti Ghosh { 4252302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 4262302aa1bSDebojyoti Ghosh 4272302aa1bSDebojyoti Ghosh PetscFunctionBegin; 4282302aa1bSDebojyoti Ghosh TSGLEEPackageInitialized = PETSC_FALSE; 4292302aa1bSDebojyoti Ghosh ierr = TSGLEERegisterDestroy();CHKERRQ(ierr); 4302302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 4312302aa1bSDebojyoti Ghosh } 4322302aa1bSDebojyoti Ghosh 4332302aa1bSDebojyoti Ghosh #undef __FUNCT__ 4342302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERegister" 4352302aa1bSDebojyoti Ghosh /*@C 4362302aa1bSDebojyoti Ghosh TSGLEERegister - register an GLEE scheme by providing the entries in the Butcher tableau 4372302aa1bSDebojyoti Ghosh 4382302aa1bSDebojyoti Ghosh Not Collective, but the same schemes should be registered on all processes on which they will be used 4392302aa1bSDebojyoti Ghosh 4402302aa1bSDebojyoti Ghosh Input Parameters: 4412302aa1bSDebojyoti Ghosh + name - identifier for method 4422302aa1bSDebojyoti Ghosh . order - order of method 4432302aa1bSDebojyoti Ghosh . s - number of stages 4442302aa1bSDebojyoti Ghosh . r - number of steps 4452302aa1bSDebojyoti Ghosh . gamma - LTE ratio 4462302aa1bSDebojyoti Ghosh . A - stage coefficients (dimension s*s, row-major) 4472302aa1bSDebojyoti Ghosh . B - step completion coefficients (dimension r*s, row-major) 4482302aa1bSDebojyoti Ghosh . U - method coefficients (dimension s*r, row-major) 4492302aa1bSDebojyoti Ghosh . V - method coefficients (dimension r*r, row-major) 4502302aa1bSDebojyoti Ghosh . S - starting coefficients 4512302aa1bSDebojyoti Ghosh . F - finishing coefficients 4522302aa1bSDebojyoti Ghosh . c - abscissa (dimension s; NULL to use row sums of A) 4532302aa1bSDebojyoti Ghosh . Fembed - step completion coefficients for embedded method 454fd96ebc2SDebojyoti Ghosh . Ferror - error computation coefficients 45557df6a1bSDebojyoti Ghosh . Serror - error initialization coefficients 4562302aa1bSDebojyoti Ghosh . pinterp - order of interpolation (0 if unavailable) 4572302aa1bSDebojyoti Ghosh - binterp - array of interpolation coefficients (NULL if unavailable) 4582302aa1bSDebojyoti Ghosh 4592302aa1bSDebojyoti Ghosh Notes: 4602302aa1bSDebojyoti Ghosh Several GLEE methods are provided, this function is only needed to create new methods. 4612302aa1bSDebojyoti Ghosh 4622302aa1bSDebojyoti Ghosh Level: advanced 4632302aa1bSDebojyoti Ghosh 4642302aa1bSDebojyoti Ghosh .keywords: TS, register 4652302aa1bSDebojyoti Ghosh 4662302aa1bSDebojyoti Ghosh .seealso: TSGLEE 4672302aa1bSDebojyoti Ghosh @*/ 4682302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEERegister(TSGLEEType name,PetscInt order,PetscInt s, PetscInt r, 4692302aa1bSDebojyoti Ghosh PetscReal gamma, 4702302aa1bSDebojyoti Ghosh const PetscReal A[],const PetscReal B[], 4712302aa1bSDebojyoti Ghosh const PetscReal U[],const PetscReal V[], 4722302aa1bSDebojyoti Ghosh const PetscReal S[],const PetscReal F[], 473fd96ebc2SDebojyoti Ghosh const PetscReal c[], 474fd96ebc2SDebojyoti Ghosh const PetscReal Fembed[],const PetscReal Ferror[], 47557df6a1bSDebojyoti Ghosh const PetscReal Serror[], 4762302aa1bSDebojyoti Ghosh PetscInt pinterp, const PetscReal binterp[]) 4772302aa1bSDebojyoti Ghosh { 4782302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 4792302aa1bSDebojyoti Ghosh GLEETableauLink link; 4802302aa1bSDebojyoti Ghosh GLEETableau t; 4812302aa1bSDebojyoti Ghosh PetscInt i,j; 4822302aa1bSDebojyoti Ghosh 4832302aa1bSDebojyoti Ghosh PetscFunctionBegin; 4842302aa1bSDebojyoti Ghosh ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 4852302aa1bSDebojyoti Ghosh ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 4862302aa1bSDebojyoti Ghosh t = &link->tab; 4872302aa1bSDebojyoti Ghosh ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 4882302aa1bSDebojyoti Ghosh t->order = order; 4892302aa1bSDebojyoti Ghosh t->s = s; 4902302aa1bSDebojyoti Ghosh t->r = r; 4912302aa1bSDebojyoti Ghosh t->gamma = gamma; 4922302aa1bSDebojyoti Ghosh ierr = PetscMalloc5(s*s,&t->A,r*r,&t->V,s,&t->c,r*s,&t->B,s*r,&t->U); CHKERRQ(ierr); 4932302aa1bSDebojyoti Ghosh ierr = PetscMalloc2(r,&t->S,r,&t->F); CHKERRQ(ierr); 4942302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 4952302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->B,B,r*s*sizeof(B[0]));CHKERRQ(ierr); 4962302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->U,U,s*r*sizeof(U[0]));CHKERRQ(ierr); 4972302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->V,V,r*r*sizeof(V[0]));CHKERRQ(ierr); 4982302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->S,S,r *sizeof(S[0]));CHKERRQ(ierr); 4992302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->F,F,r *sizeof(F[0]));CHKERRQ(ierr); 5002302aa1bSDebojyoti Ghosh if (c) { 5012302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->c,c,s*sizeof(c[0])); CHKERRQ(ierr); 5022302aa1bSDebojyoti Ghosh } else { 5032302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 5042302aa1bSDebojyoti Ghosh } 5052302aa1bSDebojyoti Ghosh ierr = PetscMalloc1(r,&t->Fembed);CHKERRQ(ierr); 506fd96ebc2SDebojyoti Ghosh ierr = PetscMalloc1(r,&t->Ferror);CHKERRQ(ierr); 50757df6a1bSDebojyoti Ghosh ierr = PetscMalloc1(r,&t->Serror);CHKERRQ(ierr); 5082302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->Fembed,Fembed,r*sizeof(Fembed[0]));CHKERRQ(ierr); 509fd96ebc2SDebojyoti Ghosh ierr = PetscMemcpy(t->Ferror,Ferror,r*sizeof(Ferror[0]));CHKERRQ(ierr); 51057df6a1bSDebojyoti Ghosh ierr = PetscMemcpy(t->Serror,Serror,r*sizeof(Serror[0]));CHKERRQ(ierr); 5112302aa1bSDebojyoti Ghosh t->pinterp = pinterp; 5122302aa1bSDebojyoti Ghosh ierr = PetscMalloc1(s*pinterp,&t->binterp);CHKERRQ(ierr); 5132302aa1bSDebojyoti Ghosh ierr = PetscMemcpy(t->binterp,binterp,s*pinterp*sizeof(binterp[0]));CHKERRQ(ierr); 5142302aa1bSDebojyoti Ghosh 5152302aa1bSDebojyoti Ghosh link->next = GLEETableauList; 5162302aa1bSDebojyoti Ghosh GLEETableauList = link; 5172302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5182302aa1bSDebojyoti Ghosh } 5192302aa1bSDebojyoti Ghosh 5202302aa1bSDebojyoti Ghosh #undef __FUNCT__ 5212302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSEvaluateStep_GLEE" 5222302aa1bSDebojyoti Ghosh static PetscErrorCode TSEvaluateStep_GLEE(TS ts,PetscInt order,Vec X,PetscBool *done) 5232302aa1bSDebojyoti Ghosh { 5242302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*) ts->data; 5252302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 5262302aa1bSDebojyoti Ghosh PetscReal h, *B = tab->B, *V = tab->V, 527fd96ebc2SDebojyoti Ghosh *F = tab->F, 528fd96ebc2SDebojyoti Ghosh *Fembed = tab->Fembed; 5292302aa1bSDebojyoti Ghosh PetscInt s = tab->s, r = tab->r, i, j; 5302302aa1bSDebojyoti Ghosh Vec *Y = glee->Y, *YdotStage = glee->YdotStage; 5312302aa1bSDebojyoti Ghosh PetscScalar *w = glee->work; 5322302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 5332302aa1bSDebojyoti Ghosh 5342302aa1bSDebojyoti Ghosh PetscFunctionBegin; 5352302aa1bSDebojyoti Ghosh 5362302aa1bSDebojyoti Ghosh switch (glee->status) { 5372302aa1bSDebojyoti Ghosh case TS_STEP_INCOMPLETE: 5382302aa1bSDebojyoti Ghosh case TS_STEP_PENDING: 5392302aa1bSDebojyoti Ghosh h = ts->time_step; break; 5402302aa1bSDebojyoti Ghosh case TS_STEP_COMPLETE: 541ed450d42SEmil Constantinescu h = ts->ptime - ts->ptime_prev; break; 5422302aa1bSDebojyoti Ghosh default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 5432302aa1bSDebojyoti Ghosh } 5442302aa1bSDebojyoti Ghosh 5452302aa1bSDebojyoti Ghosh if (order == tab->order) { 5462302aa1bSDebojyoti Ghosh 5472302aa1bSDebojyoti Ghosh /* Note: Irrespective of whether status is TS_STEP_INCOMPLETE 548fd96ebc2SDebojyoti Ghosh or TS_STEP_COMPLETE, glee->X has the solution at the 5492302aa1bSDebojyoti Ghosh beginning of the time step. So no need to roll-back. 5502302aa1bSDebojyoti Ghosh */ 5512302aa1bSDebojyoti Ghosh if (glee->status == TS_STEP_INCOMPLETE) { 5522302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 5532302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr); 5542302aa1bSDebojyoti Ghosh ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr); 5552302aa1bSDebojyoti Ghosh for (j=0; j<s; j++) w[j] = h*B[i*s+j]; 5562302aa1bSDebojyoti Ghosh ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr); 5572302aa1bSDebojyoti Ghosh } 5582302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(X);CHKERRQ(ierr); 5592302aa1bSDebojyoti Ghosh ierr = VecMAXPY(X,r,F,Y);CHKERRQ(ierr); 5602302aa1bSDebojyoti Ghosh } else {ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);} 5612302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5622302aa1bSDebojyoti Ghosh 5632302aa1bSDebojyoti Ghosh } else if (order == tab->order-1) { 5642302aa1bSDebojyoti Ghosh 5652302aa1bSDebojyoti Ghosh /* Complete with the embedded method (Fembed) */ 5662302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 5672302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(Y[i]); CHKERRQ(ierr); 5682302aa1bSDebojyoti Ghosh ierr = VecMAXPY(Y[i],r,(V+i*r),glee->X); CHKERRQ(ierr); 5692302aa1bSDebojyoti Ghosh for (j=0; j<s; j++) w[j] = h*B[i*s+j]; 5702302aa1bSDebojyoti Ghosh ierr = VecMAXPY(Y[i],s,w,YdotStage); CHKERRQ(ierr); 5712302aa1bSDebojyoti Ghosh } 5722302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(X);CHKERRQ(ierr); 5732302aa1bSDebojyoti Ghosh ierr = VecMAXPY(X,r,Fembed,Y);CHKERRQ(ierr); 574fd96ebc2SDebojyoti Ghosh 5752302aa1bSDebojyoti Ghosh if (done) *done = PETSC_TRUE; 5762302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5772302aa1bSDebojyoti Ghosh } 5782302aa1bSDebojyoti Ghosh if (done) *done = PETSC_FALSE; 5792302aa1bSDebojyoti Ghosh else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"GLEE '%s' of order %D cannot evaluate step at order %D",tab->name,tab->order,order); 5802302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 5812302aa1bSDebojyoti Ghosh } 5822302aa1bSDebojyoti Ghosh 5832302aa1bSDebojyoti Ghosh #undef __FUNCT__ 5842302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStep_GLEE" 5852302aa1bSDebojyoti Ghosh static PetscErrorCode TSStep_GLEE(TS ts) 5862302aa1bSDebojyoti Ghosh { 5872302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 5882302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 5892302aa1bSDebojyoti Ghosh const PetscInt s = tab->s, r = tab->r; 5902302aa1bSDebojyoti Ghosh PetscReal *A = tab->A, *U = tab->U, 5912302aa1bSDebojyoti Ghosh *F = tab->F, 5922302aa1bSDebojyoti Ghosh *c = tab->c; 5932302aa1bSDebojyoti Ghosh Vec *Y = glee->Y, *X = glee->X, 5942302aa1bSDebojyoti Ghosh *YStage = glee->YStage, 5952302aa1bSDebojyoti Ghosh *YdotStage = glee->YdotStage, 5962302aa1bSDebojyoti Ghosh W = glee->W; 5972302aa1bSDebojyoti Ghosh SNES snes; 5982302aa1bSDebojyoti Ghosh PetscScalar *w = glee->work; 5992302aa1bSDebojyoti Ghosh TSAdapt adapt; 6002302aa1bSDebojyoti Ghosh PetscInt i,j,reject,next_scheme,its,lits; 6012302aa1bSDebojyoti Ghosh PetscReal next_time_step; 6022302aa1bSDebojyoti Ghosh PetscReal t; 6032302aa1bSDebojyoti Ghosh PetscBool accept; 6042302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 6052302aa1bSDebojyoti Ghosh 6062302aa1bSDebojyoti Ghosh PetscFunctionBegin; 607*642f3702SEmil Constantinescu ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 608*642f3702SEmil Constantinescu 6092302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { ierr = VecCopy(Y[i],X[i]); CHKERRQ(ierr); } 6102302aa1bSDebojyoti Ghosh 6112302aa1bSDebojyoti Ghosh ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 6122302aa1bSDebojyoti Ghosh next_time_step = ts->time_step; 6132302aa1bSDebojyoti Ghosh t = ts->ptime; 6142302aa1bSDebojyoti Ghosh accept = PETSC_TRUE; 6152302aa1bSDebojyoti Ghosh glee->status = TS_STEP_INCOMPLETE; 6162302aa1bSDebojyoti Ghosh 6172302aa1bSDebojyoti Ghosh for (reject=0; reject<ts->max_reject && !ts->reason; reject++,ts->reject++) { 6182302aa1bSDebojyoti Ghosh 6192302aa1bSDebojyoti Ghosh PetscReal h = ts->time_step; 6202302aa1bSDebojyoti Ghosh ierr = TSPreStep(ts);CHKERRQ(ierr); 6212302aa1bSDebojyoti Ghosh 6222302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) { 6232302aa1bSDebojyoti Ghosh 6242302aa1bSDebojyoti Ghosh glee->stage_time = t + h*c[i]; 6252302aa1bSDebojyoti Ghosh ierr = TSPreStage(ts,glee->stage_time); CHKERRQ(ierr); 6262302aa1bSDebojyoti Ghosh 6272302aa1bSDebojyoti Ghosh if (A[i*s+i] == 0) { /* Explicit stage */ 6282302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(YStage[i]);CHKERRQ(ierr); 6292302aa1bSDebojyoti Ghosh ierr = VecMAXPY(YStage[i],r,(U+i*r),X);CHKERRQ(ierr); 6302302aa1bSDebojyoti Ghosh for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 6312302aa1bSDebojyoti Ghosh ierr = VecMAXPY(YStage[i],i,w,YdotStage);CHKERRQ(ierr); 6322302aa1bSDebojyoti Ghosh } else { /* Implicit stage */ 6332302aa1bSDebojyoti Ghosh glee->scoeff = 1.0/A[i*s+i]; 6342302aa1bSDebojyoti Ghosh /* compute right-hand-side */ 6352302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(W);CHKERRQ(ierr); 6362302aa1bSDebojyoti Ghosh ierr = VecMAXPY(W,r,(U+i*r),X);CHKERRQ(ierr); 6372302aa1bSDebojyoti Ghosh for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 6382302aa1bSDebojyoti Ghosh ierr = VecMAXPY(W,i,w,YdotStage);CHKERRQ(ierr); 6392302aa1bSDebojyoti Ghosh ierr = VecScale(W,glee->scoeff/h);CHKERRQ(ierr); 6402302aa1bSDebojyoti Ghosh /* set initial guess */ 6412302aa1bSDebojyoti Ghosh ierr = VecCopy(i>0 ? YStage[i-1] : ts->vec_sol,YStage[i]);CHKERRQ(ierr); 6422302aa1bSDebojyoti Ghosh /* solve for this stage */ 6432302aa1bSDebojyoti Ghosh ierr = SNESSolve(snes,W,YStage[i]);CHKERRQ(ierr); 6442302aa1bSDebojyoti Ghosh ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 6452302aa1bSDebojyoti Ghosh ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 6462302aa1bSDebojyoti Ghosh ts->snes_its += its; ts->ksp_its += lits; 6472302aa1bSDebojyoti Ghosh } 6482302aa1bSDebojyoti Ghosh ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 6496b234c65SEmil Constantinescu ierr = TSAdaptCheckStage(adapt,ts,glee->stage_time,Y[i],&accept);CHKERRQ(ierr); 6502302aa1bSDebojyoti Ghosh if (!accept) goto reject_step; 6512302aa1bSDebojyoti Ghosh ierr = TSPostStage(ts,glee->stage_time,i,YStage); CHKERRQ(ierr); 6522302aa1bSDebojyoti Ghosh ierr = TSComputeRHSFunction(ts,t+h*c[i],YStage[i],YdotStage[i]);CHKERRQ(ierr); 6532302aa1bSDebojyoti Ghosh } 6542302aa1bSDebojyoti Ghosh ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 6552302aa1bSDebojyoti Ghosh glee->status = TS_STEP_PENDING; 6562302aa1bSDebojyoti Ghosh 6572302aa1bSDebojyoti Ghosh /* Register only the current method as a candidate because we're not supporting multiple candidates yet. */ 6582302aa1bSDebojyoti Ghosh ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 6592302aa1bSDebojyoti Ghosh ierr = TSAdaptCandidatesClear(adapt);CHKERRQ(ierr); 6602302aa1bSDebojyoti Ghosh ierr = TSAdaptCandidateAdd(adapt,tab->name,tab->order,1,tab->ccfl,1.*tab->s,PETSC_TRUE);CHKERRQ(ierr); 6612302aa1bSDebojyoti Ghosh ierr = TSAdaptChoose(adapt,ts,ts->time_step,&next_scheme,&next_time_step,&accept);CHKERRQ(ierr); 6622302aa1bSDebojyoti Ghosh if (accept) { 6632302aa1bSDebojyoti Ghosh /* ignore next_scheme for now */ 6642302aa1bSDebojyoti Ghosh ts->ptime += ts->time_step; 6652302aa1bSDebojyoti Ghosh ts->time_step = next_time_step; 6662302aa1bSDebojyoti Ghosh ts->steps++; 6672302aa1bSDebojyoti Ghosh glee->status = TS_STEP_COMPLETE; 6682302aa1bSDebojyoti Ghosh ierr = PetscObjectComposedDataSetReal((PetscObject)ts->vec_sol,explicit_stage_time_id,ts->ptime);CHKERRQ(ierr); 6692302aa1bSDebojyoti Ghosh break; 6702302aa1bSDebojyoti Ghosh } else { /* Roll back the current step */ 6712302aa1bSDebojyoti Ghosh ierr = VecMAXPY(ts->vec_sol,r,F,X);CHKERRQ(ierr); 6722302aa1bSDebojyoti Ghosh ts->time_step = next_time_step; 6732302aa1bSDebojyoti Ghosh glee->status = TS_STEP_INCOMPLETE; 6742302aa1bSDebojyoti Ghosh } 6752302aa1bSDebojyoti Ghosh reject_step: continue; 6762302aa1bSDebojyoti Ghosh } 6772302aa1bSDebojyoti Ghosh if (glee->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 6782302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 6792302aa1bSDebojyoti Ghosh } 6802302aa1bSDebojyoti Ghosh 6812302aa1bSDebojyoti Ghosh #undef __FUNCT__ 6822302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSInterpolate_GLEE" 6832302aa1bSDebojyoti Ghosh static PetscErrorCode TSInterpolate_GLEE(TS ts,PetscReal itime,Vec X) 6842302aa1bSDebojyoti Ghosh { 6852302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 6862302aa1bSDebojyoti Ghosh PetscInt s=glee->tableau->s, pinterp=glee->tableau->pinterp,i,j; 6872302aa1bSDebojyoti Ghosh PetscReal h,tt,t; 6882302aa1bSDebojyoti Ghosh PetscScalar *b; 6892302aa1bSDebojyoti Ghosh const PetscReal *B = glee->tableau->binterp; 6902302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 6912302aa1bSDebojyoti Ghosh 6922302aa1bSDebojyoti Ghosh PetscFunctionBegin; 6932302aa1bSDebojyoti Ghosh if (!B) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSGLEE %s does not have an interpolation formula",glee->tableau->name); 6942302aa1bSDebojyoti Ghosh switch (glee->status) { 6952302aa1bSDebojyoti Ghosh case TS_STEP_INCOMPLETE: 6962302aa1bSDebojyoti Ghosh case TS_STEP_PENDING: 6972302aa1bSDebojyoti Ghosh h = ts->time_step; 6982302aa1bSDebojyoti Ghosh t = (itime - ts->ptime)/h; 6992302aa1bSDebojyoti Ghosh break; 7002302aa1bSDebojyoti Ghosh case TS_STEP_COMPLETE: 701ed450d42SEmil Constantinescu h = ts->ptime - ts->ptime_prev; 7022302aa1bSDebojyoti Ghosh t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 7032302aa1bSDebojyoti Ghosh break; 7042302aa1bSDebojyoti Ghosh default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 7052302aa1bSDebojyoti Ghosh } 7062302aa1bSDebojyoti Ghosh ierr = PetscMalloc1(s,&b);CHKERRQ(ierr); 7072302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) b[i] = 0; 7082302aa1bSDebojyoti Ghosh for (j=0,tt=t; j<pinterp; j++,tt*=t) { 7092302aa1bSDebojyoti Ghosh for (i=0; i<s; i++) { 7102302aa1bSDebojyoti Ghosh b[i] += h * B[i*pinterp+j] * tt; 7112302aa1bSDebojyoti Ghosh } 7122302aa1bSDebojyoti Ghosh } 7132302aa1bSDebojyoti Ghosh ierr = VecCopy(glee->YStage[0],X);CHKERRQ(ierr); 7142302aa1bSDebojyoti Ghosh ierr = VecMAXPY(X,s,b,glee->YdotStage);CHKERRQ(ierr); 7152302aa1bSDebojyoti Ghosh ierr = PetscFree(b);CHKERRQ(ierr); 7162302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7172302aa1bSDebojyoti Ghosh } 7182302aa1bSDebojyoti Ghosh 7192302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/ 7202302aa1bSDebojyoti Ghosh #undef __FUNCT__ 7212302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSReset_GLEE" 7222302aa1bSDebojyoti Ghosh static PetscErrorCode TSReset_GLEE(TS ts) 7232302aa1bSDebojyoti Ghosh { 7242302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7252302aa1bSDebojyoti Ghosh PetscInt s, r; 7262302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 7272302aa1bSDebojyoti Ghosh 7282302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7292302aa1bSDebojyoti Ghosh if (!glee->tableau) PetscFunctionReturn(0); 7302302aa1bSDebojyoti Ghosh s = glee->tableau->s; 7312302aa1bSDebojyoti Ghosh r = glee->tableau->r; 7322302aa1bSDebojyoti Ghosh ierr = VecDestroyVecs(r,&glee->Y);CHKERRQ(ierr); 7332302aa1bSDebojyoti Ghosh ierr = VecDestroyVecs(r,&glee->X);CHKERRQ(ierr); 7342302aa1bSDebojyoti Ghosh ierr = VecDestroyVecs(s,&glee->YStage);CHKERRQ(ierr); 7352302aa1bSDebojyoti Ghosh ierr = VecDestroyVecs(s,&glee->YdotStage);CHKERRQ(ierr); 7362302aa1bSDebojyoti Ghosh ierr = VecDestroy(&glee->Ydot);CHKERRQ(ierr); 7372302aa1bSDebojyoti Ghosh ierr = VecDestroy(&glee->W);CHKERRQ(ierr); 7382302aa1bSDebojyoti Ghosh ierr = PetscFree(glee->work);CHKERRQ(ierr); 7392302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7402302aa1bSDebojyoti Ghosh } 7412302aa1bSDebojyoti Ghosh 7422302aa1bSDebojyoti Ghosh #undef __FUNCT__ 7432302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSDestroy_GLEE" 7442302aa1bSDebojyoti Ghosh static PetscErrorCode TSDestroy_GLEE(TS ts) 7452302aa1bSDebojyoti Ghosh { 7462302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 7472302aa1bSDebojyoti Ghosh 7482302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7492302aa1bSDebojyoti Ghosh ierr = TSReset_GLEE(ts);CHKERRQ(ierr); 7502302aa1bSDebojyoti Ghosh ierr = PetscFree(ts->data);CHKERRQ(ierr); 7512302aa1bSDebojyoti Ghosh ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",NULL);CHKERRQ(ierr); 7522302aa1bSDebojyoti Ghosh ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",NULL);CHKERRQ(ierr); 7532302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7542302aa1bSDebojyoti Ghosh } 7552302aa1bSDebojyoti Ghosh 7562302aa1bSDebojyoti Ghosh #undef __FUNCT__ 7572302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetVecs" 7582302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEEGetVecs(TS ts,DM dm,Vec *Ydot) 7592302aa1bSDebojyoti Ghosh { 7602302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7612302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 7622302aa1bSDebojyoti Ghosh 7632302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7642302aa1bSDebojyoti Ghosh if (Ydot) { 7652302aa1bSDebojyoti Ghosh if (dm && dm != ts->dm) { 7662302aa1bSDebojyoti Ghosh ierr = DMGetNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr); 7672302aa1bSDebojyoti Ghosh } else *Ydot = glee->Ydot; 7682302aa1bSDebojyoti Ghosh } 7692302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7702302aa1bSDebojyoti Ghosh } 7712302aa1bSDebojyoti Ghosh 7722302aa1bSDebojyoti Ghosh 7732302aa1bSDebojyoti Ghosh #undef __FUNCT__ 7742302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEERestoreVecs" 7752302aa1bSDebojyoti Ghosh static PetscErrorCode TSGLEERestoreVecs(TS ts,DM dm,Vec *Ydot) 7762302aa1bSDebojyoti Ghosh { 7772302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 7782302aa1bSDebojyoti Ghosh 7792302aa1bSDebojyoti Ghosh PetscFunctionBegin; 7802302aa1bSDebojyoti Ghosh if (Ydot) { 7812302aa1bSDebojyoti Ghosh if (dm && dm != ts->dm) { 7822302aa1bSDebojyoti Ghosh ierr = DMRestoreNamedGlobalVector(dm,"TSGLEE_Ydot",Ydot);CHKERRQ(ierr); 7832302aa1bSDebojyoti Ghosh } 7842302aa1bSDebojyoti Ghosh } 7852302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 7862302aa1bSDebojyoti Ghosh } 7872302aa1bSDebojyoti Ghosh 7882302aa1bSDebojyoti Ghosh /* 7892302aa1bSDebojyoti Ghosh This defines the nonlinear equation that is to be solved with SNES 7902302aa1bSDebojyoti Ghosh */ 7912302aa1bSDebojyoti Ghosh #undef __FUNCT__ 7922302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormFunction_GLEE" 7932302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormFunction_GLEE(SNES snes,Vec X,Vec F,TS ts) 7942302aa1bSDebojyoti Ghosh { 7952302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 7962302aa1bSDebojyoti Ghosh DM dm,dmsave; 7972302aa1bSDebojyoti Ghosh Vec Ydot; 7982302aa1bSDebojyoti Ghosh PetscReal shift = glee->scoeff / ts->time_step; 7992302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 8002302aa1bSDebojyoti Ghosh 8012302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8022302aa1bSDebojyoti Ghosh ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 8032302aa1bSDebojyoti Ghosh ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr); 8042302aa1bSDebojyoti Ghosh /* Set Ydot = shift*X */ 8052302aa1bSDebojyoti Ghosh ierr = VecCopy(X,Ydot);CHKERRQ(ierr); 8062302aa1bSDebojyoti Ghosh ierr = VecScale(Ydot,shift);CHKERRQ(ierr); 8072302aa1bSDebojyoti Ghosh dmsave = ts->dm; 8082302aa1bSDebojyoti Ghosh ts->dm = dm; 8092302aa1bSDebojyoti Ghosh 8102302aa1bSDebojyoti Ghosh ierr = TSComputeIFunction(ts,glee->stage_time,X,Ydot,F,PETSC_FALSE);CHKERRQ(ierr); 8112302aa1bSDebojyoti Ghosh 8122302aa1bSDebojyoti Ghosh ts->dm = dmsave; 8132302aa1bSDebojyoti Ghosh ierr = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr); 8142302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8152302aa1bSDebojyoti Ghosh } 8162302aa1bSDebojyoti Ghosh 8172302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8182302aa1bSDebojyoti Ghosh #define __FUNCT__ "SNESTSFormJacobian_GLEE" 8192302aa1bSDebojyoti Ghosh static PetscErrorCode SNESTSFormJacobian_GLEE(SNES snes,Vec X,Mat A,Mat B,TS ts) 8202302aa1bSDebojyoti Ghosh { 8212302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 8222302aa1bSDebojyoti Ghosh DM dm,dmsave; 8232302aa1bSDebojyoti Ghosh Vec Ydot; 8242302aa1bSDebojyoti Ghosh PetscReal shift = glee->scoeff / ts->time_step; 8252302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 8262302aa1bSDebojyoti Ghosh 8272302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8282302aa1bSDebojyoti Ghosh ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 8292302aa1bSDebojyoti Ghosh ierr = TSGLEEGetVecs(ts,dm,&Ydot);CHKERRQ(ierr); 8302302aa1bSDebojyoti Ghosh /* glee->Ydot has already been computed in SNESTSFormFunction_GLEE (SNES guarantees this) */ 8312302aa1bSDebojyoti Ghosh dmsave = ts->dm; 8322302aa1bSDebojyoti Ghosh ts->dm = dm; 8332302aa1bSDebojyoti Ghosh 8342302aa1bSDebojyoti Ghosh ierr = TSComputeIJacobian(ts,glee->stage_time,X,Ydot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 8352302aa1bSDebojyoti Ghosh 8362302aa1bSDebojyoti Ghosh ts->dm = dmsave; 8372302aa1bSDebojyoti Ghosh ierr = TSGLEERestoreVecs(ts,dm,&Ydot);CHKERRQ(ierr); 8382302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8392302aa1bSDebojyoti Ghosh } 8402302aa1bSDebojyoti Ghosh 8412302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8422302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMCoarsenHook_TSGLEE" 8432302aa1bSDebojyoti Ghosh static PetscErrorCode DMCoarsenHook_TSGLEE(DM fine,DM coarse,void *ctx) 8442302aa1bSDebojyoti Ghosh { 8452302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8462302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8472302aa1bSDebojyoti Ghosh } 8482302aa1bSDebojyoti Ghosh 8492302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8502302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMRestrictHook_TSGLEE" 8512302aa1bSDebojyoti Ghosh static PetscErrorCode DMRestrictHook_TSGLEE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 8522302aa1bSDebojyoti Ghosh { 8532302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8542302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8552302aa1bSDebojyoti Ghosh } 8562302aa1bSDebojyoti Ghosh 8572302aa1bSDebojyoti Ghosh 8582302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8592302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainHook_TSGLEE" 8602302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainHook_TSGLEE(DM dm,DM subdm,void *ctx) 8612302aa1bSDebojyoti Ghosh { 8622302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8632302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8642302aa1bSDebojyoti Ghosh } 8652302aa1bSDebojyoti Ghosh 8662302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8672302aa1bSDebojyoti Ghosh #define __FUNCT__ "DMSubDomainRestrictHook_TSGLEE" 8682302aa1bSDebojyoti Ghosh static PetscErrorCode DMSubDomainRestrictHook_TSGLEE(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 8692302aa1bSDebojyoti Ghosh { 8702302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8712302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 8722302aa1bSDebojyoti Ghosh } 8732302aa1bSDebojyoti Ghosh 8742302aa1bSDebojyoti Ghosh #undef __FUNCT__ 8752302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetUp_GLEE" 8762302aa1bSDebojyoti Ghosh static PetscErrorCode TSSetUp_GLEE(TS ts) 8772302aa1bSDebojyoti Ghosh { 8782302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 8792302aa1bSDebojyoti Ghosh GLEETableau tab; 8802302aa1bSDebojyoti Ghosh PetscInt s,r; 8812302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 8822302aa1bSDebojyoti Ghosh DM dm; 8832302aa1bSDebojyoti Ghosh 8842302aa1bSDebojyoti Ghosh PetscFunctionBegin; 8852302aa1bSDebojyoti Ghosh if (!glee->tableau) { 8862302aa1bSDebojyoti Ghosh ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr); 8872302aa1bSDebojyoti Ghosh } 8882302aa1bSDebojyoti Ghosh tab = glee->tableau; 8892302aa1bSDebojyoti Ghosh s = tab->s; 8902302aa1bSDebojyoti Ghosh r = tab->r; 8912302aa1bSDebojyoti Ghosh ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->Y);CHKERRQ(ierr); 8922302aa1bSDebojyoti Ghosh ierr = VecDuplicateVecs(ts->vec_sol,r,&glee->X);CHKERRQ(ierr); 8932302aa1bSDebojyoti Ghosh ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YStage);CHKERRQ(ierr); 8942302aa1bSDebojyoti Ghosh ierr = VecDuplicateVecs(ts->vec_sol,s,&glee->YdotStage);CHKERRQ(ierr); 8952302aa1bSDebojyoti Ghosh ierr = VecDuplicate(ts->vec_sol,&glee->Ydot);CHKERRQ(ierr); 8962302aa1bSDebojyoti Ghosh ierr = VecDuplicate(ts->vec_sol,&glee->W);CHKERRQ(ierr); 8972302aa1bSDebojyoti Ghosh ierr = PetscMalloc1(s,&glee->work);CHKERRQ(ierr); 8982302aa1bSDebojyoti Ghosh ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 8992302aa1bSDebojyoti Ghosh if (dm) { 9002302aa1bSDebojyoti Ghosh ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLEE,DMRestrictHook_TSGLEE,ts);CHKERRQ(ierr); 9012302aa1bSDebojyoti Ghosh ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLEE,DMSubDomainRestrictHook_TSGLEE,ts);CHKERRQ(ierr); 9022302aa1bSDebojyoti Ghosh } 9032302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9042302aa1bSDebojyoti Ghosh } 9052302aa1bSDebojyoti Ghosh 9062302aa1bSDebojyoti Ghosh #undef __FUNCT__ 9072302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSStartingMethod_GLEE" 9082302aa1bSDebojyoti Ghosh PetscErrorCode TSStartingMethod_GLEE(TS ts) 9092302aa1bSDebojyoti Ghosh { 9102302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 9116e2e232bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 9126e2e232bSDebojyoti Ghosh PetscInt r=tab->r,i; 9136e2e232bSDebojyoti Ghosh PetscReal *S=tab->S; 9142302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 9152302aa1bSDebojyoti Ghosh 9162302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9172302aa1bSDebojyoti Ghosh for (i=0; i<r; i++) { 9182302aa1bSDebojyoti Ghosh ierr = VecZeroEntries(glee->Y[i]);CHKERRQ(ierr); 9192302aa1bSDebojyoti Ghosh ierr = VecAXPY(glee->Y[i],S[i],ts->vec_sol);CHKERRQ(ierr); 9202302aa1bSDebojyoti Ghosh } 9212302aa1bSDebojyoti Ghosh 9222302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9232302aa1bSDebojyoti Ghosh } 9242302aa1bSDebojyoti Ghosh 9252302aa1bSDebojyoti Ghosh 9262302aa1bSDebojyoti Ghosh /*------------------------------------------------------------*/ 9272302aa1bSDebojyoti Ghosh 9282302aa1bSDebojyoti Ghosh #undef __FUNCT__ 9292302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSSetFromOptions_GLEE" 9306b234c65SEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TS ts) 9312302aa1bSDebojyoti Ghosh { 9322302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 9332302aa1bSDebojyoti Ghosh char gleetype[256]; 9342302aa1bSDebojyoti Ghosh 9352302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9362302aa1bSDebojyoti Ghosh ierr = PetscOptionsHead(PetscOptionsObject,"GLEE ODE solver options");CHKERRQ(ierr); 9372302aa1bSDebojyoti Ghosh { 9382302aa1bSDebojyoti Ghosh GLEETableauLink link; 9392302aa1bSDebojyoti Ghosh PetscInt count,choice; 9402302aa1bSDebojyoti Ghosh PetscBool flg; 9412302aa1bSDebojyoti Ghosh const char **namelist; 9422302aa1bSDebojyoti Ghosh 9432302aa1bSDebojyoti Ghosh ierr = PetscStrncpy(gleetype,TSGLEEDefaultType,sizeof(gleetype));CHKERRQ(ierr); 9442302aa1bSDebojyoti Ghosh for (link=GLEETableauList,count=0; link; link=link->next,count++) ; 9452302aa1bSDebojyoti Ghosh ierr = PetscMalloc1(count,&namelist);CHKERRQ(ierr); 9462302aa1bSDebojyoti Ghosh for (link=GLEETableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 9472302aa1bSDebojyoti Ghosh ierr = PetscOptionsEList("-ts_glee_type","Family of GLEE method","TSGLEESetType",(const char*const*)namelist,count,gleetype,&choice,&flg);CHKERRQ(ierr); 9482302aa1bSDebojyoti Ghosh ierr = TSGLEESetType(ts,flg ? namelist[choice] : gleetype);CHKERRQ(ierr); 9492302aa1bSDebojyoti Ghosh ierr = PetscFree(namelist);CHKERRQ(ierr); 9502302aa1bSDebojyoti Ghosh } 9512302aa1bSDebojyoti Ghosh ierr = PetscOptionsTail();CHKERRQ(ierr); 9522302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9532302aa1bSDebojyoti Ghosh } 9542302aa1bSDebojyoti Ghosh 9552302aa1bSDebojyoti Ghosh #undef __FUNCT__ 9562302aa1bSDebojyoti Ghosh #define __FUNCT__ "PetscFormatRealArray" 9572302aa1bSDebojyoti Ghosh static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 9582302aa1bSDebojyoti Ghosh { 9592302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 9602302aa1bSDebojyoti Ghosh PetscInt i; 9612302aa1bSDebojyoti Ghosh size_t left,count; 9622302aa1bSDebojyoti Ghosh char *p; 9632302aa1bSDebojyoti Ghosh 9642302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9652302aa1bSDebojyoti Ghosh for (i=0,p=buf,left=len; i<n; i++) { 9662302aa1bSDebojyoti Ghosh ierr = PetscSNPrintfCount(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 9672302aa1bSDebojyoti Ghosh if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 9682302aa1bSDebojyoti Ghosh left -= count; 9692302aa1bSDebojyoti Ghosh p += count; 9702302aa1bSDebojyoti Ghosh *p++ = ' '; 9712302aa1bSDebojyoti Ghosh } 9722302aa1bSDebojyoti Ghosh p[i ? 0 : -1] = 0; 9732302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 9742302aa1bSDebojyoti Ghosh } 9752302aa1bSDebojyoti Ghosh 9762302aa1bSDebojyoti Ghosh #undef __FUNCT__ 9772302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSView_GLEE" 9782302aa1bSDebojyoti Ghosh static PetscErrorCode TSView_GLEE(TS ts,PetscViewer viewer) 9792302aa1bSDebojyoti Ghosh { 9802302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 9812302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 9822302aa1bSDebojyoti Ghosh PetscBool iascii; 9832302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 9842302aa1bSDebojyoti Ghosh TSAdapt adapt; 9852302aa1bSDebojyoti Ghosh 9862302aa1bSDebojyoti Ghosh PetscFunctionBegin; 9872302aa1bSDebojyoti Ghosh ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9882302aa1bSDebojyoti Ghosh if (iascii) { 9892302aa1bSDebojyoti Ghosh TSGLEEType gleetype; 9902302aa1bSDebojyoti Ghosh char buf[512]; 9912302aa1bSDebojyoti Ghosh ierr = TSGLEEGetType(ts,&gleetype);CHKERRQ(ierr); 9922302aa1bSDebojyoti Ghosh ierr = PetscViewerASCIIPrintf(viewer," GLEE %s\n",gleetype);CHKERRQ(ierr); 9932302aa1bSDebojyoti Ghosh ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 9942302aa1bSDebojyoti Ghosh ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 9952302aa1bSDebojyoti Ghosh /* Note: print out r as well */ 9962302aa1bSDebojyoti Ghosh } 9972302aa1bSDebojyoti Ghosh ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 9982302aa1bSDebojyoti Ghosh ierr = TSAdaptView(adapt,viewer);CHKERRQ(ierr); 9992302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10002302aa1bSDebojyoti Ghosh } 10012302aa1bSDebojyoti Ghosh 10022302aa1bSDebojyoti Ghosh #undef __FUNCT__ 10032302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSLoad_GLEE" 10042302aa1bSDebojyoti Ghosh static PetscErrorCode TSLoad_GLEE(TS ts,PetscViewer viewer) 10052302aa1bSDebojyoti Ghosh { 10062302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 10072302aa1bSDebojyoti Ghosh SNES snes; 10082302aa1bSDebojyoti Ghosh TSAdapt tsadapt; 10092302aa1bSDebojyoti Ghosh 10102302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10112302aa1bSDebojyoti Ghosh ierr = TSGetAdapt(ts,&tsadapt);CHKERRQ(ierr); 10122302aa1bSDebojyoti Ghosh ierr = TSAdaptLoad(tsadapt,viewer);CHKERRQ(ierr); 10132302aa1bSDebojyoti Ghosh ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 10142302aa1bSDebojyoti Ghosh ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 10152302aa1bSDebojyoti Ghosh /* function and Jacobian context for SNES when used with TS is always ts object */ 10162302aa1bSDebojyoti Ghosh ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 10172302aa1bSDebojyoti Ghosh ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 10182302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10192302aa1bSDebojyoti Ghosh } 10202302aa1bSDebojyoti Ghosh 10212302aa1bSDebojyoti Ghosh #undef __FUNCT__ 10222302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType" 10232302aa1bSDebojyoti Ghosh /*@C 10242302aa1bSDebojyoti Ghosh TSGLEESetType - Set the type of GLEE scheme 10252302aa1bSDebojyoti Ghosh 10262302aa1bSDebojyoti Ghosh Logically collective 10272302aa1bSDebojyoti Ghosh 10282302aa1bSDebojyoti Ghosh Input Parameter: 10292302aa1bSDebojyoti Ghosh + ts - timestepping context 10302302aa1bSDebojyoti Ghosh - gleetype - type of GLEE-scheme 10312302aa1bSDebojyoti Ghosh 10322302aa1bSDebojyoti Ghosh Level: intermediate 10332302aa1bSDebojyoti Ghosh 10342302aa1bSDebojyoti Ghosh .seealso: TSGLEEGetType(), TSGLEE 10352302aa1bSDebojyoti Ghosh @*/ 10362302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType(TS ts,TSGLEEType gleetype) 10372302aa1bSDebojyoti Ghosh { 10382302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 10392302aa1bSDebojyoti Ghosh 10402302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10412302aa1bSDebojyoti Ghosh PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10422302aa1bSDebojyoti Ghosh ierr = PetscTryMethod(ts,"TSGLEESetType_C",(TS,TSGLEEType),(ts,gleetype));CHKERRQ(ierr); 10432302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10442302aa1bSDebojyoti Ghosh } 10452302aa1bSDebojyoti Ghosh 10462302aa1bSDebojyoti Ghosh #undef __FUNCT__ 10472302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType" 10482302aa1bSDebojyoti Ghosh /*@C 10492302aa1bSDebojyoti Ghosh TSGLEEGetType - Get the type of GLEE scheme 10502302aa1bSDebojyoti Ghosh 10512302aa1bSDebojyoti Ghosh Logically collective 10522302aa1bSDebojyoti Ghosh 10532302aa1bSDebojyoti Ghosh Input Parameter: 10542302aa1bSDebojyoti Ghosh . ts - timestepping context 10552302aa1bSDebojyoti Ghosh 10562302aa1bSDebojyoti Ghosh Output Parameter: 10572302aa1bSDebojyoti Ghosh . gleetype - type of GLEE-scheme 10582302aa1bSDebojyoti Ghosh 10592302aa1bSDebojyoti Ghosh Level: intermediate 10602302aa1bSDebojyoti Ghosh 10612302aa1bSDebojyoti Ghosh .seealso: TSGLEESetType() 10622302aa1bSDebojyoti Ghosh @*/ 10632302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType *gleetype) 10642302aa1bSDebojyoti Ghosh { 10652302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 10662302aa1bSDebojyoti Ghosh 10672302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10682302aa1bSDebojyoti Ghosh PetscValidHeaderSpecific(ts,TS_CLASSID,1); 10692302aa1bSDebojyoti Ghosh ierr = PetscUseMethod(ts,"TSGLEEGetType_C",(TS,TSGLEEType*),(ts,gleetype));CHKERRQ(ierr); 10702302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10712302aa1bSDebojyoti Ghosh } 10722302aa1bSDebojyoti Ghosh 10732302aa1bSDebojyoti Ghosh #undef __FUNCT__ 10742302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEEGetType_GLEE" 10752302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEEGetType_GLEE(TS ts,TSGLEEType *gleetype) 10762302aa1bSDebojyoti Ghosh { 10772302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10782302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 10792302aa1bSDebojyoti Ghosh 10802302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10812302aa1bSDebojyoti Ghosh if (!glee->tableau) { 10822302aa1bSDebojyoti Ghosh ierr = TSGLEESetType(ts,TSGLEEDefaultType);CHKERRQ(ierr); 10832302aa1bSDebojyoti Ghosh } 10842302aa1bSDebojyoti Ghosh *gleetype = glee->tableau->name; 10852302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 10862302aa1bSDebojyoti Ghosh } 10872302aa1bSDebojyoti Ghosh #undef __FUNCT__ 10882302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGLEESetType_GLEE" 10892302aa1bSDebojyoti Ghosh PetscErrorCode TSGLEESetType_GLEE(TS ts,TSGLEEType gleetype) 10902302aa1bSDebojyoti Ghosh { 10912302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 10922302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 10932302aa1bSDebojyoti Ghosh PetscBool match; 10942302aa1bSDebojyoti Ghosh GLEETableauLink link; 10952302aa1bSDebojyoti Ghosh 10962302aa1bSDebojyoti Ghosh PetscFunctionBegin; 10972302aa1bSDebojyoti Ghosh if (glee->tableau) { 10982302aa1bSDebojyoti Ghosh ierr = PetscStrcmp(glee->tableau->name,gleetype,&match);CHKERRQ(ierr); 10992302aa1bSDebojyoti Ghosh if (match) PetscFunctionReturn(0); 11002302aa1bSDebojyoti Ghosh } 11012302aa1bSDebojyoti Ghosh for (link = GLEETableauList; link; link=link->next) { 11022302aa1bSDebojyoti Ghosh ierr = PetscStrcmp(link->tab.name,gleetype,&match);CHKERRQ(ierr); 11032302aa1bSDebojyoti Ghosh if (match) { 11042302aa1bSDebojyoti Ghosh ierr = TSReset_GLEE(ts);CHKERRQ(ierr); 11052302aa1bSDebojyoti Ghosh glee->tableau = &link->tab; 11062302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 11072302aa1bSDebojyoti Ghosh } 11082302aa1bSDebojyoti Ghosh } 11092302aa1bSDebojyoti Ghosh SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",gleetype); 11102302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 11112302aa1bSDebojyoti Ghosh } 11122302aa1bSDebojyoti Ghosh 11132302aa1bSDebojyoti Ghosh #undef __FUNCT__ 11142302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSGetStages_GLEE" 11152302aa1bSDebojyoti Ghosh static PetscErrorCode TSGetStages_GLEE(TS ts,PetscInt *ns,Vec **Y) 11162302aa1bSDebojyoti Ghosh { 11172302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 11182302aa1bSDebojyoti Ghosh 11192302aa1bSDebojyoti Ghosh PetscFunctionBegin; 11202302aa1bSDebojyoti Ghosh *ns = glee->tableau->s; 11212302aa1bSDebojyoti Ghosh if(Y) *Y = glee->Y; 11222302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 11232302aa1bSDebojyoti Ghosh } 11242302aa1bSDebojyoti Ghosh 11252302aa1bSDebojyoti Ghosh #undef __FUNCT__ 1126b2bf4f3aSDebojyoti Ghosh #define __FUNCT__ "TSGetSolutionComponents_GLEE" 1127b2bf4f3aSDebojyoti Ghosh PetscErrorCode TSGetSolutionComponents_GLEE(TS ts,PetscInt *n,Vec *Y) 11282302aa1bSDebojyoti Ghosh { 11292302aa1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 11302302aa1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 11316e2e232bSDebojyoti Ghosh PetscErrorCode ierr; 11322302aa1bSDebojyoti Ghosh 11332302aa1bSDebojyoti Ghosh PetscFunctionBegin; 11344cdd57e5SDebojyoti Ghosh if (!Y) *n = tab->r; 11352302aa1bSDebojyoti Ghosh else { 11364cdd57e5SDebojyoti Ghosh if ((*n >= 0) && (*n < tab->r)) { 11374cdd57e5SDebojyoti Ghosh ierr = VecCopy(glee->Y[*n],*Y); CHKERRQ(ierr); 11389657682dSDebojyoti Ghosh } else SETERRQ3(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Second argument (%d) out of range[%d,%d].",*n,0,tab->r-1); 11392302aa1bSDebojyoti Ghosh } 11402302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 11412302aa1bSDebojyoti Ghosh } 11422302aa1bSDebojyoti Ghosh 11434cdd57e5SDebojyoti Ghosh #undef __FUNCT__ 11444cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetAuxSolution_GLEE" 11454cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetAuxSolution_GLEE(TS ts,Vec *X) 11464cdd57e5SDebojyoti Ghosh { 11474cdd57e5SDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 11484cdd57e5SDebojyoti Ghosh GLEETableau tab = glee->tableau; 11494cdd57e5SDebojyoti Ghosh PetscReal *F = tab->Fembed; 11504cdd57e5SDebojyoti Ghosh PetscInt r = tab->r; 11514cdd57e5SDebojyoti Ghosh Vec *Y = glee->Y; 11524cdd57e5SDebojyoti Ghosh PetscErrorCode ierr; 11534cdd57e5SDebojyoti Ghosh 11544cdd57e5SDebojyoti Ghosh PetscFunctionBegin; 11554cdd57e5SDebojyoti Ghosh ierr = VecZeroEntries(*X);CHKERRQ(ierr); 11564cdd57e5SDebojyoti Ghosh ierr = VecMAXPY((*X),r,F,Y);CHKERRQ(ierr); 11574cdd57e5SDebojyoti Ghosh PetscFunctionReturn(0); 11584cdd57e5SDebojyoti Ghosh } 11594cdd57e5SDebojyoti Ghosh 11604cdd57e5SDebojyoti Ghosh #undef __FUNCT__ 11614cdd57e5SDebojyoti Ghosh #define __FUNCT__ "TSGetTimeError_GLEE" 11624cdd57e5SDebojyoti Ghosh PetscErrorCode TSGetTimeError_GLEE(TS ts,Vec *X) 11634cdd57e5SDebojyoti Ghosh { 11644cdd57e5SDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 11654cdd57e5SDebojyoti Ghosh GLEETableau tab = glee->tableau; 11664cdd57e5SDebojyoti Ghosh PetscReal *F = tab->Ferror; 11674cdd57e5SDebojyoti Ghosh PetscInt r = tab->r; 11684cdd57e5SDebojyoti Ghosh Vec *Y = glee->Y; 11694cdd57e5SDebojyoti Ghosh PetscErrorCode ierr; 11704cdd57e5SDebojyoti Ghosh 11714cdd57e5SDebojyoti Ghosh PetscFunctionBegin; 11724cdd57e5SDebojyoti Ghosh ierr = VecZeroEntries(*X);CHKERRQ(ierr); 11734cdd57e5SDebojyoti Ghosh ierr = VecMAXPY((*X),r,F,Y);CHKERRQ(ierr); 11744cdd57e5SDebojyoti Ghosh PetscFunctionReturn(0); 11754cdd57e5SDebojyoti Ghosh } 11764cdd57e5SDebojyoti Ghosh 117757df6a1bSDebojyoti Ghosh #undef __FUNCT__ 117857df6a1bSDebojyoti Ghosh #define __FUNCT__ "TSSetTimeError_GLEE" 117957df6a1bSDebojyoti Ghosh PetscErrorCode TSSetTimeError_GLEE(TS ts,Vec X) 118057df6a1bSDebojyoti Ghosh { 118157df6a1bSDebojyoti Ghosh TS_GLEE *glee = (TS_GLEE*)ts->data; 118257df6a1bSDebojyoti Ghosh GLEETableau tab = glee->tableau; 118357df6a1bSDebojyoti Ghosh PetscReal *S = tab->Serror; 118457df6a1bSDebojyoti Ghosh PetscInt r = tab->r,i; 118557df6a1bSDebojyoti Ghosh Vec *Y = glee->Y; 118657df6a1bSDebojyoti Ghosh PetscErrorCode ierr; 118757df6a1bSDebojyoti Ghosh 118857df6a1bSDebojyoti Ghosh PetscFunctionBegin; 11899657682dSDebojyoti Ghosh if (r != 2) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSSetTimeError_GLEE not supported for '%s' with r=%D.",tab->name,tab->r); 119057df6a1bSDebojyoti Ghosh for (i=1; i<r; i++) { 119157df6a1bSDebojyoti Ghosh ierr = VecCopy(ts->vec_sol,Y[i]); CHKERRQ(ierr); 119257df6a1bSDebojyoti Ghosh ierr = VecAXPBY(Y[i],S[0],S[1],X); CHKERRQ(ierr); 119357df6a1bSDebojyoti Ghosh } 119457df6a1bSDebojyoti Ghosh PetscFunctionReturn(0); 119557df6a1bSDebojyoti Ghosh } 119657df6a1bSDebojyoti Ghosh 11972302aa1bSDebojyoti Ghosh /* ------------------------------------------------------------ */ 11982302aa1bSDebojyoti Ghosh /*MC 11992302aa1bSDebojyoti Ghosh TSGLEE - ODE and DAE solver using General Linear with Error Estimation schemes 12002302aa1bSDebojyoti Ghosh 12012302aa1bSDebojyoti Ghosh The user should provide the right hand side of the equation 12022302aa1bSDebojyoti Ghosh using TSSetRHSFunction(). 12032302aa1bSDebojyoti Ghosh 12042302aa1bSDebojyoti Ghosh Notes: 12052302aa1bSDebojyoti Ghosh The default is TSGLEE35, it can be changed with TSGLEESetType() or -ts_glee_type 12062302aa1bSDebojyoti Ghosh 12072302aa1bSDebojyoti Ghosh Level: beginner 12082302aa1bSDebojyoti Ghosh 12092302aa1bSDebojyoti Ghosh .seealso: TSCreate(), TS, TSSetType(), TSGLEESetType(), TSGLEEGetType(), 12102302aa1bSDebojyoti Ghosh TSGLEE23, TTSGLEE24, TSGLEE35, TSGLEE25I, TSGLEEEXRK2A, 12112302aa1bSDebojyoti Ghosh TSGLEERK32G1, TSGLEERK285EX, TSGLEEType, TSGLEERegister() 12122302aa1bSDebojyoti Ghosh 12132302aa1bSDebojyoti Ghosh M*/ 12142302aa1bSDebojyoti Ghosh #undef __FUNCT__ 12152302aa1bSDebojyoti Ghosh #define __FUNCT__ "TSCreate_GLEE" 12162302aa1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSCreate_GLEE(TS ts) 12172302aa1bSDebojyoti Ghosh { 12182302aa1bSDebojyoti Ghosh TS_GLEE *th; 12192302aa1bSDebojyoti Ghosh PetscErrorCode ierr; 12202302aa1bSDebojyoti Ghosh 12212302aa1bSDebojyoti Ghosh PetscFunctionBegin; 12222302aa1bSDebojyoti Ghosh #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 12232302aa1bSDebojyoti Ghosh ierr = TSGLEEInitializePackage();CHKERRQ(ierr); 12242302aa1bSDebojyoti Ghosh #endif 12252302aa1bSDebojyoti Ghosh 12262302aa1bSDebojyoti Ghosh ts->ops->reset = TSReset_GLEE; 12272302aa1bSDebojyoti Ghosh ts->ops->destroy = TSDestroy_GLEE; 12282302aa1bSDebojyoti Ghosh ts->ops->view = TSView_GLEE; 12292302aa1bSDebojyoti Ghosh ts->ops->load = TSLoad_GLEE; 12302302aa1bSDebojyoti Ghosh ts->ops->setup = TSSetUp_GLEE; 12312302aa1bSDebojyoti Ghosh ts->ops->step = TSStep_GLEE; 12322302aa1bSDebojyoti Ghosh ts->ops->interpolate = TSInterpolate_GLEE; 12332302aa1bSDebojyoti Ghosh ts->ops->evaluatestep = TSEvaluateStep_GLEE; 12342302aa1bSDebojyoti Ghosh ts->ops->setfromoptions = TSSetFromOptions_GLEE; 12352302aa1bSDebojyoti Ghosh ts->ops->getstages = TSGetStages_GLEE; 12362302aa1bSDebojyoti Ghosh ts->ops->snesfunction = SNESTSFormFunction_GLEE; 12372302aa1bSDebojyoti Ghosh ts->ops->snesjacobian = SNESTSFormJacobian_GLEE; 1238b2bf4f3aSDebojyoti Ghosh ts->ops->getsolutioncomponents = TSGetSolutionComponents_GLEE; 12394cdd57e5SDebojyoti Ghosh ts->ops->getauxsolution = TSGetAuxSolution_GLEE; 12404cdd57e5SDebojyoti Ghosh ts->ops->gettimeerror = TSGetTimeError_GLEE; 124157df6a1bSDebojyoti Ghosh ts->ops->settimeerror = TSSetTimeError_GLEE; 12422302aa1bSDebojyoti Ghosh ts->ops->startingmethod = TSStartingMethod_GLEE; 12432302aa1bSDebojyoti Ghosh 12442302aa1bSDebojyoti Ghosh ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 12452302aa1bSDebojyoti Ghosh ts->data = (void*)th; 12462302aa1bSDebojyoti Ghosh 12472302aa1bSDebojyoti Ghosh ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEEGetType_C",TSGLEEGetType_GLEE);CHKERRQ(ierr); 12482302aa1bSDebojyoti Ghosh ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLEESetType_C",TSGLEESetType_GLEE);CHKERRQ(ierr); 12492302aa1bSDebojyoti Ghosh PetscFunctionReturn(0); 12502302aa1bSDebojyoti Ghosh } 1251