8.3
general documentation
cs_navsto_system.h File Reference
#include "cs_advection_field.h"
#include "cs_cdo_turbulence.h"
#include "cs_equation.h"
#include "cs_field.h"
#include "cs_param_types.h"
#include "cs_property.h"
#include "cs_maxwell.h"
#include "cs_mesh.h"
#include "cs_navsto_param.h"
#include "cs_time_plot.h"
#include "cs_time_step.h"
#include "cs_thermal_system.h"
#include "cs_xdef.h"
+ Include dependency graph for cs_navsto_system.h:

Go to the source code of this file.

Data Structures

struct  cs_navsto_system_t
 Structure managing the Navier-Stokes system. More...
 

Typedefs

typedef void *() cs_navsto_init_scheme_context_t(const cs_navsto_param_t *nsp, cs_adv_field_t *adv_field, cs_real_t *mflux, cs_real_t *mflux_pre, cs_boundary_type_t *fb_type, void *nscc)
 Allocate and initialize the context structure related to a given discretization scheme for the resolution of the Navier-Stokes system with a specified coupling algorithm. More...
 
typedef void *() cs_navsto_free_scheme_context_t(void *scheme_context)
 Free the context structure related to a given discretization scheme for the resolution of the Navier-Stokes system with a specified coupling algorithm. More...
 
typedef void() cs_navsto_init_values_t(const cs_navsto_param_t *nsp, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts, cs_field_t *field)
 According to the model, coupling algorithm and the space discretization, initialize the field values which are not associated to a cs_equation_t structure (which manages the initialization) More...
 
typedef void() cs_navsto_compute_t(const cs_mesh_t *mesh, const cs_navsto_param_t *nsp, void *scheme_context)
 Compute for the current time step the new state for the Navier-Stokes system. This means that equations are built and then solved. More...
 

Functions

bool cs_navsto_system_is_activated (void)
 Check if the resolution of the Navier-Stokes system has been activated. More...
 
void cs_navsto_system_update_model (bool with_thermal)
 Update the flag associated to the modelling options. More...
 
cs_navsto_system_tcs_navsto_system_activate (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)
 Allocate and initialize the Navier-Stokes (NS) system. More...
 
void cs_navsto_system_destroy (void)
 Free the main structure related to the Navier-Stokes system. More...
 
cs_navsto_param_tcs_navsto_system_get_param (void)
 Retrieve the structure storing the parameters for the Navier–Stokes system. More...
 
cs_equation_tcs_navsto_system_get_momentum_eq (void)
 Retrieve a pointer to the equation related to the momentum equation. More...
 
cs_adv_field_tcs_navsto_get_adv_field (void)
 Retrieve the advection field structure (the mass flux) related to the Navier-Stokes system. More...
 
cs_real_tcs_navsto_get_mass_flux (bool previous)
 Retrieve the mass flux array related to the Navier-Stokes system. More...
 
void cs_navsto_system_init_setup (void)
 Start setting-up the Navier-Stokes system At this stage, numerical settings should be completely determined but connectivity and geometrical information is not yet available. More...
 
