7.1
general documentation
cs_navsto_param.h File Reference
#include "cs_boundary.h"
#include "cs_cdo_turbulence.h"
#include "cs_equation_param.h"
#include "cs_iter_algo.h"
#include "cs_math.h"
#include "cs_param_sles.h"
#include "cs_physical_constants.h"
+ Include dependency graph for cs_navsto_param.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_boussinesq_t
 Structure storing the parameters related to the Boussinesq source term in the momentum equation. 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
 

Enumerations

enum  cs_navsto_param_model_t { CS_NAVSTO_MODEL_STOKES, CS_NAVSTO_MODEL_OSEEN, CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES, CS_NAVSTO_N_MODELS }
 Describe the system of equations related to the Navier-Stokes to be solved. More...
 
enum  cs_navsto_param_model_bit_t {
  CS_NAVSTO_MODEL_STEADY = 1<<0, CS_NAVSTO_MODEL_GRAVITY_EFFECTS = 1<<1, CS_NAVSTO_MODEL_CORIOLIS_EFFECTS = 1<<2, CS_NAVSTO_MODEL_PASSIVE_THERMAL_TRACER = 1<<3,
  CS_NAVSTO_MODEL_BOUSSINESQ = 1<<4, CS_NAVSTO_MODEL_WITH_SOLIDIFICATION = 1<<5
}
 Bit values for additional physical modelling related to the Navier-Stokes system of equations. More...
 
enum  cs_navsto_param_post_bit_t {
  CS_NAVSTO_POST_VELOCITY_DIVERGENCE = 1<< 0, CS_NAVSTO_POST_KINETIC_ENERGY = 1<< 1, CS_NAVSTO_POST_VORTICITY = 1<< 2, CS_NAVSTO_POST_VELOCITY_GRADIENT = 1<< 3,
  CS_NAVSTO_POST_STREAM_FUNCTION = 1<< 4, CS_NAVSTO_POST_HELICITY = 1<< 5, CS_NAVSTO_POST_ENSTROPHY = 1<< 6, CS_NAVSTO_POST_MASS_DENSITY = 1<< 7,
  CS_NAVSTO_POST_CELL_MASS_FLUX_BALANCE = 1<< 8, CS_NAVSTO_POST_PRESSURE_GRADIENT = 1<< 9
}
 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. More...
 
enum  cs_navsto_sles_t {
  CS_NAVSTO_SLES_ADDITIVE_GMRES_BY_BLOCK, CS_NAVSTO_SLES_BLOCK_MULTIGRID_CG, CS_NAVSTO_SLES_BY_BLOCKS, CS_NAVSTO_SLES_DIAG_SCHUR_GCR,
  CS_NAVSTO_SLES_DIAG_SCHUR_GMRES, CS_NAVSTO_SLES_DIAG_SCHUR_MINRES, CS_NAVSTO_SLES_EQ_WITHOUT_BLOCK, CS_NAVSTO_SLES_GCR,
  CS_NAVSTO_SLES_GKB_PETSC, CS_NAVSTO_SLES_GKB_GMRES, CS_NAVSTO_SLES_GKB_SATURNE, CS_NAVSTO_SLES_LOWER_SCHUR_GCR,
  CS_NAVSTO_SLES_MINRES, CS_NAVSTO_SLES_MULTIPLICATIVE_GMRES_BY_BLOCK, CS_NAVSTO_SLES_MUMPS, CS_NAVSTO_SLES_NOTAY_TRANSFORM,
  CS_NAVSTO_SLES_SGS_SCHUR_GCR, CS_NAVSTO_SLES_UPPER_SCHUR_GCR, CS_NAVSTO_SLES_UPPER_SCHUR_GMRES, CS_NAVSTO_SLES_USER,
  CS_NAVSTO_SLES_UZAWA_SCHUR_GCR, CS_NAVSTO_SLES_UZAWA_AL, CS_NAVSTO_SLES_UZAWA_CG, CS_NAVSTO_SLES_N_TYPES
}
 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). More...
 
enum  cs_navsto_param_coupling_t { CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY, CS_NAVSTO_COUPLING_MONOLITHIC, CS_NAVSTO_COUPLING_PROJECTION, CS_NAVSTO_N_COUPLINGS }
 Choice of algorithm for solving the system. More...
 
