My Project
programmer's documentation
|
Go to the source code of this file.
Data Structures | |
struct | cs_navsto_param_t |
Structure storing the parameters related to the resolution of the Navier-Stokes system. More... | |
Functions | |
static bool | cs_navsto_param_is_steady (cs_navsto_param_t *nsp) |
Ask cs_navsto_param_t structure if the settings correspond to a steady computation. More... | |
cs_navsto_param_t * | cs_navsto_param_create (const cs_boundary_t *boundaries, cs_navsto_param_model_t model, cs_navsto_param_time_state_t time_state, cs_navsto_param_coupling_t algo_coupling) |
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... | |
const char * | cs_navsto_param_get_coupling_name (cs_navsto_param_coupling_t coupling) |
Retrieve the name of the coupling algorithm. 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 unkowns. 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 unkowns. 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... | |
void | cs_navsto_set_pressure_bc_by_value (cs_navsto_param_t *nsp, const char *z_name, cs_real_t *values) |
Define the pressure field on a boundary using a uniform value. More... | |
void | 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... | |
void | 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... | |
void | 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_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... | |
enum cs_navsto_key_t |
List of available keys for setting the parameters of the Navier-Stokes system.
Enumerator | |
---|---|
CS_NSKEY_ADVECTION_FORMULATION | Set the type of formulation for the advection term, for example in the Oseen problem . cf. cs_param.h |
CS_NSKEY_ADVECTION_SCHEME | Set the type of scheme for the advection term, for example in the Oseen problem . cf. cs_param.h |
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_MAX_ALGO_ITER | Set the maximal number of iteration of the coupling algorithm. Not useful for a monolithic approach. In this case, only the maximal number of iterations for the iterative solver is taken into account |
CS_NSKEY_QUADRATURE | Set the type to use in all routines involving quadrature (similar to CS_EQKEY_BC_QUADRATURE) |
CS_NSKEY_RESIDUAL_TOLERANCE | Tolerance at which the Navier–Stokes is resolved (apply to the residual of the coupling algorithm chosen to solve the Navier–Stokes system) |
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_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 |
Choice of algorithm for solving the system.
Modelling related to the Navier-Stokes system of equations.
High-level information about the way of settings the SLES for solving the Navier-Stokes system. When a 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 preconditionner may have entries for the approximation of the inverse of the Schur complement).
Enumerator | |
---|---|
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 preconditionner or the type of residual normalization |
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 or CS_NAVSTO_COUPLING_UZAWA is set as coupling algorithm). This option is only available with the support to the PETSc library up to now. |
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} preconditionned by one multigrid iteration and the block 11 is set to the identity. This option is only available with the support to the PETSc library up to now. |
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 preconditionned 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_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 preconditionned 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_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 preconditionner 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 | Associated keyword: "gkb" 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 conjuguate gradient algorithm preconditionned 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_N_TYPES |
Status of the time for the Navier-Stokes system of equations.
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 unkowns. 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 unkowns. 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_time_state_t | time_state, | ||
cs_navsto_param_coupling_t | algo_coupling | ||
) |
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 | model related to the NS system to solve |
[in] | time_state | state of the time for the NS equations |
[in] | algo_coupling | algorithm used for solving the NS system |
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 |
|
inlinestatic |
Ask cs_navsto_param_t structure if the settings correspond to a steady computation.
[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 |
void cs_navsto_set_pressure_bc_by_value | ( | cs_navsto_param_t * | nsp, |
const char * | z_name, | ||
cs_real_t * | values | ||
) |
Define 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_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 |
void 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 |
void 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 |
void 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 |