My Project
programmer's documentation
Functions | Variables
cs_equation_param.c File Reference
#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_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"
Include dependency graph for cs_equation_param.c:

Functions

static void _set_key (const char *label, 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_equation_param_tcs_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_param_update_from (const cs_equation_param_t *ref, cs_equation_param_t *dst)
 Copy the settings from one cs_equation_param_t structure to another one. More...
 
cs_equation_param_tcs_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, int field_id)
 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_summary_param (const cs_equation_param_t *eqp)
 Summary of a cs_equation_param_t structure. More...
 
cs_xdef_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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_tcs_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...
 
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)
 Define and initialize a new structure to store parameters related to a diffusion term. More...
 
void cs_equation_add_time (cs_equation_param_t *eqp, cs_property_t *property)
 Define and initialize a new structure to store parameters related to an unsteady term. More...
 
void cs_equation_add_advection (cs_equation_param_t *eqp, cs_adv_field_t *adv_field)
 Define and initialize a new structure to store parameters related to an advection term. More...
 
int cs_equation_add_reaction (cs_equation_param_t *eqp, cs_property_t *property)
 Define and initialize a new structure to store parameters related to a reaction term. More...
 
cs_xdef_tcs_equation_add_source_term_by_val (cs_equation_param_t *eqp, const char *z_name, cs_real_t *val)
 Define a new source term structure and initialize it by value. More...
 
cs_xdef_tcs_equation_add_source_term_by_analytic (cs_equation_param_t *eqp, const char *z_name, cs_analytic_func_t *ana, void *input)
 Define a new source term structure and initialize it by an analytical function. More...
 
cs_xdef_tcs_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)
 Define a new source term defined by an array. 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 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. We assume that values are interlaced (if eqp->dim > 1) More...
 

Variables

static const cs_real_t _weak_pena_bc_coef_by_default = 100.
 
static const cs_real_t _strong_pena_bc_coef_by_default = 1e12
 
static const char _err_empty_eqp []
 

Function Documentation

◆ _set_key()

static void _set_key ( const char *  label,
cs_equation_param_t eqp,
cs_equation_key_t  key,
const char *  keyval 
)
static

Set a parameter attached to a keyname in a cs_equation_param_t structure.

Parameters
[in]labellabel to identify the error message
[in,out]eqppointer to a cs_equation_param_t structure
[in]keykey related to the member of eq to set
[in]keyvalaccessor to the value to set

◆ cs_equation_add_advection()

void cs_equation_add_advection ( cs_equation_param_t eqp,
cs_adv_field_t adv_field 
)

Define and initialize a new structure to store parameters related to an advection term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]adv_fieldpointer to a cs_adv_field_t structure

◆ cs_equation_add_bc_by_analytic()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function defining the value
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_bc_by_array()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the related boundary zone
[in]locinformation to know where are located values
[in]arraypointer to an array
[in]is_ownertransfer the lifecycle to the cs_xdef_t structure (true or false)
[in]indexoptional pointer to the array index
Returns
a pointer to the new allocated cs_xdef_t structure

◆ cs_equation_add_bc_by_value()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]bc_typetype of boundary condition to add
[in]z_namename of the related boundary zone
[in]valuespointer to a array storing the values
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_diffusion()

void cs_equation_add_diffusion ( cs_equation_param_t eqp,
cs_property_t property 
)

Define and initialize a new structure to store parameters related to a diffusion term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_ic_by_analytic()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]analyticpointer to an analytic function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_ic_by_qov()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]quantityquantity to distribute over the mesh location
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_ic_by_value()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_reaction()

int cs_equation_add_reaction ( cs_equation_param_t eqp,
cs_property_t property 
)

Define and initialize a new structure to store parameters related to a reaction term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure
Returns
the id related to the reaction term

◆ cs_equation_add_sliding_condition()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the related boundary zone

◆ cs_equation_add_source_term_by_analytic()

cs_xdef_t* cs_equation_add_source_term_by_analytic ( cs_equation_param_t eqp,
const char *  z_name,
cs_analytic_func_t ana,
void *  input 
)

