Gerris
Functions

fluid.c File Reference

Low-level discrete operators. More...

Functions

FttCellFace gfs_cell_face (FttCell *cell, FttDirection d)
gdouble gfs_neighbor_value (const FttCellFace *face, guint v, gdouble *x)
gdouble gfs_center_gradient (FttCell *cell, FttComponent c, guint v)
void gfs_center_gradient_stencil (FttCell *cell, FttComponent c, guint v)
gdouble gfs_center_van_leer_gradient (FttCell *cell, FttComponent c, guint v)
gdouble gfs_center_minmod_gradient (FttCell *cell, FttComponent c, guint v)
gdouble gfs_center_superbee_gradient (FttCell *cell, FttComponent c, guint v)
gdouble gfs_center_sweby_gradient (FttCell *cell, FttComponent c, guint v)
gdouble gfs_center_regular_gradient (FttCell *cell, FttComponent c, GfsVariable *v)
gdouble gfs_center_regular_2nd_derivative (FttCell *cell, FttComponent c, GfsVariable *v)
void gfs_face_gradient (const FttCellFace *face, GfsGradient *g, guint v, gint max_level)
void gfs_face_weighted_gradient (const FttCellFace *face, GfsGradient *g, guint v, gint max_level)
void gfs_face_cm_gradient (const FttCellFace *face, GfsGradient *g, guint v, gint max_level)
void gfs_face_cm_weighted_gradient (const FttCellFace *face, GfsGradient *g, guint v, gint max_level)
void gfs_cell_dirichlet_gradient (FttCell *cell, guint v, gint max_level, gdouble v0, FttVector *grad)
void gfs_mixed_cell_gradient (FttCell *cell, GfsVariable *v, FttVector *g)
gdouble gfs_mixed_cell_interpolate (FttCell *cell, FttVector p, GfsVariable *v)
gdouble gfs_cell_dirichlet_gradient_flux (FttCell *cell, guint v, gint max_level, gdouble v0)
gdouble gfs_cell_dirichlet_gradient_flux_stencil (FttCell *cell, gint max_level, gdouble v0, GfsLinearProblem *lp, GfsStencil *stencil)
gdouble gfs_cell_dirichlet_value (FttCell *cell, GfsVariable *v, gint max_level)
void gfs_shear_strain_rate_tensor (FttCell *cell, GfsVariable **u, gdouble t[FTT_DIMENSION][FTT_DIMENSION])
gdouble gfs_2nd_principal_invariant (FttCell *cell, GfsVariable **u)
void gfs_get_from_below_intensive (FttCell *cell, const GfsVariable *v)
void gfs_cell_coarse_fine (FttCell *parent, GfsVariable *v)
void gfs_cell_cleanup (FttCell *cell, GfsDomain *domain)
void gfs_cell_reset (FttCell *cell, GfsVariable *v)
void gfs_face_reset (FttCellFace *face, GfsVariable *v)
GtsRange gfs_stats_variable (FttCell *root, GfsVariable *v, FttTraverseFlags flags, gint max_depth)
GfsNorm gfs_norm_variable (FttCell *root, GfsVariable *v, FttTraverseFlags flags, gint max_depth)
void gfs_norm_init (GfsNorm *n)
void gfs_norm_reset (GfsNorm *n)
void gfs_norm_add (GfsNorm *n, gdouble val, gdouble weight)
void gfs_norm_update (GfsNorm *n)
gdouble gfs_face_interpolated_value (const FttCellFace *face, guint v)
gdouble gfs_face_interpolated_value_generic (const FttCellFace *face, const GfsVariable *v)
gdouble gfs_face_weighted_interpolated_value (const FttCellFace *face, guint v)
void gfs_normal_divergence (FttCell *cell, GfsVariable *v)
void gfs_normal_divergence_2D (FttCell *cell, GfsVariable *v)
gdouble gfs_divergence (FttCell *cell, GfsVariable **v)
gdouble gfs_vorticity (FttCell *cell, GfsVariable **v)
gdouble gfs_vector_norm2 (FttCell *cell, GfsVariable **v)
gdouble gfs_vector_norm (FttCell *cell, GfsVariable **v)
gdouble gfs_vector_lambda2 (FttCell *cell, GfsVariable **v)
void gfs_pressure_force (FttCell *cell, GfsVariable *p, FttVector *f)
void gfs_cell_traverse_mixed (FttCell *root, FttTraverseType order, FttTraverseFlags flags, FttCellTraverseFunc func, gpointer data)
void gfs_cell_corner_values (FttCell *cell, GfsVariable *v, gint max_level, gdouble f[4 *(FTT_DIMENSION-1)+1])
gdouble gfs_interpolate_from_corners (FttCell *cell, FttVector p, gdouble *f)
gdouble gfs_interpolate (FttCell *cell, FttVector p, GfsVariable *v)
void gfs_interpolate_stencil (FttCell *cell, GfsVariable *v)
gdouble gfs_center_curvature (FttCell *cell, FttComponent c, guint v)
gdouble gfs_streamline_curvature (FttCell *cell, GfsVariable **v)
void gfs_cell_corner_interpolator (FttCell *cell, FttDirection d[FTT_DIMENSION], gint max_level, gboolean centered, GfsInterpolator *inter)
gdouble gfs_cell_corner_value (FttCell *cell, FttDirection d[FTT_DIMENSION], GfsVariable *v, gint max_level)
GfsStencil * gfs_stencil_new (FttCell *cell, GfsLinearProblem *lp, gdouble coeff)
void gfs_stencil_add_element (GfsStencil *stencil, FttCell *cell, GfsLinearProblem *lp, gdouble coeff)
void gfs_stencil_destroy (GfsStencil *stencil)
void gfs_face_weighted_gradient_stencil (const FttCellFace *face, GfsGradient *g, gint max_level, GfsLinearProblem *lp, GfsStencil *stencil)
void gfs_face_cm_weighted_gradient_stencil (const FttCellFace *face, GfsGradient *g, gint max_level, GfsLinearProblem *lp, GfsStencil *stencil)
void gfs_cm_gradient (FttCell *cell, GfsVariable *v, FttVector *g)

