1*e27a552bSJed Brown /* 2*e27a552bSJed Brown Code for timestepping with additive Runge-Kutta IMEX method 3*e27a552bSJed Brown 4*e27a552bSJed Brown Notes: 5*e27a552bSJed Brown The general system is written as 6*e27a552bSJed Brown 7*e27a552bSJed Brown F(t,X,Xdot) = G(t,X) 8*e27a552bSJed Brown 9*e27a552bSJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 10*e27a552bSJed Brown 11*e27a552bSJed Brown */ 12*e27a552bSJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 13*e27a552bSJed Brown 14*e27a552bSJed Brown static const TSRosWType TSRosWDefault = TSROSW2E; 15*e27a552bSJed Brown static PetscBool TSRosWRegisterAllCalled; 16*e27a552bSJed Brown static PetscBool TSRosWPackageInitialized; 17*e27a552bSJed Brown 18*e27a552bSJed Brown typedef struct _ARKTableau *ARKTableau; 19*e27a552bSJed Brown struct _ARKTableau { 20*e27a552bSJed Brown char *name; 21*e27a552bSJed Brown PetscInt order; /* Classical approximation order of the method */ 22*e27a552bSJed Brown PetscInt s; /* Number of stages */ 23*e27a552bSJed Brown PetscInt pinterp; /* Interpolation order */ 24*e27a552bSJed Brown PetscReal *At,*bt,*ct; /* Stiff tableau */ 25*e27a552bSJed Brown PetscReal *A,*b,*c; /* Non-stiff tableau */ 26*e27a552bSJed Brown PetscReal *binterpt,*binterp; /* Dense output formula */ 27*e27a552bSJed Brown }; 28*e27a552bSJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 29*e27a552bSJed Brown struct _ARKTableauLink { 30*e27a552bSJed Brown struct _ARKTableau tab; 31*e27a552bSJed Brown ARKTableauLink next; 32*e27a552bSJed Brown }; 33*e27a552bSJed Brown static ARKTableauLink ARKTableauList; 34*e27a552bSJed Brown 35*e27a552bSJed Brown typedef struct { 36*e27a552bSJed Brown ARKTableau tableau; 37*e27a552bSJed Brown Vec *Y; /* States computed during the step */ 38*e27a552bSJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 39*e27a552bSJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 40*e27a552bSJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 41*e27a552bSJed Brown Vec Work; /* Generic work vector */ 42*e27a552bSJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 43*e27a552bSJed Brown PetscScalar *work; /* Scalar work */ 44*e27a552bSJed Brown PetscReal shift; 45*e27a552bSJed Brown PetscReal stage_time; 46*e27a552bSJed Brown PetscBool imex; 47*e27a552bSJed Brown } TS_RosW; 48*e27a552bSJed Brown 49*e27a552bSJed Brown #undef __FUNCT__ 50*e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterAll" 51*e27a552bSJed Brown /*@C 52*e27a552bSJed Brown TSRosWRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSRosW 53*e27a552bSJed Brown 54*e27a552bSJed Brown Not Collective, but should be called by all processes which will need the schemes to be registered 55*e27a552bSJed Brown 56*e27a552bSJed Brown Level: advanced 57*e27a552bSJed Brown 58*e27a552bSJed Brown .keywords: TS, TSRosW, register, all 59*e27a552bSJed Brown 60*e27a552bSJed Brown .seealso: TSRosWRegisterDestroy() 61*e27a552bSJed Brown @*/ 62*e27a552bSJed Brown PetscErrorCode TSRosWRegisterAll(void) 63*e27a552bSJed Brown { 64*e27a552bSJed Brown PetscErrorCode ierr; 65*e27a552bSJed Brown 66*e27a552bSJed Brown PetscFunctionBegin; 67*e27a552bSJed Brown if (TSRosWRegisterAllCalled) PetscFunctionReturn(0); 68*e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_TRUE; 69*e27a552bSJed Brown 70*e27a552bSJed Brown { 71*e27a552bSJed Brown const PetscReal 72*e27a552bSJed Brown A[3][3] = {{0,0,0}, 73*e27a552bSJed Brown {0.41421356237309504880,0,0}, 74*e27a552bSJed Brown {0.75,0.25,0}}, 75*e27a552bSJed Brown At[3][3] = {{0,0,0}, 76*e27a552bSJed Brown {0.12132034355964257320,0.29289321881345247560,0}, 77*e27a552bSJed Brown {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}, 78*e27a552bSJed Brown binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 79*e27a552bSJed Brown ierr = TSRosWRegister(TSROSW2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 80*e27a552bSJed Brown } 81*e27a552bSJed Brown { /* Optimal for linear implicit part */ 82*e27a552bSJed Brown const PetscReal s2 = PetscSqrtReal((PetscReal)2.0), 83*e27a552bSJed Brown A[3][3] = {{0,0,0}, 84*e27a552bSJed Brown {2-s2,0,0}, 85*e27a552bSJed Brown {(3-2*s2)/6,(3+2*s2)/6,0}}, 86*e27a552bSJed Brown At[3][3] = {{0,0,0}, 87*e27a552bSJed Brown {1-1/s2,1-1/s2,0}, 88*e27a552bSJed Brown {1/(2*s2),1/(2*s2),1-1/s2}}, 89*e27a552bSJed Brown binterpt[3][2] = {{1,-0.5},{0,0},{0,0.5}}; 90*e27a552bSJed Brown ierr = TSRosWRegister(TSROSW2E,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 91*e27a552bSJed Brown } 92*e27a552bSJed Brown { 93*e27a552bSJed Brown const PetscReal 94*e27a552bSJed Brown A[4][4] = {{0,0,0,0}, 95*e27a552bSJed Brown {1767732205903./2027836641118.,0,0,0}, 96*e27a552bSJed Brown {5535828885825./10492691773637.,788022342437./10882634858940.,0,0}, 97*e27a552bSJed Brown {6485989280629./16251701735622.,-4246266847089./9704473918619.,10755448449292./10357097424841.,0}}, 98*e27a552bSJed Brown At[4][4] = {{0,0,0,0}, 99*e27a552bSJed Brown {1767732205903./4055673282236.,1767732205903./4055673282236.,0,0}, 100*e27a552bSJed Brown {2746238789719./10658868560708.,-640167445237./6845629431997.,1767732205903./4055673282236.,0}, 101*e27a552bSJed Brown {1471266399579./7840856788654.,-4482444167858./7529755066697.,11266239266428./11593286722821.,1767732205903./4055673282236.}}, 102*e27a552bSJed Brown binterpt[4][2] = {{4655552711362./22874653954995., -215264564351./13552729205753.}, 103*e27a552bSJed Brown {-18682724506714./9892148508045.,17870216137069./13817060693119.}, 104*e27a552bSJed Brown {34259539580243./13192909600954.,-28141676662227./17317692491321.}, 105*e27a552bSJed Brown {584795268549./6622622206610., 2508943948391./7218656332882.}}; 106*e27a552bSJed Brown ierr = TSRosWRegister(TSROSW3,3,4,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,2,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 107*e27a552bSJed Brown } 108*e27a552bSJed Brown { 109*e27a552bSJed Brown const PetscReal 110*e27a552bSJed Brown A[6][6] = {{0,0,0,0,0,0}, 111*e27a552bSJed Brown {1./2,0,0,0,0,0}, 112*e27a552bSJed Brown {13861./62500.,6889./62500.,0,0,0,0}, 113*e27a552bSJed Brown {-116923316275./2393684061468.,-2731218467317./15368042101831.,9408046702089./11113171139209.,0,0,0}, 114*e27a552bSJed Brown {-451086348788./2902428689909.,-2682348792572./7519795681897.,12662868775082./11960479115383.,3355817975965./11060851509271.,0,0}, 115*e27a552bSJed Brown {647845179188./3216320057751.,73281519250./8382639484533.,552539513391./3454668386233.,3354512671639./8306763924573.,4040./17871.,0}}, 116*e27a552bSJed Brown At[6][6] = {{0,0,0,0,0,0}, 117*e27a552bSJed Brown {1./4,1./4,0,0,0,0}, 118*e27a552bSJed Brown {8611./62500.,-1743./31250.,1./4,0,0,0}, 119*e27a552bSJed Brown {5012029./34652500.,-654441./2922500.,174375./388108.,1./4,0,0}, 120*e27a552bSJed Brown {15267082809./155376265600.,-71443401./120774400.,730878875./902184768.,2285395./8070912.,1./4,0}, 121*e27a552bSJed Brown {82889./524892.,0,15625./83664.,69875./102672.,-2260./8211,1./4}}, 122*e27a552bSJed Brown binterpt[6][3] = {{6943876665148./7220017795957.,-54480133./30881146.,6818779379841./7100303317025.}, 123*e27a552bSJed Brown {0,0,0}, 124*e27a552bSJed Brown {7640104374378./9702883013639.,-11436875./14766696.,2173542590792./12501825683035.}, 125*e27a552bSJed Brown {-20649996744609./7521556579894.,174696575./18121608.,-31592104683404./5083833661969.}, 126*e27a552bSJed Brown {8854892464581./2390941311638.,-12120380./966161.,61146701046299./7138195549469.}, 127*e27a552bSJed Brown {-11397109935349./6675773540249.,3843./706.,-17219254887155./4939391667607.}}; 128*e27a552bSJed Brown ierr = TSRosWRegister(TSROSW4,4,6,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 129*e27a552bSJed Brown } 130*e27a552bSJed Brown { 131*e27a552bSJed Brown const PetscReal 132*e27a552bSJed Brown A[8][8] = {{0,0,0,0,0,0,0,0}, 133*e27a552bSJed Brown {41./100,0,0,0,0,0,0,0}, 134*e27a552bSJed Brown {367902744464./2072280473677.,677623207551./8224143866563.,0,0,0,0,0,0}, 135*e27a552bSJed Brown {1268023523408./10340822734521.,0,1029933939417./13636558850479.,0,0,0,0,0}, 136*e27a552bSJed Brown {14463281900351./6315353703477.,0,66114435211212./5879490589093.,-54053170152839./4284798021562.,0,0,0,0}, 137*e27a552bSJed Brown {14090043504691./34967701212078.,0,15191511035443./11219624916014.,-18461159152457./12425892160975.,-281667163811./9011619295870.,0,0,0}, 138*e27a552bSJed Brown {19230459214898./13134317526959.,0,21275331358303./2942455364971.,-38145345988419./4862620318723.,-1./8,-1./8,0,0}, 139*e27a552bSJed Brown {-19977161125411./11928030595625.,0,-40795976796054./6384907823539.,177454434618887./12078138498510.,782672205425./8267701900261.,-69563011059811./9646580694205.,7356628210526./4942186776405.,0}}, 140*e27a552bSJed Brown At[8][8] = {{0,0,0,0,0,0,0,0}, 141*e27a552bSJed Brown {41./200.,41./200.,0,0,0,0,0,0}, 142*e27a552bSJed Brown {41./400.,-567603406766./11931857230679.,41./200.,0,0,0,0,0}, 143*e27a552bSJed Brown {683785636431./9252920307686.,0,-110385047103./1367015193373.,41./200.,0,0,0,0}, 144*e27a552bSJed Brown {3016520224154./10081342136671.,0,30586259806659./12414158314087.,-22760509404356./11113319521817.,41./200.,0,0,0}, 145*e27a552bSJed Brown {218866479029./1489978393911.,0,638256894668./5436446318841.,-1179710474555./5321154724896.,-60928119172./8023461067671.,41./200.,0,0}, 146*e27a552bSJed Brown {1020004230633./5715676835656.,0,25762820946817./25263940353407.,-2161375909145./9755907335909.,-211217309593./5846859502534.,-4269925059573./7827059040749.,41./200,0}, 147*e27a552bSJed Brown {-872700587467./9133579230613.,0,0,22348218063261./9555858737531.,-1143369518992./8141816002931.,-39379526789629./19018526304540.,32727382324388./42900044865799.,41./200.}}, 148*e27a552bSJed Brown binterpt[8][3] = {{-17674230611817./10670229744614. , 43486358583215./12773830924787. , -9257016797708./5021505065439.}, 149*e27a552bSJed Brown {0 , 0 , 0 }, 150*e27a552bSJed Brown {0 , 0 , 0 }, 151*e27a552bSJed Brown {65168852399939./7868540260826. , -91478233927265./11067650958493., 26096422576131./11239449250142.}, 152*e27a552bSJed Brown {15494834004392./5936557850923. , -79368583304911./10890268929626., 92396832856987./20362823103730.}, 153*e27a552bSJed Brown {-99329723586156./26959484932159., -12239297817655./9152339842473. , 30029262896817./10175596800299.}, 154*e27a552bSJed Brown {-19024464361622./5461577185407. , 115839755401235./10719374521269., -26136350496073./3983972220547.}, 155*e27a552bSJed Brown {-6511271360970./6095937251113. , 5843115559534./2180450260947. , -5289405421727./3760307252460. }}; 156*e27a552bSJed Brown ierr = TSRosWRegister(TSROSW5,5,8,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL,3,binterpt[0],PETSC_NULL);CHKERRQ(ierr); 157*e27a552bSJed Brown } 158*e27a552bSJed Brown 159*e27a552bSJed Brown PetscFunctionReturn(0); 160*e27a552bSJed Brown } 161*e27a552bSJed Brown 162*e27a552bSJed Brown #undef __FUNCT__ 163*e27a552bSJed Brown #define __FUNCT__ "TSRosWRegisterDestroy" 164*e27a552bSJed Brown /*@C 165*e27a552bSJed Brown TSRosWRegisterDestroy - Frees the list of schemes that were registered by TSRosWRegister(). 166*e27a552bSJed Brown 167*e27a552bSJed Brown Not Collective 168*e27a552bSJed Brown 169*e27a552bSJed Brown Level: advanced 170*e27a552bSJed Brown 171*e27a552bSJed Brown .keywords: TSRosW, register, destroy 172*e27a552bSJed Brown .seealso: TSRosWRegister(), TSRosWRegisterAll(), TSRosWRegisterDynamic() 173*e27a552bSJed Brown @*/ 174*e27a552bSJed Brown PetscErrorCode TSRosWRegisterDestroy(void) 175*e27a552bSJed Brown { 176*e27a552bSJed Brown PetscErrorCode ierr; 177*e27a552bSJed Brown ARKTableauLink link; 178*e27a552bSJed Brown 179*e27a552bSJed Brown PetscFunctionBegin; 180*e27a552bSJed Brown while ((link = ARKTableauList)) { 181*e27a552bSJed Brown ARKTableau t = &link->tab; 182*e27a552bSJed Brown ARKTableauList = link->next; 183*e27a552bSJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 184*e27a552bSJed Brown ierr = PetscFree2(t->binterpt,t->binterp);CHKERRQ(ierr); 185*e27a552bSJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 186*e27a552bSJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 187*e27a552bSJed Brown } 188*e27a552bSJed Brown TSRosWRegisterAllCalled = PETSC_FALSE; 189*e27a552bSJed Brown PetscFunctionReturn(0); 190*e27a552bSJed Brown } 191*e27a552bSJed Brown 192*e27a552bSJed Brown #undef __FUNCT__ 193*e27a552bSJed Brown #define __FUNCT__ "TSRosWInitializePackage" 194*e27a552bSJed Brown /*@C 195*e27a552bSJed Brown TSRosWInitializePackage - This function initializes everything in the TSRosW package. It is called 196*e27a552bSJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_RosW() 197*e27a552bSJed Brown when using static libraries. 198*e27a552bSJed Brown 199*e27a552bSJed Brown Input Parameter: 200*e27a552bSJed Brown path - The dynamic library path, or PETSC_NULL 201*e27a552bSJed Brown 202*e27a552bSJed Brown Level: developer 203*e27a552bSJed Brown 204*e27a552bSJed Brown .keywords: TS, TSRosW, initialize, package 205*e27a552bSJed Brown .seealso: PetscInitialize() 206*e27a552bSJed Brown @*/ 207*e27a552bSJed Brown PetscErrorCode TSRosWInitializePackage(const char path[]) 208*e27a552bSJed Brown { 209*e27a552bSJed Brown PetscErrorCode ierr; 210*e27a552bSJed Brown 211*e27a552bSJed Brown PetscFunctionBegin; 212*e27a552bSJed Brown if (TSRosWPackageInitialized) PetscFunctionReturn(0); 213*e27a552bSJed Brown TSRosWPackageInitialized = PETSC_TRUE; 214*e27a552bSJed Brown ierr = TSRosWRegisterAll();CHKERRQ(ierr); 215*e27a552bSJed Brown ierr = PetscRegisterFinalize(TSRosWFinalizePackage);CHKERRQ(ierr); 216*e27a552bSJed Brown PetscFunctionReturn(0); 217*e27a552bSJed Brown } 218*e27a552bSJed Brown 219*e27a552bSJed Brown #undef __FUNCT__ 220*e27a552bSJed Brown #define __FUNCT__ "TSRosWFinalizePackage" 221*e27a552bSJed Brown /*@C 222*e27a552bSJed Brown TSRosWFinalizePackage - This function destroys everything in the TSRosW package. It is 223*e27a552bSJed Brown called from PetscFinalize(). 224*e27a552bSJed Brown 225*e27a552bSJed Brown Level: developer 226*e27a552bSJed Brown 227*e27a552bSJed Brown .keywords: Petsc, destroy, package 228*e27a552bSJed Brown .seealso: PetscFinalize() 229*e27a552bSJed Brown @*/ 230*e27a552bSJed Brown PetscErrorCode TSRosWFinalizePackage(void) 231*e27a552bSJed Brown { 232*e27a552bSJed Brown PetscErrorCode ierr; 233*e27a552bSJed Brown 234*e27a552bSJed Brown PetscFunctionBegin; 235*e27a552bSJed Brown TSRosWPackageInitialized = PETSC_FALSE; 236*e27a552bSJed Brown ierr = TSRosWRegisterDestroy();CHKERRQ(ierr); 237*e27a552bSJed Brown PetscFunctionReturn(0); 238*e27a552bSJed Brown } 239*e27a552bSJed Brown 240*e27a552bSJed Brown #undef __FUNCT__ 241*e27a552bSJed Brown #define __FUNCT__ "TSRosWRegister" 242*e27a552bSJed Brown /*@C 243*e27a552bSJed Brown TSRosWRegister - register an ARK IMEX scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation 244*e27a552bSJed Brown 245*e27a552bSJed Brown Not Collective, but the same schemes should be registered on all processes on which they will be used 246*e27a552bSJed Brown 247*e27a552bSJed Brown Input Parameters: 248*e27a552bSJed Brown + name - identifier for method 249*e27a552bSJed Brown . order - approximation order of method 250*e27a552bSJed Brown . s - number of stages, this is the dimension of the matrices below 251*e27a552bSJed Brown . At - Butcher table of stage coefficients for stiff part (dimension s*s, row-major) 252*e27a552bSJed Brown . bt - Butcher table for completing the stiff part of the step (dimension s; PETSC_NULL to use the last row of At) 253*e27a552bSJed Brown . ct - Abscissa of each stiff stage (dimension s, PETSC_NULL to use row sums of At) 254*e27a552bSJed Brown . A - Non-stiff stage coefficients (dimension s*s, row-major) 255*e27a552bSJed Brown . b - Non-stiff step completion table (dimension s; PETSC_NULL to use last row of At) 256*e27a552bSJed Brown . c - Non-stiff abscissa (dimension s; PETSC_NULL to use row sums of A) 257*e27a552bSJed Brown . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt and binterp 258*e27a552bSJed Brown . binterpt - Coefficients of the interpolation formula for the stiff part (dimension s*pinterp) 259*e27a552bSJed Brown - binterp - Coefficients of the interpolation formula for the non-stiff part (dimension s*pinterp; PETSC_NULL to reuse binterpt) 260*e27a552bSJed Brown 261*e27a552bSJed Brown Notes: 262*e27a552bSJed Brown Several ARK IMEX methods are provided, this function is only needed to create new methods. 263*e27a552bSJed Brown 264*e27a552bSJed Brown Level: advanced 265*e27a552bSJed Brown 266*e27a552bSJed Brown .keywords: TS, register 267*e27a552bSJed Brown 268*e27a552bSJed Brown .seealso: TSRosW 269*e27a552bSJed Brown @*/ 270*e27a552bSJed Brown PetscErrorCode TSRosWRegister(const TSRosWType name,PetscInt order,PetscInt s, 271*e27a552bSJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 272*e27a552bSJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[], 273*e27a552bSJed Brown PetscInt pinterp,const PetscReal binterpt[],const PetscReal binterp[]) 274*e27a552bSJed Brown { 275*e27a552bSJed Brown PetscErrorCode ierr; 276*e27a552bSJed Brown ARKTableauLink link; 277*e27a552bSJed Brown ARKTableau t; 278*e27a552bSJed Brown PetscInt i,j; 279*e27a552bSJed Brown 280*e27a552bSJed Brown PetscFunctionBegin; 281*e27a552bSJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 282*e27a552bSJed Brown ierr = PetscMemzero(link,sizeof(*link));CHKERRQ(ierr); 283*e27a552bSJed Brown t = &link->tab; 284*e27a552bSJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 285*e27a552bSJed Brown t->order = order; 286*e27a552bSJed Brown t->s = s; 287*e27a552bSJed 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); 288*e27a552bSJed Brown ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 289*e27a552bSJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 290*e27a552bSJed Brown if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);} 291*e27a552bSJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 292*e27a552bSJed Brown if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);} 293*e27a552bSJed Brown else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i]; 294*e27a552bSJed Brown if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);} 295*e27a552bSJed 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]; 296*e27a552bSJed Brown if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);} 297*e27a552bSJed 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]; 298*e27a552bSJed Brown t->pinterp = pinterp; 299*e27a552bSJed Brown ierr = PetscMalloc2(s*pinterp,PetscReal,&t->binterpt,s*pinterp,PetscReal,&t->binterp);CHKERRQ(ierr); 300*e27a552bSJed Brown ierr = PetscMemcpy(t->binterpt,binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 301*e27a552bSJed Brown ierr = PetscMemcpy(t->binterp,binterp?binterp:binterpt,s*pinterp*sizeof(binterpt[0]));CHKERRQ(ierr); 302*e27a552bSJed Brown link->next = ARKTableauList; 303*e27a552bSJed Brown ARKTableauList = link; 304*e27a552bSJed Brown PetscFunctionReturn(0); 305*e27a552bSJed Brown } 306*e27a552bSJed Brown 307*e27a552bSJed Brown #undef __FUNCT__ 308*e27a552bSJed Brown #define __FUNCT__ "TSStep_RosW" 309*e27a552bSJed Brown static PetscErrorCode TSStep_RosW(TS ts) 310*e27a552bSJed Brown { 311*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 312*e27a552bSJed Brown ARKTableau tab = ark->tableau; 313*e27a552bSJed Brown const PetscInt s = tab->s; 314*e27a552bSJed Brown const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c; 315*e27a552bSJed Brown PetscScalar *w = ark->work; 316*e27a552bSJed Brown Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z; 317*e27a552bSJed Brown SNES snes; 318*e27a552bSJed Brown PetscInt i,j,its,lits; 319*e27a552bSJed Brown PetscReal h,t; 320*e27a552bSJed Brown PetscErrorCode ierr; 321*e27a552bSJed Brown 322*e27a552bSJed Brown PetscFunctionBegin; 323*e27a552bSJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 324*e27a552bSJed Brown h = ts->time_step = ts->next_time_step; 325*e27a552bSJed Brown t = ts->ptime; 326*e27a552bSJed Brown 327*e27a552bSJed Brown for (i=0; i<s; i++) { 328*e27a552bSJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 329*e27a552bSJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 330*e27a552bSJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 331*e27a552bSJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 332*e27a552bSJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 333*e27a552bSJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 334*e27a552bSJed Brown } else { 335*e27a552bSJed Brown ark->stage_time = t + h*ct[i]; 336*e27a552bSJed Brown ark->shift = 1./(h*At[i*s+i]); 337*e27a552bSJed Brown /* Affine part */ 338*e27a552bSJed Brown ierr = VecZeroEntries(W);CHKERRQ(ierr); 339*e27a552bSJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 340*e27a552bSJed Brown ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 341*e27a552bSJed Brown /* Ydot = shift*(Y-Z) */ 342*e27a552bSJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 343*e27a552bSJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 344*e27a552bSJed Brown ierr = VecMAXPY(Z,i,w,YdotI);CHKERRQ(ierr); 345*e27a552bSJed Brown /* Initial guess taken from last stage */ 346*e27a552bSJed Brown ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr); 347*e27a552bSJed Brown ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 348*e27a552bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 349*e27a552bSJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 350*e27a552bSJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 351*e27a552bSJed Brown } 352*e27a552bSJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 353*e27a552bSJed Brown ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],ark->imex);CHKERRQ(ierr); 354*e27a552bSJed Brown if (ark->imex) { 355*e27a552bSJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 356*e27a552bSJed Brown } else { 357*e27a552bSJed Brown ierr = VecZeroEntries(YdotRHS[i]);CHKERRQ(ierr); 358*e27a552bSJed Brown } 359*e27a552bSJed Brown } 360*e27a552bSJed Brown for (j=0; j<s; j++) w[j] = -h*bt[j]; 361*e27a552bSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 362*e27a552bSJed Brown for (j=0; j<s; j++) w[j] = h*b[j]; 363*e27a552bSJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 364*e27a552bSJed Brown 365*e27a552bSJed Brown ts->ptime += ts->time_step; 366*e27a552bSJed Brown ts->next_time_step = ts->time_step; 367*e27a552bSJed Brown ts->steps++; 368*e27a552bSJed Brown PetscFunctionReturn(0); 369*e27a552bSJed Brown } 370*e27a552bSJed Brown 371*e27a552bSJed Brown #undef __FUNCT__ 372*e27a552bSJed Brown #define __FUNCT__ "TSInterpolate_RosW" 373*e27a552bSJed Brown static PetscErrorCode TSInterpolate_RosW(TS ts,PetscReal itime,Vec X) 374*e27a552bSJed Brown { 375*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 376*e27a552bSJed Brown PetscInt s = ark->tableau->s,pinterp = ark->tableau->pinterp,i,j; 377*e27a552bSJed Brown PetscReal tt,t = (itime - ts->ptime)/ts->time_step + 1; /* In the interval [0,1] */ 378*e27a552bSJed Brown PetscScalar *bt,*b; 379*e27a552bSJed Brown const PetscReal *Bt = ark->tableau->binterpt,*B = ark->tableau->binterp; 380*e27a552bSJed Brown PetscErrorCode ierr; 381*e27a552bSJed Brown 382*e27a552bSJed Brown PetscFunctionBegin; 383*e27a552bSJed Brown if (!Bt || !B) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSRosW %s does not have an interpolation formula",ark->tableau->name); 384*e27a552bSJed Brown ierr = PetscMalloc2(s,PetscScalar,&bt,s,PetscScalar,&b);CHKERRQ(ierr); 385*e27a552bSJed Brown for (i=0; i<s; i++) bt[i] = b[i] = 0; 386*e27a552bSJed Brown for (j=0,tt=t; j<pinterp; j++,tt*=t) { 387*e27a552bSJed Brown for (i=0; i<s; i++) { 388*e27a552bSJed Brown bt[i] += ts->time_step * Bt[i*pinterp+j] * tt * -1.0; 389*e27a552bSJed Brown b[i] += ts->time_step * B[i*pinterp+j] * tt; 390*e27a552bSJed Brown } 391*e27a552bSJed Brown } 392*e27a552bSJed 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"); 393*e27a552bSJed Brown ierr = VecCopy(ark->Y[0],X);CHKERRQ(ierr); 394*e27a552bSJed Brown ierr = VecMAXPY(X,s,bt,ark->YdotI);CHKERRQ(ierr); 395*e27a552bSJed Brown ierr = VecMAXPY(X,s,b,ark->YdotRHS);CHKERRQ(ierr); 396*e27a552bSJed Brown ierr = PetscFree2(bt,b);CHKERRQ(ierr); 397*e27a552bSJed Brown PetscFunctionReturn(0); 398*e27a552bSJed Brown } 399*e27a552bSJed Brown 400*e27a552bSJed Brown /*------------------------------------------------------------*/ 401*e27a552bSJed Brown #undef __FUNCT__ 402*e27a552bSJed Brown #define __FUNCT__ "TSReset_RosW" 403*e27a552bSJed Brown static PetscErrorCode TSReset_RosW(TS ts) 404*e27a552bSJed Brown { 405*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 406*e27a552bSJed Brown PetscInt s; 407*e27a552bSJed Brown PetscErrorCode ierr; 408*e27a552bSJed Brown 409*e27a552bSJed Brown PetscFunctionBegin; 410*e27a552bSJed Brown if (!ark->tableau) PetscFunctionReturn(0); 411*e27a552bSJed Brown s = ark->tableau->s; 412*e27a552bSJed Brown ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 413*e27a552bSJed Brown ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 414*e27a552bSJed Brown ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 415*e27a552bSJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 416*e27a552bSJed Brown ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 417*e27a552bSJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 418*e27a552bSJed Brown ierr = PetscFree(ark->work);CHKERRQ(ierr); 419*e27a552bSJed Brown PetscFunctionReturn(0); 420*e27a552bSJed Brown } 421*e27a552bSJed Brown 422*e27a552bSJed Brown #undef __FUNCT__ 423*e27a552bSJed Brown #define __FUNCT__ "TSDestroy_RosW" 424*e27a552bSJed Brown static PetscErrorCode TSDestroy_RosW(TS ts) 425*e27a552bSJed Brown { 426*e27a552bSJed Brown PetscErrorCode ierr; 427*e27a552bSJed Brown 428*e27a552bSJed Brown PetscFunctionBegin; 429*e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 430*e27a552bSJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 431*e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","",PETSC_NULL);CHKERRQ(ierr); 432*e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","",PETSC_NULL);CHKERRQ(ierr); 433*e27a552bSJed Brown PetscFunctionReturn(0); 434*e27a552bSJed Brown } 435*e27a552bSJed Brown 436*e27a552bSJed Brown /* 437*e27a552bSJed Brown This defines the nonlinear equation that is to be solved with SNES 438*e27a552bSJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 439*e27a552bSJed Brown */ 440*e27a552bSJed Brown #undef __FUNCT__ 441*e27a552bSJed Brown #define __FUNCT__ "SNESTSFormFunction_RosW" 442*e27a552bSJed Brown static PetscErrorCode SNESTSFormFunction_RosW(SNES snes,Vec X,Vec F,TS ts) 443*e27a552bSJed Brown { 444*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 445*e27a552bSJed Brown PetscErrorCode ierr; 446*e27a552bSJed Brown 447*e27a552bSJed Brown PetscFunctionBegin; 448*e27a552bSJed Brown ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 449*e27a552bSJed Brown ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,F,ark->imex);CHKERRQ(ierr); 450*e27a552bSJed Brown PetscFunctionReturn(0); 451*e27a552bSJed Brown } 452*e27a552bSJed Brown 453*e27a552bSJed Brown #undef __FUNCT__ 454*e27a552bSJed Brown #define __FUNCT__ "SNESTSFormJacobian_RosW" 455*e27a552bSJed Brown static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 456*e27a552bSJed Brown { 457*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 458*e27a552bSJed Brown PetscErrorCode ierr; 459*e27a552bSJed Brown 460*e27a552bSJed Brown PetscFunctionBegin; 461*e27a552bSJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */ 462*e27a552bSJed Brown ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 463*e27a552bSJed Brown PetscFunctionReturn(0); 464*e27a552bSJed Brown } 465*e27a552bSJed Brown 466*e27a552bSJed Brown #undef __FUNCT__ 467*e27a552bSJed Brown #define __FUNCT__ "TSSetUp_RosW" 468*e27a552bSJed Brown static PetscErrorCode TSSetUp_RosW(TS ts) 469*e27a552bSJed Brown { 470*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 471*e27a552bSJed Brown ARKTableau tab = ark->tableau; 472*e27a552bSJed Brown PetscInt s = tab->s; 473*e27a552bSJed Brown PetscErrorCode ierr; 474*e27a552bSJed Brown 475*e27a552bSJed Brown PetscFunctionBegin; 476*e27a552bSJed Brown if (!ark->tableau) { 477*e27a552bSJed Brown ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr); 478*e27a552bSJed Brown } 479*e27a552bSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 480*e27a552bSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 481*e27a552bSJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 482*e27a552bSJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 483*e27a552bSJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 484*e27a552bSJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 485*e27a552bSJed Brown ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr); 486*e27a552bSJed Brown PetscFunctionReturn(0); 487*e27a552bSJed Brown } 488*e27a552bSJed Brown /*------------------------------------------------------------*/ 489*e27a552bSJed Brown 490*e27a552bSJed Brown #undef __FUNCT__ 491*e27a552bSJed Brown #define __FUNCT__ "TSSetFromOptions_RosW" 492*e27a552bSJed Brown static PetscErrorCode TSSetFromOptions_RosW(TS ts) 493*e27a552bSJed Brown { 494*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 495*e27a552bSJed Brown PetscErrorCode ierr; 496*e27a552bSJed Brown char arktype[256]; 497*e27a552bSJed Brown 498*e27a552bSJed Brown PetscFunctionBegin; 499*e27a552bSJed Brown ierr = PetscOptionsHead("RosW ODE solver options");CHKERRQ(ierr); 500*e27a552bSJed Brown { 501*e27a552bSJed Brown ARKTableauLink link; 502*e27a552bSJed Brown PetscInt count,choice; 503*e27a552bSJed Brown PetscBool flg; 504*e27a552bSJed Brown const char **namelist; 505*e27a552bSJed Brown ierr = PetscStrncpy(arktype,TSRosWDefault,sizeof arktype);CHKERRQ(ierr); 506*e27a552bSJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 507*e27a552bSJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 508*e27a552bSJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 509*e27a552bSJed Brown ierr = PetscOptionsEList("-ts_rosw_type","Family of ARK IMEX method","TSRosWSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 510*e27a552bSJed Brown ierr = TSRosWSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 511*e27a552bSJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 512*e27a552bSJed Brown flg = (PetscBool)!ark->imex; 513*e27a552bSJed Brown ierr = PetscOptionsBool("-ts_rosw_fully_implicit","Solve the problem fully implicitly","TSRosWSetFullyImplicit",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 514*e27a552bSJed Brown ark->imex = (PetscBool)!flg; 515*e27a552bSJed Brown ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); 516*e27a552bSJed Brown } 517*e27a552bSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 518*e27a552bSJed Brown PetscFunctionReturn(0); 519*e27a552bSJed Brown } 520*e27a552bSJed Brown 521*e27a552bSJed Brown #undef __FUNCT__ 522*e27a552bSJed Brown #define __FUNCT__ "PetscFormatRealArray" 523*e27a552bSJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 524*e27a552bSJed Brown { 525*e27a552bSJed Brown PetscErrorCode ierr; 526*e27a552bSJed Brown int i,left,count; 527*e27a552bSJed Brown char *p; 528*e27a552bSJed Brown 529*e27a552bSJed Brown PetscFunctionBegin; 530*e27a552bSJed Brown for (i=0,p=buf,left=(int)len; i<n; i++) { 531*e27a552bSJed Brown ierr = PetscSNPrintf(p,left,fmt,&count,x[i]);CHKERRQ(ierr); 532*e27a552bSJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 533*e27a552bSJed Brown left -= count; 534*e27a552bSJed Brown p += count; 535*e27a552bSJed Brown *p++ = ' '; 536*e27a552bSJed Brown } 537*e27a552bSJed Brown p[i ? 0 : -1] = 0; 538*e27a552bSJed Brown PetscFunctionReturn(0); 539*e27a552bSJed Brown } 540*e27a552bSJed Brown 541*e27a552bSJed Brown #undef __FUNCT__ 542*e27a552bSJed Brown #define __FUNCT__ "TSView_RosW" 543*e27a552bSJed Brown static PetscErrorCode TSView_RosW(TS ts,PetscViewer viewer) 544*e27a552bSJed Brown { 545*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 546*e27a552bSJed Brown ARKTableau tab = ark->tableau; 547*e27a552bSJed Brown PetscBool iascii; 548*e27a552bSJed Brown PetscErrorCode ierr; 549*e27a552bSJed Brown 550*e27a552bSJed Brown PetscFunctionBegin; 551*e27a552bSJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 552*e27a552bSJed Brown if (iascii) { 553*e27a552bSJed Brown const TSRosWType arktype; 554*e27a552bSJed Brown char buf[512]; 555*e27a552bSJed Brown ierr = TSRosWGetType(ts,&arktype);CHKERRQ(ierr); 556*e27a552bSJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 557*e27a552bSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->ct);CHKERRQ(ierr); 558*e27a552bSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa ct = %s\n",buf);CHKERRQ(ierr); 559*e27a552bSJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"% 8.6f",tab->s,tab->c);CHKERRQ(ierr); 560*e27a552bSJed Brown ierr = PetscViewerASCIIPrintf(viewer," Nonstiff abscissa c = %s\n",buf);CHKERRQ(ierr); 561*e27a552bSJed Brown } 562*e27a552bSJed Brown ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 563*e27a552bSJed Brown PetscFunctionReturn(0); 564*e27a552bSJed Brown } 565*e27a552bSJed Brown 566*e27a552bSJed Brown #undef __FUNCT__ 567*e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType" 568*e27a552bSJed Brown /*@C 569*e27a552bSJed Brown TSRosWSetType - Set the type of ARK IMEX scheme 570*e27a552bSJed Brown 571*e27a552bSJed Brown Logically collective 572*e27a552bSJed Brown 573*e27a552bSJed Brown Input Parameter: 574*e27a552bSJed Brown + ts - timestepping context 575*e27a552bSJed Brown - arktype - type of ARK-IMEX scheme 576*e27a552bSJed Brown 577*e27a552bSJed Brown Level: intermediate 578*e27a552bSJed Brown 579*e27a552bSJed Brown .seealso: TSRosWGetType() 580*e27a552bSJed Brown @*/ 581*e27a552bSJed Brown PetscErrorCode TSRosWSetType(TS ts,const TSRosWType arktype) 582*e27a552bSJed Brown { 583*e27a552bSJed Brown PetscErrorCode ierr; 584*e27a552bSJed Brown 585*e27a552bSJed Brown PetscFunctionBegin; 586*e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 587*e27a552bSJed Brown ierr = PetscTryMethod(ts,"TSRosWSetType_C",(TS,const TSRosWType),(ts,arktype));CHKERRQ(ierr); 588*e27a552bSJed Brown PetscFunctionReturn(0); 589*e27a552bSJed Brown } 590*e27a552bSJed Brown 591*e27a552bSJed Brown #undef __FUNCT__ 592*e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType" 593*e27a552bSJed Brown /*@C 594*e27a552bSJed Brown TSRosWGetType - Get the type of ARK IMEX scheme 595*e27a552bSJed Brown 596*e27a552bSJed Brown Logically collective 597*e27a552bSJed Brown 598*e27a552bSJed Brown Input Parameter: 599*e27a552bSJed Brown . ts - timestepping context 600*e27a552bSJed Brown 601*e27a552bSJed Brown Output Parameter: 602*e27a552bSJed Brown . arktype - type of ARK-IMEX scheme 603*e27a552bSJed Brown 604*e27a552bSJed Brown Level: intermediate 605*e27a552bSJed Brown 606*e27a552bSJed Brown .seealso: TSRosWGetType() 607*e27a552bSJed Brown @*/ 608*e27a552bSJed Brown PetscErrorCode TSRosWGetType(TS ts,const TSRosWType *arktype) 609*e27a552bSJed Brown { 610*e27a552bSJed Brown PetscErrorCode ierr; 611*e27a552bSJed Brown 612*e27a552bSJed Brown PetscFunctionBegin; 613*e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 614*e27a552bSJed Brown ierr = PetscUseMethod(ts,"TSRosWGetType_C",(TS,const TSRosWType*),(ts,arktype));CHKERRQ(ierr); 615*e27a552bSJed Brown PetscFunctionReturn(0); 616*e27a552bSJed Brown } 617*e27a552bSJed Brown 618*e27a552bSJed Brown #undef __FUNCT__ 619*e27a552bSJed Brown #define __FUNCT__ "TSRosWSetFullyImplicit" 620*e27a552bSJed Brown /*@C 621*e27a552bSJed Brown TSRosWSetFullyImplicit - Solve both parts of the equation implicitly 622*e27a552bSJed Brown 623*e27a552bSJed Brown Logically collective 624*e27a552bSJed Brown 625*e27a552bSJed Brown Input Parameter: 626*e27a552bSJed Brown + ts - timestepping context 627*e27a552bSJed Brown - flg - PETSC_TRUE for fully implicit 628*e27a552bSJed Brown 629*e27a552bSJed Brown Level: intermediate 630*e27a552bSJed Brown 631*e27a552bSJed Brown .seealso: TSRosWGetType() 632*e27a552bSJed Brown @*/ 633*e27a552bSJed Brown PetscErrorCode TSRosWSetFullyImplicit(TS ts,PetscBool flg) 634*e27a552bSJed Brown { 635*e27a552bSJed Brown PetscErrorCode ierr; 636*e27a552bSJed Brown 637*e27a552bSJed Brown PetscFunctionBegin; 638*e27a552bSJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 639*e27a552bSJed Brown ierr = PetscTryMethod(ts,"TSRosWSetFullyImplicit_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 640*e27a552bSJed Brown PetscFunctionReturn(0); 641*e27a552bSJed Brown } 642*e27a552bSJed Brown 643*e27a552bSJed Brown EXTERN_C_BEGIN 644*e27a552bSJed Brown #undef __FUNCT__ 645*e27a552bSJed Brown #define __FUNCT__ "TSRosWGetType_RosW" 646*e27a552bSJed Brown PetscErrorCode TSRosWGetType_RosW(TS ts,const TSRosWType *arktype) 647*e27a552bSJed Brown { 648*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 649*e27a552bSJed Brown PetscErrorCode ierr; 650*e27a552bSJed Brown 651*e27a552bSJed Brown PetscFunctionBegin; 652*e27a552bSJed Brown if (!ark->tableau) {ierr = TSRosWSetType(ts,TSRosWDefault);CHKERRQ(ierr);} 653*e27a552bSJed Brown *arktype = ark->tableau->name; 654*e27a552bSJed Brown PetscFunctionReturn(0); 655*e27a552bSJed Brown } 656*e27a552bSJed Brown #undef __FUNCT__ 657*e27a552bSJed Brown #define __FUNCT__ "TSRosWSetType_RosW" 658*e27a552bSJed Brown PetscErrorCode TSRosWSetType_RosW(TS ts,const TSRosWType arktype) 659*e27a552bSJed Brown { 660*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 661*e27a552bSJed Brown PetscErrorCode ierr; 662*e27a552bSJed Brown PetscBool match; 663*e27a552bSJed Brown ARKTableauLink link; 664*e27a552bSJed Brown 665*e27a552bSJed Brown PetscFunctionBegin; 666*e27a552bSJed Brown if (ark->tableau) { 667*e27a552bSJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 668*e27a552bSJed Brown if (match) PetscFunctionReturn(0); 669*e27a552bSJed Brown } 670*e27a552bSJed Brown for (link = ARKTableauList; link; link=link->next) { 671*e27a552bSJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 672*e27a552bSJed Brown if (match) { 673*e27a552bSJed Brown ierr = TSReset_RosW(ts);CHKERRQ(ierr); 674*e27a552bSJed Brown ark->tableau = &link->tab; 675*e27a552bSJed Brown PetscFunctionReturn(0); 676*e27a552bSJed Brown } 677*e27a552bSJed Brown } 678*e27a552bSJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 679*e27a552bSJed Brown PetscFunctionReturn(0); 680*e27a552bSJed Brown } 681*e27a552bSJed Brown #undef __FUNCT__ 682*e27a552bSJed Brown #define __FUNCT__ "TSRosWSetFullyImplicit_RosW" 683*e27a552bSJed Brown PetscErrorCode TSRosWSetFullyImplicit_RosW(TS ts,PetscBool flg) 684*e27a552bSJed Brown { 685*e27a552bSJed Brown TS_RosW *ark = (TS_RosW*)ts->data; 686*e27a552bSJed Brown 687*e27a552bSJed Brown PetscFunctionBegin; 688*e27a552bSJed Brown ark->imex = (PetscBool)!flg; 689*e27a552bSJed Brown PetscFunctionReturn(0); 690*e27a552bSJed Brown } 691*e27a552bSJed Brown EXTERN_C_END 692*e27a552bSJed Brown 693*e27a552bSJed Brown /* ------------------------------------------------------------ */ 694*e27a552bSJed Brown /*MC 695*e27a552bSJed Brown TSRosW - ODE solver using Rosenbrock-W schemes 696*e27a552bSJed Brown 697*e27a552bSJed Brown These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly 698*e27a552bSJed Brown nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part 699*e27a552bSJed Brown of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 700*e27a552bSJed Brown 701*e27a552bSJed Brown Notes: 702*e27a552bSJed Brown This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 703*e27a552bSJed Brown 704*e27a552bSJed Brown Level: beginner 705*e27a552bSJed Brown 706*e27a552bSJed Brown .seealso: TSCreate(), TS, TSSetType(), TSRosWRegister() 707*e27a552bSJed Brown 708*e27a552bSJed Brown M*/ 709*e27a552bSJed Brown EXTERN_C_BEGIN 710*e27a552bSJed Brown #undef __FUNCT__ 711*e27a552bSJed Brown #define __FUNCT__ "TSCreate_RosW" 712*e27a552bSJed Brown PetscErrorCode TSCreate_RosW(TS ts) 713*e27a552bSJed Brown { 714*e27a552bSJed Brown TS_RosW *th; 715*e27a552bSJed Brown PetscErrorCode ierr; 716*e27a552bSJed Brown 717*e27a552bSJed Brown PetscFunctionBegin; 718*e27a552bSJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 719*e27a552bSJed Brown ierr = TSRosWInitializePackage(PETSC_NULL);CHKERRQ(ierr); 720*e27a552bSJed Brown #endif 721*e27a552bSJed Brown 722*e27a552bSJed Brown ts->ops->reset = TSReset_RosW; 723*e27a552bSJed Brown ts->ops->destroy = TSDestroy_RosW; 724*e27a552bSJed Brown ts->ops->view = TSView_RosW; 725*e27a552bSJed Brown ts->ops->setup = TSSetUp_RosW; 726*e27a552bSJed Brown ts->ops->step = TSStep_RosW; 727*e27a552bSJed Brown ts->ops->interpolate = TSInterpolate_RosW; 728*e27a552bSJed Brown ts->ops->setfromoptions = TSSetFromOptions_RosW; 729*e27a552bSJed Brown ts->ops->snesfunction = SNESTSFormFunction_RosW; 730*e27a552bSJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_RosW; 731*e27a552bSJed Brown 732*e27a552bSJed Brown ierr = PetscNewLog(ts,TS_RosW,&th);CHKERRQ(ierr); 733*e27a552bSJed Brown ts->data = (void*)th; 734*e27a552bSJed Brown th->imex = PETSC_TRUE; 735*e27a552bSJed Brown 736*e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWGetType_C","TSRosWGetType_RosW",TSRosWGetType_RosW);CHKERRQ(ierr); 737*e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetType_C","TSRosWSetType_RosW",TSRosWSetType_RosW);CHKERRQ(ierr); 738*e27a552bSJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSRosWSetFullyImplicit_C","TSRosWSetFullyImplicit_RosW",TSRosWSetFullyImplicit_RosW);CHKERRQ(ierr); 739*e27a552bSJed Brown PetscFunctionReturn(0); 740*e27a552bSJed Brown } 741*e27a552bSJed Brown EXTERN_C_END 742