Define a new source term structure and initialize it by an analytical function.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]anapointer to an analytical function
[in]inputNULL or pointer to a structure cast on-the-fly
Returns
a pointer to the new cs_source_term_t structure

◆ cs_equation_add_source_term_by_array()

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 
)

Define a new source term defined by an array.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" if all cells are considered)
[in]locinformation to know where are located values
[in]arraypointer to an array
[in]is_ownertransfer the lifecycle to the cs_xdef_t structure (true or false)
[in]indexoptional pointer to the array index
Returns
a pointer to the new cs_source_term_t structure

◆ cs_equation_add_source_term_by_val()

cs_xdef_t* cs_equation_add_source_term_by_val ( cs_equation_param_t eqp,
const char *  z_name,
cs_real_t val 
)

Define a new source term structure and initialize it by value.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]z_namename of the associated zone (if NULL or "" all cells are considered)
[in]valpointer to the value
Returns
a pointer to the new cs_xdef_t structure

◆ cs_equation_add_time()

void cs_equation_add_time ( cs_equation_param_t eqp,
cs_property_t property 
)

Define and initialize a new structure to store parameters related to an unsteady term.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]propertypointer to a cs_property_t structure

◆ cs_equation_add_xdef_bc()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]xdefpointer to the cs_xdef_t structure to transfer

◆ cs_equation_create_param()

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.

Parameters
[in]namename of the equation
[in]typetype of equation
[in]dimdim of the variable associated to this equation
[in]default_bctype of boundary condition set by default
Returns
a pointer to a new allocated cs_equation_param_t structure

◆ cs_equation_enforce_vertex_dofs()

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  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. We assume that values are interlaced (if eqp->dim > 1)

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]n_eltsnumber of vertices to enforce
[in]elt_idslist of vertices
[in]elt_valueslist of associated values

◆ cs_equation_free_param()

cs_equation_param_t* cs_equation_free_param ( cs_equation_param_t eqp)

Free a cs_equation_param_t.

Parameters
[in,out]eqppointer to a cs_equation_param_t
Returns
a NULL pointer

◆ cs_equation_param_last_stage()

void cs_equation_param_last_stage ( cs_equation_param_t eqp)

Last modification of the cs_equation_param_t structure before launching the computation.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure

◆ cs_equation_param_set_sles()

void cs_equation_param_set_sles ( cs_equation_param_t eqp,
int  field_id 
)

Set parameters for initializing SLES structures used for the resolution of the linear system. Settings are related to this equation.

Parameters
[in]eqppointer to a cs_equation_param_t struct.
[in]field_idid of the cs_field_t struct. for this equation

◆ cs_equation_param_update_from()

void cs_equation_param_update_from ( const cs_equation_param_t ref,
cs_equation_param_t dst 
)

Copy the settings from one cs_equation_param_t structure to another one.

Parameters
[in]refpointer to the reference cs_equation_param_t
[in,out]dstpointer to the cs_equation_param_t to update

◆ cs_equation_set_param()

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.

Parameters
[in,out]eqppointer to a cs_equation_param_t structure
[in]keykey related to the member of eq to set
[in]keyvalaccessor to the value to set

◆ cs_equation_summary_param()

void cs_equation_summary_param ( const cs_equation_param_t eqp)

Summary of a cs_equation_param_t structure.

Parameters
[in]eqppointer to a cs_equation_param_t structure

Variable Documentation

◆ _err_empty_eqp

const char _err_empty_eqp[]
static
Initial value:
=
N_(" Stop setting an empty cs_equation_param_t structure.\n"
" Please check your settings.\n")

◆ _strong_pena_bc_coef_by_default

const cs_real_t _strong_pena_bc_coef_by_default = 1e12
static

◆ _weak_pena_bc_coef_by_default

const cs_real_t _weak_pena_bc_coef_by_default = 100.
static
N_
#define N_(String)
Definition: cs_defs.h:56