Detailed Description

Low-level discrete operators.


Function Documentation

gdouble gfs_2nd_principal_invariant ( FttCell *  cell,
GfsVariable **  u 
)
Parameters:
cella #FttCell.
uthe velocity.
Returns:
the second principal invariant of the shear strain rate tensor of cell: sqrt(D:D).

Here is the call graph for this function:

void gfs_cell_cleanup ( FttCell *  cell,
GfsDomain *  domain 
)
Parameters:
cella #FttCell.
domaina #GfsDomain.

Frees the memory allocated for extra data associated with cell.

This function must be used as "cleanup function" when using ftt_cell_destroy().

Here is the caller graph for this function:

void gfs_cell_coarse_fine ( FttCell *  parent,
GfsVariable *  v 
)
Parameters:
parenta #FttCell.
va #GfsVariable.

Initializes v on the children of parent using interpolation.

First-order interpolation (straight injection) is used for boundary cells and second-order interpolation for the other cells.

Here is the call graph for this function:

void gfs_cell_corner_interpolator ( FttCell *  cell,
FttDirection  d[FTT_DIMENSION],
gint  max_level,
gboolean  centered,
GfsInterpolator *  inter 
)
Parameters:
cella #FttCell.
da set of perpendicular directions.
max_levelthe maximum cell level to consider (-1 means no restriction).
centeredTRUE if the interpolator is cell-centered.
intera #GfsInterpolator.

Fills inter with the interpolator for the corner of cell defined by d.

Here is the call graph for this function:

Here is the caller graph for this function:

gdouble gfs_cell_corner_value ( FttCell *  cell,
FttDirection  d[FTT_DIMENSION],
GfsVariable *  v,
gint  max_level 
)
Parameters:
cella #FttCell.
da set of perpendicular directions.
va #GfsVariable.
max_levelthe maximum cell level to consider (-1 means no restriction).
Returns:
the value of variable v interpolated at the corner of cell defined by d.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_cell_corner_values ( FttCell *  cell,
GfsVariable *  v,
gint  max_level,
gdouble  f[4 *(FTT_DIMENSION-1)+1] 
)
Parameters:
cella #FttCell containing location p.
va #GfsVariable.
max_levelthe maximum cell level to consider (-1 means no restriction).
fan array with the correct size (4*(FTT_DIMENSION - 1) + 1).

