My Project
programmer's documentation
Data Structures | Enumerations | Functions
cs_navsto_param.h File Reference
#include "cs_boundary.h"
#include "cs_equation_param.h"
Include dependency graph for cs_navsto_param.h:

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...
 

Enumerations

enum  cs_navsto_param_model_t {
  CS_NAVSTO_MODEL_STOKES, CS_NAVSTO_MODEL_OSEEN, CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES, CS_NAVSTO_MODEL_BOUSSINESQ_NAVIER_STOKES,
  CS_NAVSTO_N_MODELS
}
 Modelling related to the Navier-Stokes system of equations. More...
 
enum  cs_navsto_param_sles_t {
  CS_NAVSTO_SLES_EQ_WITHOUT_BLOCK, CS_NAVSTO_SLES_BLOCK_MULTIGRID_CG, CS_NAVSTO_SLES_ADDITIVE_GMRES_BY_BLOCK, CS_NAVSTO_SLES_DIAG_SCHUR_GMRES,
  CS_NAVSTO_SLES_UPPER_SCHUR_GMRES, CS_NAVSTO_SLES_GKB_GMRES, CS_NAVSTO_SLES_GKB, CS_NAVSTO_SLES_N_TYPES
}
 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). More...
 
enum  cs_navsto_param_time_state_t { CS_NAVSTO_TIME_STATE_FULL_STEADY, CS_NAVSTO_TIME_STATE_LIMIT_STEADY, CS_NAVSTO_TIME_STATE_UNSTEADY, CS_NAVSTO_N_TIME_STATES }
 Status of the time for the Navier-Stokes system of equations. More...
 
enum  cs_navsto_param_coupling_t {
  CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY, CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY_VPP, CS_NAVSTO_COUPLING_MONOLITHIC, CS_NAVSTO_COUPLING_PROJECTION,
  CS_NAVSTO_COUPLING_UZAWA, CS_NAVSTO_N_COUPLINGS
}
 Choice of algorithm for solving the system. More...
 
enum  cs_navsto_key_t {
  CS_NSKEY_ADVECTION_FORMULATION, CS_NSKEY_ADVECTION_SCHEME, CS_NSKEY_DOF_REDUCTION, CS_NSKEY_GD_SCALE_COEF,
  CS_NSKEY_MAX_ALGO_ITER, CS_NSKEY_QUADRATURE, CS_NSKEY_RESIDUAL_TOLERANCE, CS_NSKEY_SLES_STRATEGY,
  CS_NSKEY_SPACE_SCHEME, CS_NSKEY_TIME_SCHEME, CS_NSKEY_TIME_THETA, CS_NSKEY_VERBOSITY,
  CS_NSKEY_N_KEYS
}
 List of available keys for setting the parameters 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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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...
 

Enumeration Type Documentation

◆ 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:

  • "cdo_fb" for CDO face-based scheme
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"

  • Example: "0.75" (keyval must be between 0 and 1)
CS_NSKEY_VERBOSITY 

Set the level of details for the specific part related to the Navier-Stokes system

CS_NSKEY_N_KEYS 

◆ cs_navsto_param_coupling_t

Choice of algorithm for solving the system.

Enumerator
CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY 

The system is solved using an artificial compressibility algorithm. One vectorial equation is solved followed by a pressure update.

CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY_VPP 

The system is solved using an artificial compressibility algorithm with a Vector Penalty Projection splitting. Two vectorial equations are solved: a momentum-like one and another one involving a grad-div operator.

CS_NAVSTO_COUPLING_MONOLITHIC 

The system is treated as a "monolithic" matrix

CS_NAVSTO_COUPLING_PROJECTION 

The system is solved using an incremental projection algorithm

CS_NAVSTO_COUPLING_UZAWA 

The system is solved without decoupling the equations using a Uzawa algorithm and an Augmented Lagrangian approach inside each sub-iteration.

CS_NAVSTO_N_COUPLINGS 

◆ cs_navsto_param_model_t

Modelling related to the Navier-Stokes system of equations.

Enumerator
CS_NAVSTO_MODEL_STOKES 

