18a381b04SJed Brown /* 28a381b04SJed Brown Code for timestepping with additive Runge-Kutta IMEX method 38a381b04SJed Brown 48a381b04SJed Brown Notes: 58a381b04SJed Brown The general system is written as 68a381b04SJed Brown 78a381b04SJed Brown F(t,X,Xdot) = G(t,X) 88a381b04SJed Brown 98a381b04SJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 108a381b04SJed Brown 118a381b04SJed Brown */ 128a381b04SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 138a381b04SJed Brown 144f385281SJed Brown static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2E; 158a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled; 168a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized; 178a381b04SJed Brown 188a381b04SJed Brown typedef struct _ARKTableau *ARKTableau; 198a381b04SJed Brown struct _ARKTableau { 208a381b04SJed Brown char *name; 214f385281SJed Brown PetscInt order; /* Classical approximation order of the method */ 224f385281SJed Brown PetscInt s; /* Number of stages */ 234f385281SJed Brown PetscInt pinterp; /* Interpolation order */ 244f385281SJed Brown PetscReal *At,*bt,*ct; /* Stiff tableau */ 258a381b04SJed Brown PetscReal *A,*b,*c; /* Non-stiff tableau */ 26cd652676SJed Brown PetscReal *binterpt,*binterp; /* Dense output formula */ 278a381b04SJed Brown }; 288a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 298a381b04SJed Brown struct _ARKTableauLink { 308a381b04SJed Brown struct _ARKTableau tab; 318a381b04SJed Brown ARKTableauLink next; 328a381b04SJed Brown }; 338a381b04SJed Brown static ARKTableauLink ARKTableauList; 348a381b04SJed Brown 358a381b04SJed Brown typedef struct { 368a381b04SJed Brown ARKTableau tableau; 378a381b04SJed Brown Vec *Y; /* States computed during the step */ 388a381b04SJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 398a381b04SJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 408a381b04SJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 418a381b04SJed Brown Vec Work; /* Generic work vector */ 428a381b04SJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 438a381b04SJed Brown PetscScalar *work; /* Scalar work */ 448a381b04SJed Brown PetscReal shift; 458a381b04SJed Brown PetscReal stage_time; 464cc180ffSJed Brown PetscBool imex; 478a381b04SJed Brown } TS_ARKIMEX; 488a381b04SJed Brown 4964f491ddSJed Brown /*MC 5064f491ddSJed Brown TSARKIMEX2D - Second order ARK IMEX scheme with L-stable implicit part. 5164f491ddSJed Brown 5264f491ddSJed Brown This method has one explicit stage and two implicit stages. 5364f491ddSJed Brown 54*b330ce4dSSatish Balay Level: advanced 55*b330ce4dSSatish Balay 5664f491ddSJed Brown .seealso: TSARKIMEX 5764f491ddSJed Brown M*/ 5864f491ddSJed Brown /*MC 5964f491ddSJed Brown TSARKIMEX2E - Second order ARK IMEX scheme with L-stable implicit part. 6064f491ddSJed Brown 6164f491ddSJed Brown This method has one explicit stage and two implicit stages. It is is an optimal method developed by Emil Constantinescu. 6264f491ddSJed Brown 63*b330ce4dSSatish Balay Level: advanced 64*b330ce4dSSatish Balay 6564f491ddSJed Brown .seealso: TSARKIMEX 6664f491ddSJed Brown M*/ 6764f491ddSJed Brown /*MC 6864f491ddSJed Brown TSARKIMEX3 - Third order ARK IMEX scheme with L-stable implicit part. 6964f491ddSJed Brown 7064f491ddSJed Brown This method has one explicit stage and three implicit stages. 7164f491ddSJed Brown 7264f491ddSJed Brown References: 7364f491ddSJed Brown Kennedy and Carpenter 2003. 7464f491ddSJed Brown 75*b330ce4dSSatish Balay Level: advanced 76*b330ce4dSSatish Balay 7764f491ddSJed Brown .seealso: TSARKIMEX 7864f491ddSJed Brown M*/ 7964f491ddSJed Brown /*MC 8064f491ddSJed Brown TSARKIMEX4 - Fourth order ARK IMEX scheme with L-stable implicit part. 8164f491ddSJed Brown 8264f491ddSJed Brown This method has one explicit stage and four implicit stages. 8364f491ddSJed Brown 8464f491ddSJed Brown References: 8564f491ddSJed Brown Kennedy and Carpenter 2003. 8664f491ddSJed Brown 87*b330ce4dSSatish Balay Level: advanced 88*b330ce4dSSatish Balay 8964f491ddSJed Brown .seealso: TSARKIMEX 9064f491ddSJed Brown M*/ 9164f491ddSJed Brown /*MC 9264f491ddSJed Brown TSARKIMEX5 - Fifth order ARK IMEX scheme with L-stable implicit part. 9364f491ddSJed Brown 9464f491ddSJed Brown This method has one explicit stage and five implicit stages. 9564f491ddSJed Brown 9664f491ddSJed Brown References: 9764f491ddSJed Brown Kennedy and Carpenter 2003. 9864f491ddSJed Brown 99*b330ce4dSSatish Balay Level: advanced 100*b330ce4dSSatish Balay 10164f491ddSJed Brown .seealso: TSARKIMEX 10264f491ddSJed Brown M*/ 10364f491ddSJed Brown 1048a381b04SJed Brown #undef __FUNCT__ 1058a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll" 1068a381b04SJed Brown /*@C 1078a381b04SJed Brown TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 1088a381b04SJed Brown 109fca742c7SJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 1108a381b04SJed Brown 1118a381b04SJed Brown Level: advanced 1128a381b04SJed Brown 1138a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all 1148a381b04SJed Brown 1158a381b04SJed Brown .seealso: TSARKIMEXRegisterDestroy() 1168a381b04SJed Brown @*/ 1178a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 1188a381b04SJed Brown { 1198a381b04SJed Brown PetscErrorCode ierr; 1208a381b04SJed Brown 1218a381b04SJed Brown PetscFunctionBegin; 1228a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 1238a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 1248a381b04SJed Brown 1258a381b04SJed Brown { 1268a381b04SJed Brown const PetscReal 1278a381b04SJed Brown A[3][3] = {{0,0,0}, 1288a381b04SJed Brown {0.41421356237309504880,0,0}, 1298a381b04SJed Brown {0.75,0.25,0}}, 1308a381b04SJed Brown At[3][3] = {{0,0,0}, 1318a381b04SJed Brown {0.12132034355964257320,0.29289321881345247560,0}, 13206db7b1cSJed Brown {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, 13306db7b1cSJed Brown binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 13406db7b1cSJed Brown ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 1358a381b04SJed Brown } 13606db7b1cSJed Brown { /* Optimal for linear implicit part */ 13703403c7fSJed Brown const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 138a3a57f36SJed Brown A[3][3] = {{0,0,0}, 139a3a57f36SJed Brown {2-s2,0,0}, 140a3a57f36SJed Brown {(3-2*s2)/6,(3+2*s2)/6,0}}, 141a3a57f36SJed Brown At[3][3] = {{0,0,0}, 142a3a57f36SJed Brown {1-1/s2,1-1/s2,0}, 143cd652676SJed Brown {1/(2*s2),1/(2*s2),1-1/s2}}, 144cd652676SJed Brown binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 145cd652676SJed Brown ierr = TSARKIMEXRegister(TSARKIMEX2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 146a3a57f36SJed Brown } 147a3a57f36SJed Brown { 148a3a57f36SJed Brown const PetscReal 149a3a57f36SJed Brown A[4][4] = {{0,0,0,0}, 1504040e9f2SJed Brown {1767732205903./2027836641118.,0,0,0}, 1514040e9f2SJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 1524040e9f2SJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 153a3a57f36SJed Brown At[4][4] = {{0,0,0,0}, 1544040e9f2SJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 1554040e9f2SJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 1564040e9f2SJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 1574040e9f2SJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 1584040e9f2SJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 1594040e9f2SJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 1604040e9f2SJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 161cd652676SJed Brown ierr = TSARKIMEXRegister(TSARKIMEX3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 162a3a57f36SJed Brown } 163a3a57f36SJed Brown { 164a3a57f36SJed Brown const PetscReal 165a3a57f36SJed Brown A[6][6] = {{0,0,0,0,0,0}, 166a3a57f36SJed Brown {1./2,0,0,0,0,0}, 1674040e9f2SJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 1684040e9f2SJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 1694040e9f2SJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 1704040e9f2SJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 171a3a57f36SJed Brown At[6][6] = {{0,0,0,0,0,0}, 172a3a57f36SJed Brown {1./4,1./4,0,0,0,0}, 1734040e9f2SJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 1744040e9f2SJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 1754040e9f2SJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 1764040e9f2SJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 1774040e9f2SJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 178cd652676SJed Brown {0,0,0}, 1794040e9f2SJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 1804040e9f2SJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 1814040e9f2SJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 1824040e9f2SJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 183cd652676SJed Brown ierr = TSARKIMEXRegister(TSARKIMEX4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 184a3a57f36SJed Brown } 185a3a57f36SJed Brown { 186a3a57f36SJed Brown const PetscReal 187a3a57f36SJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 188a3a57f36SJed Brown {41./100,0,0,0,0,0,0,0}, 1894040e9f2SJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 1904040e9f2SJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 1914040e9f2SJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 1924040e9f2SJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 1934040e9f2SJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 1944040e9f2SJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 195a3a57f36SJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 1964040e9f2SJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 1974040e9f2SJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 1984040e9f2SJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 1994040e9f2SJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 2004040e9f2SJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 2014040e9f2SJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 2024040e9f2SJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 2034040e9f2SJed Brown binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.}, 204cd652676SJed Brown {0 , 0 , 0 }, 205cd652676SJed Brown {0 , 0 , 0 }, 2064040e9f2SJed Brown {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.}, 2074040e9f2SJed Brown {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.}, 2084040e9f2SJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.}, 2094040e9f2SJed Brown {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.}, 2104040e9f2SJed Brown {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }}; 211cd652676SJed Brown ierr = TSARKIMEXRegister(TSARKIMEX5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 212a3a57f36SJed Brown } 213a3a57f36SJed Brown 2148a381b04SJed Brown PetscFunctionReturn(0); 2158a381b04SJed Brown } 2168a381b04SJed Brown 2178a381b04SJed Brown #undef __FUNCT__ 2188a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy" 2198a381b04SJed Brown /*@C 2208a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 2218a381b04SJed Brown 2228a381b04SJed Brown Not Collective 2238a381b04SJed Brown 2248a381b04SJed Brown Level: advanced 2258a381b04SJed Brown 2268a381b04SJed Brown .keywords: TSARKIMEX, register, destroy 2278a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic() 2288a381b04SJed Brown @*/ 2298a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 2308a381b04SJed Brown { 2318a381b04SJed Brown PetscErrorCode ierr; 2328a381b04SJed Brown ARKTableauLink link; 2338a381b04SJed Brown 2348a381b04SJed Brown PetscFunctionBegin; 2358a381b04SJed Brown while ((link = ARKTableauList)) { 2368a381b04SJed Brown ARKTableau t = &link->tab; 2378a381b04SJed Brown ARKTableauList = link->next; 2388a381b04SJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 239cd652676SJed Brown ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 2408a381b04SJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 2418a381b04SJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 2428a381b04SJed Brown } 2438a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 2448a381b04SJed Brown PetscFunctionReturn(0); 2458a381b04SJed Brown } 2468a381b04SJed Brown 2478a381b04SJed Brown #undef __FUNCT__ 2488a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage" 2498a381b04SJed Brown /*@C 2508a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 2518a381b04SJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 2528a381b04SJed Brown when using static libraries. 2538a381b04SJed Brown 2548a381b04SJed Brown Input Parameter: 2558a381b04SJed Brown path - The dynamic library path, or PETSC_NULL 2568a381b04SJed Brown 2578a381b04SJed Brown Level: developer 2588a381b04SJed Brown 2598a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package 2608a381b04SJed Brown .seealso: PetscInitialize() 2618a381b04SJed Brown @*/ 2628a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[]) 2638a381b04SJed Brown { 2648a381b04SJed Brown PetscErrorCode ierr; 2658a381b04SJed Brown 2668a381b04SJed Brown PetscFunctionBegin; 2678a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 2688a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 2698a381b04SJed Brown ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 2708a381b04SJed Brown ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 2718a381b04SJed Brown PetscFunctionReturn(0); 2728a381b04SJed Brown } 2738a381b04SJed Brown 2748a381b04SJed Brown #undef __FUNCT__ 2758a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage" 2768a381b04SJed Brown /*@C 2778a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 2788a381b04SJed Brown called from PetscFinalize(). 2798a381b04SJed Brown 2808a381b04SJed Brown Level: developer 2818a381b04SJed Brown 2828a381b04SJed Brown .keywords: Petsc, destroy, package 2838a381b04SJed Brown .seealso: PetscFinalize() 2848a381b04SJed Brown @*/ 2858a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 2868a381b04SJed Brown { 2878a381b04SJed Brown PetscErrorCode ierr; 2888a381b04SJed Brown 2898a381b04SJed Brown PetscFunctionBegin; 2908a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 2918a381b04SJed Brown ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 2928a381b04SJed Brown PetscFunctionReturn(0); 2938a381b04SJed Brown } 2948a381b04SJed Brown 2958a381b04SJed Brown #undef __FUNCT__ 2968a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister" 297cd652676SJed Brown /*@C 298cd652676SJed Brown TSARKIMEXRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 299cd652676SJed Brown 300cd652676SJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 301cd652676SJed Brown 302cd652676SJed Brown Input Parameters: 303cd652676SJed Brown + name - identifier for method 304cd652676SJed Brown . order - approximation order of method 305cd652676SJed Brown . s - number of stages, this is the dimension of the matrices below 306cd652676SJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 307cd652676SJed Brown . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At) 308cd652676SJed Brown . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At) 309cd652676SJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 310cd652676SJed Brown . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At) 311cd652676SJed Brown . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A) 312cd652676SJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 313cd652676SJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 314cd652676SJed Brown - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt) 315cd652676SJed Brown 316cd652676SJed Brown Notes: 317cd652676SJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 318cd652676SJed Brown 319cd652676SJed Brown Level: advanced 320cd652676SJed Brown 321cd652676SJed Brown .keywords: TS, register 322cd652676SJed Brown 323cd652676SJed Brown .seealso: TSARKIMEX 324cd652676SJed Brown @*/ 3258a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s, 3268a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 327cd652676SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 328cd652676SJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 3298a381b04SJed Brown { 3308a381b04SJed Brown PetscErrorCode ierr; 3318a381b04SJed Brown ARKTableauLink link; 3328a381b04SJed Brown ARKTableau t; 3338a381b04SJed Brown PetscInt i,j; 3348a381b04SJed Brown 3358a381b04SJed Brown PetscFunctionBegin; 3368a381b04SJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 337cd652676SJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 3388a381b04SJed Brown t = &link->tab; 3398a381b04SJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 3408a381b04SJed Brown t->order = order; 3418a381b04SJed Brown t->s = s; 3428a381b04SJed Brown ierr = PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 3438a381b04SJed Brown ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 3448a381b04SJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 3458a381b04SJed Brown if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);} 3468a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 3478a381b04SJed Brown if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);} 3488a381b04SJed Brown else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i]; 3498a381b04SJed Brown if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);} 3508a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 3518a381b04SJed Brown if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);} 3528a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 3534f385281SJed Brown t->pinterp = pinterp; 354cd652676SJed Brown ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr); 355cd652676SJed Brown ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 356cd652676SJed Brown ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 3578a381b04SJed Brown link->next = ARKTableauList; 3588a381b04SJed Brown ARKTableauList = link; 3598a381b04SJed Brown PetscFunctionReturn(0); 3608a381b04SJed Brown } 3618a381b04SJed Brown 3628a381b04SJed Brown #undef __FUNCT__ 3638a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX" 3648a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 3658a381b04SJed Brown { 3668a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 3678a381b04SJed Brown ARKTableau tab = ark->tableau; 3688a381b04SJed Brown const PetscInt s = tab->s; 3698a381b04SJed Brown const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c; 370406d0ec2SJed Brown PetscScalar *w = ark->work; 3718a381b04SJed Brown Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z; 3728a381b04SJed Brown SNES snes; 3738a381b04SJed Brown PetscInt i,j,its,lits; 374cdbf8f93SLisandro Dalcin PetscReal next_time_step; 3758a381b04SJed Brown PetscReal h,t; 3768a381b04SJed Brown PetscErrorCode ierr; 3778a381b04SJed Brown 3788a381b04SJed Brown PetscFunctionBegin; 3798a381b04SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 380cdbf8f93SLisandro Dalcin next_time_step = ts->time_step; 381cdbf8f93SLisandro Dalcin h = ts->time_step; 3828a381b04SJed Brown t = ts->ptime; 3838a381b04SJed Brown 3848a381b04SJed Brown for (i=0; i<s; i++) { 3858a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 3868a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 3878a381b04SJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 3888a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 3898a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 3908a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 3918a381b04SJed Brown } else { 3928a381b04SJed Brown ark->stage_time = t + h*ct[i]; 3938a381b04SJed Brown ark->shift = 1./(h*At[i*s+i]); 3948a381b04SJed Brown /* Affine part */ 3958a381b04SJed Brown ierr = VecZeroEntries(W);CHKERRQ(ierr); 3968a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 3978a381b04SJed Brown ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 3988a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 3998a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 4004f385281SJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 4014f385281SJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 4028a381b04SJed Brown /* Initial guess taken from last stage */ 4038a381b04SJed Brown ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr); 4048a381b04SJed Brown ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 4058a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 4068a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 4078a381b04SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 4088a381b04SJed Brown } 4098a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 4104cc180ffSJed Brown ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr); 4114cc180ffSJed Brown if (ark->imex) { 4128a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 4134cc180ffSJed Brown } else { 4144cc180ffSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 4154cc180ffSJed Brown } 4168a381b04SJed Brown } 417a17b68daSJed Brown for (j=0; j<s; j++) w[j] = -h*bt[j]; 4188a381b04SJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 4198a381b04SJed Brown for (j=0; j<s; j++) w[j] = h*b[j]; 4208a381b04SJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 4218a381b04SJed Brown 4228a381b04SJed Brown ts->ptime += ts->time_step; 423cdbf8f93SLisandro Dalcin ts->time_step = next_time_step; 4248a381b04SJed Brown ts->steps++; 4258a381b04SJed Brown PetscFunctionReturn(0); 4268a381b04SJed Brown } 4278a381b04SJed Brown 428cd652676SJed Brown #undef __FUNCT__ 429cd652676SJed Brown #define __FUNCT__ "TSInterpolate_ARKIMEX" 430cd652676SJed Brown static PetscErrorCode TSInterpolate_ARKIMEX(TS ts,PetscReal itime,Vec X) 431cd652676SJed Brown { 432cd652676SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 4334f385281SJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 4344f385281SJed Brown PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */ 435cd652676SJed Brown PetscScalar *bt,*b; 436cd652676SJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 437cd652676SJed Brown PetscErrorCode ierr; 438cd652676SJed Brown 439cd652676SJed Brown PetscFunctionBegin; 440cd652676SJed Brown if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSARKIMEX %s does not have an interpolation formula",ark->tableau->name); 441cd652676SJed Brown ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr); 442cd652676SJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 4434f385281SJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 444cd652676SJed Brown for (i=0; i<s; i++) { 4454f385281SJed Brown bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0; 4464f385281SJed Brown b[i] += ts->time_step * B[i*pinterp+j] * tt; 447cd652676SJed Brown } 448cd652676SJed Brown } 449cd652676SJed Brown if (ark->tableau->At[0*s+0] != 0.0) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_SUP,"First stage not explicit so starting stage not saved"); 450cd652676SJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 451cd652676SJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 452cd652676SJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 453cd652676SJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 454cd652676SJed Brown PetscFunctionReturn(0); 455cd652676SJed Brown } 456cd652676SJed Brown 4578a381b04SJed Brown /*------------------------------------------------------------*/ 4588a381b04SJed Brown #undef __FUNCT__ 4598a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX" 4608a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 4618a381b04SJed Brown { 4628a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 4638a381b04SJed Brown PetscInt s; 4648a381b04SJed Brown PetscErrorCode ierr; 4658a381b04SJed Brown 4668a381b04SJed Brown PetscFunctionBegin; 4678a381b04SJed Brown if (!ark->tableau) PetscFunctionReturn(0); 4688a381b04SJed Brown s = ark->tableau->s; 4698a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 4708a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 4718a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 4728a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 4738a381b04SJed Brown ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 4748a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 4758a381b04SJed Brown ierr = PetscFree(ark->work);CHKERRQ(ierr); 4768a381b04SJed Brown PetscFunctionReturn(0); 4778a381b04SJed Brown } 4788a381b04SJed Brown 4798a381b04SJed Brown #undef __FUNCT__ 4808a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX" 4818a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 4828a381b04SJed Brown { 4838a381b04SJed Brown PetscErrorCode ierr; 4848a381b04SJed Brown 4858a381b04SJed Brown PetscFunctionBegin; 4868a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 4878a381b04SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 4888a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr); 4898a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr); 4908a381b04SJed Brown PetscFunctionReturn(0); 4918a381b04SJed Brown } 4928a381b04SJed Brown 4938a381b04SJed Brown /* 4948a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 4958a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 4968a381b04SJed Brown */ 4978a381b04SJed Brown #undef __FUNCT__ 4988a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 4998a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 5008a381b04SJed Brown { 5018a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 5028a381b04SJed Brown PetscErrorCode ierr; 5038a381b04SJed Brown 5048a381b04SJed Brown PetscFunctionBegin; 5058a381b04SJed Brown ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 5064cc180ffSJed Brown ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr); 5078a381b04SJed Brown PetscFunctionReturn(0); 5088a381b04SJed Brown } 5098a381b04SJed Brown 5108a381b04SJed Brown #undef __FUNCT__ 5118a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 5128a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 5138a381b04SJed Brown { 5148a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 5158a381b04SJed Brown PetscErrorCode ierr; 5168a381b04SJed Brown 5178a381b04SJed Brown PetscFunctionBegin; 5188a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 5198a381b04SJed Brown ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 5208a381b04SJed Brown PetscFunctionReturn(0); 5218a381b04SJed Brown } 5228a381b04SJed Brown 5238a381b04SJed Brown #undef __FUNCT__ 5248a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX" 5258a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 5268a381b04SJed Brown { 5278a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 5288a381b04SJed Brown ARKTableau tab = ark->tableau; 5298a381b04SJed Brown PetscInt s = tab->s; 5308a381b04SJed Brown PetscErrorCode ierr; 5318a381b04SJed Brown 5328a381b04SJed Brown PetscFunctionBegin; 5338a381b04SJed Brown if (!ark->tableau) { 534e24355feSJed Brown ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 5358a381b04SJed Brown } 5368a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 5378a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 5388a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 5398a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 5408a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 5418a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 5428a381b04SJed Brown ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr); 5438a381b04SJed Brown PetscFunctionReturn(0); 5448a381b04SJed Brown } 5458a381b04SJed Brown /*------------------------------------------------------------*/ 5468a381b04SJed Brown 5478a381b04SJed Brown #undef __FUNCT__ 5488a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 5498a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts) 5508a381b04SJed Brown { 5514cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 5528a381b04SJed Brown PetscErrorCode ierr; 5538a381b04SJed Brown char arktype[256]; 5548a381b04SJed Brown 5558a381b04SJed Brown PetscFunctionBegin; 5568a381b04SJed Brown ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr); 5578a381b04SJed Brown { 5588a381b04SJed Brown ARKTableauLink link; 5598a381b04SJed Brown PetscInt count,choice; 5608a381b04SJed Brown PetscBool flg; 5618a381b04SJed Brown const char **namelist; 5628a381b04SJed Brown ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr); 5638a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 5648a381b04SJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 5658a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 5668a381b04SJed Brown ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 5678a381b04SJed Brown ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 5688a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 5694cc180ffSJed Brown flg = (PetscBool)!ark->imex; 5704cc180ffSJed Brown ierr = PetscOptionsBool("-ts_arkimex_fully_implicit","Solve the problem fully implicitly","TSARKIMEXSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 5714cc180ffSJed Brown ark->imex = (PetscBool)!flg; 572d52bd9f3SBarry Smith ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 5738a381b04SJed Brown } 5748a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 5758a381b04SJed Brown PetscFunctionReturn(0); 5768a381b04SJed Brown } 5778a381b04SJed Brown 5788a381b04SJed Brown #undef __FUNCT__ 5798a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray" 5808a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 5818a381b04SJed Brown { 582257d2499SJed Brown PetscErrorCode ierr; 5838a381b04SJed Brown int i,left,count; 5848a381b04SJed Brown char *p; 5858a381b04SJed Brown 5868a381b04SJed Brown PetscFunctionBegin; 5878a381b04SJed Brown for (i=0,p=buf,left=(int)len; i<n; i++) { 588257d2499SJed Brown ierr = PetscSNPrintf(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 5898a381b04SJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 5908a381b04SJed Brown left -= count; 5918a381b04SJed Brown p += count; 5928a381b04SJed Brown *p++ = ' '; 5938a381b04SJed Brown } 5948a381b04SJed Brown p[i ? 0 : -1] = 0; 5958a381b04SJed Brown PetscFunctionReturn(0); 5968a381b04SJed Brown } 5978a381b04SJed Brown 5988a381b04SJed Brown #undef __FUNCT__ 5998a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX" 6008a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 6018a381b04SJed Brown { 6028a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 6038a381b04SJed Brown ARKTableau tab = ark->tableau; 6048a381b04SJed Brown PetscBool iascii; 6058a381b04SJed Brown PetscErrorCode ierr; 6068a381b04SJed Brown 6078a381b04SJed Brown PetscFunctionBegin; 6088a381b04SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6098a381b04SJed Brown if (iascii) { 6108a381b04SJed Brown const TSARKIMEXType arktype; 6118a381b04SJed Brown char buf[512]; 6128a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 6138a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 6148a381b04SJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 61531f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 61631f6fcc0SJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 61731f6fcc0SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 6188a381b04SJed Brown } 619d52bd9f3SBarry Smith ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 6208a381b04SJed Brown PetscFunctionReturn(0); 6218a381b04SJed Brown } 6228a381b04SJed Brown 6238a381b04SJed Brown #undef __FUNCT__ 6248a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType" 6258a381b04SJed Brown /*@C 6268a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 6278a381b04SJed Brown 6288a381b04SJed Brown Logically collective 6298a381b04SJed Brown 6308a381b04SJed Brown Input Parameter: 6318a381b04SJed Brown + ts - timestepping context 6328a381b04SJed Brown - arktype - type of ARK-IMEX scheme 6338a381b04SJed Brown 6348a381b04SJed Brown Level: intermediate 6358a381b04SJed Brown 6368a381b04SJed Brown .seealso: TSARKIMEXGetType() 6378a381b04SJed Brown @*/ 6388a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype) 6398a381b04SJed Brown { 6408a381b04SJed Brown PetscErrorCode ierr; 6418a381b04SJed Brown 6428a381b04SJed Brown PetscFunctionBegin; 6438a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6448a381b04SJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 6458a381b04SJed Brown PetscFunctionReturn(0); 6468a381b04SJed Brown } 6478a381b04SJed Brown 6488a381b04SJed Brown #undef __FUNCT__ 6498a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType" 6508a381b04SJed Brown /*@C 6518a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 6528a381b04SJed Brown 6538a381b04SJed Brown Logically collective 6548a381b04SJed Brown 6558a381b04SJed Brown Input Parameter: 6568a381b04SJed Brown . ts - timestepping context 6578a381b04SJed Brown 6588a381b04SJed Brown Output Parameter: 6598a381b04SJed Brown . arktype - type of ARK-IMEX scheme 6608a381b04SJed Brown 6618a381b04SJed Brown Level: intermediate 6628a381b04SJed Brown 6638a381b04SJed Brown .seealso: TSARKIMEXGetType() 6648a381b04SJed Brown @*/ 6658a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype) 6668a381b04SJed Brown { 6678a381b04SJed Brown PetscErrorCode ierr; 6688a381b04SJed Brown 6698a381b04SJed Brown PetscFunctionBegin; 6708a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6718a381b04SJed Brown ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 6728a381b04SJed Brown PetscFunctionReturn(0); 6738a381b04SJed Brown } 6748a381b04SJed Brown 6754cc180ffSJed Brown #undef __FUNCT__ 6764cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit" 6774cc180ffSJed Brown /*@C 6784cc180ffSJed Brown TSARKIMEXSetFullyImplicit - Solve both parts of the equation implicitly 6794cc180ffSJed Brown 6804cc180ffSJed Brown Logically collective 6814cc180ffSJed Brown 6824cc180ffSJed Brown Input Parameter: 6834cc180ffSJed Brown + ts - timestepping context 6844cc180ffSJed Brown - flg - PETSC_TRUE for fully implicit 6854cc180ffSJed Brown 6864cc180ffSJed Brown Level: intermediate 6874cc180ffSJed Brown 6884cc180ffSJed Brown .seealso: TSARKIMEXGetType() 6894cc180ffSJed Brown @*/ 6904cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit(TS ts,PetscBool flg) 6914cc180ffSJed Brown { 6924cc180ffSJed Brown PetscErrorCode ierr; 6934cc180ffSJed Brown 6944cc180ffSJed Brown PetscFunctionBegin; 6954cc180ffSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 6964cc180ffSJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 6974cc180ffSJed Brown PetscFunctionReturn(0); 6984cc180ffSJed Brown } 6994cc180ffSJed Brown 7008a381b04SJed Brown EXTERN_C_BEGIN 7018a381b04SJed Brown #undef __FUNCT__ 7028a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 7038a381b04SJed Brown PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype) 7048a381b04SJed Brown { 7058a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7068a381b04SJed Brown PetscErrorCode ierr; 7078a381b04SJed Brown 7088a381b04SJed Brown PetscFunctionBegin; 7098a381b04SJed Brown if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);} 7108a381b04SJed Brown *arktype = ark->tableau->name; 7118a381b04SJed Brown PetscFunctionReturn(0); 7128a381b04SJed Brown } 7138a381b04SJed Brown #undef __FUNCT__ 7148a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 7158a381b04SJed Brown PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype) 7168a381b04SJed Brown { 7178a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7188a381b04SJed Brown PetscErrorCode ierr; 7198a381b04SJed Brown PetscBool match; 7208a381b04SJed Brown ARKTableauLink link; 7218a381b04SJed Brown 7228a381b04SJed Brown PetscFunctionBegin; 7238a381b04SJed Brown if (ark->tableau) { 7248a381b04SJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 7258a381b04SJed Brown if (match) PetscFunctionReturn(0); 7268a381b04SJed Brown } 7278a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 7288a381b04SJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 7298a381b04SJed Brown if (match) { 7308a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 7318a381b04SJed Brown ark->tableau = &link->tab; 7328a381b04SJed Brown PetscFunctionReturn(0); 7338a381b04SJed Brown } 7348a381b04SJed Brown } 7358a381b04SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 7368a381b04SJed Brown PetscFunctionReturn(0); 7378a381b04SJed Brown } 7384cc180ffSJed Brown #undef __FUNCT__ 7394cc180ffSJed Brown #define __FUNCT__ "TSARKIMEXSetFullyImplicit_ARKIMEX" 7404cc180ffSJed Brown PetscErrorCode TSARKIMEXSetFullyImplicit_ARKIMEX(TS ts,PetscBool flg) 7414cc180ffSJed Brown { 7424cc180ffSJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 7434cc180ffSJed Brown 7444cc180ffSJed Brown PetscFunctionBegin; 7454cc180ffSJed Brown ark->imex = (PetscBool)!flg; 7464cc180ffSJed Brown PetscFunctionReturn(0); 7474cc180ffSJed Brown } 7488a381b04SJed Brown EXTERN_C_END 7498a381b04SJed Brown 7508a381b04SJed Brown /* ------------------------------------------------------------ */ 7518a381b04SJed Brown /*MC 7528a381b04SJed Brown TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes 7538a381b04SJed Brown 754fca742c7SJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 755fca742c7SJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 756fca742c7SJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 757fca742c7SJed Brown 758fca742c7SJed Brown Notes: 759c8058688SBarry Smith The default is TSARKIMEX2E, it can be changed with TSARKIMEXSetType() or -ts_arkimex_type 760c8058688SBarry Smith 761fca742c7SJed Brown This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 762fca742c7SJed Brown 7638a381b04SJed Brown Level: beginner 7648a381b04SJed Brown 765c8058688SBarry Smith .seealso: TSCreate(), TS, TSSetType(), TSARKIMEXSetType(), TSARKIMEXGetType(), TSARKIMEXSetFullyImplicit(), TSARKIMEX2D, TTSARKIMEX2E, TSARKIMEX3, 76626885f59SBarry Smith TSARKIMEX4, TSARKIMEX5, TSARKIMEXType, TSARKIMEXRegister() 7678a381b04SJed Brown 7688a381b04SJed Brown M*/ 7698a381b04SJed Brown EXTERN_C_BEGIN 7708a381b04SJed Brown #undef __FUNCT__ 7718a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX" 7728a381b04SJed Brown PetscErrorCode TSCreate_ARKIMEX(TS ts) 7738a381b04SJed Brown { 7748a381b04SJed Brown TS_ARKIMEX *th; 7758a381b04SJed Brown PetscErrorCode ierr; 7768a381b04SJed Brown 7778a381b04SJed Brown PetscFunctionBegin; 7788a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 7798a381b04SJed Brown ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr); 7808a381b04SJed Brown #endif 7818a381b04SJed Brown 7828a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 7838a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 7848a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 7858a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 7868a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 787cd652676SJed Brown ts->ops->interpolate = TSInterpolate_ARKIMEX; 7888a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 7898a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 7908a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 7918a381b04SJed Brown 7928a381b04SJed Brown ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr); 7938a381b04SJed Brown ts->data = (void*)th; 7944cc180ffSJed Brown th->imex = PETSC_TRUE; 7958a381b04SJed Brown 7968a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 7978a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 7984cc180ffSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetFullyImplicit_C","TSARKIMEXSetFullyImplicit_ARKIMEX",TSARKIMEXSetFullyImplicit_ARKIMEX);CHKERRQ(ierr); 7998a381b04SJed Brown PetscFunctionReturn(0); 8008a381b04SJed Brown } 8018a381b04SJed Brown EXTERN_C_END 802