Fills f with the values of v interpolated at the corners of cell.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_cell_dirichlet_gradient ( FttCell *  cell,
guint  v,
gint  max_level,
gdouble  v0,
FttVector *  grad 
)
Parameters:
cella #FttCell.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).
v0the Dirichlet value on the boundary.
grada #FttVector.

Fills grad with components of the gradient of variable v interpolated at the center of area of the solid boundary contained in cell. The gradient is scaled by the size of the cell.

Here is the call graph for this function:

Here is the caller graph for this function:

gdouble gfs_cell_dirichlet_gradient_flux ( FttCell *  cell,
guint  v,
gint  max_level,
gdouble  v0 
)
Parameters:
cella #FttCell.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).
v0the Dirichlet value on the boundary.
Returns:
the flux of the gradient of variable v through the solid boundary contained in cell.

Here is the call graph for this function:

gdouble gfs_cell_dirichlet_gradient_flux_stencil ( FttCell *  cell,
gint  max_level,
gdouble  v0,
GfsLinearProblem *  lp,
GfsStencil *  stencil 
)
Parameters:
cella #FttCell.
max_levelthe maximum cell level to consider (-1 means no restriction).
v0the Dirichlet value on the boundary.
lpthe #GfsLinearProblem.
stencila GfsStencil.

Returns the stencil accounting for the flux of the gradient of a given variable through the solid boundary contained in cell. Returns the stencil equivalent to the action of gfs_cell_dirichlet_gradient_flux().

Returns:
the flux of the gradient of variable v through the solid boundary contained in cell.
gdouble gfs_cell_dirichlet_value ( FttCell *  cell,
GfsVariable *  v,
gint  max_level 
)
Parameters:
cella #FttCell.
va #GfsVariable.
max_levelthe maximum cell level to consider (-1 means no restriction).
Returns:
the value of variable v interpolated at the center of area of the solid boundary contained in cell.

Here is the call graph for this function:

Here is the caller graph for this function:

FttCellFace gfs_cell_face ( FttCell *  cell,
FttDirection  d 
)
Parameters:
cella #FttCell.
da direction.

This function is different from ftt_cell_face() because it takes into account the solid fractions.

Returns:
the face of cell in direction d.

Here is the caller graph for this function:

