• Main Page
  • Classes
  • Files
  • File List
  • File Members

athena.h

Go to the documentation of this file.
00001 #ifndef ATHENA_H
00002 #define ATHENA_H 
00003 /*============================================================================*/
00004 /*! \file athena.h
00005  *  \brief Contains definitions of many data types and structures.
00006  *
00007  * PURPOSE: Contains definitions of the following data types and structures:
00008  * - Real    - either float or double, depending on configure option
00009  * - ConsS   - cell-centered conserved variables
00010  * - PrimS   - cell-centered primitive variables
00011  * - Cons1DS - conserved variables in 1D: same as ConsS minus Bx
00012  * - Prim1DS - primitive variables in 1D: same as PrimS minus Bx
00013  * - GrainS  - basic properties of particles
00014  * - GridS   - everything in a single Grid
00015  * - DomainS - everything in a single Domain (potentially many Grids)
00016  * - MeshS   - everything across whole Mesh (potentially many Domains)
00017  * - OutputS - everything associated with an individual output                */
00018 /*============================================================================*/
00019 #include "defs.h"
00020 
00021 #ifdef MPI_PARALLEL
00022 #include "mpi.h"
00023 #endif
00024 
00025 /*! \typedef Real
00026  *  \brief Variable precision float, depends on macro set by configure.
00027  */ 
00028 #if defined(SINGLE_PREC)
00029 #ifdef MPI_PARALLEL
00030 #error: MPI requires double precision
00031 #endif /*MPI_PARALLEL */
00032 typedef float  Real;
00033 #elif defined(DOUBLE_PREC)
00034 typedef double Real;
00035 #else
00036 # error "Not a valid precision flag"
00037 #endif
00038 
00039 /*! \struct Real3Vect
00040  *  \brief General 3-vectors of Reals.
00041  */
00042 typedef struct Real3Vect_s{
00043   Real x, y, z;
00044 }Real3Vect;
00045 /*! \struct Int3Vect
00046  *  \brief General 3-vectors of ints.
00047  */
00048 typedef struct Int3Vect_s{
00049   int i, j, k;
00050 }Int3Vect;
00051 
00052 /*! \struct SideS
00053  *  \brief Sides of a cube, used to find overlaps between Grids 
00054  *   at different levels.
00055  */
00056 typedef struct Side_s{
00057   int ijkl[3];    /*!< indices of left-sides  in each dir [0,1,2]=[i,j,k] */ 
00058   int ijkr[3];    /*!< indices of right-sides in each dir [0,1,2]=[i,j,k] */ 
00059 }SideS;
00060 
00061 /*! \struct GridsDataS
00062  *  \brief Number of zones in, and identifying information about, a Grid.
00063  */
00064 typedef struct GridsData_s{
00065   int Nx[3];       /*!< number of zones in each dir [0,1,2]=[x1,x2,x3] */
00066   int Disp[3];   /*!< i,j,k displacements from origin of root [0,1,2]=[i,j,k] */
00067   int ID_Comm_world;      /*!< ID of process for this Grid in MPI_COMM_WORLD */
00068   int ID_Comm_Domain;     /*!< ID of process for this Grid in Comm_Domain    */
00069 #ifdef STATIC_MESH_REFINEMENT
00070   int ID_Comm_Children;     /*!< ID updating this Grid in Comm_Domain    */
00071   int ID_Comm_Parent;     /*!< ID updating this Grid in Comm_Domain    */
00072 #endif
00073 }GridsDataS;
00074 
00075 /*----------------------------------------------------------------------------*/
00076 /*! \struct ConsS
00077  *  \brief Conserved variables.
00078  *  IMPORTANT!! The order of the elements in ConsS CANNOT be changed.
00079  */
00080 typedef struct Cons_s{
00081   Real d;                       /*!< density */
00082   Real M1;                      /*!< momentum density in 1-direction*/
00083   Real M2;                      /*!< momentum density in 2-direction*/
00084   Real M3;                      /*!< momentum density in 3-direction*/
00085 #ifndef BAROTROPIC
00086   Real E;                       /*!< total energy density */
00087 #endif /* BAROTROPIC */
00088 #ifdef MHD
00089   Real B1c;                     /*!< cell centered magnetic fields in 1-dir*/
00090   Real B2c;                     /*!< cell centered magnetic fields in 2-dir*/
00091   Real B3c;                     /*!< cell centered magnetic fields in 3-dir*/
00092 #endif /* MHD */
00093 #if (NSCALARS > 0)
00094   Real s[NSCALARS];             /*!< passively advected scalars */
00095 #endif
00096 #ifdef CYLINDRICAL
00097   Real Pflux;                   /*!< pressure component of flux */
00098 #endif
00099 }ConsS;
00100 
00101 /*----------------------------------------------------------------------------*/
00102 /*! \struct PrimS
00103  *  \brief Primitive variables.
00104  *  IMPORTANT!! The order of the elements in PrimS CANNOT be changed.
00105  */
00106 typedef struct Prim_s{
00107   Real d;                       /*!< density */
00108   Real V1;                      /*!< velocity in 1-direction */
00109   Real V2;                      /*!< velocity in 2-direction */
00110   Real V3;                      /*!< velocity in 3-direction */
00111 #ifndef BAROTROPIC
00112   Real P;                       /*!< pressure */
00113 #endif /* BAROTROPIC */
00114 #ifdef MHD
00115   Real B1c;                     /*!< cell centered magnetic fields in 1-dir */
00116   Real B2c;                     /*!< cell centered magnetic fields in 2-dir */
00117   Real B3c;                     /*!< cell centered magnetic fields in 3-dir */
00118 #endif /* MHD */
00119 #if (NSCALARS > 0)
00120   Real r[NSCALARS];             /*!< density-normalized advected scalars */
00121 #endif
00122 }PrimS;
00123 
00124 /*----------------------------------------------------------------------------*/
00125 /*! \struct Cons1DS
00126  *  \brief Conserved variables in 1D (does not contain Bx).
00127  *  IMPORTANT!! The order of the elements in Cons1DS CANNOT be changed.
00128  */
00129 typedef struct Cons1D_s{
00130   Real d;                       /*!< density */
00131   Real Mx;                      /*!< momentum density in X,Y,Z; where X is    */
00132   Real My;                      /*!< direction longitudinal to 1D slice; which*/
00133   Real Mz;                      /*!< can be in any dimension: 1,2,or 3        */
00134 #ifndef BAROTROPIC
00135   Real E;                       /*!< total energy density */
00136 #endif /* BAROTROPIC */
00137 #ifdef MHD
00138   Real By;                      /*!< cell centered magnetic fields in Y */
00139   Real Bz;                      /*!< cell centered magnetic fields in Z */
00140 #endif /* MHD */
00141 #if (NSCALARS > 0)
00142   Real s[NSCALARS];             /*!< passively advected scalars */
00143 #endif
00144 #ifdef CYLINDRICAL
00145   Real Pflux;                   /*!< pressure component of flux */
00146 #endif
00147 }Cons1DS;
00148 
00149 /*----------------------------------------------------------------------------*/
00150 /*! \struct Prim1DS
00151  *  \brief Primitive variables in 1D (does not contain Bx).
00152  *  IMPORTANT!! The order of the elements in Prim1DS CANNOT be changed.
00153  */
00154 typedef struct Prim1D_s{
00155   Real d;                       /*!< density */
00156   Real Vx;                      /*!< velocity in X-direction */
00157   Real Vy;                      /*!< velocity in Y-direction */
00158   Real Vz;                      /*!< velocity in Z-direction */
00159 #ifndef BAROTROPIC
00160   Real P;                       /*!< pressure */
00161 #endif /* BAROTROPIC */
00162 #ifdef MHD
00163   Real By;                      /*!< cell centered magnetic fields in Y-dir */
00164   Real Bz;                      /*!< cell centered magnetic fields in Z-dir */
00165 #endif /* MHD */
00166 #if (NSCALARS > 0)
00167   Real r[NSCALARS];             /*!< density-normalized advected scalars */
00168 #endif
00169 }Prim1DS;
00170 
00171 /*----------------------------------------------------------------------------*/
00172 /*! \struct GrainS
00173  *  \brief Basic quantities for one pseudo-particle.
00174  */
00175 #ifdef PARTICLES
00176 
00177 /* Physical quantities of a dust particle */
00178 typedef struct Grain_s{
00179   Real x1,x2,x3;        /*!< coordinate in X,Y,Z */
00180   Real v1,v2,v3;        /*!< velocity in X,Y,Z */
00181   int property;         /*!< index of grain properties */
00182   short pos;            /*!< position: 0: ghost; 1: grid; >=10: cross out/in; */
00183   long my_id;           /*!< particle id */
00184 #ifdef MPI_PARALLEL
00185   int init_id;          /*!< particle's initial host processor id */
00186 #endif
00187 #ifdef FARGO
00188   Real shift;           /*!< amount of shift in x2 direction */
00189 #endif
00190 }GrainS;
00191 
00192 /*! \struct GrainAux
00193  *  \brief Auxilary quantities for a dust particle. */
00194 typedef struct GrainAux_s{
00195   Real dpar;            /*!< local particle density */
00196 #ifdef FARGO
00197   Real shift;           /*!< amount of shift in x2 direction */
00198 #endif
00199 }GrainAux;
00200 
00201 /*! \struct Grain_Property
00202  *  \brief List of physical grain properties. */
00203 typedef struct Grain_Property_s{
00204 #ifdef FEEDBACK
00205   Real m;               /*!< mass of this type of particle */
00206 #endif
00207   Real rad;             /*!< radius of this type of particle (cm) */
00208   Real rho;             /*!< solid density of this type of particle (g/cm^3) */
00209   long num;             /*!< number of particles with this property */
00210   short integrator;     /*!< integrator type: exp (1), semi (2) or full (3) */
00211 }Grain_Property;
00212 
00213 /*! \struct GPCouple
00214  *  \brief Grid elements for gas-particle coupling. */
00215 typedef struct GPCouple_s{
00216   Real grid_d;          /*!< gas density (at 1/2 step) */
00217   Real grid_v1;         /*!< gas 1-velocity (at 1/2 step) */
00218   Real grid_v2;         /*!< gas 2-velocity (at 1/2 step) */
00219   Real grid_v3;         /*!< gas 3-velocity (at 1/2 step) */
00220 #ifndef BAROTROPIC
00221   Real grid_cs;         /*!< gas sound speed */
00222 #endif
00223 #ifdef FEEDBACK
00224   Real fb1;             /*!< 1-momentum feedback to the grid */
00225   Real fb2;             /*!< 2-momentum feedback to the grid */
00226   Real fb3;             /*!< 3-momentum feedback to the grid */
00227   Real FBstiff;         /*!< stiffness of the feedback term */
00228   Real Eloss;           /*!< energy dissipation */
00229 #endif
00230 }GPCouple;
00231 
00232 #endif /* PARTICLES */
00233 
00234 /*----------------------------------------------------------------------------*/
00235 /*! \struct GridOvrlpS
00236  *  \brief Contains information about Grid overlaps, used for SMR.
00237  */
00238 #ifdef STATIC_MESH_REFINEMENT
00239 typedef struct GridOvrlp_s{
00240   int ijks[3];         /*!< start ijk on this Grid of overlap [0,1,2]=[i,j,k] */
00241   int ijke[3];         /*!< end   ijk on this Grid of overlap [0,1,2]=[i,j,k] */
00242   int ID, DomN;        /*!< processor ID, and Domain #, of OVERLAP Grid */
00243   int nWordsRC, nWordsP; /*!< # of words communicated for Rest/Corr and Prol */
00244   ConsS **myFlx[6];   /*!< fluxes of conserved variables at 6 boundaries */
00245 #ifdef MHD
00246   Real **myEMF1[6];      /*!< fluxes of magnetic field (EMF1) at 6 boundaries */
00247   Real **myEMF2[6];      /*!< fluxes of magnetic field (EMF2) at 6 boundaries */
00248   Real **myEMF3[6];      /*!< fluxes of magnetic field (EMF3) at 6 boundaries */
00249 #endif
00250 }GridOvrlpS;
00251 #endif /* STATIC_MESH_REFINEMENT */
00252 
00253 /*----------------------------------------------------------------------------*/
00254 /*! \struct GridS
00255  *  \brief 3D arrays of dependent variables, plus grid data, plus particle data,
00256  *   plus data about child and parent Grids, plus MPI rank information for a
00257  *   Grid.
00258  *
00259  *   Remember a Grid is defined as the region of a Domain at some
00260  *   refinement level being updated by a single processor.  Uses an array of
00261  *   ConsS, rather than arrays of each variable, to increase locality of data
00262  *   for a given cell in memory.  */
00263 
00264 typedef struct Grid_s{
00265   ConsS ***U;                /*!< conserved variables */
00266 #ifdef MHD
00267   Real ***B1i,***B2i,***B3i;    /*!< interface magnetic fields */
00268 #ifdef RESISTIVITY
00269   Real ***eta_Ohm,***eta_Hall,***eta_AD; /*!< magnetic diffusivities */ 
00270 #endif
00271 #endif /* MHD */
00272 #ifdef SELF_GRAVITY
00273   Real ***Phi, ***Phi_old;      /*!< gravitational potential */
00274   Real ***x1MassFlux;           /*!< x1 mass flux for source term correction */
00275   Real ***x2MassFlux;           /*!< x2 mass flux for source term correction */
00276   Real ***x3MassFlux;           /*!< x3 mass flux for source term correction */
00277 #endif /* GRAVITY */
00278   Real MinX[3];       /*!< min(x) in each dir on this Grid [0,1,2]=[x1,x2,x3] */
00279   Real MaxX[3];       /*!< max(x) in each dir on this Grid [0,1,2]=[x1,x2,x3] */
00280   Real dx1,dx2,dx3;   /*!< cell size on this Grid */
00281   Real time, dt;           /*!< current time and timestep  */
00282   int is,ie;               /*!< start/end cell index in x1 direction */
00283   int js,je;               /*!< start/end cell index in x2 direction */
00284   int ks,ke;               /*!< start/end cell index in x3 direction */
00285   int Nx[3];     /*!< # of zones in each dir on Grid [0,1,2]=[x1,x2,x3] */
00286   int Disp[3];   /*!< i,j,k displacements of Grid from origin [0,1,2]=[i,j,k] */
00287 
00288   int rx1_id, lx1_id;  /*!< ID of Grid to R/L in x1-dir (default=-1; no Grid) */
00289   int rx2_id, lx2_id;  /*!< ID of Grid to R/L in x2-dir (default=-1; no Grid) */
00290   int rx3_id, lx3_id;  /*!< ID of Grid to R/L in x3-dir (default=-1; no Grid) */
00291 
00292 #ifdef PARTICLES
00293   int partypes;              /*!< number of particle types */
00294   Grain_Property *grproperty;/*!< array of particle properties of all types */
00295   long nparticle;            /*!< number of particles */
00296   long arrsize;              /*!< size of the particle array */
00297   Grain *particle;           /*!< array of all particles */
00298   GrainAux *parsub;          /*!< supplemental particle information */
00299   GPCouple ***Coup;          /*!< array of gas-particle coupling */
00300 #endif /* PARTICLES */
00301 
00302 #ifdef STATIC_MESH_REFINEMENT
00303   int NCGrid;        /*!< # of child  Grids that overlap this Grid */
00304   int NPGrid;        /*!< # of parent Grids that this Grid overlaps */
00305   int NmyCGrid;      /*!< # of child  Grids on same processor as this Grid */
00306   int NmyPGrid;      /*!< # of parent Grids on same processor (either 0 or 1) */
00307 
00308   GridOvrlpS *CGrid;  /*!< 1D array of data for NCGrid child  overlap regions */
00309   GridOvrlpS *PGrid;  /*!< 1D array of data for NPGrid parent overlap regions */
00310 /* NB: The first NmyCGrid[NmyPGrid] elements of these arrays contain overlap
00311  * regions being updated by the same processor as this Grid */
00312 #endif /* STATIC_MESH_REFINEMENT */
00313 
00314 #ifdef CYLINDRICAL
00315   Real *r,*ri;                  /*!< cylindrical scaling factors */ 
00316 #endif /* CYLINDRICAL */
00317 
00318 }GridS;
00319 
00320 /*! \fn void (*VGFun_t)(GridS *pG)
00321  *  \brief Generic void function of Grid. */
00322 typedef void (*VGFun_t)(GridS *pG);    /* generic void function of Grid */
00323 
00324 /*----------------------------------------------------------------------------*/
00325 /*! \struct DomainS
00326  *  \brief Information about one region of Mesh at some particular level.
00327  *
00328  * Contains pointer to a single Grid, even though the Domain may contain many
00329  * Grids, because for any general parallelization mode, no more than one Grid
00330  * can exist per Domain per processor.
00331  *
00332  * The i,j,k displacements are measured in units of grid cells on this Domain
00333  */
00334 typedef struct Domain_s{
00335   Real RootMinX[3]; /*!< min(x) in each dir on root Domain [0,1,2]=[x1,x2,x3] */
00336   Real RootMaxX[3]; /*!< max(x) in each dir on root Domain [0,1,2]=[x1,x2,x3] */
00337   Real MinX[3];     /*!< min(x) in each dir on this Domain [0,1,2]=[x1,x2,x3] */
00338   Real MaxX[3];     /*!< max(x) in each dir on this Domain [0,1,2]=[x1,x2,x3] */
00339   Real dx[3];                /*!< cell size in this Domain [0,1,2]=[x1,x2,x3] */
00340   int Nx[3];    /*!< # of zones in each dir in this Domain [0,1,2]=[x1,x2,x3] */
00341   int NGrid[3]; /*!< # of Grids in each dir in this Domain [0,1,2]=[x1,x2,x3] */
00342   int Disp[3]; /*!< i,j,k displacements of Domain from origin [0,1,2]=[i,j,k] */
00343   int Level,DomNumber;   /*!< level and ID number of this Domain */
00344   int InputBlock;      /*!< # of <domain> block in input file for this Domain */
00345   GridS *Grid;     /*!< pointer to Grid in this Dom updated on this processor */
00346 
00347   GridsDataS ***GData;/*!< size,location,& processor IDs of Grids in this Dom */
00348 
00349   VGFun_t ix1_BCFun, ox1_BCFun;/*!< ix1/ox1 BC function pointers for this Dom */
00350   VGFun_t ix2_BCFun, ox2_BCFun;/*!< ix1/ox1 BC function pointers for this Dom */
00351   VGFun_t ix3_BCFun, ox3_BCFun;/*!< ix1/ox1 BC function pointers for this Dom */
00352 
00353 #ifdef MPI_PARALLEL
00354   MPI_Comm Comm_Domain;      /*!< MPI communicator between Grids on this Dom */
00355   MPI_Group Group_Domain;    /*!< MPI group for Domain communicator */
00356 #ifdef STATIC_MESH_REFINEMENT
00357   MPI_Comm Comm_Parent;      /*!< MPI communicator to Grids in parent Domain  */
00358   MPI_Comm Comm_Children;    /*!< MPI communicator to Grids in  child Domains */
00359   MPI_Group Group_Children;  /*!< MPI group for Children communicator */
00360 #endif /* STATIC_MESH_REFINEMENT */
00361 #endif /* MPI_PARALLEL */
00362 }DomainS;
00363 
00364 /*! \fn void (*VDFun_t)(DomainS *pD)
00365  *  \brief Generic void function of Domain. */
00366 typedef void (*VDFun_t)(DomainS *pD);  /* generic void function of Domain */
00367 
00368 /*----------------------------------------------------------------------------*/
00369 /*! \struct MeshS
00370  *  \brief Information about entire mesh hierarchy, including array of Domains.
00371  */
00372 
00373 typedef struct Mesh_s{
00374   Real RootMinX[3]; /*!< min(x) in each dir on root Domain [0,1,2]=[x1,x2,x3] */
00375   Real RootMaxX[3]; /*!< max(x) in each dir on root Domain [0,1,2]=[x1,x2,x3] */
00376   Real dx[3];     /*!< cell size on root Domain [0,1,2]=[x1,x2,x3] */
00377   Real time, dt;  /*!< current time and timestep for entire Mesh */
00378   int Nx[3];    /*!< # of zones in each dir on root Domain [0,1,2]=[x1,x2,x3] */
00379   int nstep;                 /*!< number of integration steps taken */
00380   int BCFlag_ix1, BCFlag_ox1;  /*!< BC flag on root domain for inner/outer x1 */
00381   int BCFlag_ix2, BCFlag_ox2;  /*!< BC flag on root domain for inner/outer x2 */
00382   int BCFlag_ix3, BCFlag_ox3;  /*!< BC flag on root domain for inner/outer x3 */
00383   int NLevels;               /*!< overall number of refinement levels in mesh */
00384   int *DomainsPerLevel;      /*!< number of Domains per level (DPL) */
00385   DomainS **Domain;        /*!< array of Domains, indexed over levels and DPL */
00386   char *outfilename;         /*!< basename for output files containing -id#  */
00387 }MeshS;
00388 
00389 /*----------------------------------------------------------------------------*/
00390 /* OutputS: Everything for outputs. */
00391 
00392 struct Output_s;
00393 /*! \fn void (*VOutFun_t)(MeshS *pM, struct Output_s *pout)
00394  *  \brief Output function pointer. */
00395 typedef void (*VOutFun_t)(MeshS *pM, struct Output_s *pout);
00396 /*! \fn void (*VResFun_t)(MeshS *pM, struct Output_s *pout)
00397  *  \brief Restart function pointer. */
00398 typedef void (*VResFun_t)(MeshS *pM, struct Output_s *pout);
00399 /*! \fn Real (*ConsFun_t)(const GridS *pG, const int i,const int j,const int k) 
00400  *  \brief Pointer to expression that computes quant for output.*/
00401 typedef Real (*ConsFun_t)(const GridS *pG, const int i,const int j,const int k);
00402 #ifdef PARTICLES
00403 /*! \fn int (*PropFun_t)(const Grain *gr, const GrainAux *grsub)
00404  *  \brief Particle property selection function */
00405 typedef int (*PropFun_t)(const Grain *gr, const GrainAux *grsub);
00406 typedef Real (*Parfun_t)(const Grid *pG, const Grain *gr);
00407 typedef Real (*Parfun_t)(const Grid *pG, const Grain *gr);
00408 #endif
00409 
00410 /*! \struct OutputS
00411  *  \brief Everything for outputs. */
00412 typedef struct Output_s{
00413   int n;          /*!< the N from the <outputN> block of this output */
00414   Real dt;        /*!< time interval between outputs  */
00415   Real t;         /*!< next time to output */
00416   int num;        /*!< dump number (0=first) */
00417   char *out;      /*!< variable (or user fun) to be output */
00418   char *id;       /*!< filename is of the form <basename>[.idump][.id].<ext> */
00419 #ifdef PARTICLES
00420   int out_pargrid;    /*!< output grid binned particles (=1) or not (=0) */
00421   PropFun_t par_prop; /*!< particle property selection function */
00422 #endif
00423 
00424 /* level and domain number of output (default = [-1,-1] = output all levels) */
00425 
00426   int nlevel, ndomain;
00427 
00428 /* variables which describe data min/max */
00429   Real dmin,dmax;   /*!< user defined min/max for scaling data */
00430   Real gmin,gmax;   /*!< computed global min/max (over all output data) */
00431   int sdmin,sdmax;  /*!< 0 = auto scale, otherwise use dmin/dmax */
00432 
00433 /* variables which describe coordinates of output data volume */
00434   int ndim;       /*!< 3=cube 2=slice 1=vector 0=scalar */
00435   int reduce_x1;  /*!< flag to denote reduction in x1 (0=no reduction) */
00436   int reduce_x2;  /*!< flag to denote reduction in x2 (0=no reduction) */
00437   int reduce_x3;  /*!< flag to denote reduction in x3 (0=no reduction) */
00438   Real x1l, x1u;  /*!< lower/upper x1 range for data slice  */
00439   Real x2l, x2u;  /*!< lower/upper x2 range for data slice  */
00440   Real x3l, x3u;  /*!< lower/upper x3 range for data slice  */
00441 
00442 /* variables which describe output format */
00443   char *out_fmt;  /*!< output format = {bin, tab, hdf, hst, pgm, ppm, ...} */
00444   char *dat_fmt;  /*!< format string for tabular type output, e.g. "%10.5e" */
00445   char *palette;  /*!< name of palette for RGB conversions */
00446   float *rgb;     /*!< array of RGB[256*3] values derived from palette */
00447   float *der;     /*!< helper array of derivatives for interpolation into RGB */
00448 
00449 /* pointers to output functions; data expressions */
00450   VOutFun_t out_fun; /*!< output function pointer */
00451   VResFun_t res_fun; /*!< restart function pointer */
00452   ConsFun_t expr;   /*!< pointer to expression that computes quant for output */
00453 
00454 }OutputS;
00455 
00456 
00457 /*----------------------------------------------------------------------------*/
00458 /* typedefs for functions:
00459  */
00460 /* for static gravitational potential and cooling, set in problem generator,
00461  * and used by integrators */
00462 
00463 /*! \fn Real (*GravPotFun_t)(const Real x1, const Real x2, const Real x3)
00464  *  \brief Gravitational potential function. */
00465 typedef Real (*GravPotFun_t)(const Real x1, const Real x2, const Real x3);
00466 #ifdef CYLINDRICAL
00467 /*! \fn Real (*StaticGravAcc_t)(const Real x1, const Real x2, const Real x3)
00468  *  \brief Static gravitational acceleration. */
00469 typedef Real (*StaticGravAcc_t)(const Real x1, const Real x2, const Real x3);
00470 #ifdef FARGO
00471 /*! \fn Real (*OrbitalFun_t)(const Real x1)
00472  *  \brief Orbital function for FARGO */
00473 typedef Real (*OrbitalFun_t)(const Real x1);
00474 /*! \fn Real (*ShearFun_t)(const Real x1)
00475  *  \brief Shear function for FARGO */
00476 typedef Real (*ShearFun_t)(const Real x1);
00477 #endif
00478 #endif /* Cylindrical */
00479 /*! \fn Real (*CoolingFun_t)(const Real d, const Real p, const Real dt);
00480  *  \brief Cooling function. */
00481 typedef Real (*CoolingFun_t)(const Real d, const Real p, const Real dt);
00482 #ifdef RESISTIVITY
00483 /*! \fn void (*EtaFun_t)(GridS *pG, int i, int j, int k,
00484                          Real *eta_O, Real *eta_H, Real *eta_A)
00485  *  \brief Resistivity Eta Function. */
00486 typedef void (*EtaFun_t)(GridS *pG, int i, int j, int k,
00487                          Real *eta_O, Real *eta_H, Real *eta_A);
00488 #endif /* RESISTIVITY */
00489 
00490 #ifdef PARTICLES
00491 /* function types for interpolation schemes and stopping time */
00492 /*! \fn void (*WeightFun_t)(GridS *pG, Real x1, Real x2, Real x3,
00493   Real3Vector cell1, Real weight[3][3][3], int *is, int *js, int *ks);
00494  *  \brief Interpolation scheme for particles. */
00495 typedef void (*WeightFun_t)(GridS *pG, Real x1, Real x2, Real x3,
00496   Real3Vector cell1, Real weight[3][3][3], int *is, int *js, int *ks);
00497 /*! \fn Real (*TSFun_t)(GridS *pG, int type, Real rho, Real cs, Real vd)
00498  *  \brief Stopping time function for particles. */
00499 typedef Real (*TSFun_t)(GridS *pG, int type, Real rho, Real cs, Real vd);
00500 #endif /* PARTICLES */
00501 
00502 /*----------------------------------------------------------------------------*/
00503 /*! \enum BCDirection
00504  *  \brief Directions for the set_bvals_fun() function */
00505 enum BCDirection {left_x1, right_x1, left_x2, right_x2, left_x3, right_x3};
00506 
00507 #endif /* ATHENA_H */

Generated on Mon Sep 27 2010 23:03:06 for Athena by  doxygen 1.7.1