#include "cs_defs.h"
#include <assert.h>
#include <ctype.h>
#include <stdlib.h>
#include <string.h>
#include <bft_error.h>
#include <bft_mem.h>
#include "cs_boundary_zone.h"
#include "cs_cdo_bc.h"
#include "cs_fp_exception.h"
#include "cs_hodge.h"
#include "cs_log.h"
#include "cs_mesh_location.h"
#include "cs_multigrid.h"
#include "cs_sles.h"
#include "cs_source_term.h"
#include "cs_volume_zone.h"
#include "cs_equation_param.h"
Functions | |
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 with the given parameters. The remaining parameters are set with default values;. More... | |
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_clear_param (cs_equation_param_t *eqp) |
Free the contents of a cs_equation_param_t. More... | |
cs_equation_param_t * | cs_equation_free_param (cs_equation_param_t *eqp) |
Free a cs_equation_param_t. More... | |
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... | |
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_last_stage (cs_equation_param_t *eqp) |
Last modification of the cs_equation_param_t structure before launching the computation. 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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here a constant value is set to all the entities belonging to the given mesh location. 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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the value related to all the entities belonging to the given mesh location is such that the integral over these cells 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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function. More... | |
void | cs_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, cs_lnum_t *index) |
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_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_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_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, cs_lnum_t *index) |
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... | |
void | cs_equation_enforce_vertex_dofs (cs_equation_param_t *eqp, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_value[], const cs_real_t elt_values[]) |
Add an enforcement of the value of degrees of freedom located at mesh vertices. The spatial discretization scheme for the given equation has to be CDO-Vertex based or CDO-Vertex+Cell-based schemes. More... | |
void | cs_equation_enforce_value_on_cell_selection (cs_equation_param_t *eqp, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], const cs_real_t ref_value[], const cs_real_t elt_values[]) |
Add an enforcement of the value related to the degrees of freedom associated to the list of selected cells. 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.
[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 |
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, | ||
cs_lnum_t * | index | ||
) |
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 structure (true or false) |
[in] | index | optional pointer to the array index |
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 |
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 |
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 | > 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 |
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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the initial value is set according to an analytical function.
[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_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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here the value related to all the entities belonging to the given mesh location is such that the integral over these cells 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 can be done on a specified mesh location. By default, the unknown is set to zero everywhere. Here a constant value is set to all the entities belonging to the given mesh location.
[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 |
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 |
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, | ||
cs_lnum_t * | index | ||
) |
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 structure (true or false) |
[in] | index | optional pointer to the array index |
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 |
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 |
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_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 |
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 |
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 |
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 |
void cs_equation_clear_param | ( | 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_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.
[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 |
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 with the given parameters. The remaining parameters are set with default values;.
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_enforce_value_on_cell_selection | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_elts, | ||
const cs_lnum_t | elt_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | elt_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_elts | number of selected cells |
[in] | elt_ids | list of cell ids |
[in] | ref_value | ignored if NULL |
[in] | elt_values | list of associated values, ignored if NULL |
void cs_equation_enforce_vertex_dofs | ( | cs_equation_param_t * | eqp, |
cs_lnum_t | n_elts, | ||
const cs_lnum_t | elt_ids[], | ||
const cs_real_t | ref_value[], | ||
const cs_real_t | elt_values[] | ||
) |
Add an enforcement of the value of degrees of freedom located at 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_elts | number of vertices to enforce |
[in] | elt_ids | list of vertices |
[in] | ref_value | ignored if NULL |
[in] | elt_values | list of associated values, ignored if NULL |
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) |
cs_equation_param_t* cs_equation_free_param | ( | cs_equation_param_t * | eqp | ) |
Free a cs_equation_param_t.
[in,out] | 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 |
void cs_equation_param_last_stage | ( | cs_equation_param_t * | eqp | ) |
Last modification of the cs_equation_param_t structure before launching the computation.
[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_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_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) |
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.
[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 |