#include "cs_advection_field.h"
#include "cs_enforcement.h"
#include "cs_hodge.h"
#include "cs_iter_algo.h"
#include "cs_param_cdo.h"
#include "cs_param_saddle.h"
#include "cs_param_sles.h"
#include "cs_property.h"
#include "cs_xdef.h"
Go to the source code of this file.
Data Structures | |
struct | cs_equation_param_t |
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources. More... | |
Macros | |
Flags specifying the behavior of an equation (activated terms like the | |
convection, diffusion, etc.). | |
#define | CS_EQUATION_LOCKED (1 << 0) /* 1 */ |
Parameters for setting an equation are not modifiable anymore. More... | |
#define | CS_EQUATION_UNSTEADY (1 << 1) /* 2 */ |
Unsteady term is needed. More... | |
#define | CS_EQUATION_CONVECTION (1 << 2) /* 4 */ |
Convection term is needed. More... | |
#define | CS_EQUATION_DIFFUSION (1 << 3) /* 8 */ |
Diffusion term is needed. A scalar-/vector-valued Laplacian with div .grad. More... | |
#define | CS_EQUATION_CURLCURL (1 << 4) /* 16 */ |
The term related to the curl-curl operator is needed. More... | |
#define | CS_EQUATION_GRADDIV (1 << 5) /* 32 */ |
The term related to the grad-div operator is needed. More... | |
#define | CS_EQUATION_REACTION (1 << 6) /* 64 */ |
Reaction term is needed. More... | |
#define | CS_EQUATION_FORCE_VALUES (1 << 7) /* 128 */ |
Add an algebraic manipulation to set the value of a given set of interior degrees of freedom. More... | |
#define | CS_EQUATION_INSIDE_SYSTEM (1 << 8) /* 256 */ |
The current equation settings belong to a system of equation. More... | |
#define | CS_EQUATION_BUILD_HOOK (1 << 9) /* 512 */ |
Activate a build hook function to get a fine control of the discretization process during the cellwise building of the linear system. Need a function matching the cs_equation_build_hook_t prototype. More... | |
#define | CS_EQUATION_USER_TRIGGERED (1 << 10) /* 1024 */ |
The stage at which the equation is solved is driven by the user. So, this user-defined equation is not solved at the same time as the other user-defined equations. More... | |
Flags specifying which extra operation related to an equation is | |
needed. | |
#define | CS_EQUATION_POST_BALANCE (1 << 0) /* 1 */ |
Compute and postprocess the equation balance. More... | |
#define | CS_EQUATION_POST_PECLET (1 << 1) /* 2 */ |
Compute and postprocess the Peclet number. More... | |
#define | CS_EQUATION_POST_UPWIND_COEF (1 << 2) /* 4 */ |
Postprocess the value of the upwinding coefficient. Only useful if variable coefficient controls the amount of upwinding. For instant, in schemes such as Schaffeter-Gummel, the upwind portion is driven by the local Peclet number. More... | |
#define | CS_EQUATION_POST_NORMAL_FLUX (1 << 3) /* 8 */ |
Postprocess the value of the normal flux across the boundary faces. More... | |
Functions | |
static void | cs_equation_param_set_flag (cs_equation_param_t *eqp, cs_flag_t flag) |
Update the flag related to a cs_equation_param_t structure. More... | |
static bool | cs_equation_param_has_diffusion (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a diffusion term. More... | |
static bool | cs_equation_param_has_curlcurl (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a curl-curl term. More... | |
static bool | cs_equation_param_has_graddiv (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a grad-div term. More... | |
static bool | cs_equation_param_has_convection (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a convection term. More... | |
static bool | cs_equation_param_has_reaction (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a reaction term. More... | |
static bool | cs_equation_param_has_time (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs an unsteady term. More... | |
static bool | cs_equation_param_has_sourceterm (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation needs a source term. More... | |
static bool | cs_equation_param_has_internal_enforcement (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation has an internal enforcement of the degrees of freedom. More... | |
static bool | cs_equation_param_has_implicit_advection (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation induces an implicit treatment of the advection term. More... | |
static bool | cs_equation_param_has_user_hook (const cs_equation_param_t *eqp) |
Ask if the parameters of the equation has activated a user hook to get a fine tuning of the cellwise system building. More... | |
static bool | cs_equation_param_has_name (cs_equation_param_t *eqp, const char *name) |
Check if a cs_equation_param_t structure has its name member equal to the given name. More... | |
cs_equation_param_t * | cs_equation_param_create (const char *name, cs_equation_type_t type, int dim, cs_param_bc_type_t default_bc) |
Create a cs_equation_param_t structure. More... | |
static cs_equation_param_t * | cs_equation_create_param (const char *name, cs_equation_type_t type, int dim, cs_param_bc_type_t default_bc) |
Create a cs_equation_param_t structure. More... | |
void | cs_equation_param_copy_from (const cs_equation_param_t *ref, cs_equation_param_t *dst, bool copy_fld_id) |
Copy the settings from one cs_equation_param_t structure to another one. The name is not copied. More... | |
void | cs_equation_param_copy_bc (const cs_equation_param_t *ref, cs_equation_param_t *dst) |
Copy only the part dedicated to the boundary conditions and the DoF (degrees of freedom) enforcement from one cs_equation_param_t structure to another one. More... | |
static void | cs_equation_copy_param_from (const cs_equation_param_t *ref, cs_equation_param_t *dst, bool copy_fid) |
Copy the settings from one cs_equation_param_t structure to another one. The name is not copied. More... | |
void | cs_equation_param_clear (cs_equation_param_t *eqp) |
Free the contents of a cs_equation_param_t. More... | |
cs_equation_param_t * | cs_equation_param_free (cs_equation_param_t *eqp) |
Free a cs_equation_param_t. More... | |
void | cs_equation_param_set (cs_equation_param_t *eqp, cs_equation_key_t key, const char *keyval) |
Set a parameter attached to a keyname in a cs_equation_param_t structure. More... | |
static void | cs_equation_set_param (cs_equation_param_t *eqp, cs_equation_key_t key, const char *keyval) |
Set a parameter attached to a keyname in a cs_equation_param_t structure. More... | |
cs_param_sles_t * | cs_equation_param_get_sles_param (cs_equation_param_t *eqp) |
Get the pointer to the set of parameters to handle the SLES solver associated to this set of equation parameters. More... | |
cs_param_saddle_t * | cs_equation_param_get_saddle_param (cs_equation_param_t *eqp) |
Get the pointer to the set of parameters to handle the SLES solver associated to this set of equation parameters. More... | |
void | cs_equation_param_set_sles (cs_equation_param_t *eqp) |
Set parameters for initializing SLES structures used for the resolution of the linear system. Settings are related to this equation. More... | |
void | cs_equation_param_set_quadrature_to_all (cs_equation_param_t *eqp, cs_quadrature_type_t qtype) |
Apply the given quadrature rule to all existing definitions under the cs_equation_param_t structure. To get a more detailed control of the quadrature rule, please consider the function cs_xdef_set_quadrature. More... | |
void | cs_equation_param_lock_settings (cs_equation_param_t *eqp) |
Lock settings to prevent from unwanted modifications. More... | |
void | cs_equation_param_unlock_settings (cs_equation_param_t *eqp) |
Unlock settings. Be sure that is really wanted (inconsistency between the setup logging and what is used may appear) More... | |
void | cs_equation_param_ensure_consistent_settings (cs_equation_param_t *eqp) |
At this stage, the numerical settings should not be modified anymore by the user. One makes a last set of modifications if needed to ensure a consistent numerical settings. More... | |
void | cs_equation_param_log (const cs_equation_param_t *eqp) |
Print the detail of a cs_equation_param_t structure. More... | |
bool | cs_equation_param_has_robin_bc (const cs_equation_param_t *eqp) |
Ask if the parameter settings of the equation has requested the treatment of Robin boundary conditions. More... | |
cs_xdef_t * | cs_equation_add_ic_by_value (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here a constant value is set to all the unknows belonging to the given zone with name z_name. More... | |
cs_xdef_t * | cs_equation_add_ic_by_qov (cs_equation_param_t *eqp, const char *z_name, double quantity) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here the value set to each unknown belonging to the given zone with name z_name is such that the integral over the cells of the zone returns the requested quantity. More... | |
cs_xdef_t * | cs_equation_add_ic_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *analytic, void *input) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here the initial value for each unknown associated to the zone with name z_name is set according to an analytical function. More... | |
cs_xdef_t * | cs_equation_add_ic_by_dof_func (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Case of a definition by a DoF function. More... | |
void | cs_equation_add_xdef_bc (cs_equation_param_t *eqp, cs_xdef_t *xdef) |
Set a boundary condition from an existing cs_xdef_t structure The lifecycle of the cs_xdef_t structure is now managed by the current cs_equation_param_t structure. More... | |
cs_xdef_t * | cs_equation_add_bc_by_value (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_real_t *values) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More... | |
cs_xdef_t * | cs_equation_add_bc_by_array (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, bool full_length) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More... | |
cs_xdef_t * | cs_equation_add_bc_by_field (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_field_t *field) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More... | |
cs_xdef_t * | cs_equation_add_bc_by_analytic (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_analytic_func_t *analytic, void *input) |
Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t. More... | |
cs_xdef_t * | cs_equation_add_bc_by_time_func (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_time_func_t *t_func, void *input) |
Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t Definition relying on a cs_time_func_t function pointer. More... | |
cs_xdef_t * | cs_equation_add_bc_by_dof_func (cs_equation_param_t *eqp, const cs_param_bc_type_t bc_type, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input) |
Define and initialize a new structure to set a boundary condition related to the given cs_equation_param_t structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t. More... | |
cs_xdef_t * | cs_equation_add_volume_mass_injection_by_array (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc_flag, cs_real_t *array, bool is_owner, bool full_length) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an array. More... | |
cs_xdef_t * | cs_equation_find_bc (cs_equation_param_t *eqp, const char *z_name) |
Return pointer to existing boundary condition definition structure for the given equation param structure and zone. More... | |
void | cs_equation_remove_bc (cs_equation_param_t *eqp, const char *z_name) |
Remove boundary condition from the given equation param structure for a given zone. More... | |
void | cs_equation_remove_ic (cs_equation_param_t *eqp, const char *z_name) |
Remove initial condition from the given equation param structure for a given zone. More... | |
void | cs_equation_add_sliding_condition (cs_equation_param_t *eqp, const char *z_name) |
Define and initialize a new structure to set a sliding boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t. More... | |
void | cs_equation_add_diffusion (cs_equation_param_t *eqp, cs_property_t *property) |
Associate a new term related to the Laplacian operator for the equation associated to the given cs_equation_param_t structure Laplacian means div-grad (either for vector-valued or scalar-valued equations) More... | |
void | cs_equation_add_curlcurl (cs_equation_param_t *eqp, cs_property_t *property, int inversion) |
Associate a new term related to the curl-curl operator for the equation associated to the given cs_equation_param_t structure. More... | |
void | cs_equation_add_graddiv (cs_equation_param_t *eqp, cs_property_t *property) |
Associate a new term related to the grad-div operator for the equation associated to the given cs_equation_param_t structure. More... | |
void | cs_equation_add_time (cs_equation_param_t *eqp, cs_property_t *property) |
Associate a new term related to the time derivative operator for the equation associated to the given cs_equation_param_t structure. More... | |
void | cs_equation_add_advection (cs_equation_param_t *eqp, cs_adv_field_t *adv_field) |
Associate a new term related to the advection operator for the equation associated to the given cs_equation_param_t structure. More... | |
void | cs_equation_add_advection_scaling_property (cs_equation_param_t *eqp, cs_property_t *property) |
Associate a scaling property to the advection. More... | |
int | cs_equation_add_reaction (cs_equation_param_t *eqp, cs_property_t *property) |
Associate a new term related to the reaction operator for the equation associated to the given cs_equation_param_t structure. More... | |
cs_xdef_t * | cs_equation_add_source_term_by_val (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a constant value. More... | |
cs_xdef_t * | cs_equation_add_source_term_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *func, void *input) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an analytical function. More... | |
cs_xdef_t * | cs_equation_add_source_term_by_dof_func (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a DoF function. More... | |
cs_xdef_t * | cs_equation_add_source_term_by_array (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc, cs_real_t *array, bool is_owner, bool full_length) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an array. More... | |
cs_xdef_t * | cs_equation_add_volume_mass_injection_by_value (cs_equation_param_t *eqp, const char *z_name, double *val) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant value. More... | |
cs_xdef_t * | cs_equation_add_volume_mass_injection_by_qov (cs_equation_param_t *eqp, const char *z_name, double *quantity) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant quantity distributed over the associated zone's volume. More... | |
cs_xdef_t * | cs_equation_add_volume_mass_injection_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *func, void *input) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an analytical function. More... | |
cs_xdef_t * | cs_equation_add_volume_mass_injection_by_dof_func (cs_equation_param_t *eqp, const char *z_name, cs_flag_t loc_flag, cs_dof_func_t *func, void *input) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a DoF function. More... | |
cs_enforcement_param_t * | cs_equation_add_vertex_dof_enforcement (cs_equation_param_t *eqp, cs_lnum_t n_vertices, const cs_lnum_t vertex_ids[], const cs_real_t ref_value[], const cs_real_t vtx_values[]) |
Add an enforcement of the value of degrees of freedom located at the mesh vertices. The spatial discretization scheme for the given equation has to be CDO vertex-based or CDO vertex+cell-based schemes. More... | |
cs_enforcement_param_t * | cs_equation_add_edge_dof_enforcement (cs_equation_param_t *eqp, cs_lnum_t n_edges, const cs_lnum_t edge_ids[], const cs_real_t ref_value[], const cs_real_t edge_values[]) |
Add an enforcement of the value of degrees of freedom located at the mesh edges. The spatial discretization scheme for the given equation has to be CDO edge-based schemes. More... | |
cs_enforcement_param_t * | cs_equation_add_face_dof_enforcement (cs_equation_param_t *eqp, cs_lnum_t n_faces, const cs_lnum_t face_ids[], const cs_real_t ref_value[], const cs_real_t face_values[]) |
Add an enforcement of the value of degrees of freedom located at the mesh faces. The spatial discretization scheme for the given equation has to be CDO face-based schemes. More... | |
cs_enforcement_param_t * | cs_equation_add_cell_enforcement (cs_equation_param_t *eqp, cs_lnum_t n_cells, const cs_lnum_t elt_ids[], const cs_real_t ref_value[], const cs_real_t cell_values[]) |
Add an enforcement of the value related to the degrees of freedom associated to the list of selected cells. More... | |
cs_enforcement_param_t * | cs_equation_add_or_replace_cell_enforcement (cs_equation_param_t *eqp, int enforcement_id, cs_lnum_t n_cells, const cs_lnum_t cell_ids[], const cs_real_t ref_value[], const cs_real_t cell_values[]) |
Add a new enforcement if enforcement_id does not exist or replace it otherwise. Enforcement of the value related to the degrees of freedom associated to the list of selected cells. More... | |
#define CS_EQUATION_BUILD_HOOK (1 << 9) /* 512 */ |
Activate a build hook function to get a fine control of the discretization process during the cellwise building of the linear system. Need a function matching the cs_equation_build_hook_t prototype.
#define CS_EQUATION_CONVECTION (1 << 2) /* 4 */ |
Convection term is needed.
#define CS_EQUATION_CURLCURL (1 << 4) /* 16 */ |
The term related to the curl-curl operator is needed.
#define CS_EQUATION_DIFFUSION (1 << 3) /* 8 */ |
Diffusion term is needed. A scalar-/vector-valued Laplacian with div .grad.
#define CS_EQUATION_FORCE_VALUES (1 << 7) /* 128 */ |
Add an algebraic manipulation to set the value of a given set of interior degrees of freedom.
#define CS_EQUATION_GRADDIV (1 << 5) /* 32 */ |
The term related to the grad-div operator is needed.
#define CS_EQUATION_INSIDE_SYSTEM (1 << 8) /* 256 */ |
The current equation settings belong to a system of equation.
#define CS_EQUATION_LOCKED (1 << 0) /* 1 */ |
Parameters for setting an equation are not modifiable anymore.
#define CS_EQUATION_POST_BALANCE (1 << 0) /* 1 */ |
Compute and postprocess the equation balance.
#define CS_EQUATION_POST_NORMAL_FLUX (1 << 3) /* 8 */ |
Postprocess the value of the normal flux across the boundary faces.
#define CS_EQUATION_POST_PECLET (1 << 1) /* 2 */ |
Compute and postprocess the Peclet number.
#define CS_EQUATION_POST_UPWIND_COEF (1 << 2) /* 4 */ |
Postprocess the value of the upwinding coefficient. Only useful if variable coefficient controls the amount of upwinding. For instant, in schemes such as Schaffeter-Gummel, the upwind portion is driven by the local Peclet number.
#define CS_EQUATION_REACTION (1 << 6) /* 64 */ |
Reaction term is needed.
#define CS_EQUATION_UNSTEADY (1 << 1) /* 2 */ |
Unsteady term is needed.
#define CS_EQUATION_USER_TRIGGERED (1 << 10) /* 1024 */ |
The stage at which the equation is solved is driven by the user. So, this user-defined equation is not solved at the same time as the other user-defined equations.
enum cs_equation_key_t |
List of available keys for setting the parameters of an equation.
Enumerator | |
---|---|
CS_EQKEY_ADV_CIP_COEF | Set the value of the stabilization scaling coefficient when a CIP advection scheme is used. This value should be > 0 |
CS_EQKEY_ADV_EXTRAPOL | Choice in the way to extrapolate the advection field when building the advection operator
|
CS_EQKEY_ADV_FORMULATION | Kind of formulation of the advective term. There are two possible choices:
|
CS_EQKEY_ADV_SCHEME | Specify the type of numerical scheme for the advective term. The available choices depend on the space discretization scheme. Please refer to the section Set the advection scheme of the user guide for more details |
CS_EQKEY_ADV_STRATEGY | Strategy used to handle the advection term
|
CS_EQKEY_ADV_UPWIND_PORTION | Value between 0 and 1 specifying the portion of upwind added to a centered discretization. |
CS_EQKEY_AMG_TYPE | Specify which type of algebraic multigrid (AMG) to choose. Available choices are:
|
CS_EQKEY_BC_ENFORCEMENT | Set the type of enforcement of the boundary conditions. Available choices are:
For HHO and CDO-Face based schemes, only the "penalization" and "algebraic" technique is available up to now. |
CS_EQKEY_BC_STRONG_PENA_COEFF | Set the value of the penalization coefficient when "penalization" is activated The default value is 1e12. cf. CS_PARAM_BC_ENFORCE_PENALIZED |
CS_EQKEY_BC_WEAK_PENA_COEFF | Set the value of the penalization coefficient when "weak" or "weak_sym" is activated. The default value is 100. cf. CS_PARAM_BC_ENFORCE_WEAK_NITSCHE or CS_PARAM_BC_ENFORCE_WEAK_SYM |
CS_EQKEY_DO_LUMPING | |
CS_EQKEY_DOF_REDUCTION | Set how is defined each degree of freedom (DoF).
|
CS_EQKEY_EXTRA_OP | Set the additional post-processing to perform. Available choices are:
|
CS_EQKEY_HODGE_DIFF_ALGO | Set the algorithm used for building the discrete Hodge operator used in the diffusion term. Available choices are:
|
CS_EQKEY_HODGE_DIFF_COEF | This key is only useful if the key CS_EQKEY_HODGE_DIFF_ALGO is set to "cost" or "ocs". keyval is either a name or a value:
|
CS_EQKEY_HODGE_TIME_ALGO | Set the algorithm used for building the discrete Hodge operator used in the unsteady term. Available choices are:
|
CS_EQKEY_HODGE_REAC_ALGO | Set the algorithm used for building the discrete Hodge operator used in the reaction term. Available choices are:
|
CS_EQKEY_ITSOL | --> deprecated key (use CS_EQKEY_SOLVER instead) |
CS_EQKEY_ITSOL_ATOL | --> deprecated key (use CS_EQKEY_SOLVER_ATOL instead) |
CS_EQKEY_ITSOL_DTOL | --> deprecated key (use CS_EQKEY_SOLVER_DTOL instead) |
CS_EQKEY_ITSOL_EPS | --> deprecated key (use CS_EQKEY_SOLVER_RTOL instead) |
CS_EQKEY_ITSOL_MAX_ITER | --> deprecated key (use CS_EQKEY_SOLVER_MAX_ITER instead) |
CS_EQKEY_ITSOL_RESNORM_TYPE | --> deprecated key (use CS_EQKEY_SOLVER_RESNORM_TYPE instead) |
CS_EQKEY_ITSOL_RESTART | --> deprecated key (use CS_EQKEY_SOLVER_RESTART instead) |
CS_EQKEY_ITSOL_RTOL | --> deprecated key (use CS_EQKEY_SOLVER_RTOL instead) |
CS_EQKEY_SOLVER | Specify the solver for the resolution of the linear system related to an equation. Please refer to the section Set the linear solver of the user guide for more details |
CS_EQKEY_SOLVER_ATOL | Absolute tolerance factor for stopping the iterative process during the iterative resolution of a linear system related to an equation. Most iterative solver are not using this tolerance (the relative tolerance is always used). PETSc solvers use this information for instance. Please refer to CS_EQKEY_SOLVER_RTOL
|
CS_EQKEY_SOLVER_DTOL | Divergence tolerance factor for stopping the iterative process during the iterative resolution of a linear system related to an equation. Most iterative solver are not using this tolerance (the relative tolerance is always used). PETSc solvers use this information for instance. Please refer to CS_EQKEY_SOLVER_RTOL
|
CS_EQKEY_SOLVER_MAX_ITER | Maximum number of iterations for solving the linear system
|
CS_EQKEY_SOLVER_NO_OP | Set if the solver can be skipt if some specific situation
|
CS_EQKEY_SOLVER_RESNORM_TYPE | Normalized or not the residual before testing if one continues iterating for solving the linear system. This normalization is performed before applying the boundary conditions to avoid handling the penalization of Dirichlet boundary conditions. If the RHS norm is equal to zero, then the "vol_tot" option is used for rescaling the residual. Available choices are: "false" or "none" (default) "rhs" "weighted_rhs" or "weighted" "filtered_rhs" or "filtered_rhs" |
CS_EQKEY_SOLVER_RESTART | Maximum number of iterations before restarting a Krylov solver Only useful with GMRES, flexible GMRES or GCR solvers. This value has an impact on the memory footprint since one has to store a buffer of double with a size equal to restart*sizeof(solution array)
|
CS_EQKEY_SOLVER_RTOL | Relative tolerance factor for stopping the iterative process during the iterative resolution of a linear system related to an equation.
|
CS_EQKEY_OMP_ASSEMBLY_STRATEGY | |
CS_EQKEY_PRECOND | Specify the preconditioner associated to an iterative solver. Please refer to the section Set the preconditioner of the user guide for more details. Be careful since some options are only available with a given family of solvers. If this is the case, be sure that your installation has been installed with the appropriate library. |
CS_EQKEY_PRECOND_BLOCK_TYPE | Specify the type of block preconditioner associated to a preconditioner. Useful for a vector-valued equation for instance, where each component can be associated to a block. Warning: Most of these options are available only with PETSc/HYPRE.
|
CS_EQKEY_SADDLE_ATOL | Absolute tolerance under which the iterative process stops when solving a saddle-point problem related to this equation. Read the description of the structure cs_param_saddle_t for more details.
|
CS_EQKEY_SADDLE_AUGMENT_SCALING | Value of the scaling coefficient for the augmentation when solving a saddle-point problem. The initial system can namely be augmented with an operator for instance the grad-div operator when using an Augmented Lagrangian approach or the GKB solver.
|
CS_EQKEY_SADDLE_DTOL | Divergence tolerance above which the iterative process stops when solving a saddle-point problem related to this equation. This triggers an error for divergence of the solver when detected. Read the description of the structure cs_param_saddle_t for more details.
|
CS_EQKEY_SADDLE_MAX_ITER | Max. number of iterations before stopping the "saddle-point" solver when solving a saddle-point problem related to this equation. This is different from CS_EQKEY_SOLVER_MAX_ITER This latter key is dedicated to the preconditioning of the (1,1)-block when a saddle-point system has to be solved. Read the description of the structure cs_param_saddle_t for more details.
|
CS_EQKEY_SADDLE_PRECOND | Block preconditioner used to solve a saddle-point system.
|
CS_EQKEY_SADDLE_RTOL | Relative tolerance under which the iterative process stops when solving a saddle-point problem related to this equation. Read the description of the structure cs_param_saddle_t for more details.
|
CS_EQKEY_SADDLE_SCHUR_APPROX | Choice of the approximation of the Schur complement when solving a saddle-point problem. Not useful when there is no saddle-point problem to be solved. Some solvers or choices of preconditioner/solver may not need an approximation of the schur complement (ALU or GKB for instance). Please refer to cs_param_saddle_schur_approx_t for more details. |
CS_EQKEY_SADDLE_SOLVER | Set the strategy/solver used to solve a saddle-point system. |
CS_EQKEY_SADDLE_SOLVER_CLASS | Strategy/solver used to solve a saddle-point system. |
CS_EQKEY_SADDLE_SOLVER_RESTART | Number of iterations before restarting a Krylov solver used as the main solver for a saddle-point problem. This key is useful when the main solver is a GMRES, a flexible GMRES or a GCR. Read the description of the structure cs_param_saddle_t for more details.
|
CS_EQKEY_SADDLE_VERBOSITY | Level of details displayed for the resolution of a saddle-point system
|
CS_EQKEY_SLES_VERBOSITY | Level of details written by the code for the resolution of the linear system
|
CS_EQKEY_SOLVER_FAMILY | |
CS_EQKEY_SPACE_SCHEME | Set the space discretization scheme. Available choices are:
|
CS_EQKEY_TIME_SCHEME | Set the scheme for the temporal discretization. Available choices are:
|
CS_EQKEY_TIME_THETA | Set a value betwwen 0 and 1 for the theta parameter when a "theta_scheme" is set for the time discretization. Only useful if CS_EQKEY_TIME_SCHEME is set to "theta_scheme".
|
CS_EQKEY_VERBOSITY | Set the level of details written by the code for an equation. The higher the more detailed information is displayed.
|
CS_EQKEY_N_KEYS |
enum cs_equation_type_t |
Type of equations managed by the solver.
void cs_equation_add_advection | ( | cs_equation_param_t * | eqp, |
cs_adv_field_t * | adv_field | ||
) |
Associate a new term related to the advection operator for the equation associated to the given cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | adv_field | pointer to a cs_adv_field_t structure |
void cs_equation_add_advection_scaling_property | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property | ||
) |
Associate a scaling property to the advection.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
cs_xdef_t * cs_equation_add_bc_by_analytic | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_analytic_func_t * | analytic, | ||
void * | input | ||
) |
Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | analytic | pointer to an analytic function defining the value |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | analytic | pointer to an analytic function defining the value |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_bc_by_array | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_flag_t | loc, | ||
cs_real_t * | array, | ||
bool | is_owner, | ||
bool | full_length | ||
) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the related boundary zone |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t struct. (true or false) |
[in] | full_length | if true, size of "array" should be allocated to the total numbers of entities related to the given location. If false, a new list is allocated and filled with the related subset indirection. |
cs_xdef_t * cs_equation_add_bc_by_dof_func | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_flag_t | loc_flag, | ||
cs_dof_func_t * | func, | ||
void * | input | ||
) |
Define and initialize a new structure to set a boundary condition related to the given cs_equation_param_t structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | loc_flag | location where values are computed |
[in] | func | pointer to cs_dof_func_t function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | loc_flag | location where values are computed |
[in] | func | pointer to cs_dof_func_t function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_bc_by_field | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_field_t * | field | ||
) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the related boundary zone |
[in] | field | pointer to a cs_field_t structure |
cs_xdef_t * cs_equation_add_bc_by_time_func | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_time_func_t * | t_func, | ||
void * | input | ||
) |
Define and initialize a new structure to set a boundary condition related to the given equation param structure ml_name corresponds to the name of a pre-existing cs_mesh_location_t Definition relying on a cs_time_func_t function pointer.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | t_func | pointer to an analytic function defining the value |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | t_func | pointer to an analytic function defining the value |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_bc_by_value | ( | cs_equation_param_t * | eqp, |
const cs_param_bc_type_t | bc_type, | ||
const char * | z_name, | ||
cs_real_t * | values | ||
) |
Define and initialize a new structure to set a boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the related boundary zone |
[in] | values | pointer to a array storing the values |
cs_enforcement_param_t * cs_equation_add_cell_enforcement | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_cells, | ||
const cs_lnum_t | cell_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | cell_values[] | ||
) |
Add an enforcement of the value related to the degrees of freedom associated to the list of selected cells.
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_cells | number of selected cells |
[in] | cell_ids | list of cell ids |
[in] | ref_value | ignored if null |
[in] | cell_values | list of associated values, ignored if null |
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_cells | number of selected cells |
[in] | cell_ids | list of cell ids |
[in] | ref_value | ignored if nullptr |
[in] | cell_values | list of associated values, ignored if nullptr |
void cs_equation_add_curlcurl | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property, | ||
int | inversion | ||
) |
Associate a new term related to the curl-curl operator for the equation associated to the given cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
[in] | inversion | 1: true, false otherwise |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
[in] | inversion | > 0: true, false otherwise |
void cs_equation_add_diffusion | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property | ||
) |
Associate a new term related to the Laplacian operator for the equation associated to the given cs_equation_param_t structure Laplacian means div-grad (either for vector-valued or scalar-valued equations)
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
cs_enforcement_param_t * cs_equation_add_edge_dof_enforcement | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_edges, | ||
const cs_lnum_t | edge_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | edge_values[] | ||
) |
Add an enforcement of the value of degrees of freedom located at the mesh edges. The spatial discretization scheme for the given equation has to be CDO edge-based schemes.
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_edges | number of edges to enforce |
[in] | edge_ids | list of edges |
[in] | ref_value | default values or ignored (may be null) |
[in] | edge_values | list of associated values, ignored if null |
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_edges | number of edges to enforce |
[in] | edge_ids | list of edges |
[in] | ref_value | default values or ignored (may be nullptr) |
[in] | edge_values | list of associated values, ignored if nullptr |
cs_enforcement_param_t * cs_equation_add_face_dof_enforcement | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_faces, | ||
const cs_lnum_t | face_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | face_values[] | ||
) |
Add an enforcement of the value of degrees of freedom located at the mesh faces. The spatial discretization scheme for the given equation has to be CDO face-based schemes.
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_faces | number of faces to enforce |
[in] | face_ids | list of faces |
[in] | ref_value | default values or ignored (may be null) |
[in] | face_values | list of associated values, ignored if null |
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_faces | number of faces to enforce |
[in] | face_ids | list of faces |
[in] | ref_value | default values or ignored (may be nullptr) |
[in] | face_values | list of associated values, ignored if nullptr |
void cs_equation_add_graddiv | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property | ||
) |
Associate a new term related to the grad-div operator for the equation associated to the given cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
cs_xdef_t * cs_equation_add_ic_by_analytic | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_analytic_func_t * | analytic, | ||
void * | input | ||
) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here the initial value for each unknown associated to the zone with name z_name is set according to an analytical function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | analytic | pointer to an analytic function |
[in] | input | null or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_ic_by_dof_func | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_flag_t | loc_flag, | ||
cs_dof_func_t * | func, | ||
void * | input | ||
) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Case of a definition by a DoF function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | loc_flag | where information is computed |
[in] | func | pointer to a DoF function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | loc_flag | where information is computed |
[in] | func | pointer to a DoF function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_ic_by_qov | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
double | quantity | ||
) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here the value set to each unknown belonging to the given zone with name z_name is such that the integral over the cells of the zone returns the requested quantity.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" all cells are considered) |
[in] | quantity | quantity to distribute over the mesh location |
cs_xdef_t * cs_equation_add_ic_by_value | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_real_t * | val | ||
) |
Define the initial condition for the unknown related to this equation. This definition applies to a volume zone. By default, the unknown is set to zero everywhere. Here a constant value is set to all the unknows belonging to the given zone with name z_name.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" all cells are considered) |
[in] | val | pointer to the value |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" all cells are considered) |
[in] | val | pointer to the value |
cs_enforcement_param_t * cs_equation_add_or_replace_cell_enforcement | ( | cs_equation_param_t * | eqp, |
int | enforcement_id, | ||
cs_lnum_t | n_cells, | ||
const cs_lnum_t | cell_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | cell_values[] | ||
) |
Add a new enforcement if enforcement_id does not exist or replace it otherwise. Enforcement of the value related to the degrees of freedom associated to the list of selected cells.
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | enforcement_id | id of the enforcement to handle |
[in] | n_cells | number of selected cells |
[in] | cell_ids | list of cell ids |
[in] | ref_value | ignored if null |
[in] | cell_values | list of associated values, ignored if null |
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | enforcement_id | id of the enforcement to handle |
[in] | n_cells | number of selected cells |
[in] | cell_ids | list of cell ids |
[in] | ref_value | ignored if nullptr |
[in] | cell_values | list of associated values, ignored if nullptr |
int cs_equation_add_reaction | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property | ||
) |
Associate a new term related to the reaction operator for the equation associated to the given cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
void cs_equation_add_sliding_condition | ( | cs_equation_param_t * | eqp, |
const char * | z_name | ||
) |
Define and initialize a new structure to set a sliding boundary condition related to the given equation structure z_name corresponds to the name of a pre-existing cs_zone_t.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the related boundary zone |
cs_xdef_t * cs_equation_add_source_term_by_analytic | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_analytic_func_t * | func, | ||
void * | input | ||
) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an analytical function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | func | pointer to an analytical function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | func | pointer to an analytical function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_source_term_by_array | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_flag_t | loc, | ||
cs_real_t * | array, | ||
bool | is_owner, | ||
bool | full_length | ||
) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by an array.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t struct. (true or false) |
[in] | full_length | if true, the size of "array" should be allocated to the total numbers of entities related to the given location. If false, a new list is allocated and filled with the related subset indirection. |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t struct. (true or false) |
[in] | full_length | if true, the size of "array" should be allocated to the total numbers of entities related to the given location. If false, a new list is allocated and filled with the related subset indirection. |
cs_xdef_t * cs_equation_add_source_term_by_dof_func | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_flag_t | loc_flag, | ||
cs_dof_func_t * | func, | ||
void * | input | ||
) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a DoF function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | loc_flag | location of element ids given as parameter |
[in] | func | pointer to a DoF function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | loc_flag | location of element ids given as parameter |
[in] | func | pointer to a DoF function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_source_term_by_val | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_real_t * | val | ||
) |
Add a new source term by initializing a cs_xdef_t structure. Case of a definition by a constant value.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" all cells are considered) |
[in] | val | pointer to the value |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" all cells are considered) |
[in] | val | pointer to the value |
void cs_equation_add_time | ( | cs_equation_param_t * | eqp, |
cs_property_t * | property | ||
) |
Associate a new term related to the time derivative operator for the equation associated to the given cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | property | pointer to a cs_property_t structure |
cs_enforcement_param_t * cs_equation_add_vertex_dof_enforcement | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_vertices, | ||
const cs_lnum_t | vertex_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | vtx_values[] | ||
) |
Add an enforcement of the value of degrees of freedom located at the mesh vertices. The spatial discretization scheme for the given equation has to be CDO vertex-based or CDO vertex+cell-based schemes.
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_vertices | number of vertices to enforce |
[in] | vertex_ids | list of vertices |
[in] | ref_value | default values or ignored (may be null) |
[in] | vtx_values | list of associated values, ignored if null |
One assumes that values are interlaced if eqp->dim > 1 ref_value or elt_values has to be defined. If both parameters are defined, one keeps the definition in elt_values
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | n_vertices | number of vertices to enforce |
[in] | vertex_ids | list of vertices |
[in] | ref_value | default values or ignored (may be nullptr) |
[in] | vtx_values | list of associated values, ignored if nullptr |
cs_xdef_t * cs_equation_add_volume_mass_injection_by_analytic | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_analytic_func_t * | func, | ||
void * | input | ||
) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an analytical function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | func | pointer to an analytical function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | func | pointer to an analytical function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_volume_mass_injection_by_array | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_flag_t | loc_flag, | ||
cs_real_t * | array, | ||
bool | is_owner, | ||
bool | full_length | ||
) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using an array.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | bc_type | type of boundary condition to add |
[in] | z_name | name of the related boundary zone |
[in] | loc | information to know where are located values |
[in] | array | pointer to an array |
[in] | is_owner | transfer the lifecycle to the cs_xdef_t struct. (true or false) |
[in] | full_length | if true, size of "array" should be allocated to the total numbers of entities related to the given location. If false, a new list is allocated and filled with the related subset indirection. |
cs_xdef_t * cs_equation_add_volume_mass_injection_by_dof_func | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
cs_flag_t | loc_flag, | ||
cs_dof_func_t * | func, | ||
void * | input | ||
) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a DoF function.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | loc_flag | where information is computed |
[in] | func | pointer to an analytical function |
[in] | input | null or pointer to a structure cast on-the-fly |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | loc_flag | where information is computed |
[in] | func | pointer to an analytical function |
[in] | input | nullptr or pointer to a structure cast on-the-fly |
cs_xdef_t * cs_equation_add_volume_mass_injection_by_qov | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
double * | quantity | ||
) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant quantity distributed over the associated zone's volume.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | quantity | pointer to quantity to distribute over the zone |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | quantity | pointer to quantity to distribute over the zone |
cs_xdef_t * cs_equation_add_volume_mass_injection_by_value | ( | cs_equation_param_t * | eqp, |
const char * | z_name, | ||
double * | val | ||
) |
Add a new volume mass injection definition source term by initializing a cs_xdef_t structure, using a constant value.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in] | val | pointer to the value |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
[in] | val | pointer to the value |
void cs_equation_add_xdef_bc | ( | cs_equation_param_t * | eqp, |
cs_xdef_t * | xdef | ||
) |
Set a boundary condition from an existing cs_xdef_t structure The lifecycle of the cs_xdef_t structure is now managed by the current cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | xdef | pointer to the cs_xdef_t structure to transfer |
|
inlinestatic |
Copy the settings from one cs_equation_param_t structure to another one. The name is not copied.
[in] | ref | pointer to the reference cs_equation_param_t |
[in,out] | dst | pointer to the cs_equation_param_t to update |
[in] | copy_fid | copy also the field id or not |
|
inlinestatic |
Create a cs_equation_param_t structure.
[in] | name | name of the equation |
[in] | type | type of equation |
[in] | dim | dim of the variable associated to this equation |
[in] | default_bc | type of boundary condition set by default |
cs_xdef_t * cs_equation_find_bc | ( | cs_equation_param_t * | eqp, |
const char * | z_name | ||
) |
Return pointer to existing boundary condition definition structure for the given equation param structure and zone.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
void cs_equation_param_clear | ( | cs_equation_param_t * | eqp | ) |
Free the contents of a cs_equation_param_t.
The cs_equation_param_t structure itself is not freed, but all its sub-structures are freed.
This is useful for equation parameters which are accessed through field keywords.
[in,out] | eqp | pointer to a cs_equation_param_t |
void cs_equation_param_copy_bc | ( | const cs_equation_param_t * | ref, |
cs_equation_param_t * | dst | ||
) |
Copy only the part dedicated to the boundary conditions and the DoF (degrees of freedom) enforcement from one cs_equation_param_t structure to another one.
[in] | ref | pointer to the reference cs_equation_param_t |
[in,out] | dst | pointer to the cs_equation_param_t to update |
void cs_equation_param_copy_from | ( | const cs_equation_param_t * | ref, |
cs_equation_param_t * | dst, | ||
bool | copy_fld_id | ||
) |
Copy the settings from one cs_equation_param_t structure to another one. The name is not copied.
[in] | ref | pointer to the reference cs_equation_param_t |
[in,out] | dst | pointer to the cs_equation_param_t to update |
[in] | copy_fld_id | copy also the field id or not |
cs_equation_param_t * cs_equation_param_create | ( | const char * | name, |
cs_equation_type_t | type, | ||
int | dim, | ||
cs_param_bc_type_t | default_bc | ||
) |
Create a cs_equation_param_t structure.
[in] | name | name of the equation |
[in] | type | type of equation |
[in] | dim | dim of the variable associated to this equation |
[in] | default_bc | type of boundary condition set by default |
Create a cs_equation_param_t structure.
[in] | name | name of the equation |
[in] | type | type of equation |
[in] | dim | dim of the variable associated to this equation |
[in] | default_bc | type of boundary condition set by default |
void cs_equation_param_ensure_consistent_settings | ( | cs_equation_param_t * | eqp | ) |
At this stage, the numerical settings should not be modified anymore by the user. One makes a last set of modifications if needed to ensure a consistent numerical settings.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
cs_equation_param_t * cs_equation_param_free | ( | cs_equation_param_t * | eqp | ) |
Free a cs_equation_param_t.
[in] | eqp | pointer to a cs_equation_param_t |
[in,out] | eqp | pointer to a cs_equation_param_t |
cs_param_saddle_t * cs_equation_param_get_saddle_param | ( | cs_equation_param_t * | eqp | ) |
Get the pointer to the set of parameters to handle the SLES solver associated to this set of equation parameters.
[in] | eqp | pointer to a cs_equation_param_t structure |
Get the pointer to the set of parameters to handle the SLES solver associated to this set of equation parameters.
[in] | eqp | pointer to a cs_equation_param_t structure |
cs_param_sles_t * cs_equation_param_get_sles_param | ( | cs_equation_param_t * | eqp | ) |
Get the pointer to the set of parameters to handle the SLES solver associated to this set of equation parameters.
[in] | eqp | pointer to a cs_equation_param_t structure |
|
inlinestatic |
Ask if the parameters of the equation needs a convection term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation needs a curl-curl term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation needs a diffusion term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation needs a grad-div term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation induces an implicit treatment of the advection term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation has an internal enforcement of the degrees of freedom.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Check if a cs_equation_param_t structure has its name member equal to the given name.
[in] | eqp | pointer to a cs_equation_param_t structure |
[in] | name | name of the equation |
|
inlinestatic |
Ask if the parameters of the equation needs a reaction term.
[in] | eqp | pointer to a cs_equation_param_t |
bool cs_equation_param_has_robin_bc | ( | const cs_equation_param_t * | eqp | ) |
Ask if the parameter settings of the equation has requested the treatment of Robin boundary conditions.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation needs a source term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation needs an unsteady term.
[in] | eqp | pointer to a cs_equation_param_t |
|
inlinestatic |
Ask if the parameters of the equation has activated a user hook to get a fine tuning of the cellwise system building.
[in] | eqp | pointer to a cs_equation_param_t |
void cs_equation_param_lock_settings | ( | cs_equation_param_t * | eqp | ) |
Lock settings to prevent from unwanted modifications.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
void cs_equation_param_log | ( | const cs_equation_param_t * | eqp | ) |
Print the detail of a cs_equation_param_t structure.
[in] | eqp | pointer to a cs_equation_param_t structure |
void cs_equation_param_set | ( | cs_equation_param_t * | eqp, |
cs_equation_key_t | key, | ||
const char * | keyval | ||
) |
Set a parameter attached to a keyname in a cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | key | key related to the member of eq to set |
[in] | keyval | accessor to the value to set |
|
inlinestatic |
Update the flag related to a cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t |
[in] | flag | flag to add |
void cs_equation_param_set_quadrature_to_all | ( | cs_equation_param_t * | eqp, |
cs_quadrature_type_t | qtype | ||
) |
Apply the given quadrature rule to all existing definitions under the cs_equation_param_t structure. To get a more detailed control of the quadrature rule, please consider the function cs_xdef_set_quadrature.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | qtype | type of quadrature to apply |
Apply the given quadrature rule to all existing definitions under the cs_equation_param_t structure. To get a more detailed control of the quadrature rule, please consider the function cs_xdef_set_quadrature.
If the default quadrature has been modified by the code, this function set the value of the quadrture to the given parameter whatever was the previous value.
To get a more detailed control of the quadrature rule, please consider the function cs_xdef_set_quadrature
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | qtype | type of quadrature to apply |
void cs_equation_param_set_sles | ( | cs_equation_param_t * | eqp | ) |
Set parameters for initializing SLES structures used for the resolution of the linear system. Settings are related to this equation.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
void cs_equation_param_unlock_settings | ( | cs_equation_param_t * | eqp | ) |
Unlock settings. Be sure that is really wanted (inconsistency between the setup logging and what is used may appear)
[in,out] | eqp | pointer to a cs_equation_param_t structure |
void cs_equation_remove_bc | ( | cs_equation_param_t * | eqp, |
const char * | z_name | ||
) |
Remove boundary condition from the given equation param structure for a given zone.
If no matching boundary condition is found, the function returns silently.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
If no matching boundary condition is found, the function returns silently.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
void cs_equation_remove_ic | ( | cs_equation_param_t * | eqp, |
const char * | z_name | ||
) |
Remove initial condition from the given equation param structure for a given zone.
If no matching boundary condition is found, the function returns silently.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if null or "" if all cells are considered) |
If no matching boundary condition is found, the function returns silently.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | z_name | name of the associated zone (if nullptr or "" if all cells are considered) |
|
inlinestatic |
Set a parameter attached to a keyname in a cs_equation_param_t structure.
[in,out] | eqp | pointer to a cs_equation_param_t structure |
[in] | key | key related to the member of eq to set |
[in] | keyval | accessor to the value to set |