void cs_navsto_system_finalize_setup (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Last step of the setup of the Navier-Stokes system. More...
 
void cs_navsto_system_define_context (const cs_mesh_t *mesh)
 Initialize the scheme context structure used to build the algebraic system. This is done after the setup step. More...
 
void cs_navsto_system_init_values (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Set an initial value for the velocity and pressure fields as well as mass fluxes and tubulent quantities if needed. More...
 
void cs_navsto_system_set_solid_cells (cs_lnum_t n_solid_cells, cs_lnum_t solid_cell_ids[])
 Set a solid zone related to the Navier-Stokes equations. More...
 
void cs_navsto_system_update (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Update variables and related quantities when a new state of the Navier-Stokes system has been computed. More...
 
void cs_navsto_system_compute_steady_state (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Build, solve and update the Navier-Stokes system in case of a steady-state approach. More...
 
void cs_navsto_system_compute (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Build, solve and update the Navier-Stokes system. More...
 
void cs_navsto_system_extra_op (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Predefined extra-operations for the Navier-Stokes system. More...
 
void cs_navsto_system_extra_post (void *input, int mesh_id, int cat_id, int ent_flag[5], cs_lnum_t n_cells, cs_lnum_t n_i_faces, cs_lnum_t n_b_faces, const cs_lnum_t cell_ids[], const cs_lnum_t i_face_ids[], const cs_lnum_t b_face_ids[], const cs_time_step_t *time_step)
 Predefined post-processing output for the Navier-Stokes system. The prototype of this function is fixed since it is a function pointer defined in cs_post.h (cs_post_time_mesh_dep_output_t) More...
 
void cs_navsto_system_log_setup (void)
 Summary of the main cs_navsto_system_t structure. More...
 

Typedef Documentation

◆ cs_navsto_compute_t

typedef void() cs_navsto_compute_t(const cs_mesh_t *mesh, const cs_navsto_param_t *nsp, void *scheme_context)

Compute for the current time step the new state for the Navier-Stokes system. This means that equations are built and then solved.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]nsppointer to a cs_navsto_param_t structure
[in,out]scheme_contextpointer to a structure cast on-the-fly

◆ cs_navsto_free_scheme_context_t

typedef void *() cs_navsto_free_scheme_context_t(void *scheme_context)

Free the context structure related to a given discretization scheme for the resolution of the Navier-Stokes system with a specified coupling algorithm.

Parameters
[in,out]scheme_contextpointer to a structure cast on-the-fly
Returns
a null pointer

◆ cs_navsto_init_scheme_context_t

typedef void *() cs_navsto_init_scheme_context_t(const cs_navsto_param_t *nsp, cs_adv_field_t *adv_field, cs_real_t *mflux, cs_real_t *mflux_pre, cs_boundary_type_t *fb_type, void *nscc)

Allocate and initialize the context structure related to a given discretization scheme for the resolution of the Navier-Stokes system with a specified coupling algorithm.

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]adv_fieldpointer to cs_adv_field_t structure
[in]mfluxcurrent values of the mass flux
[in]mflux_precurrent values of the mass flux
[in]fb_typetype of boundary for each boundary face
[in,out]nsccNavier-Stokes coupling context: pointer to a structure cast on-the-fly
Returns
a pointer to a new allocated scheme context structure

◆ cs_navsto_init_values_t

typedef void() cs_navsto_init_values_t(const cs_navsto_param_t *nsp, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts, cs_field_t *field)

According to the model, coupling algorithm and the space discretization, initialize the field values which are not associated to a cs_equation_t structure (which manages the initialization)

Parameters
[in]nsppointer to a cs_navsto_param_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]tspointer to a cs_time_step_t structure
[in,out]fieldpointer to a cs_field_t structure

Function Documentation

◆ cs_navsto_get_adv_field()

cs_adv_field_t * cs_navsto_get_adv_field ( void  )

Retrieve the advection field structure (the mass flux) related to the Navier-Stokes system.

Returns
a pointer to the advection field structure

◆ cs_navsto_get_mass_flux()

cs_real_t * cs_navsto_get_mass_flux ( bool  previous)

Retrieve the mass flux array related to the Navier-Stokes system.

Parameters
[in]previousif true return the previous state otherwise the current state.
Returns
a pointer to an array of cs_real_t

◆ cs_navsto_system_activate()

cs_navsto_system_t * cs_navsto_system_activate ( 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 
)

Allocate and initialize the Navier-Stokes (NS) system.

Parameters
[in]boundariespointer to the domain boundaries
[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 cs_navsto_system_t structure

◆ cs_navsto_system_compute()

void cs_navsto_system_compute ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Build, solve and update the Navier-Stokes system.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_stepstructure managing the time stepping

◆ cs_navsto_system_compute_steady_state()

void cs_navsto_system_compute_steady_state ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Build, solve and update the Navier-Stokes system in case of a steady-state approach.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_stepstructure managing the time stepping

◆ cs_navsto_system_define_context()

void cs_navsto_system_define_context ( const cs_mesh_t mesh)

Initialize the scheme context structure used to build the algebraic system. This is done after the setup step.

Parameters
[in]meshpointer to a cs_mesh_t structure

◆ cs_navsto_system_destroy()

void cs_navsto_system_destroy ( void  )

Free the main structure related to the Navier-Stokes system.

◆ cs_navsto_system_extra_op()

void cs_navsto_system_extra_op ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Predefined extra-operations for the Navier-Stokes system.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_steppointer to a cs_time_step_t structure

◆ cs_navsto_system_extra_post()

void cs_navsto_system_extra_post ( void *  input,
int  mesh_id,
int  cat_id,
int  ent_flag[5],
cs_lnum_t  n_cells,
cs_lnum_t  n_i_faces,
cs_lnum_t  n_b_faces,
const cs_lnum_t  cell_ids[],
const cs_lnum_t  i_face_ids[],
const cs_lnum_t  b_face_ids[],
const cs_time_step_t time_step 
)

Predefined post-processing output for the Navier-Stokes system. The prototype of this function is fixed since it is a function pointer defined in cs_post.h (cs_post_time_mesh_dep_output_t)

Parameters
[in,out]inputpointer to a optional structure (here a cs_navsto_system_t structure)
[in]mesh_idid of the output mesh for the current call
[in]cat_idcategory id of the output mesh for this call
[in]ent_flagindicate global presence of cells (ent_flag[0]), interior faces (ent_flag[1]), boundary faces (ent_flag[2]), particles (ent_flag[3]) or probes (ent_flag[4])
[in]n_cellslocal number of cells of post_mesh
[in]n_i_faceslocal number of interior faces of post_mesh
[in]n_b_faceslocal number of boundary faces of post_mesh
[in]cell_idslist of cells (0 to n-1)
[in]i_face_idslist of interior faces (0 to n-1)
[in]b_face_idslist of boundary faces (0 to n-1)
[in]time_steppointer to a cs_time_step_t struct.

◆ cs_navsto_system_finalize_setup()

void cs_navsto_system_finalize_setup ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Last step of the setup of the Navier-Stokes system.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_steppointer to a cs_time_step_t structure

◆ cs_navsto_system_get_momentum_eq()

cs_equation_t * cs_navsto_system_get_momentum_eq ( void  )

Retrieve a pointer to the equation related to the momentum equation.

Returns
nullptr or the pointer
null or the pointer

◆ cs_navsto_system_get_param()

cs_navsto_param_t * cs_navsto_system_get_param ( void  )

Retrieve the structure storing the parameters for the Navier–Stokes system.

Returns
nullptr or the pointer to a cs_navsto_param_t structure
null or the pointer to a cs_navsto_param_t structure

◆ cs_navsto_system_init_setup()

void cs_navsto_system_init_setup ( void  )

Start setting-up the Navier-Stokes system At this stage, numerical settings should be completely determined but connectivity and geometrical information is not yet available.

◆ cs_navsto_system_init_values()

void cs_navsto_system_init_values ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Set an initial value for the velocity and pressure fields as well as mass fluxes and tubulent quantities if needed.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_steppointer to a cs_time_step_t structure

◆ cs_navsto_system_is_activated()

bool cs_navsto_system_is_activated ( void  )

Check if the resolution of the Navier-Stokes system has been activated.

Returns
true or false

◆ cs_navsto_system_log_setup()

void cs_navsto_system_log_setup ( void  )

Summary of the main cs_navsto_system_t structure.

◆ cs_navsto_system_set_solid_cells()

void cs_navsto_system_set_solid_cells ( cs_lnum_t  n_solid_cells,
cs_lnum_t  solid_cell_ids[] 
)

Set a solid zone related to the Navier-Stokes equations.

Parameters
[in]n_solid_cellsnumber of solid cells
[in]solid_cell_idslist of cell ids (local numbering)

◆ cs_navsto_system_update()

void cs_navsto_system_update ( const cs_mesh_t mesh,
const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t time_step 
)

Update variables and related quantities when a new state of the Navier-Stokes system has been computed.

Parameters
[in]meshpointer to a cs_mesh_t structure
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]time_stepstructure managing the time stepping

◆ cs_navsto_system_update_model()

void cs_navsto_system_update_model ( bool  with_thermal)

Update the flag associated to the modelling options.

Parameters
[in]with_thermaltrue or false