7.2
general documentation
cs_solidification.h File Reference
#include "cs_base.h"
#include "cs_boundary.h"
#include "cs_equation.h"
#include "cs_iter_algo.h"
#include "cs_navsto_param.h"
#include "cs_time_plot.h"
#include "cs_time_step.h"
+ Include dependency graph for cs_solidification.h:

Go to the source code of this file.

Data Structures

struct  cs_solidification_stefan_t
 
struct  cs_solidification_voller_t
 
struct  cs_solidification_binary_alloy_t
 
struct  cs_solidification_t
 

Macros

Flags specifying automatic post-processing for the solidification

module

#define CS_SOLIDIFICATION_POST_CELL_STATE   (1 << 0) /* = 1 */
 State related to each cell between (solid, mushy, liquid or eutectic) More...
 
#define CS_SOLIDIFICATION_POST_ENTHALPY   (1 << 1) /* = 2 */
 Enthalpy in each cell. By default, only the temperature is post-processed. More...
 
#define CS_SOLIDIFICATION_POST_CBULK_ADIM   (1 << 2) /* = 4 */
 Compute and post-process (C_bulk - C_0)/C_0 Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated C_0 is the reference concentration. More...
 
#define CS_SOLIDIFICATION_POST_CLIQ   (1 << 3) /* = 8 */
 Post-process Cliq the liquid solute distribution (wt %) Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated. More...
 
#define CS_SOLIDIFICATION_POST_LIQUIDUS_TEMPERATURE   (1 << 4) /* = 16 */
 Activate the (volumic) post-processing of the liquidus temperature in each cell. More...
 
#define CS_SOLIDIFICATION_POST_SEGREGATION_INDEX   (1 << 5) /* = 32 */
 Activate the computation and output in the file solidification.dat for each time step of the segregation index defined by sqrt( 1/|Domaine| * {Domain} ((C_bulk - C_0)/C_0)**2 ) Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated. More...
 
#define CS_SOLIDIFICATION_POST_SOLIDIFICATION_RATE   (1 << 6) /* = 64 */
 Activate the computation and output in the file solidification.dat for each time step of the integral over the computational domain of the solid fraction divided by the volume of the domain. More...
 
#define CS_SOLIDIFICATION_ADVANCED_ANALYSIS   (1 << 7) /* = 128 */
 Activate a set of post-processing (Advanced usage. Only for the understanding of the solidification process) More...
 
Flags specifying options specific to the solidification module
#define CS_SOLIDIFICATION_USE_ENTHALPY_VARIABLE   (1 << 0) /*= 1 */
 The dynamic system of equations is associated with an energy equation solved using the enthalpy as variable (not fully available). More...
 
#define CS_SOLIDIFICATION_NO_VELOCITY_FIELD   (1 << 1) /*= 2 */
 The system of equations does not involve the Navier-Stokes equations. No velocity is taken into account. More...
 
#define CS_SOLIDIFICATION_WITH_SOLUTE_SOURCE_TERM   (1 << 2) /*= 4 */
 The solute equation related to the transport of the bulk concentration is treated with a source term related to an explicit advection of the quantity (C - Cl). The default behavior is to add a weighting coefficient to the (implicit) advection term related to the liquid fraction This option is related to the LEGACY strategy. More...
 
#define CS_SOLIDIFICATION_USE_EXTRAPOLATION   (1 << 3) /*= 8 */
 Use an extrapolation during the computation of different terms according to the strategy. This extrapolation of variable at time step n+1 uses values at n and n-1: 2*u^n - u^{n-1}. More...
 
#define CS_SOLIDIFICATION_WITH_PENALIZED_EUTECTIC   (1 << 4) /*= 16 */
 Option related to the PATH strategy. Introduce a reaction term and a source term in order to remain on the eutectic plateau. More...
 
#define CS_SOLIDIFICATION_BINARY_ALLOY_M_FUNC   (1 << 7) /*= 128 */
 