enum  cs_navsto_key_t {
  CS_NSKEY_ADVECTION_EXTRAPOL, CS_NSKEY_ADVECTION_FORMULATION, CS_NSKEY_ADVECTION_SCHEME, CS_NSKEY_ADVECTION_STRATEGY,
  CS_NSKEY_DOF_REDUCTION, CS_NSKEY_GD_SCALE_COEF, CS_NSKEY_IL_ALGO_ATOL, CS_NSKEY_IL_ALGO_DTOL,
  CS_NSKEY_IL_ALGO_RTOL, CS_NSKEY_IL_ALGO_RESTART, CS_NSKEY_IL_ALGO_VERBOSITY, CS_NSKEY_MAX_IL_ALGO_ITER,
  CS_NSKEY_MAX_NL_ALGO_ITER, CS_NSKEY_MAX_OUTER_ITER, CS_NSKEY_NL_ALGO, CS_NSKEY_NL_ALGO_ATOL,
  CS_NSKEY_NL_ALGO_DTOL, CS_NSKEY_NL_ALGO_RTOL, CS_NSKEY_NL_ALGO_VERBOSITY, CS_NSKEY_QUADRATURE,
  CS_NSKEY_SCHUR_STRATEGY, CS_NSKEY_SLES_STRATEGY, CS_NSKEY_SPACE_SCHEME, CS_NSKEY_THERMAL_TOLERANCE,
  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

double cs_navsto_param_get_notay_scaling (void)
 Retrieve the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha. More...
 
void cs_navsto_param_set_notay_scaling (double scaling_coef)
 Set the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha. More...
 
cs_navsto_param_tcs_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_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...
 
cs_navsto_param_boussinesq_tcs_navsto_param_add_boussinesq_term (cs_navsto_param_t *nsp, cs_real_t dilatation_coef, cs_real_t reference_value)
 Add a new Boussinesq term (source term for the momemtum equation) More...
 
void cs_navsto_param_set_boussinesq_array (cs_navsto_param_boussinesq_t *bp, const cs_real_t *var)
 Set the array of values to consider in the Boussinesq term. More...
 
cs_navsto_param_sles_tcs_navsto_param_get_sles_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...
 
cs_equation_param_tcs_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_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 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_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 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_tcs_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_tcs_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_tcs_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_tcs_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_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_tcs_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_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...
 

Macro Definition Documentation

◆ CS_NAVSTO_STREAM_EQNAME

#define CS_NAVSTO_STREAM_EQNAME   "streamfunction_eq"

Typedef Documentation

◆ cs_navsto_param_model_flag_t

◆ cs_navsto_param_post_flag_t

Enumeration Type Documentation

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

Number of iterations before restarting a Krylov solver as the main solver (useful if the strategy implied a GMRES, flexible GMRES or GCR)

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 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 (Picard or Anderson)

CS_NSKEY_NL_ALGO_ATOL 

Absolute tolerance at which the non-linearity arising from the advection term is resolved

CS_NSKEY_NL_ALGO_DTOL 

Divergence tolerance at which the non-linearity arising from the advection term is detected

CS_NSKEY_NL_ALGO_RTOL 

Relative tolerance at which the non-linearity arising from the advection term is resolved

CS_NSKEY_NL_ALGO_VERBOSITY 

Level of verbosity related to the non-linear algorithm

CS_NSKEY_QUADRATURE 

Set the type to use in all functions involving quadrature (similar to CS_EQKEY_BC_QUADRATURE)

CS_NSKEY_SCHUR_STRATEGY 

Set the way to define the Schur complement approximation (cf. cs_param_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:

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

  • 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_MONOLITHIC 

The system is treated as a "monolithic" matrix

CS_NAVSTO_COUPLING_PROJECTION 

The system is solved using an incremental projection algorithm

CS_NAVSTO_N_COUPLINGS 

◆ cs_navsto_param_model_bit_t

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. This option is not compatible with the option CS_NAVSTO_MODEL_BOUSSINESQ

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_WITH_SOLIDIFICATION 

The Navier-Stokes is modified to take into account solidification process. A boussinesq term is added as well as a head loss term derived from a Darcy approximation. Please see the cs_solidification_model_t structure and related structures for more details

◆ cs_navsto_param_model_t

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 

◆ cs_navsto_param_post_bit_t

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.

Enumerator
CS_NAVSTO_POST_VELOCITY_DIVERGENCE 

Compute div(u) and associate a field to this quantity

CS_NAVSTO_POST_KINETIC_ENERGY 

Compute rho/2 u u and associate a field to this quantity

CS_NAVSTO_POST_VORTICITY 

Compute w = curl(u) and associate a field to this quantity

CS_NAVSTO_POST_VELOCITY_GRADIENT 

Compute the tensor grad(u) and associate a field to this quantity

CS_NAVSTO_POST_STREAM_FUNCTION 

Add an equation to compute the stream function and associate a field to this quantity. This equation is a scalar valued diffusion equation: -Lap(Psi) = w_z This is only useful for a 2D computation (assuming that the z-azis is the extruded one, i.e. the flow is in the x-y plane

CS_NAVSTO_POST_HELICITY 

Compute h = u w and associate a field to this quantity

CS_NAVSTO_POST_ENSTROPHY 

Compute w w and associate a field to this quantity

CS_NAVSTO_POST_MASS_DENSITY 

Compute the mass density in each cell and monitor the evolution of the integral of the mass in the full computational domain. If a Boussinesq approximation is set, then the mass density used in the buoyancy term is considered (not the reference value assumed to be constant).

CS_NAVSTO_POST_CELL_MASS_FLUX_BALANCE 

Compute the balance of the mass flux for each cell. When the mass density is constant, this quantity is only a scaling (by the mass density) of the quantity computed with the flag CS_NAVSTO_POST_VELOCITY_DIVERGENCE.

CS_NAVSTO_POST_PRESSURE_GRADIENT 

Compute the presure gradient in each cell.

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

Associated keyword: "diag_schur_gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm with a block diagonal preconditioner using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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_DIAG_SCHUR_MINRES 

Associated keyword: "diag_schur_minres"

The Stokes or Navier-Stokes system with an explicit advection is solved using a MINRES algorithm with a block diagonal preconditioner using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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_GCR 

Associated keyword: "gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm without preconditioning. The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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_LOWER_SCHUR_GCR 

Associated keyword: "lower_schur_gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm with an lower triangular block preconditioner using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

CS_NAVSTO_SLES_MINRES 

Associated keyword: "minres"

The Stokes or Navier-Stokes system with an explicit advection is solved using a MINRES algorithm without preconditioning. The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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 the full (saddle-point) system arising from the discretization of the Navier-Stokes equations

CS_NAVSTO_SLES_NOTAY_TRANSFORM 

Associated keyword: "notay"

Transform the saddle-point problem into an equivalent system without a zero block for the (2,2) block. This transformation allows one to consider standard preconditionner and iterative Krylow solver.

CS_NAVSTO_SLES_SGS_SCHUR_GCR 

Associated keyword: "sgs_schur_gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm with a symmetric Gauss-Seidel block preconditioner using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

CS_NAVSTO_SLES_UPPER_SCHUR_GCR 

Associated keyword: "upper_schur_gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm with an upper triangular block preconditioner using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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_USER 

Associated keyword: "user"

Resolution using a user-defined function. This function fulfills a pre-defined prototype

CS_NAVSTO_SLES_UZAWA_SCHUR_GCR 

Associated keyword: "uza_schur_gcr"

The Stokes or Navier-Stokes system is solved using a GCR algorithm with an Uzawa algorithm tuned for block preconditioning and using a Schur approximation for the pressure block (the block 22). The system is stored using a hybrid assembled/unassembled blocks. The velocity block is assembled (with potentially sub-blocks for each component) and the velocity divergence/pressure gradient operators are unassembled.

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 

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

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

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_add_boussinesq_term()

cs_navsto_param_boussinesq_t* cs_navsto_param_add_boussinesq_term ( cs_navsto_param_t nsp,
cs_real_t  dilatation_coef,
cs_real_t  reference_value 
)

Add a new Boussinesq term (source term for the momemtum equation)

Parameters
[in,out]nsppointer to a cs_navsto_param_t structure
Returns
a pointer to the newly added 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_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.

Parameters
[in]boundariespointer to a cs_boundary_t structure
[in]modeltype of model related to the NS system
[in]model_flagadditional high-level model options
[in]algo_couplingalgorithm used for solving the NS system
[in]post_flagpredefined post-processings options
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_get_model_name()

const char* cs_navsto_param_get_model_name ( cs_navsto_param_model_t  model)

Retrieve the name of the model system of equations.

Parameters
[in]modela cs_navsto_param_model_t
Returns
the corresponding name

◆ cs_navsto_param_get_notay_scaling()

double cs_navsto_param_get_notay_scaling ( void  )

Retrieve the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha.

Returns
the value of the scaling coefficient

◆ cs_navsto_param_get_sles_param()

cs_navsto_param_sles_t* cs_navsto_param_get_sles_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)

Parameters
[in]nsppointer to a cs_navsto_param_t structure
Returns
a pointer to the set of SLES parameters

◆ cs_navsto_param_get_velocity_param()

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)

Parameters
[in]nsppointer to a cs_navsto_param_t structure
Returns
a pointer to the set of parameters related to the momentum equation

◆ 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_set_boussinesq_array()

void cs_navsto_param_set_boussinesq_array ( cs_navsto_param_boussinesq_t bp,
const cs_real_t var 
)

Set the array of values to consider in the Boussinesq term.

Parameters
[in,out]bppointer to a cs_navsto_param_boussinesq_t structure
[in]varshared pointer to the array of values to consider

◆ cs_navsto_param_set_notay_scaling()

void cs_navsto_param_set_notay_scaling ( double  scaling_coef)

Set the scaling coefficient used in the Notay's transformation devised in "Algebraic multigrid for Stokes equations" SIAM J. Sci. Comput. Vol. 39 (5), 2017 In this article, this scaling is denoted by alpha.

Parameters
[in]scaling_coefvalued of the scaling coefficient

◆ 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()

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.

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
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_set_reference_pressure()

void cs_navsto_set_reference_pressure ( cs_navsto_param_t nsp,
cs_real_t  pref 
)

Set the value to consider for the reference pressure.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]prefvalue of the reference pressure

◆ 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()

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.

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
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_set_velocity_inlet_by_array()

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.

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]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_set_velocity_inlet_by_dof_func()

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.

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]funcpointer to a cs_dof_function_t
[in]func_inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_set_velocity_inlet_by_value()

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.

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
Returns
a pointer to the new cs_xdef_t structure

◆ cs_navsto_set_velocity_wall_by_value()

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.

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
Returns
a pointer to the new cs_xdef_t structure