Gerris
|
Poisson and diffusion solvers. More...
Functions | |
void | gfs_multilevel_params_write (GfsMultilevelParams *par, FILE *fp) |
void | gfs_multilevel_params_stats_write (GfsMultilevelParams *par, FILE *fp) |
GfsLinearProblem * | gfs_linear_problem_new (GfsDomain *domain) |
void | gfs_linear_problem_add_stencil (GfsLinearProblem *lp, GfsStencil *stencil) |
void | gfs_linear_problem_destroy (GfsLinearProblem *lp) |
GfsLinearProblem * | gfs_get_poisson_problem (GfsDomain *domain, GfsVariable *rhs, GfsVariable *lhs, GfsVariable *dia, gint maxlevel, GfsVariable *v) |
void | gfs_relax (GfsDomain *domain, guint d, gint max_depth, gdouble omega, GfsVariable *u, GfsVariable *rhs, GfsVariable *dia) |
void | gfs_residual (GfsDomain *domain, guint d, FttTraverseFlags flags, gint max_depth, GfsVariable *u, GfsVariable *rhs, GfsVariable *dia, GfsVariable *res) |
void | gfs_poisson_coefficients (GfsDomain *domain, GfsFunction *alpha, gboolean positive, gboolean centered, gboolean reset) |
void | gfs_source_tension_coefficients (GfsSourceTension *s, GfsDomain *domain, GfsFunction *alpha) |
void | gfs_poisson_cycle (GfsDomain *domain, GfsMultilevelParams *p, GfsVariable *u, GfsVariable *rhs, GfsVariable *dia, GfsVariable *res) |
gdouble | gfs_poisson_compatibility (GfsDomain *domain, GfsVariable *lhs, GfsVariable *rhs, gdouble dt) |
void | gfs_poisson_solve (GfsDomain *domain, GfsMultilevelParams *par, GfsVariable *lhs, GfsVariable *rhs, GfsVariable *res, GfsVariable *dia, gdouble dt) |
void | gfs_diffusion_coefficients (GfsDomain *domain, GfsSourceDiffusion *d, gdouble dt, GfsVariable *rhoc, GfsVariable *metric, GfsFunction *alpha, gdouble beta) |
void | gfs_diffusion_rhs (GfsDomain *domain, GfsVariable *v, GfsVariable *rhs, GfsVariable *rhoc, GfsVariable *metric, gdouble beta) |
void | gfs_diffusion_residual (GfsDomain *domain, GfsVariable *u, GfsVariable *rhs, GfsVariable *rhoc, GfsVariable *metric, GfsVariable *res) |
void | gfs_diffusion_cycle (GfsDomain *domain, guint levelmin, guint depth, guint nrelax, GfsVariable *u, GfsVariable *rhs, GfsVariable *rhoc, GfsVariable *metric, GfsVariable *res) |
GfsLinearProblem * | gfs_get_diffusion_problem (GfsDomain *domain, GfsVariable *rhs, GfsVariable *lhs, GfsVariable *rhoc, GfsVariable *metric, gint maxlevel, GfsVariable *v) |
Poisson and diffusion solvers.
void gfs_diffusion_coefficients | ( | GfsDomain * | domain, |
GfsSourceDiffusion * | d, | ||
gdouble | dt, | ||
GfsVariable * | rhoc, | ||
GfsVariable * | metric, | ||
GfsFunction * | alpha, | ||
gdouble | beta | ||
) |
domain | a #GfsDomain. |
d | a #GfsSourceDiffusion. |
dt | the time-step. |
rhoc | where to store the mass. |
metric | where to store the implicit metric term (or NULL). |
alpha | the inverse of density or NULL. |
beta | the implicitness parameter (0.5 Crank-Nicholson, 1. backward Euler). |
Initializes the face coefficients for the diffusion equation.
void gfs_diffusion_cycle | ( | GfsDomain * | domain, |
guint | levelmin, | ||
guint | depth, | ||
guint | nrelax, | ||
GfsVariable * | u, | ||
GfsVariable * | rhs, | ||
GfsVariable * | rhoc, | ||
GfsVariable * | metric, | ||
GfsVariable * | res | ||
) |
domain | the domain on which to solve the diffusion equation. |
levelmin | the top level of the multigrid hierarchy. |
depth | the total depth of the domain. |
nrelax | the number of relaxations to apply at each level. |
u | the variable to use as left-hand side. |
rhs | the right-hand side. |
rhoc | the mass. |
metric | the metric term. |
res | the residual. |
Apply one multigrid iteration to the diffusion equation for u.
The initial value of res on the leaves of root must be set to the residual of the diffusion equation using gfs_diffusion_residual().
The diffusion coefficients must be set using gfs_diffusion_coefficients().
The values of u on the leaf cells are updated as well as the values of res (i.e. the cell tree is ready for another iteration).
void gfs_diffusion_residual | ( | GfsDomain * | domain, |
GfsVariable * | u, | ||
GfsVariable * | rhs, | ||
GfsVariable * | rhoc, | ||
GfsVariable * | metric, | ||
GfsVariable * | res | ||
) |
domain | a #GfsDomain. |
u | the variable to use as left-hand side. |
rhs | the right-hand side. |
rhoc | the mass. |
metric | the metric term. |
res | the residual. |
Sets the res variable of each leaf cell of domain to the residual of the diffusion equation for v.
The diffusion coefficients must have been set using gfs_diffusion_coefficients() and the right-hand side using gfs_diffusion_rhs().
void gfs_diffusion_rhs | ( | GfsDomain * | domain, |
GfsVariable * | v, | ||
GfsVariable * | rhs, | ||
GfsVariable * | rhoc, | ||
GfsVariable * | metric, | ||
gdouble | beta | ||
) |
domain | a #GfsDomain. |
v | a #GfsVariable. |
rhs | a #GfsVariable. |
rhoc | the mass. |
metric | the metric term. |
beta | the implicitness parameter (0.5 Crank-Nicholson, 1. backward Euler). |
Adds to the rhs variable of cell the right-hand side of the diffusion equation for variable v.
The diffusion coefficients must have been already set using gfs_diffusion_coefficients().
GfsLinearProblem* gfs_get_diffusion_problem | ( | GfsDomain * | domain, |
GfsVariable * | rhs, | ||
GfsVariable * | lhs, | ||
GfsVariable * | rhoc, | ||
GfsVariable * | metric, | ||
gint | maxlevel, | ||
GfsVariable * | v | ||
) |
domain | the domain over which the poisson problem is defined |
rhs | the variable to use as right-hand side |
lhs | the variable to use as left-hand side |
rhoc | the mass. |
metric | the metric term (or NULL). |
maxlevel | the maximum level to consider (or -1). |
v | a #GfsVariable of which lhs is an homogeneous version. |
Extracts the diffusion problem associated with lhs and rhs.
GfsLinearProblem* gfs_get_poisson_problem | ( | GfsDomain * | domain, |
GfsVariable * | rhs, | ||
GfsVariable * | lhs, | ||
GfsVariable * | dia, | ||
gint | maxlevel, | ||
GfsVariable * | v | ||
) |
domain | the domain over which the poisson problem is defined |
rhs | the variable to use as right-hand side |
lhs | the variable to use as left-hand side |
dia | the diagonal weight |
maxlevel | the maximum level to consider (or -1). |
v | a #GfsVariable of which lhs is an homogeneous version. |
Extracts the poisson problem associated with lhs and rhs.
void gfs_linear_problem_add_stencil | ( | GfsLinearProblem * | lp, |
GfsStencil * | stencil | ||
) |
lp | a #GfsLinearProblem. |
stencil | a #GfsStencil. |
Adds a stencil to the linear problem.
void gfs_linear_problem_destroy | ( | GfsLinearProblem * | lp | ) |
lp | a #GfsLinearProblem. |
Destroys a #GfsLinearProblem.
GfsLinearProblem* gfs_linear_problem_new | ( | GfsDomain * | domain | ) |
domain | a #GfsDomain. |
void gfs_multilevel_params_stats_write | ( | GfsMultilevelParams * | par, |
FILE * | fp | ||
) |
par | the multilevel parameters. |
fp | a file pointer. |
Writes in fp the statistics contained in p.
void gfs_multilevel_params_write | ( | GfsMultilevelParams * | par, |
FILE * | fp | ||
) |
par | the multilevel parameters. |
fp | a file pointer. |
Writes in fp a text representation of the multilevel parameters par.
void gfs_poisson_coefficients | ( | GfsDomain * | domain, |
GfsFunction * | alpha, | ||
gboolean | positive, | ||
gboolean | centered, | ||
gboolean | reset | ||
) |
domain | a #GfsDomain. |
alpha | the inverse of density or NULL. |
positive | if TRUE, alpha must be strictly positive. |
centered | TRUE if solving for a centered variable. |
reset | TRUE if resetting previous coefficients. |
Initializes the face coefficients for the Poisson equation .
If alpha is NULL, it is taken to be unity.
gdouble gfs_poisson_compatibility | ( | GfsDomain * | domain, |
GfsVariable * | lhs, | ||
GfsVariable * | rhs, | ||
gdouble | dt | ||
) |
domain | the domain over which the poisson problem is solved. |
lhs | the variable to use as left-hand side. |
rhs | the variable to use as right-hand side. |
dt | the timestep. |
void gfs_poisson_cycle | ( | GfsDomain * | domain, |
GfsMultilevelParams * | p, | ||
GfsVariable * | u, | ||
GfsVariable * | rhs, | ||
GfsVariable * | dia, | ||
GfsVariable * | res | ||
) |
domain | the domain on which to solve the Poisson equation. |
p | the #GfsMultilevelParams. |
u | the variable to use as left-hand side. |
rhs | the variable to use as right-hand side. |
dia | the diagonal weight. |
res | the residual. |
Apply one multigrid iteration to the Poisson equation defined by u and rhs.
The initial value of res on the leaves of root must be set to the residual of the Poisson equation (using gfs_residual()).
The face coefficients must be set using gfs_poisson_coefficients().
The values of u on the leaf cells are updated as well as the values of res (i.e. the cell tree is ready for another iteration).
void gfs_poisson_solve | ( | GfsDomain * | domain, |
GfsMultilevelParams * | par, | ||
GfsVariable * | lhs, | ||
GfsVariable * | rhs, | ||
GfsVariable * | res, | ||
GfsVariable * | dia, | ||
gdouble | dt | ||
) |
domain | the domain over which the poisson problem is solved. |
par | the parameters of the poisson problem. |
lhs | the variable to use as left-hand side. |
rhs | the variable to use as right-hand side. |
res | the variable in which to store the residual |
dia | the diagonal weight. |
dt | the length of the time-step. |
Solves the poisson problem over domain using Gerris' native multigrid poisson solver.
void gfs_relax | ( | GfsDomain * | domain, |
guint | d, | ||
gint | max_depth, | ||
gdouble | omega, | ||
GfsVariable * | u, | ||
GfsVariable * | rhs, | ||
GfsVariable * | dia | ||
) |
domain | the domain to relax. |
d | number of dimensions (2 or 3). |
max_depth | the maximum depth of the domain to relax. |
omega | the over-relaxation parameter. |
u | the variable to use as left-hand side. |
rhs | the variable to use as right-hand side. |
dia | the diagonal weight. |
Apply one pass of a Jacobi relaxation to all the leaf cells of domain with a level inferior or equal to max_depth and to all the cells at level max_depth. The relaxation should converge (if the right-hand-side rhs verifies the solvability conditions) toward the solution of a Poisson equation for u at the maximum depth.
void gfs_residual | ( | GfsDomain * | domain, |
guint | d, | ||
FttTraverseFlags | flags, | ||
gint | max_depth, | ||
GfsVariable * | u, | ||
GfsVariable * | rhs, | ||
GfsVariable * | dia, | ||
GfsVariable * | res | ||
) |
domain | a domain. |
d | number of dimensions (2 or 3). |
flags | which types of cells are to be visited. |
max_depth | maximum depth of the traversal. |
u | the variable to use as left-hand side. |
rhs | the variable to use as right-hand side. |
dia | the diagonal weight. |
res | the variable to use to store the residual. |
For each cell of domain, computes the sum of the residual over the volume of the cell for a Poisson equation with u as left-hand-side and rhs as right-hand-side. Stores the result in res.
void gfs_source_tension_coefficients | ( | GfsSourceTension * | s, |
GfsDomain * | domain, | ||
GfsFunction * | alpha | ||
) |
s | a #GfsSourceTension. |
domain | a #GfsDomain. |
alpha | the inverse of density or NULL. |
Initializes the face coefficients with the surface tension term (interface curvature times surface tension coefficient).
If alpha is NULL, it is taken to be unity.