void gfs_cell_reset ( FttCell *  cell,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
va #GfsVariable to reset.

Sets the value of the variable v of cell to zero.

Here is the caller graph for this function:

void gfs_cell_traverse_mixed ( FttCell *  root,
FttTraverseType  order,
FttTraverseFlags  flags,
FttCellTraverseFunc  func,
gpointer  data 
)
Parameters:
rootthe root #FttCell of the tree to traverse.
orderthe order in which the cells are visited - FTT_PRE_ORDER, FTT_POST_ORDER.
flagswhich types of children are to be visited.
functhe function to call for each visited #FttCell.
datauser data to pass to func.

Traverses a cell tree starting at the given root #FttCell. Calls the given function for each mixed cell.

gdouble gfs_center_curvature ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The curvature is normalized by the square of the size of the cell.

Returns:
the value of the c component of the curvature of variable v at the center of the cell.

Here is the call graph for this function:

gdouble gfs_center_gradient ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The gradient is normalized by the size of the cell.

Returns:
the value of the c component of the gradient of variable v at the center of the cell.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_center_gradient_stencil ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

Sets to 1. the v variable of all the cells which would be used if gfs_center_gradient() was called with identical arguments.

gdouble gfs_center_minmod_gradient ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The gradient is normalized by the size of the cell and is limited using a minmod limiter.

Returns:
the value of the c component of the gradient of variable v at the center of the cell.
gdouble gfs_center_regular_2nd_derivative ( FttCell *  cell,
FttComponent  c,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable.

The derivative is normalized by the size of the cell squared. Only regular Cartesian stencils are used to compute the derivative.

Returns:
the value of the c component of the 2nd derivative of variable v at the center of the cell.
gdouble gfs_center_regular_gradient ( FttCell *  cell,
FttComponent  c,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable.

The gradient is normalized by the size of the cell. Only regular Cartesian stencils are used to compute the gradient.

Returns:
the value of the c component of the gradient of variable v at the center of the cell.
gdouble gfs_center_superbee_gradient ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The gradient is normalized by the size of the cell and is limited using a superbee limiter.

Returns:
the value of the c component of the gradient of variable v at the center of the cell.
gdouble gfs_center_sweby_gradient ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The gradient is normalized by the size of the cell and is limited using a Sweby limiter (beta = 1.5).

Returns:
the value of the c component of the gradient of variable v at the center of the cell.
gdouble gfs_center_van_leer_gradient ( FttCell *  cell,
FttComponent  c,
guint  v 
)
Parameters:
cella #FttCell.
ca component.
va #GfsVariable index.

The gradient is normalized by the size of the cell and is limited using van Leer's limiter.

Returns:
the value of the c component of the gradient of variable v at the center of the cell.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_cm_gradient ( FttCell *  cell,
GfsVariable *  v,
FttVector *  g 
)
Parameters:
cella #FttCell.
va #GfsVariable.
ga #FttVector.

Fills g with the components of the gradient of v at the center of mass of cell.

The gradient is normalized by the size of the cell.

Here is the call graph for this function:

gdouble gfs_divergence ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.
Returns:
the divergence of the (centered) vector field v in cell.

Here is the call graph for this function:

void gfs_face_cm_gradient ( const FttCellFace *  face,
GfsGradient *  g,
guint  v,
gint  max_level 
)
Parameters:
facea #FttCellFace.
gthe #GfsGradient.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).

Set the value of g as the gradient of variable v on the face. Variable v is defined at the center of mass of its cell. Linear interpolation is used to evaluate the gradient in the vicinity of cut cells.

Here is the caller graph for this function:

void gfs_face_cm_weighted_gradient ( const FttCellFace *  face,
GfsGradient *  g,
guint  v,
gint  max_level 
)
Parameters:
facea #FttCellFace.
gthe #GfsGradient.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).

Set the value of g as the gradient of variable v on the face weighted by the value of the v field of the face state vector of the corresponding cell. Variable v is defined at the center of mass of its cell. Linear interpolation is used to evaluate the gradient in the vicinity of cut cells.

void gfs_face_cm_weighted_gradient_stencil ( const FttCellFace *  face,
GfsGradient *  g,
gint  max_level,
GfsLinearProblem *  lp,
GfsStencil *  stencil 
)
Parameters:
facea #FttCellFace.
ga GfsGradient.
max_levelthe maximum cell level to consider (-1 means no restriction).
lpthe linear problem.
stencila stencil.

Fills stencil with the stencil corresponding to the weighted gradient on face. The weights are defined at the center of mass of the cell. Linear interpolation is used to evaluate the gradient in the vicinity of cut cells. This function returns the stencil equivalent to the action of gfs_face_cm_weighted_gradient().

void gfs_face_gradient ( const FttCellFace *  face,
GfsGradient *  g,
guint  v,
gint  max_level 
)
Parameters:
facea #FttCellFace.
gthe #GfsGradient.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).

Set the value of g as the gradient of variable v on the face. The value returned is second order accurate in space and conservative, in the sense that values at a coarse/fine cell boundary are consistent.

The value of the gradient (normalised by the size of face->cell) is given by: g->b - g->a*GFS_VALUE (face->cell, v).

Here is the caller graph for this function:

gdouble gfs_face_interpolated_value ( const FttCellFace *  face,
guint  v 
)
Parameters:
facea #FttFace.
va #GfsVariable index.

Computes the value of variable v on the face using second-order interpolation from the cell-centered values.

Note that this function cannot be called for coarse -> fine faces (use gfs_face_interpolated_value_generic() instead).

Returns:
the value of variable v on the face.

Here is the call graph for this function:

Here is the caller graph for this function:

gdouble gfs_face_interpolated_value_generic ( const FttCellFace *  face,
const GfsVariable *  v 
)
Parameters:
facea #FttFace.
va #GfsVariable.

Same as gfs_face_interpolated_value() but also works for coarse -> fine faces.

Returns:
the value of variable v on the face.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_face_reset ( FttCellFace *  face,
GfsVariable *  v 
)
Parameters:
facea #FttCellFace.
va #GfsVariable to reset.

Sets the value of the variable v of face to zero (on both sides).

void gfs_face_weighted_gradient ( const FttCellFace *  face,
GfsGradient *  g,
guint  v,
gint  max_level 
)
Parameters:
facea #FttCellFace.
gthe #GfsGradient.
va #GfsVariable index.
max_levelthe maximum cell level to consider (-1 means no restriction).

Set the value of g as the gradient of variable v on the face weighted by the value of the v field of the face state vector of the corresponding cell. The value returned is second order accurate in space and conservative, in the sense that values at a coarse/fine cell boundary are consistent.

void gfs_face_weighted_gradient_stencil ( const FttCellFace *  face,
GfsGradient *  g,
gint  max_level,
GfsLinearProblem *  lp,
GfsStencil *  stencil 
)
Parameters:
facea FttCellFace
ga GfsGradient
max_levelan integer
lpthe linear problem
stencila stencil

Fills stencil with the stencil corresponding to the weighted gradient on face.

gdouble gfs_face_weighted_interpolated_value ( const FttCellFace *  face,
guint  v 
)
Parameters:
facea #FttFace.
va #GfsVariable index.

Computes the value of variable v on the face weighted by the value of the v field of the face state vector using interpolation from the cell-centered values. The value returned is second order accurate in space and conservative, in the sense that values at a coarse/fine cell boundary are consistent.

Returns:
the weighted value of variable v on the face.

Here is the call graph for this function:

void gfs_get_from_below_intensive ( FttCell *  cell,
const GfsVariable *  v 
)
Parameters:
cella #FttCell.
va #GfsVariable to "get from below".

Sets the value of the "intensive" variable v of cell by taking the volume weighted average of the values of its children cells.

This functions fails if cell is a leaf of the cell tree.

Here is the caller graph for this function:

gdouble gfs_interpolate ( FttCell *  cell,
FttVector  p,
GfsVariable *  v 
)
Parameters:
cella #FttCell containing location p.
pthe location at which to interpolate.
va #GfsVariable.

Interpolates the v variable of cell, at location p. Linear interpolation is used and the boundaries of the domain are treated as planes of symmetry for all variables.

Returns:
the interpolated value of variable v at location p.

Here is the call graph for this function:

Here is the caller graph for this function:

gdouble gfs_interpolate_from_corners ( FttCell *  cell,
FttVector  p,
gdouble *  f 
)
Parameters:
cella #FttCell containing location p.
pthe location at which to interpolate.
fvalues at the corners of cell.
Returns:
the interpolated value at location p.

Here is the call graph for this function:

Here is the caller graph for this function:

void gfs_interpolate_stencil ( FttCell *  cell,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
va #GfsVariable.

Sets to 1. the v variable of all the cells which would be used by a call to gfs_interpolate().

Here is the call graph for this function:

void gfs_mixed_cell_gradient ( FttCell *  cell,
GfsVariable *  v,
FttVector *  g 
)
Parameters:
cella mixed #FttCell.
va #GfsVariable.
gthe gradient.

Fills g with the components of the gradient of v at the center of mass of cell.

The gradient is normalized by the size of the cell.

Here is the call graph for this function:

Here is the caller graph for this function:

gdouble gfs_mixed_cell_interpolate ( FttCell *  cell,
FttVector  p,
GfsVariable *  v 
)
Parameters:
cella mixed #FttCell.
pa #FttVector.
va #GfsVariable.
Returns:
the value of variable v interpolated at position p within cell.

Here is the call graph for this function:

gdouble gfs_neighbor_value ( const FttCellFace *  face,
guint  v,
gdouble *  x 
)
Parameters:
facea #FttCellFace.
va variable index.
xrelative position of the neigboring value.
Returns:
the (averaged) neighboring value of variable v at face and its position x.

Here is the caller graph for this function:

void gfs_norm_add ( GfsNorm *  n,
gdouble  val,
gdouble  weight 
)
Parameters:
na #GfsNorm.
vala value to add to n.
weightweight of val.

Adds val to n.

void gfs_norm_init ( GfsNorm *  n)
Parameters:
na #GfsNorm.

Initializes a #GfsNorm.

Here is the caller graph for this function:

void gfs_norm_reset ( GfsNorm *  n)
Parameters:
na #GfsNorm.

Sets all the fields of n to 0.

void gfs_norm_update ( GfsNorm *  n)
Parameters:
na #GfsNorm.

Updates the fields of n.

Here is the caller graph for this function:

GfsNorm gfs_norm_variable ( FttCell *  root,
GfsVariable *  v,
FttTraverseFlags  flags,
gint  max_depth 
)
Parameters:
rootthe root #FttCell of the tree to obtain norm from.
vthe variable to consider for norm statistics.
flagswhich types of cells are to be visited.
max_depthmaximum depth of the traversal.

Traverses the cell tree defined by root using ftt_cell_traverse() and gathers norm statistics about variable v.

Returns:
a #GfsNorm containing the norm statistics about v.

Here is the call graph for this function:

void gfs_normal_divergence ( FttCell *  cell,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
va #GfsVariable.

Adds to variable v of cell the integral of the divergence of the (MAC) velocity field in this cell.

void gfs_normal_divergence_2D ( FttCell *  cell,
GfsVariable *  v 
)
Parameters:
cella #FttCell.
va #GfsVariable.

Fills variable v of cell with the integral of the 2D divergence of the (MAC) velocity field in this cell.

void gfs_pressure_force ( FttCell *  cell,
GfsVariable *  p,
FttVector *  f 
)
Parameters:
cella #FttCell.
pa #GfsVariable.
fa #FttVector.

Fills f with the pressure p component of the force exerted by the fluid on the fraction of embedded boundary contained in cell.

Here is the call graph for this function:

void gfs_shear_strain_rate_tensor ( FttCell *  cell,
GfsVariable **  u,
gdouble  t[FTT_DIMENSION][FTT_DIMENSION] 
)
Parameters:
cella #FttCell.
uthe velocity.
tthe shear strain rate tensor t[i][j] = (d_iu_j+d_ju_i)/2.

Fills t with the shear strain rate tensor at the center of mass of cell, normalised by the size of the cell.

Here is the call graph for this function:

Here is the caller graph for this function:

GtsRange gfs_stats_variable ( FttCell *  root,
GfsVariable *  v,
FttTraverseFlags  flags,
gint  max_depth 
)
Parameters:
rootthe root #FttCell of the tree to obtain statistics from.
vthe variable to consider for statistics.
flagswhich types of cells are to be visited.
max_depthmaximum depth of the traversal.

Traverses the cell tree defined by root using ftt_cell_traverse() and gathers statistics about variable v.

Returns:
a #GtsRange containing the statistics about v.

Here is the call graph for this function:

void gfs_stencil_add_element ( GfsStencil *  stencil,
FttCell *  cell,
GfsLinearProblem *  lp,
gdouble  coeff 
)
Parameters:
stencilthe stencil
cella #FttCell.
lpthe linear problem
coeffthe cell coefficient.

Adds the contribution of a cell to stencil.

Here is the caller graph for this function:

void gfs_stencil_destroy ( GfsStencil *  stencil)
Parameters:
stencila stencil

Destroys a stencil.

Here is the caller graph for this function:

GfsStencil* gfs_stencil_new ( FttCell *  cell,
GfsLinearProblem *  lp,
gdouble  coeff 
)
Parameters:
cellthe #FttCell associated with this stencil.
lpthe linear problem.
coeffthe cell coefficient.
Returns:
a new #GfsStencil.

Here is the call graph for this function:

gdouble gfs_streamline_curvature ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.

The curvature is normalized by the size of the cell.

Returns:
the value of the curvature of the streamline defined by v passing through the center of the cell.

Here is the call graph for this function:

gdouble gfs_vector_lambda2 ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.
Returns:
The value of the lambda2 eigenvalue used by Jeong and Hussain as vortex criterion (JFM 285, 69-94, 1995), normalized by the square of the size of the cell.

Here is the call graph for this function:

gdouble gfs_vector_norm ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.
Returns:
the norm of the vector field v in cell.

Here is the call graph for this function:

gdouble gfs_vector_norm2 ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.
Returns:
the squared norm of the vector field v in cell.

Here is the caller graph for this function:

gdouble gfs_vorticity ( FttCell *  cell,
GfsVariable **  v 
)
Parameters:
cella #FttCell.
vthe components of the vector.
Returns:
the vorticity (norm of the vorticity vector in 3D) of the vector field v in cell.

Here is the call graph for this function:

 All Data Structures Files Functions Variables