#define CS_SOLIDIFICATION_BINARY_ALLOY_C_FUNC   (1 << 8) /*= 256 */
 
#define CS_SOLIDIFICATION_BINARY_ALLOY_G_FUNC   (1 << 9) /*= 512 */
 
#define CS_SOLIDIFICATION_BINARY_ALLOY_T_FUNC   (1 <<10) /*= 1024 */
 

Typedefs

typedef void() cs_solidification_func_t(const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts)
 Function pointer associated to a solidification model aiming at updating/initializing the solidification variables/properties dedicated to the model. More...
 

Enumerations

enum  cs_solidification_model_t {
  CS_SOLIDIFICATION_MODEL_STEFAN, CS_SOLIDIFICATION_MODEL_VOLLER_PRAKASH_87, CS_SOLIDIFICATION_MODEL_VOLLER_NL, CS_SOLIDIFICATION_MODEL_BINARY_ALLOY,
  CS_SOLIDIFICATION_N_MODELS
}
 Type of physical model used to simulate the solidifcation/fusion process. More...
 
enum  cs_solidification_state_t {
  CS_SOLIDIFICATION_STATE_SOLID = 0, CS_SOLIDIFICATION_STATE_MUSHY = 1, CS_SOLIDIFICATION_STATE_LIQUID = 2, CS_SOLIDIFICATION_STATE_EUTECTIC = 3,
  CS_SOLIDIFICATION_N_STATES = 4
}
 Kind of state in which a cell or an entity is. More...
 
enum  cs_solidification_strategy_t { CS_SOLIDIFICATION_STRATEGY_LEGACY, CS_SOLIDIFICATION_STRATEGY_TAYLOR, CS_SOLIDIFICATION_STRATEGY_PATH, CS_SOLIDIFICATION_N_STRATEGIES }
 Kind of strategy to use to model the segregation/solidification process. This implies a setting of functions to update the liquid fraction, the thermal source terms, the liquid concentration and its related quantities. More...
 

Functions

bool cs_solidification_is_activated (void)
 Test if solidification module is activated. More...
 
cs_solidification_tcs_solidification_get_structure (void)
 Retrieve the main structure to deal with solidification process. More...
 
void cs_solidification_set_verbosity (int verbosity)
 Set the level of verbosity for the solidification module. More...
 
void cs_solidification_set_forcing_eps (cs_real_t forcing_eps)
 Set the value of the epsilon parameter used in the forcing term of the momentum equation. More...
 
cs_solidification_tcs_solidification_activate (cs_solidification_model_t model, cs_flag_t options, cs_flag_t post_flag, const cs_boundary_t *boundaries, cs_navsto_param_model_t ns_model, cs_navsto_param_model_flag_t ns_model_flag, cs_navsto_param_coupling_t algo_coupling, cs_navsto_param_post_flag_t ns_post_flag)
 Activate the solidification module. More...
 
cs_solidification_stefan_tcs_solidification_get_stefan_struct (void)
 Get the structure defining the Stefan model. More...
 
cs_solidification_stefan_tcs_solidification_check_stefan_model (void)
 Sanity checks on the consistency of the Stefan's model settings. More...
 
void cs_solidification_set_stefan_model (cs_real_t t_change, cs_real_t latent_heat)
 Set the main physical parameters which describe the Stefan model. More...
 
cs_solidification_voller_tcs_solidification_get_voller_struct (void)
 Get the structure defining the Voller model. More...
 
cs_solidification_voller_tcs_solidification_check_voller_model (void)
 Sanity checks on the consistency of the Voller's model settings. More...
 
void cs_solidification_set_voller_model (cs_real_t beta, cs_real_t t_ref, cs_real_t t_solidus, cs_real_t t_liquidus, cs_real_t latent_heat, cs_real_t s_das)
 Set the main physical parameters which describe the Voller and Prakash modelling. More...
 
