My Project
programmer's documentation
|
Routines to handle the cs_navsto_system_t structure. More...
#include "cs_advection_field.h"
#include "cs_equation.h"
#include "cs_field.h"
#include "cs_param.h"
#include "cs_property.h"
#include "cs_mesh.h"
#include "cs_navsto_param.h"
#include "cs_time_step.h"
#include "cs_xdef.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_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... | |
cs_navsto_system_t * | cs_navsto_system_activate (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) |
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_t * | cs_navsto_system_get_param (void) |
Retrieve the structure storing the parameters for the Navier–Stokes system. More... | |
cs_equation_t * | cs_navsto_system_get_momentum_eq (void) |
Retrieve a pointer to the equation related to the momentum equation. 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_set_sles (void) |
Define the settings for SLES related to the Navier-Stokes system. More... | |
void | cs_navsto_system_initialize (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts) |
Initialize the context structure used to build the algebraic system This is done after the setup step. More... | |
void | cs_navsto_system_compute_steady_state (const cs_mesh_t *mesh, 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_time_step_t *time_step) |
Build, solve and update the Navier-Stokes system. More... | |
void | cs_navsto_system_extra_op (const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *cdoq) |
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... | |
Routines to handle the cs_navsto_system_t structure.
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.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in,out] | scheme_context | pointer to a structure cast on-the-fly |
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.
[in,out] | scheme_context | pointer to a structure cast on-the-fly |
typedef void*() cs_navsto_init_scheme_context_t(const cs_navsto_param_t *nsp, 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.
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | fb_type | type of boundary for each boundary face |
[in,out] | nscc | Navier-Stokes coupling context: pointer to a structure cast on-the-fly |
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)
[in] | nsp | pointer to a cs_navsto_param_t structure |
[in] | quant | pointer to a cs_cdo_quantities_t structure |
[in] | ts | pointer to a cs_time_step_t structure |
[in,out] | field | pointer to a cs_field_t structure |
cs_navsto_system_t* cs_navsto_system_activate | ( | 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 | ||
) |
Allocate and initialize the Navier-Stokes (NS) system.
[in] | boundaries | pointer to the domain boundaries |
[in] | model | type of model related to the NS system |
[in] | time_state | state of the time for the NS equations |
[in] | algo_coupling | algorithm used for solving the NS system |
void cs_navsto_system_compute | ( | const cs_mesh_t * | mesh, |
const cs_time_step_t * | time_step | ||
) |
Build, solve and update the Navier-Stokes system.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | time_step | structure managing the time stepping |
void cs_navsto_system_compute_steady_state | ( | const cs_mesh_t * | mesh, |
const cs_time_step_t * | time_step | ||
) |
Build, solve and update the Navier-Stokes system in case of a steady-state approach.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | time_step | structure managing the time stepping |
void cs_navsto_system_destroy | ( | void | ) |
Free the main structure related to the Navier-Stokes system.
void cs_navsto_system_extra_op | ( | const cs_cdo_connect_t * | connect, |
const cs_cdo_quantities_t * | cdoq | ||
) |
Predefined extra-operations for the Navier-Stokes system.
[in] | connect | pointer to a cs_cdo_connect_t structure |
[in] | cdoq | pointer to a cs_cdo_quantities_t structure |
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)
[in,out] | input | pointer to a optional structure (here a cs_gwf_t structure) |
[in] | mesh_id | id of the output mesh for the current call |
[in] | cat_id | category id of the output mesh for this call |
[in] | ent_flag | indicate 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_cells | local number of cells of post_mesh |
[in] | n_i_faces | local number of interior faces of post_mesh |
[in] | n_b_faces | local number of boundary faces of post_mesh |
[in] | cell_ids | list of cells (0 to n-1) |
[in] | i_face_ids | list of interior faces (0 to n-1) |
[in] | b_face_ids | list of boundary faces (0 to n-1) |
[in] | time_step | pointer to a cs_time_step_t struct. |
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.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | connect | pointer to a cs_cdo_connect_t structure |
[in] | quant | pointer to a cs_cdo_quantities_t structure |
[in] | time_step | pointer to a cs_time_step_t structure |
cs_equation_t* cs_navsto_system_get_momentum_eq | ( | void | ) |
Retrieve a pointer to the equation related to the momentum equation.
cs_navsto_param_t* cs_navsto_system_get_param | ( | void | ) |
Retrieve the structure storing the parameters for the Navier–Stokes system.
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.
void cs_navsto_system_initialize | ( | const cs_mesh_t * | mesh, |
const cs_cdo_connect_t * | connect, | ||
const cs_cdo_quantities_t * | quant, | ||
const cs_time_step_t * | ts | ||
) |
Initialize the context structure used to build the algebraic system This is done after the setup step.
[in] | mesh | pointer to a cs_mesh_t structure |
[in] | connect | pointer to a cs_cdo_connect_t structure |
[in] | quant | pointer to a cs_cdo_quantities_t structure |
[in] | ts | pointer to a cs_time_step_t structure |
bool cs_navsto_system_is_activated | ( | void | ) |
Check if the resolution of the Navier-Stokes system has been activated.
void cs_navsto_system_log_setup | ( | void | ) |
Summary of the main cs_navsto_system_t structure.
void cs_navsto_system_set_sles | ( | void | ) |
Define the settings for SLES related to the Navier-Stokes system.