Stokes equations (mass and momentum) with the classical choice of variables i.e. velocity and pressure

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

CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES 

Navier-Stokes equations: mass and momentum with a constant mass density

CS_NAVSTO_MODEL_BOUSSINESQ_NAVIER_STOKES 

Navier-Stokes equations: mass and momentum with a constant mass density 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 equation related to the energy is considered.

CS_NAVSTO_N_MODELS 

◆ cs_navsto_param_sles_t

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 

◆ cs_navsto_param_time_state_t

Status of the time for the Navier-Stokes system of equations.

Enumerator
CS_NAVSTO_TIME_STATE_FULL_STEADY 

The Navier-Stokes system of equations is solved without taking into account the time effect

CS_NAVSTO_TIME_STATE_LIMIT_STEADY 

The Navier-Stokes system of equations is solved as a limit of a unsteady process

CS_NAVSTO_TIME_STATE_UNSTEADY 

The Navier-Stokes system of equations is time-dependent

CS_NAVSTO_N_TIME_STATES 

Function Documentation

◆ cs_navsto_add_oseen_field()

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.

Parameters
[in,out]nsppointer to a cs_navsto_param_t
[in,out]adv_fldpointer to a cs_adv_field_t

◆ cs_navsto_add_pressure_ic_by_analytic()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_pressure_ic_by_value()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_source_term_by_analytic()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]anapointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_source_term_by_array()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]locinformation to know where are located values
[in]arraypointer to an array
[in]is_ownertransfer the lifecycle to the cs_xdef_t structure (true or false)
[in]indexoptional pointer to the array index
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_source_term_by_val()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]valpointer to the value to set
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_velocity_ic_by_analytic()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_add_velocity_ic_by_value()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_param_create()

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.

Parameters
[in]boundariespointer to a cs_boundary_t structure
[in]modelmodel related to the NS system to solve
[in]time_statestate of the time for the NS equations
[in]algo_couplingalgorithm used for solving the NS system
Returns
a pointer to a new allocated structure

◆ cs_navsto_param_free()

cs_navsto_param_t* cs_navsto_param_free ( cs_navsto_param_t param)

Free a cs_navsto_param_t structure.

Parameters
[in,out]parampointer to a cs_navsto_param_t structure
Returns
a NULL pointer

◆ cs_navsto_param_get_coupling_name()

const char* cs_navsto_param_get_coupling_name ( cs_navsto_param_coupling_t  coupling)

Retrieve the name of the coupling algorithm.

Parameters
[in]couplinga cs_navsto_param_coupling_t
Returns
the name

◆ cs_navsto_param_is_steady()

static bool cs_navsto_param_is_steady ( cs_navsto_param_t nsp)
inlinestatic

Ask cs_navsto_param_t structure if the settings correspond to a steady computation.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
Returns
true or false

◆ cs_navsto_param_log()

void cs_navsto_param_log ( const cs_navsto_param_t nsp)

Summary of the main cs_navsto_param_t structure.

Parameters
[in]nsppointer to a cs_navsto_param_t structure

◆ cs_navsto_param_set()

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.

Parameters
[in,out]nsppointer to a cs_navsto_param_t structure to set
[in]keykey related to the member of eq to set
[in]keyvalaccessor to the value to set

◆ cs_navsto_param_transfer()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in,out]eqppointer to a cs_equation_param_t structure

◆ cs_navsto_set_fixed_walls()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure

◆ cs_navsto_set_outlets()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure

◆ cs_navsto_set_pressure_bc_by_value()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all boundary faces are considered)
[in]valuevalue to set

◆ cs_navsto_set_symmetries()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure

◆ cs_navsto_set_velocity_inlet_by_analytic()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all boundary faces are considered)
[in]anapointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly

◆ cs_navsto_set_velocity_inlet_by_value()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all boundary faces are considered)
[in]valuesarray of three real values

◆ cs_navsto_set_velocity_wall_by_value()

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.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]z_namename of the associated zone (if NULL or "" all boundary faces are considered)
[in]valuesarray of three real values