Values, macros and structures in Gerris
From Gerris
Values defined in Gerris
N_CELLS = FTT_CELLS = <math>2^{dimension}</math>
FTT_NEIGHBORS_2D = (FTT_BOTTOM) + 1
GFS_STATUS_UNDEFINED = 0
GFS_STATUS_SOLID = 1
GFS_STATUS_FLUID = 2
Gerris Macros
Ftt Macros
Name of the macro | Macro definition |
---|---|
FTT_CELLS_DIRECTION(d) | (FTT_CELLS/2) |
ftt_cell_dz(c) | (1.) |
FTT_CELL_ID(c) | ((c)->flags & FTT_FLAG_ID) |
FTT_CELL_IS_DESTROYED(c) | (((c)->flags & FTT_FLAG_DESTROYED) != 0) |
FTT_CELL_IS_FLAGGED_LEAF(cell) | (((cell)->flags & FTT_FLAG_LEAF) != 0) |
FTT_CELL_IS_LEAF(c) | ((c)->children == NULL) |
FTT_CELL_IS_ROOT(c) | ((c)->parent == NULL) |
ftt_cell_level(c) | ((c)-> parent ? (c)->parent->level + 1 : ((struct _FttRootCell *) c)->level) |
ftt_cell_parent(c) | ((c)->parent ? (c)->parent->parent : NULL) |
FTT_CELLS_DIRECTION(d) | (FTT_CELLS/2) |
FTT_FACE_DIRECT(f) | ((f)->d % 2 == 0) |
FTT_FACE_REVERSE(dst, src) | ((dst)->cell = (src)->neighbor, (dst)->neighbor = (src)->cell, (dst)->d = FTT_OPPOSITE_DIRECTION((src)->d)) |
FTT_OPPOSITE_DIRECTION(d) | (ftt_opposite_direction[d]) |
FTT_ORTHOGONAL_COMPONENT(c) | (((c) + 1) % FTT_DIMENSION) |
FTT_ROOT_CELL(c) | ((struct _FttRootCell *) c) |
gint ftt_opposite_direction[FTT_NEIGHBORS] =
#if FTT_2D
{1,0,3,2};
#else /* FTT_3D*/
{1,0,3,2,5,4};
#endif /* FTT_3D */
Gfs Macros
Name of the macro | Macro definition |
---|---|
GFS_BOUNDARY_PERIODIC(obj) | GTS_OBJECT_CAST (obj, GfsBoundaryPeriodic, gfs_boundary_periodic_class()) |
GFS_CELL_IS_BOUNDARY(cell) | (((cell)->flags & GFS_FLAG_BOUNDARY) != 0) |
GFS_CELL_IS_GRADIENT_BOUNDARY(cell) | (((cell)->flags & GFS_FLAG_GRADIENT_BOUNDARY) != 0) |
GFS_CELL_IS_PERMANENT(cell) | (((cell)->flags & GFS_FLAG_PERMANENT) != 0) |
GFS_DOUBLE_TO_POINTER(d) | (*((gpointer *) &(d))) |
GFS_FACE_FRACTION(fa) | (GFS_IS_MIXED ((fa)->cell) ? GFS_STATE ((fa)->cell)->solid->s[(fa)->d] : 1.) |
GFS_FACE_FRACTION_LEFT(fa) | GFS_FACE_FRACTION(fa) |
GFS_FACE_FRACTION_RIGHT(fa) | (GFS_IS_MIXED ((fa)->neighbor) ? GFS_STATE ((fa)->neighbor)->solid->s[FTT_OPPOSITE_DIRECTION((fa)->d)] : 1.) |
GFS_FACE_NORMAL_VELOCITY(fa) | (GFS_STATE ((fa)->cell)->f[(fa)->d].un) |
GFS_FACE_NORMAL_VELOCITY_LEFT(fa) | (GFS_STATE ((fa)->cell)->f[(fa)->d].un) |
GFS_FACE_NORMAL_VELOCITY_RIGHT(fa) | (GFS_STATE ((fa)->neighbor)->f[FTT_OPPOSITE_DIRECTION((fa)->d)].un) |
GFS_GENERIC_SURFACE_CLASS(k) | GTS_OBJECT_CLASS_CAST (k, GfsGenericSurfaceClass, gfs_generic_surface_class()) |
GFS_IS_BOUNDARY_PERIODIC(obj) | (gts_object_is_from_class (obj, gfs_boundary_periodic_class())) |
GFS_IS_FLUID(cell) | ((cell) != NULL && GFS_STATE (cell)->solid == NULL) |
GFS_IS_MIXED(cell) | ((cell) != NULL && GFS_STATE (cell)->solid != NULL) |
GFS_IS_SOLID(obj) | (gts_object_is_from_class (obj, gfs_solid_class ())) |
GFS_STATE(cell) | ((GfsStateVector *) (cell)->data) |
GFS_VALUE(cell,v) | ((&GFS_STATE (cell)->place_holder)[(v)->i]) |
GFS_VALUEI(cell, index) | ((&GFS_STATE (cell)->place_holder)[index]) |
Other Macros
Name of the macro | Macro definition |
---|---|
REFINE_CORNER(cell) | {if (cell && FTTT_CELL_IS_LEAF (cell) && ftt_cell_level (cell) < level - 1) ftt_cell_refine_single (cell, init, data);} |
Gerris structures
CellCube
typedef struct {
GtsPoint p[8];
GfsSegment s[12];
}
The vertices of the cube are stored in the GtsPoint array and the edges in the GfsSegment array.
CellFace
typedef struct {
GtsPoint p[4];
GfsSegment s[4];
}
The vertices of the face are stored in the GtsPoint array and the edges in the GfsSegment array.
CompatPar
typedef struct {
GfsVariable * lhs;
gboolean dirichlet;
}
The homogeneous version of the considered GfsVariable is stored in lhs.
DiffusionCoeff
typedef struct {
GfsSourceDiffusion * d;
gdouble lambda2[FTT_DIMENSION];
gdouble dt;
GfsVariable * rhoc, * metric;
GfsFunction * alpha;
GfsDomain * domain;
}
d is a GfsSourceDiffusion, dt is the timestep, the GfsVariable rhoc is the mass, metric is the implicit metric term, the GfsFunction alpha is the inverse of the density and domain is the considered GfsDomain.
FttCell
struct {
/*<public>*/
guint flags;
gpointer data;
/*<private>*/
struct _FttOct * parent, * children;
}
children and parent are the two FttOcts associated with the FttCell.
FttCellChildren
struct {
FttCell * c[FTT_CELLS];
};
The children of the considered FttCell are stored in c (4 in 2D and 8 in 3D).
FttCellFace
See FttCellFace.
FttCellFlags
typedef enum {
FTT_FLAG_ID = 7,
FTT_FLAG_DESTROYED = 1 << 3,
FTT_FLAG_LEAF = 1 << 4, /* used only for I/O operations */
FTT_FLAG_TRAVERSED = FTT_FLAG_LEAF, /* used for face traversal */
FTT_FLAG_USER = 5 /* user flags start here */
}
FttCellNeighbors
struct {
/* right, left, top, bottom, front, back */
FttCell * c[FTT_NEIGHBORS];
}
c is the array containing the neighbors of a considered FttCell.
FttCellTraverse
struct {
FttCell ** cells;
FttCell ** curent;
}
FttComponent
typedef enum
{
FTT_X = 0,
FTT_Y,
#if (!FTT_2D)
FTT_Z,
#endif /* FTT_3D */
FTT_DIMENSION,
FTT_XY,
#if FTT_2D
FTT_XYZ = FTT_XY
#else /* FTT_3D */
FTT_XYZ
#endif /* FTT_3D */
}
FttDirection
See FttDirection.
FttFaceType
typedef enum {
FTT_BOUNDARY,
FTT_FINE_FINE,
FTT_FINE_COARSE
}
FttOct
struct {
guint level;
FttCell * parent;
FttCellNeighbors neighbors;
FttVector pos;
FttCell cells[FTT_CELLS];
}
parent is the FttCell to which the FttOct is linked, level its level (0 for the root, 1 for the children of the root, and so on).
The FttCellNeighbors neighbors stores the neighbors of parent and the FttVector pos gives the coordinates of the center of parent. The array cells contains the children of parent.
FttRootCell
struct{
FttCell cell;
FttCellNeighbors neighbors;
FttVector pos;
guint level;
gpointer parent;
}
The FttCellNeighbors neighbors stores the neighbors of the FttCell cell (which should be FttRootCells too).
The FttVector pos gives the coordinates of the center of cell and level is the level of cell (the size of cell is <math>2^{level}</math>).
FttTraverseType
typedef enum {
FTT_PRE_ORDER,
FTT_POST_ORDER
}
When the type FTT_PRE_ORDER is used, the traversed FttCell is considered before its FttCellChildren (and all functions are applied in this order). It happens the other way when the type FTT_POST_ORDER is used.
FttVector
struct {
gdouble x, y, z;
}
FttVector defines a vector of coordinates x, y, and z.
GfsBc
struct {
/*< private >*/
GtsObject parent;
GfsLinearProblem * lp;
/*< public >*/
GfsBoundary * b;
GfsVariable * v;
gboolean extra;
FttFaceTraverseFunc bc, homogeneous_bc;
FttFaceTraverseFunc homogeneous_bc_stencil;
FttFaceTraverseFunc face_bc;
}
The GfsLinearProblem lp is the linear problem being solved by Gerris.
The boundary which is concerned by the given condition is stored in the GfsBoundary b.
The boundary condition is applied to the GfsVariable v.
bc, homogeneous_bc, homogeneous_bc_stencil and face_bc are the functions to call when traversing the tree.
GfsBc is used to define a default symmetry boundary condition.
GfsBoundary
struct {
/*< private >*/
GtsObject parent;
FttCell * root;
GfsBox * box;
FttDirection d;
guint depth;
GfsBc * default_bc;
gboolean changed;
/*< public >*/
GfsVariable * v;
GfsBoundaryVariableType type;
GHashTable * bc;
}
The boundary defined by this GfsBoundary structure is the one associated with the FttCellFace in FttDirection d of FttCell root.
default_bc is the default symmetry boundary condition (see GfsBc) and box the GfsBox associated with root.
The boundary condition is applied to the GfsVariable v.
GfsBoundary is used to define one of the boundaries of a GfsBox.
GfsBoundaryPeriodic
struct {
/*> private >*/
GfsBoundary parent;
GfsBox * matching;
FttDirection d;
GArray * sndbuf, * rcvbuf;
guint sndcount, rcvcount;
gdouble rotate;
}
parent is the GfsBoundary concerned by the periodic condition.
This boundary is associated with the GfsBox matching in FttDirection d.
GfsBox
struct {
GtsGNode parent;
FttCell * root;
GtsObject * neighbor[FTT_NEIGHBORS];
guint id;
int pid;
gint size;
}
A GfsBox is attached to a FttCell root.
The neighbors are stored in the GtsObject array neighbor.
Its size is given by size, the process identifier associated to this GfsBox by pid, and its id by id.
GfsClock
struct {
gboolean started;
glong start, stop;
}
GfsDomain
struct {
GtsWGraph parent;
int pid;
GfsClock * timer;
GHashTable * timers;
GtsRange timestep;
GtsRange size;
gboolean profile_bc;
GtsRange mpi_messages;
GtsRange mpi_wait;
guint rootlevel;
FttVector refpos;
FttVector lambda;
GArray * allocated;
GSList * variables;
GSList * derived_variables;
GfsVariable * velocity[FTT_DIMENSION];
GSList * variables_io;
gboolean binary;
gint max_depth_write;
FttCellInitFunc cell_init;
gpointer cell_init_data;
gint version;
gpointer array;
gboolean overlap; /* whether to overlap MPI communications with computation */
/* coordinate metrics */
gpointer metric_data;
gdouble (* face_metric) (const GfsDomain *, const FttCellFace *);
gdouble (* cell_metric) (const GfsDomain *, const FttCell *);
gdouble (* solid_metric) (const GfsDomain *, const FttCell *);
gdouble (* scale_metric) (const GfsDomain *, const FttCell *, FttComponent);
gdouble (* face_scale_metric) (const GfsDomain *, const FttCellFace *, FttComponent);
gdouble (* viscous_metric_implicit) (const GfsDomain * domain,
FttCell * cell,
FttComponent component);
gdouble (* viscous_metric_explicit) (const GfsDomain * domain,
FttCell * cell,
GfsVariable * v,
GfsDiffusion * d);
void (* advection_metric) (const GfsDomain * domain,
FttCell * cell,
FttComponent c1,
gdouble m[2]);
/* Object hash table for (optional) objects IDs */
GHashTable * objects;
/* total number of parallel processes */
int np;
/* real time */
GTimer * clock;
GPtrArray * sorted; /**< array of sorted boxes */
gboolean dirty; /**< whether the sorted array needs updating */
}
parent is the graph defining the geometry of the considered domain. All computations are realised with pid pid.
rootlevel is the root level (default 0). The FttVector refpos gives the position of the center of the domain (default (0,0,0) ). variables is a pointer to the list of all variables computed in the domain, and derived_variables a pointer to the list of all derivatives.
velocity is a GfsVariable pointer to the table of all FTT_DIMENSION components of the velocity of the fluid.
The GfsClock timer is used to compute the CPU time, whereas the GTimer clock is used for the real time.
GfsEvent
struct {
GtsSListContainee parent;
gdouble t, start, end, step;
guint i, istart, iend, istep;
guint n;
gboolean end_event, realised, redo;
gchar * name;
}
Events are used to control any action which needs to be performed at a given time during a simulation. This includes one-off actions as well as periodically repeated actions.
name is the (optional) name of the event, start the physical starting time (default is zero). The end keyword can be used to indicate the end of the simulation. istart is the time step starting time (default is zero).
step means repeat every step physical time units (default is infinity), istep repeat every istep time steps (default is infinity).
end is the ending physical time (default is infinity; the event can stop before this time if it is never reached) and iend the ending number of time steps (default is infinity).
GfsFaceStateVector
struct {
gdouble un;
gdouble v;
}
GfsFlags
typedef enum {
GFS_FLAG_USED = 1 << FTT_FLAG_USER,
GFS_FLAG_BOUNDARY = 1 << (FTT_FLAG_USER + 1),
GFS_FLAG_DIRICHLET = 1 << (FTT_FLAG_USER + 2),
GFS_FLAG_GRADIENT_BOUNDARY = 1 << (FTT_FLAG_USER + 3),
GFS_FLAG_PERMANENT = 1 << (FTT_FLAG_USER + 4),
GFS_FLAG_THIN = 1 << (FTT_FLAG_USER + 5),
GFS_FLAG_USER = FTT_FLAG_USER + 6 /* user flags start here */
}
GfsGradient
struct {
gdouble a,b;
}
Used for interpolations : v = a*v(cell) + b
GfsInterpolator
struct {
#if FTT_2D
FttCell * c[7];
gdouble w[7];
#else /* 3D */
FttCell * c[29];
gdouble w[29];
#endif /* 3D */
guint n;
}
GfsLinearProblem
struct {
GPtrArray * LP;
GArray * rhs, * lhs;
GfsVariable * id, * neighbor, * neighborw;
gint istart;
}
GfsMultilevelParams
struct {
gdouble tolerance;
guint nrelax, erelax;
guint minlevel ;
guint nitermax, nitermin;
guint dimension;
guint niter;
guint depth;
gboolean weighted, function;
gdouble beta, omega;
GfsNorm residual_before, residual;
GfsPoissonSoolverFunc poisson_solve;
}
This object stores the different parameters used for the multigrid relaxation : the authorised tolerance for the residual, the minimum and maximum numbers of iterations nitermin and nitermax, the minimum level to use (root is 0) minlevel.
It also keeps the dimension, the number of iterations niter, the different norms of the residual in the GfsNorm residual_before for the beginning and residual for the current one.
GfsNorm
struct {
gdouble bias, first, second, infty, w;
}
GfsNorm stores the values of the different norms of the residual.
GfsSegment
struct {
GtsPoint * E, * D;
gdouble x;
guint n;
gint inside;
⇐ previous: GfsView, ⇑ up: Gerris Flow Solver Programming Course for Dummies