void cs_solidification_set_voller_model_no_velocity (cs_real_t t_solidus, cs_real_t t_liquidus, cs_real_t latent_heat)
 Set the main physical parameters which describe the Voller and Prakash modelling. More...
 
cs_solidification_binary_alloy_tcs_solidification_get_binary_alloy_struct (void)
 Get the structure defining the binary alloy model. More...
 
cs_solidification_binary_alloy_tcs_solidification_check_binary_alloy_model (void)
 Sanity checks on the consistency of the settings of the binary alloy model. More...
 
void cs_solidification_set_binary_alloy_model (const char *name, const char *varname, cs_real_t beta_t, cs_real_t temp0, cs_real_t beta_c, cs_real_t conc0, cs_real_t kp, cs_real_t mliq, cs_real_t t_eutec, cs_real_t t_melt, cs_real_t solute_diff, cs_real_t latent_heat, cs_real_t s_das)
 Set the main physical parameters which describe a solidification process with a binary alloy (with components A and B) Add a transport equation for the solute concentration to simulate the conv/diffusion of the alloy ratio between the two components of the alloy. More...
 
void cs_solidification_set_strategy (cs_solidification_strategy_t strategy)
 Set the strategy to update quantitiess (liquid fraction and the thermal source term for the two main quantities) More...
 
void cs_solidification_set_segr_functions (cs_solidification_func_t *vel_forcing, cs_solidification_func_t *cliq_update, cs_solidification_func_t *gliq_update, cs_solidification_func_t *thm_st_update)
 Set the functions to perform the update of physical properties and/or the computation of the thermal source term or quantities and/or the way to perform the coupling between the thermal equation and the bulk concentration computation. All this setting defines the way to compute the solidification process of a binary alloy. If a function is set to NULL then the automatic settings are kept. More...
 
cs_solidification_tcs_solidification_destroy_all (void)
 Free the main structure related to the solidification module. More...
 
void cs_solidification_init_setup (void)
 Setup equations/properties related to the solidification module. More...
 
void cs_solidification_finalize_setup (const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant)
 Finalize the setup stage for equations related to the solidification module. More...
 
void cs_solidification_log_setup (void)
 Summarize the solidification module in the log file dedicated to the setup. More...
 
