#include "cs_boundary.h"
#include "cs_cdo_turbulence.h"
#include "cs_equation_param.h"
#include "cs_math.h"
#include "cs_param_sles.h"
#include "cs_physical_constants.h"
Go to the source code of this file.
Data Structures | |
struct | cs_navsto_param_sles_t |
Structure storing the parameters for solving the Navier-Stokes system. More... | |
struct | cs_navsto_param_t |
Structure storing the parameters related to the resolution of the Navier-Stokes system. More... | |
Macros | |
#define | CS_NAVSTO_STREAM_EQNAME "streamfunction_eq" |
Typedefs | |
typedef cs_flag_t | cs_navsto_param_model_flag_t |
typedef cs_flag_t | cs_navsto_param_post_flag_t |
Functions | |
cs_navsto_param_t * | cs_navsto_param_create (const cs_boundary_t *boundaries, cs_navsto_param_model_t model, cs_navsto_param_model_flag_t model_flag, cs_navsto_param_coupling_t algo_coupling, cs_navsto_param_post_flag_t post_flag) |
Create a new structure to store all numerical parameters related to the resolution of the Navier-Stokes (NS) system. More... | |
cs_navsto_param_t * | cs_navsto_param_free (cs_navsto_param_t *param) |
Free a cs_navsto_param_t structure. More... | |
void | cs_navsto_param_set (cs_navsto_param_t *nsp, cs_navsto_key_t key, const char *keyval) |
Set a parameter attached to a keyname in a cs_navsto_param_t structure. More... | |
void | cs_navsto_param_transfer (const cs_navsto_param_t *nsp, cs_equation_param_t *eqp) |
Apply the numerical settings defined for the Navier-Stokes system to an equation related to this system. More... | |
void | cs_navsto_param_log (const cs_navsto_param_t *nsp) |
Summary of the main cs_navsto_param_t structure. More... | |
cs_equation_param_t * | cs_navsto_param_get_velocity_param (const cs_navsto_param_t *nsp) |
Retrieve the cs_equation_param_t structure related to the velocity equation (momentum equation in most of the cases) More... | |
const char * | cs_navsto_param_get_model_name (cs_navsto_param_model_t model) |
Retrieve the name of the model system of equations. More... | |
const char * | cs_navsto_param_get_coupling_name (cs_navsto_param_coupling_t coupling) |
Retrieve the name of the coupling algorithm. More... | |
void | cs_navsto_set_reference_pressure (cs_navsto_param_t *nsp, cs_real_t pref) |
Set the value to consider for the reference pressure. More... | |
cs_xdef_t * | cs_navsto_add_velocity_ic_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *val) |
Define the initial condition for the velocity unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set to a constant value. More... | |
cs_xdef_t * | cs_navsto_add_velocity_ic_by_analytic (cs_navsto_param_t *nsp, const char *z_name, cs_analytic_func_t *analytic, void *input) |
Define the initial condition for the velocity unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function. More... | |
cs_xdef_t * | cs_navsto_add_pressure_ic_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *val) |
Define the initial condition for the pressure unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set to a constant value. More... | |
cs_xdef_t * | cs_navsto_add_pressure_ic_by_analytic (cs_navsto_param_t *nsp, const char *z_name, cs_analytic_func_t *analytic, void *input) |
Define the initial condition for the pressure unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function. More... | |
void | cs_navsto_set_fixed_walls (cs_navsto_param_t *nsp) |
Add the definition of boundary conditions related to a fixed wall into the set of parameters for the management of the Navier-Stokes system of equations. More... | |
void | cs_navsto_set_symmetries (cs_navsto_param_t *nsp) |
Add the definition of boundary conditions related to a symmetry into the set of parameters for the management of the Navier-Stokes system of equations. More... | |
void | cs_navsto_set_outlets (cs_navsto_param_t *nsp) |
Add the definition of boundary conditions related to outlets into the set of parameters for the management of the Navier-Stokes system of equations. More... | |
cs_xdef_t * | cs_navsto_set_pressure_bc_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *values) |
Set the pressure field on a boundary using a uniform value. More... | |
cs_xdef_t * | cs_navsto_set_velocity_wall_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *values) |
Define the velocity field for a sliding wall boundary using a uniform value. More... | |
cs_xdef_t * | cs_navsto_set_velocity_inlet_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *values) |
Define the velocity field for an inlet boundary using a uniform value. More... | |
cs_xdef_t * | cs_navsto_set_velocity_inlet_by_analytic (cs_navsto_param_t *nsp, const char *z_name, cs_analytic_func_t *ana, void *input) |
Define the velocity field for an inlet boundary using an analytical function. More... | |
cs_xdef_t * | cs_navsto_set_velocity_inlet_by_array (cs_navsto_param_t *nsp, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index) |
Define the velocity field for an inlet boundary using an array of values. More... | |
cs_xdef_t * | cs_navsto_set_velocity_inlet_by_dof_func (cs_navsto_param_t *nsp, const char *z_name, cs_dof_func_t *func, void *func_input) |
Define the velocity field for an inlet boundary using a DoF function. More... | |
cs_xdef_t * | cs_navsto_add_source_term_by_analytic (cs_navsto_param_t *nsp, const char *z_name, cs_analytic_func_t *ana, void *input) |
Define a new source term structure defined by an analytical function. More... | |
cs_xdef_t * | cs_navsto_add_source_term_by_val (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *val) |
Define a new source term structure defined by a constant value. More... | |
cs_xdef_t * | cs_navsto_add_source_term_by_array (cs_navsto_param_t *nsp, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index) |
Define a new source term structure defined by an array. More... | |
void | cs_navsto_add_oseen_field (cs_navsto_param_t *nsp, cs_adv_field_t *adv_fld) |
Add a advection field for the Oseen problem. More... | |
#define CS_NAVSTO_STREAM_EQNAME "streamfunction_eq" |
typedef cs_flag_t cs_navsto_param_post_flag_t |
enum cs_navsto_key_t |
List of available keys for setting the parameters of the Navier-Stokes system.
Enumerator | |
---|---|
CS_NSKEY_ADVECTION_EXTRAPOL | Set the extrapolation to use for the estimation of the advection field (cf. cs_param_advection_extrapol_t)) |
CS_NSKEY_ADVECTION_FORMULATION | Set the type of formulation for the advection term, for example in the Oseen problem. (cf. cs_param_advection_form_t) |
CS_NSKEY_ADVECTION_SCHEME | Set the type of scheme for the advection term, for example in the Oseen problem. (cf. cs_param_advection_scheme_t) |
CS_NSKEY_ADVECTION_STRATEGY | Set the strategy to handle the advection term (cf. cs_param_advection_strategy_t) |
CS_NSKEY_DOF_REDUCTION | Set how the DoFs are defined (similar to CS_EQKEY_DOF_REDUCTION) Enable to set this type of DoFs definition for all related equations |
CS_NSKEY_GD_SCALE_COEF | Set the scaling of the grad-div term when an artificial compressibility algorithm or an Uzawa-Augmented Lagrangian method is used |
CS_NSKEY_IL_ALGO_ATOL | Absolute tolerance at which the Oseen or Stokes system is resolved. These systems corresponds to an inner linear system to solve when considering the Navier-Stokes system since one has to handle the non-linearity in addition as an outer process. |
CS_NSKEY_IL_ALGO_DTOL | Divergence tolerance at which the Oseen or Stokes system is resolved. These systems corresponds to an inner linear system to solve when considering the Navier-Stokes system since one has to handle the non-linearity in addition as an outer process. |
CS_NSKEY_IL_ALGO_RTOL | Relative tolerance at which the Oseen or Stokes system is resolved. These systems corresponds to an inner linear system to solve when considering the Navier-Stokes system since one has to handle the non-linearity in addition as an outer process. |
CS_NSKEY_IL_ALGO_VERBOSITY | Level of verbosity related to the inner linear algorithm (cf. CS_NSKEY_SLES_STRATEGY) |
CS_NSKEY_MAX_IL_ALGO_ITER | Set the maximal number of iteration for solving the inner linear system. |
CS_NSKEY_MAX_NL_ALGO_ITER | Set the maximal number of Picard iterations for solving the non-linearity arising from the advection form |
CS_NSKEY_MAX_OUTER_ITER | Set the maximal number of outer iterations for solving the full system including the turbulence modelling or the thermal system for instance |
CS_NSKEY_NL_ALGO | Type of algorithm to consider to solve the non-linearity arising from the Navier-Stokes system |
CS_NSKEY_NL_ALGO_ATOL | Absolute tolerance at which the non-linearity arising from the advection term is resolved using a Picard (fixed-point) algorithm |
CS_NSKEY_NL_ALGO_DTOL | Divergence tolerance at which the non-linearity arising from the advection term is resolved using a Picard (fixed-point) algorithm |
CS_NSKEY_NL_ALGO_RTOL | Relative tolerance at which the non-linearity arising from the advection term is resolved using a Picard (fixed-point) algorithm |
CS_NSKEY_NL_ALGO_VERBOSITY | Level of verbosity related to the non-linear algorithm (Picar for instance) |
CS_NSKEY_QUADRATURE | Set the type to use in all routines involving quadrature (similar to CS_EQKEY_BC_QUADRATURE) |
CS_NSKEY_SCHUR_STRATEGY | Set the way to define the Schur complement approximation (cf. cs_navsto_schur_approx_t) |
CS_NSKEY_SLES_STRATEGY | Strategy for solving the SLES arising from the discretization of the Navier-Stokes system |
CS_NSKEY_SPACE_SCHEME | Numerical scheme for the space discretization. Available choices are:
|
CS_NSKEY_THERMAL_TOLERANCE | Value of the tolerance criterion under which one stops iterating on the Navier-Stokes and thermal systems |
CS_NSKEY_TIME_SCHEME | Numerical scheme for the time discretization |
CS_NSKEY_TIME_THETA | Set the value of theta. Only useful if CS_NSKEY_TIME_SCHEME is set to "theta_scheme"
|
CS_NSKEY_VERBOSITY | Set the level of details for the specific part related to the Navier-Stokes system |
CS_NSKEY_N_KEYS |
enum cs_navsto_nl_algo_t |
Choice of algorithm for solving the system.
Bit values for additional physical modelling related to the Navier-Stokes system of equations.
Enumerator | |
---|---|
CS_NAVSTO_MODEL_STEADY | There is no time dependency in the model |
CS_NAVSTO_MODEL_GRAVITY_EFFECTS | Take into account the gravity effects (add a constant source term equal to rho*vect(g)) |
CS_NAVSTO_MODEL_CORIOLIS_EFFECTS | Take into account the Coriolis effects (add a source term) |
CS_NAVSTO_MODEL_PASSIVE_THERMAL_TRACER | An additional equation is created involving the thermal equation. The advection field is automatically set as the mass flux. By default, the temperature is solved in Celsius and homogeneous Neumann boundary conditions (no flux) are set. |
CS_NAVSTO_MODEL_BOUSSINESQ | An additional equation is created involving the thermal equation. The advection field is automatically set as the mass flux. By default, the temperature is solved in Celsius and homogeneous Neumann boundary conditions (no flux) are set. The variation of mass density around a reference value is related to the variation of temperature w.r.t. a reference temperature and it plays a role in the gravity effects (rho.vect(g)) The gradient of temperature is assumed to have a small norm and the mass density variates in a small range. In this case, an additional momentum source term is added. |
CS_NAVSTO_MODEL_SOLIDIFICATION_BOUSSINESQ | This option is similar to CS_NAVSTO_MODEL_BOUSSINESQ. The difference is that the variation of mass density relies not only on the temperature but also the alloy concentration(s). The gradient of temperature/alloy concentrations is assumed to have a small norm and the mass density variates in a small range. In this case, additional equations related to the temperature/alloy concetrations are considered. A momentum source term is added wich is managed by the solidification module. |
Describe the system of equations related to the Navier-Stokes to be solved.
Enumerator | |
---|---|
CS_NAVSTO_MODEL_STOKES | Stokes equations (mass and momentum) with the classical choice of variables i.e. velocity and pressure. Mass density is assumed to be constant. |
CS_NAVSTO_MODEL_OSEEN | Like the incompressible Navier-Stokes equations (mass and momentum) but with a velocity field which is given. Thus the advection term in the momentum equation is linear. Unknowns: velocity and pressure. Mass density is assumed to be constant. The advection field is set using cs_navsto_add_oseen_field |
CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES | Navier-Stokes equations: mass and momentum with a constant mass density Mass equation is equivalent to an incompressibility constraint |
CS_NAVSTO_N_MODELS |
Bit values for additional generic postprocessing related to the Navier-Stokes module. In what follows, w denotes the vorticity vector, u the velocity vector and k the kinetic energy defined by 1/2 * u u.
Strategy to build the Schur complement approximation. This appears in block preconditioning or uzawa algorithms when a monolithic (fully coupled approach) is used. | A B^t| | B 0 |.
enum cs_navsto_sles_t |
High-level information about the way of settings the SLES for solving the Navier-Stokes system. When the system is treated as a saddle-point problem (monolithic approach in what follows), then one uses these notations: A_{00} is the upper-left block and A_{11} (should be 0 but the preconditioner may have entries for the approximation of the inverse of the Schur complement).
Enumerator | |
---|---|
CS_NAVSTO_SLES_ADDITIVE_GMRES_BY_BLOCK | Associated keyword: "additive_gmres" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm) The Navier-Stokes system of equations is solved an additive preconditioner (block diagonal matrix where the block 00 is A_{00}) and the block 11 is set to the identity. Preconditioner/solver for the block 00 is set using the momentum equation. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_BLOCK_MULTIGRID_CG | Associated keyword: "block_amg_cg" The Navier-Stokes system of equations is solved using a multigrid on each diagonal block as a preconditioner and applying a conjugate gradient as solver. Use this strategy when the saddle-point problem has been reformulated into a "classical" linear system. For instance when a Uzawa or an Artificial Compressibility coupling algorithm is used. (i.e. with the parameter CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY is set as coupling algorithm). This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_BY_BLOCKS | Associated keyword: "blocks" The Navier-Stokes system is split into a 3x3 block matrix for the velocity unknows and in a non-assembly way for the divergence/pressure gradient operators. |
CS_NAVSTO_SLES_DIAG_SCHUR_GMRES | Associated keyword: "diag_schur_gmres" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The Navier-Stokes system of equations is solved using a block diagonal preconditioner where the block 00 is A_{00} preconditioned with one multigrid iteration and the block 11 is an approximation of the Schur complement preconditioned with one multigrid iteration. The main iterative solver is a flexible GMRES. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_EQ_WITHOUT_BLOCK | Associated keyword: "no_block" Use the same mechanism as for a stand-alone equation. In this case, the setting relies on the function cs_equation_set_sles and the different options for solving a linear system such as the choice of the iterative solver or the choice of the preconditioner or the type of residual normalization |
CS_NAVSTO_SLES_GKB_PETSC | Associated keyword: "gkb_petsc" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The Navier-Stokes system of equations is solved using a Golub-Kahan bi-diagonalization. One assumes that the saddle-point system is symmetric. By default, the block A_{00} may be augmented (this is not the default choice) and is solved with a conjugate gradient algorithm preconditioned with a multigrid. The residual is computed in the energy norm. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_GKB_GMRES | Associated keyword: "gkb_gmres" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The Navier-Stokes system of equations is solved using a Golub-Kahan bi-diagonalization (GKB) as preconditioner of a flexible GMRES solver. The GKB algorithm is solved with a reduced tolerance as well as the CG+Multigrid used as an inner solver in the GKB algorithm. One assumes that the saddle-point system is symmetric. The residual for the GKB part is computed in the energy norm. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_GKB_SATURNE | Associated keyword: "gkb" or "gkb_saturne" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The Navier-Stokes system of equations is solved using a Golub-Kahan bi-diagonalization. By default, the block A_{00} may be augmented (this is not the default choice) and is solved with the SLES settings given to the momentum equation A conjugate gradient algorithm preconditioned with a multigrid for Stokes for instance. The residual is computed in the energy norm. |
CS_NAVSTO_SLES_MULTIPLICATIVE_GMRES_BY_BLOCK | Associated keyword: "multiplicative_gmres" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm) The Navier-Stokes system of equations is solved a multiplicative preconditioner (block diagonal matrix where the block 00 is A_{00}) and the block 11 is set to the identity. Block 01 is also considered in the block preconditioner. Preconditioner/solver for the block 00 is set using the momentum equation. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_MUMPS | Associated keyword: "mumps" Direct solver to solve systems arising from the discretization of the Navier-Stokes equations |
CS_NAVSTO_SLES_UPPER_SCHUR_GMRES | Associated keyword: "upper_schur_gmres" Available choice when a monolithic approach is used (i.e. with the parameter CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The Navier-Stokes system of equations is solved using a upper triangular block preconditioner where the block 00 is A_{00} preconditioned with one multigrid iteration and the block 11 is an approximation of the Schur complement preconditioned with a minres. The main iterative solver is a flexible GMRES. This option is only available with the support to the PETSc library up to now. |
CS_NAVSTO_SLES_UZAWA_AL | Associated keyword: "uzawa_al" Resolution using an uzawa algorithm with an Augmented Lagrangian approach |
CS_NAVSTO_SLES_UZAWA_CG | Associated keyword: "uzawa_cg" Resolution using an uzawa algorithm optimized using a conjugate gradient reformulation. Two systems are solved at each iteration (one related to the velocity block, one related to the Schur complement approximation - size of the pressure space). |
CS_NAVSTO_SLES_N_TYPES |
void cs_navsto_add_oseen_field | ( | cs_navsto_param_t * | nsp, |
cs_adv_field_t * | adv_fld | ||
) |
Add a advection field for the Oseen problem.
[in,out] | nsp | pointer to a cs_navsto_param_t |
[in,out] | adv_fld | pointer to a cs_adv_field_t |
cs_xdef_t* cs_navsto_add_pressure_ic_by_analytic | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_analytic_func_t * | analytic, | ||
void * | input | ||
) |
Define the initial condition for the pressure unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" if all cells are considered) |
[in] | analytic | pointer to an analytic function |
[in] | input | NULL or pointer to a structure cast on-the-fly |
cs_xdef_t* cs_navsto_add_pressure_ic_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | val | ||
) |
Define the initial condition for the pressure unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set to a constant value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" if all cells are considered) |
[in] | val | pointer to the value |
cs_xdef_t* cs_navsto_add_source_term_by_analytic | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_analytic_func_t * | ana, | ||
void * | input | ||
) |
Define a new source term structure defined by an analytical function.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all cells are considered) |
[in] | ana | pointer to an analytical function |
[in] | input | NULL or pointer to a structure cast on-the-fly |
cs_xdef_t* cs_navsto_add_source_term_by_array | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_flag_t | loc, | ||
cs_real_t * | array, | ||
bool | is_owner, | ||
cs_lnum_t * | index | ||
) |
Define a new source term structure defined by an array.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all cells are considered) |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t structure (true or false) |
[in] | index | optional pointer to the array index |
cs_xdef_t* cs_navsto_add_source_term_by_val | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | val | ||
) |
Define a new source term structure defined by a constant value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all cells are considered) |
[in] | val | pointer to the value to set |
cs_xdef_t* cs_navsto_add_velocity_ic_by_analytic | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_analytic_func_t * | analytic, | ||
void * | input | ||
) |
Define the initial condition for the velocity unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" if all cells are considered) |
[in] | analytic | pointer to an analytic function |
[in] | input | NULL or pointer to a structure cast on-the-fly |
cs_xdef_t* cs_navsto_add_velocity_ic_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | val | ||
) |
Define the initial condition for the velocity unknowns. This definition can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set to a constant value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" if all cells are considered) |
[in] | val | pointer to the value |
cs_navsto_param_t* cs_navsto_param_create | ( | const cs_boundary_t * | boundaries, |
cs_navsto_param_model_t | model, | ||
cs_navsto_param_model_flag_t | model_flag, | ||
cs_navsto_param_coupling_t | algo_coupling, | ||
cs_navsto_param_post_flag_t | post_flag | ||
) |
Create a new structure to store all numerical parameters related to the resolution of the Navier-Stokes (NS) system.
[in] | boundaries | pointer to a cs_boundary_t structure |
[in] | model | type of model related to the NS system |
[in] | model_flag | additional high-level model options |
[in] | algo_coupling | algorithm used for solving the NS system |
[in] | post_flag | predefined post-processings options |
cs_navsto_param_t* cs_navsto_param_free | ( | cs_navsto_param_t * | param | ) |
Free a cs_navsto_param_t structure.
[in,out] | param | pointer to a cs_navsto_param_t structure |
const char* cs_navsto_param_get_coupling_name | ( | cs_navsto_param_coupling_t | coupling | ) |
Retrieve the name of the coupling algorithm.
[in] | coupling | a cs_navsto_param_coupling_t |
const char* cs_navsto_param_get_model_name | ( | cs_navsto_param_model_t | model | ) |
Retrieve the name of the model system of equations.
[in] | model | a cs_navsto_param_model_t |
cs_equation_param_t* cs_navsto_param_get_velocity_param | ( | const cs_navsto_param_t * | nsp | ) |
Retrieve the cs_equation_param_t structure related to the velocity equation (momentum equation in most of the cases)
[in] | nsp | pointer to a cs_navsto_param_t structure |
void cs_navsto_param_log | ( | const cs_navsto_param_t * | nsp | ) |
Summary of the main cs_navsto_param_t structure.
[in] | nsp | pointer to a cs_navsto_param_t structure |
void cs_navsto_param_set | ( | cs_navsto_param_t * | nsp, |
cs_navsto_key_t | key, | ||
const char * | keyval | ||
) |
Set a parameter attached to a keyname in a cs_navsto_param_t structure.
[in,out] | nsp | pointer to a cs_navsto_param_t structure to set |
[in] | key | key related to the member of eq to set |
[in] | keyval | accessor to the value to set |
void cs_navsto_param_transfer | ( | const cs_navsto_param_t * | nsp, |
cs_equation_param_t * | eqp | ||
) |
Apply the numerical settings defined for the Navier-Stokes system to an equation related to this system.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
void cs_navsto_set_fixed_walls | ( | cs_navsto_param_t * | nsp | ) |
Add the definition of boundary conditions related to a fixed wall into the set of parameters for the management of the Navier-Stokes system of equations.
[in] | nsp | pointer to a cs_navsto_param_t structure |
void cs_navsto_set_outlets | ( | cs_navsto_param_t * | nsp | ) |
Add the definition of boundary conditions related to outlets into the set of parameters for the management of the Navier-Stokes system of equations.
[in] | nsp | pointer to a cs_navsto_param_t structure |
cs_xdef_t* cs_navsto_set_pressure_bc_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | values | ||
) |
Set the pressure field on a boundary using a uniform value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | value | value to set |
void cs_navsto_set_reference_pressure | ( | cs_navsto_param_t * | nsp, |
cs_real_t | pref | ||
) |
Set the value to consider for the reference pressure.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | pref | value of the reference pressure |
void cs_navsto_set_symmetries | ( | cs_navsto_param_t * | nsp | ) |
Add the definition of boundary conditions related to a symmetry into the set of parameters for the management of the Navier-Stokes system of equations.
[in] | nsp | pointer to a cs_navsto_param_t structure |
cs_xdef_t* cs_navsto_set_velocity_inlet_by_analytic | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_analytic_func_t * | ana, | ||
void * | input | ||
) |
Define the velocity field for an inlet boundary using an analytical function.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | ana | pointer to an analytical function |
[in] | input | NULL or pointer to a structure cast on-the-fly |
cs_xdef_t* cs_navsto_set_velocity_inlet_by_array | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_flag_t | loc, | ||
cs_real_t * | array, | ||
bool | is_owner, | ||
cs_lnum_t * | index | ||
) |
Define the velocity field for an inlet boundary using an array of values.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t structure (true or false) |
[in] | index | optional pointer to the array index |
cs_xdef_t* cs_navsto_set_velocity_inlet_by_dof_func | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_dof_func_t * | func, | ||
void * | func_input | ||
) |
Define the velocity field for an inlet boundary using a DoF function.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | func | pointer to a cs_dof_function_t |
[in] | func_input | NULL or pointer to a structure cast on-the-fly |
cs_xdef_t* cs_navsto_set_velocity_inlet_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | values | ||
) |
Define the velocity field for an inlet boundary using a uniform value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | values | array of three real values |
cs_xdef_t* cs_navsto_set_velocity_wall_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | values | ||
) |
Define the velocity field for a sliding wall boundary using a uniform value.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | z_name | name of the associated zone (if NULL or "" all boundary faces are considered) |
[in] | values | array of three real values |