void cs_solidification_init_values (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Set an initial values for all quantities related to this module This is done after the setup step. More...
 
void cs_solidification_compute (const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *time_step)
 Solve equations related to the solidification module. More...
 
void cs_solidification_extra_op (const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts)
 Predefined extra-operations for the solidification module. More...
 
void cs_solidification_extra_post (void *input, int mesh_id, int cat_id, int ent_flag[5], cs_lnum_t n_cells, cs_lnum_t n_i_faces, cs_lnum_t n_b_faces, const cs_lnum_t cell_ids[], const cs_lnum_t i_face_ids[], const cs_lnum_t b_face_ids[], const cs_time_step_t *time_step)
 Predefined post-processing output for the solidification module. Prototype of this function is fixed since it is a function pointer defined in cs_post.h (cs_post_time_mesh_dep_output_t) More...
 

Macro Definition Documentation

◆ CS_SOLIDIFICATION_ADVANCED_ANALYSIS

#define CS_SOLIDIFICATION_ADVANCED_ANALYSIS   (1 << 7) /* = 128 */

Activate a set of post-processing (Advanced usage. Only for the understanding of the solidification process)

◆ CS_SOLIDIFICATION_BINARY_ALLOY_C_FUNC

#define CS_SOLIDIFICATION_BINARY_ALLOY_C_FUNC   (1 << 8) /*= 256 */

◆ CS_SOLIDIFICATION_BINARY_ALLOY_G_FUNC

#define CS_SOLIDIFICATION_BINARY_ALLOY_G_FUNC   (1 << 9) /*= 512 */

◆ CS_SOLIDIFICATION_BINARY_ALLOY_M_FUNC

#define CS_SOLIDIFICATION_BINARY_ALLOY_M_FUNC   (1 << 7) /*= 128 */

◆ CS_SOLIDIFICATION_BINARY_ALLOY_T_FUNC

#define CS_SOLIDIFICATION_BINARY_ALLOY_T_FUNC   (1 <<10) /*= 1024 */

◆ CS_SOLIDIFICATION_NO_VELOCITY_FIELD

#define CS_SOLIDIFICATION_NO_VELOCITY_FIELD   (1 << 1) /*= 2 */

The system of equations does not involve the Navier-Stokes equations. No velocity is taken into account.

◆ CS_SOLIDIFICATION_POST_CBULK_ADIM

#define CS_SOLIDIFICATION_POST_CBULK_ADIM   (1 << 2) /* = 4 */

Compute and post-process (C_bulk - C_0)/C_0 Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated C_0 is the reference concentration.

◆ CS_SOLIDIFICATION_POST_CELL_STATE

#define CS_SOLIDIFICATION_POST_CELL_STATE   (1 << 0) /* = 1 */

State related to each cell between (solid, mushy, liquid or eutectic)

◆ CS_SOLIDIFICATION_POST_CLIQ

#define CS_SOLIDIFICATION_POST_CLIQ   (1 << 3) /* = 8 */

Post-process Cliq the liquid solute distribution (wt %) Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated.

◆ CS_SOLIDIFICATION_POST_ENTHALPY

#define CS_SOLIDIFICATION_POST_ENTHALPY   (1 << 1) /* = 2 */

Enthalpy in each cell. By default, only the temperature is post-processed.

◆ CS_SOLIDIFICATION_POST_LIQUIDUS_TEMPERATURE

#define CS_SOLIDIFICATION_POST_LIQUIDUS_TEMPERATURE   (1 << 4) /* = 16 */

Activate the (volumic) post-processing of the liquidus temperature in each cell.

◆ CS_SOLIDIFICATION_POST_SEGREGATION_INDEX

#define CS_SOLIDIFICATION_POST_SEGREGATION_INDEX   (1 << 5) /* = 32 */

Activate the computation and output in the file solidification.dat for each time step of the segregation index defined by sqrt( 1/|Domaine| * {Domain} ((C_bulk - C_0)/C_0)**2 ) Only available if the model CS_SOLIDIFICATION_MODEL_BINARY_ALLOY is activated.

◆ CS_SOLIDIFICATION_POST_SOLIDIFICATION_RATE

#define CS_SOLIDIFICATION_POST_SOLIDIFICATION_RATE   (1 << 6) /* = 64 */

Activate the computation and output in the file solidification.dat for each time step of the integral over the computational domain of the solid fraction divided by the volume of the domain.

◆ CS_SOLIDIFICATION_USE_ENTHALPY_VARIABLE

#define CS_SOLIDIFICATION_USE_ENTHALPY_VARIABLE   (1 << 0) /*= 1 */

The dynamic system of equations is associated with an energy equation solved using the enthalpy as variable (not fully available).

◆ CS_SOLIDIFICATION_USE_EXTRAPOLATION

#define CS_SOLIDIFICATION_USE_EXTRAPOLATION   (1 << 3) /*= 8 */

Use an extrapolation during the computation of different terms according to the strategy. This extrapolation of variable at time step n+1 uses values at n and n-1: 2*u^n - u^{n-1}.

◆ CS_SOLIDIFICATION_WITH_PENALIZED_EUTECTIC

#define CS_SOLIDIFICATION_WITH_PENALIZED_EUTECTIC   (1 << 4) /*= 16 */

Option related to the PATH strategy. Introduce a reaction term and a source term in order to remain on the eutectic plateau.

◆ CS_SOLIDIFICATION_WITH_SOLUTE_SOURCE_TERM

#define CS_SOLIDIFICATION_WITH_SOLUTE_SOURCE_TERM   (1 << 2) /*= 4 */

The solute equation related to the transport of the bulk concentration is treated with a source term related to an explicit advection of the quantity (C - Cl). The default behavior is to add a weighting coefficient to the (implicit) advection term related to the liquid fraction This option is related to the LEGACY strategy.

Typedef Documentation

◆ cs_solidification_func_t

typedef void() cs_solidification_func_t(const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts)

Function pointer associated to a solidification model aiming at updating/initializing the solidification variables/properties dedicated to the model.

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

Enumeration Type Documentation

◆ cs_solidification_model_t

Type of physical model used to simulate the solidifcation/fusion process.

Enumerator
CS_SOLIDIFICATION_MODEL_STEFAN 

Phase change model without advection field. The phase change is assumed to be at a given temperature meaning that the liquid fraction is a step function w.r.t. the temperature.

CS_SOLIDIFICATION_MODEL_VOLLER_PRAKASH_87 

Modelling introduced in Voller and Prakash entitled: "A fixed grid numerical modelling methodology for convection-diffusion mushy region phase-change problems" Int. J. Heat Transfer, 30 (8), 1987. No tracer. Only physical constants describing the solidification process are used.

CS_SOLIDIFICATION_MODEL_VOLLER_NL 

Modelling based on CS_SOLIDIFICATION_MODEL_VOLLER_PRAKASH_87 but the thermal equation and the update of the liquid fraction is non-linear

CS_SOLIDIFICATION_MODEL_BINARY_ALLOY 

The basis is similar to CS_SOLIDIFICATION_MODEL_VOLLER_PRAKASH_87 A tracer equation is added corresponding to the evolution of the bulk concentration in alloy. This alloy has two chemical constituents hence the name "binary alloy".

CS_SOLIDIFICATION_N_MODELS 

◆ cs_solidification_state_t

Kind of state in which a cell or an entity is.

Enumerator
CS_SOLIDIFICATION_STATE_SOLID 

Solid state (the liquid fraction is equal to 0)

CS_SOLIDIFICATION_STATE_MUSHY 

Mushy state meaning that the liquid fraction is strictly below 1 and strictly above 0.

CS_SOLIDIFICATION_STATE_LIQUID 

Liquid state (the liquid fraction is equal to 1)

CS_SOLIDIFICATION_STATE_EUTECTIC 

Eutectic state. Only possible with a CS_SOLIDIFICATION_MODEL_BINARY_ALLOY model. A rough transition between the mushy and solid state occurs when an eutectic transition happens.

CS_SOLIDIFICATION_N_STATES 

◆ cs_solidification_strategy_t

Kind of strategy to use to model the segregation/solidification process. This implies a setting of functions to update the liquid fraction, the thermal source terms, the liquid concentration and its related quantities.

Enumerator
CS_SOLIDIFICATION_STRATEGY_LEGACY 
CS_SOLIDIFICATION_STRATEGY_TAYLOR 
CS_SOLIDIFICATION_STRATEGY_PATH 
CS_SOLIDIFICATION_N_STRATEGIES 

Function Documentation

◆ cs_solidification_activate()

cs_solidification_t* cs_solidification_activate ( cs_solidification_model_t  model,
cs_flag_t  options,
cs_flag_t  post_flag,
const cs_boundary_t boundaries,
cs_navsto_param_model_t  ns_model,
cs_navsto_param_model_flag_t  ns_model_flag,
cs_navsto_param_coupling_t  algo_coupling,
cs_navsto_param_post_flag_t  ns_post_flag 
)

Activate the solidification module.

Parameters
[in]modeltype of modelling
[in]optionsflag to handle optional parameters
[in]post_flagpredefined post-processings
[in]boundariespointer to the domain boundaries
[in]ns_modelmodel equations for the NavSto system
[in]ns_model_flagoption flag for the Navier-Stokes system
[in]algo_couplingalgorithm used for solving the NavSto system
[in]ns_post_flagpredefined post-processings for Navier-Stokes
Returns
a pointer to a new allocated solidification structure

◆ cs_solidification_check_binary_alloy_model()

cs_solidification_binary_alloy_t* cs_solidification_check_binary_alloy_model ( void  )

Sanity checks on the consistency of the settings of the binary alloy model.

Returns
a pointer to the structure

◆ cs_solidification_check_stefan_model()

cs_solidification_stefan_t* cs_solidification_check_stefan_model ( void  )

Sanity checks on the consistency of the Stefan's model settings.

Returns
a pointer to the structure

◆ cs_solidification_check_voller_model()

cs_solidification_voller_t* cs_solidification_check_voller_model ( void  )

Sanity checks on the consistency of the Voller's model settings.

Returns
a pointer to the structure

◆ cs_solidification_compute()

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

Solve equations related to the solidification module.

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

◆ cs_solidification_destroy_all()

cs_solidification_t* cs_solidification_destroy_all ( void  )

Free the main structure related to the solidification module.

Returns
a NULL pointer

◆ cs_solidification_extra_op()

void cs_solidification_extra_op ( const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant,
const cs_time_step_t ts 
)

Predefined extra-operations for the solidification module.

Parameters
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure
[in]tspointer to a cs_time_step_t structure

◆ cs_solidification_extra_post()

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

Predefined post-processing output for the solidification module. Prototype of this function is fixed since it is a function pointer defined in cs_post.h (cs_post_time_mesh_dep_output_t)

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

◆ cs_solidification_finalize_setup()

void cs_solidification_finalize_setup ( const cs_cdo_connect_t connect,
const cs_cdo_quantities_t quant 
)

Finalize the setup stage for equations related to the solidification module.

Parameters
[in]connectpointer to a cs_cdo_connect_t structure
[in]quantpointer to a cs_cdo_quantities_t structure

◆ cs_solidification_get_binary_alloy_struct()

cs_solidification_binary_alloy_t* cs_solidification_get_binary_alloy_struct ( void  )

Get the structure defining the binary alloy model.

Returns
a pointer to the structure

◆ cs_solidification_get_stefan_struct()

cs_solidification_stefan_t* cs_solidification_get_stefan_struct ( void  )

Get the structure defining the Stefan model.

Returns
a pointer to the structure

◆ cs_solidification_get_structure()

cs_solidification_t* cs_solidification_get_structure ( void  )

Retrieve the main structure to deal with solidification process.

Returns
a pointer to a new allocated solidification structure

◆ cs_solidification_get_voller_struct()

cs_solidification_voller_t* cs_solidification_get_voller_struct ( void  )

Get the structure defining the Voller model.

Returns
a pointer to the structure

◆ cs_solidification_init_setup()

void cs_solidification_init_setup ( void  )

Setup equations/properties related to the solidification module.

◆ cs_solidification_init_values()

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

Set an initial values for all quantities related to this module This is done after the setup step.

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

◆ cs_solidification_is_activated()

bool cs_solidification_is_activated ( void  )

Test if solidification module is activated.

◆ cs_solidification_log_setup()

void cs_solidification_log_setup ( void  )

Summarize the solidification module in the log file dedicated to the setup.

◆ cs_solidification_set_binary_alloy_model()

void cs_solidification_set_binary_alloy_model ( const char *  name,
const char *  varname,
cs_real_t  beta_t,
cs_real_t  temp0,
cs_real_t  beta_c,
cs_real_t  conc0,
cs_real_t  kp,
cs_real_t  mliq,
cs_real_t  t_eutec,
cs_real_t  t_melt,
cs_real_t  solute_diff,
cs_real_t  latent_heat,
cs_real_t  s_das 
)

Set the main physical parameters which describe a solidification process with a binary alloy (with components A and B) Add a transport equation for the solute concentration to simulate the conv/diffusion of the alloy ratio between the two components of the alloy.

Parameters
[in]namename of the binary alloy
[in]varnamename of the unknown related to the tracer eq.
[in]beta_tthermal dilatation coefficient
[in]temp0reference temperature (Boussinesq term)
[in]beta_csolutal dilatation coefficient
[in]conc0reference mixture concentration (Boussinesq term)
[in]kpvalue of the distribution coefficient
[in]mliqliquidus slope for the solute concentration
[in]t_eutectemperature at the eutectic point
[in]t_meltphase-change temperature for the pure material (A)
[in]solute_diffsolutal diffusion coefficient in the liquid
[in]latent_heatlatent heat
[in]s_dassecondary dendrite arm spacing

◆ cs_solidification_set_forcing_eps()

void cs_solidification_set_forcing_eps ( cs_real_t  forcing_eps)

Set the value of the epsilon parameter used in the forcing term of the momentum equation.

Parameters
[in]forcing_epsepsilon used in the penalization term to avoid a division by zero

◆ cs_solidification_set_segr_functions()

void cs_solidification_set_segr_functions ( cs_solidification_func_t vel_forcing,
cs_solidification_func_t cliq_update,
cs_solidification_func_t gliq_update,
cs_solidification_func_t thm_st_update 
)

Set the functions to perform the update of physical properties and/or the computation of the thermal source term or quantities and/or the way to perform the coupling between the thermal equation and the bulk concentration computation. All this setting defines the way to compute the solidification process of a binary alloy. If a function is set to NULL then the automatic settings are kept.

–Advanced usage– This enables to finely control the numerical or physical modelling aspects.

Parameters
[in]vel_forcingpointer to update the velocity forcing
[in]cliq_updatepointer to update the liquid concentration
[in]gliq_updatepointer to update the liquid fraction
[in]thm_st_updatepointer to update thermal source terms

◆ cs_solidification_set_stefan_model()

void cs_solidification_set_stefan_model ( cs_real_t  t_change,
cs_real_t  latent_heat 
)

Set the main physical parameters which describe the Stefan model.

Parameters
[in]t_changeliquidus/solidus temperature (in K)
[in]latent_heatlatent heat

◆ cs_solidification_set_strategy()

void cs_solidification_set_strategy ( cs_solidification_strategy_t  strategy)

Set the strategy to update quantitiess (liquid fraction and the thermal source term for the two main quantities)

Parameters
[in]strategystrategy to perform the update of quantities

◆ cs_solidification_set_verbosity()

void cs_solidification_set_verbosity ( int  verbosity)

Set the level of verbosity for the solidification module.

Parameters
[in]verbositylevel of verbosity to set

◆ cs_solidification_set_voller_model()

void cs_solidification_set_voller_model ( cs_real_t  beta,
cs_real_t  t_ref,
cs_real_t  t_solidus,
cs_real_t  t_liquidus,
cs_real_t  latent_heat,
cs_real_t  s_das 
)

Set the main physical parameters which describe the Voller and Prakash modelling.

Parameters
[in]betathermal dilatation coefficient
[in]t_refreference temperature (for the Boussinesq approx)
[in]t_solidussolidus temperature (in K)
[in]t_liquidusliquidus temperature (in K)
[in]latent_heatlatent heat
[in]s_dassecondary dendrite space arms

◆ cs_solidification_set_voller_model_no_velocity()

void cs_solidification_set_voller_model_no_velocity ( cs_real_t  t_solidus,
cs_real_t  t_liquidus,
cs_real_t  latent_heat 
)

Set the main physical parameters which describe the Voller and Prakash modelling.

Parameters
[in]betathermal dilatation coefficient
[in]t_refreference temperature (for the Boussinesq approx)
[in]t_solidussolidus temperature (in K)
[in]t_liquidusliquidus temperature (in K)
[in]latent_heatlatent heat