9.0
general documentation
Loading...
Searching...
No Matches
cs_gradient.h File Reference
Include dependency graph for cs_gradient.h:

Go to the source code of this file.

Enumerations

enum  cs_gradient_type_t {
  CS_GRADIENT_GREEN_ITER , CS_GRADIENT_LSQ , CS_GRADIENT_GREEN_LSQ , CS_GRADIENT_GREEN_VTX ,
  CS_GRADIENT_GREEN_R
}
enum  cs_gradient_limit_t { CS_GRADIENT_LIMIT_NONE = -1 , CS_GRADIENT_LIMIT_CELL = 0 , CS_GRADIENT_LIMIT_FACE = 1 }

Functions

void cs_gradient_initialize (void)
 Initialize gradient computation API.
void cs_gradient_finalize (void)
 Finalize gradient computation API.
void cs_gradient_free_quantities (void)
 Free saved gradient quantities.
void cs_gradient_scalar (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int hyd_p_flag, int w_stride, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, cs_real_3_t f_ext[], const cs_field_bc_coeffs_t *bc_coeffs, cs_real_t var[], cs_real_t *c_weight, const cs_internal_coupling_t *cpl, cs_real_t grad[][3])
void cs_gradient_vector (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, const cs_field_bc_coeffs_t *bc_coeffs_v, cs_real_t var[][3], cs_real_t *c_weight, const cs_internal_coupling_t *cpl, cs_real_t gradv[][3][3])
void cs_gradient_tensor (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, const cs_field_bc_coeffs_t *bc_coeffs_ts, cs_real_6_t *var, cs_real_63_t *grad)
void cs_gradient_scalar_synced_input (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int hyd_p_flag, int w_stride, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, cs_real_t f_ext[][3], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t var[], const cs_real_t *c_weight, const cs_internal_coupling_t *cpl, cs_real_t grad[][3])
void cs_gradient_vector_synced_input (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_t var[][3], const cs_real_t val_f[][3], const cs_real_t c_weight[], const cs_internal_coupling_t *cpl, cs_real_t grad[][3][3])
 Compute cell gradient of vector field.
void cs_gradient_tensor_synced_input (const char *var_name, cs_gradient_type_t gradient_type, cs_halo_type_t halo_type, int inc, int n_r_sweeps, int verbosity, cs_gradient_limit_t clip_mode, double epsilon, double clip_coeff, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_t var[][6], const cs_real_t val_f[][6], cs_real_63_t *grad)
 Compute cell gradient of tensor.
void cs_gradient_scalar_cell (const cs_mesh_t *m, const cs_mesh_quantities_t *fvq, cs_lnum_t c_id, cs_halo_type_t halo_type, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t var[], const cs_real_t c_weight[], cs_real_t grad[3])
 Compute the gradient of a scalar field at a given cell using least-squares reconstruction.
void cs_gradient_vector_cell (const cs_mesh_t *m, const cs_mesh_quantities_t *fvq, cs_lnum_t c_id, cs_halo_type_t halo_type, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_real_t var[][3], const cs_real_t c_weight[], cs_real_t grad[3][3])
 Compute the gradient of a vector field at a given cell using least-squares reconstruction.
void cs_gradient_tensor_cell (const cs_mesh_t *m, const cs_mesh_quantities_t *fvq, cs_lnum_t c_id, cs_halo_type_t halo_type, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_real_t var[][6], const cs_real_t c_weight[], cs_real_t grad[6][3])
 Compute the gradient of a tensor field at a given cell using least-squares reconstruction.
void cs_gradient_type_by_imrgra (int imrgra, cs_gradient_type_t *gradient_type, cs_halo_type_t *halo_type)
void cs_gradient_porosity_balance (int inc)
 compute the steady balance due to porous modelling for the pressure gradient.

Variables

const char * cs_gradient_type_name []

Enumeration Type Documentation

◆ cs_gradient_limit_t

Enumerator
CS_GRADIENT_LIMIT_NONE 

no limiter

CS_GRADIENT_LIMIT_CELL 

cell values extrapolated by cell gradient to neighbor cells can not exceed the limiting factor times the actual value

CS_GRADIENT_LIMIT_FACE 

cell values extrapolated by face gradient (mean of cell gradients) extrapolated to neighbor cells can not exceed the limiting factor times the actual value

◆ cs_gradient_type_t

Enumerator
CS_GRADIENT_GREEN_ITER 

Iterative

CS_GRADIENT_LSQ 

Least-squares

CS_GRADIENT_GREEN_LSQ 

Green-Gauss reconstruction with least squares gradient face values

CS_GRADIENT_GREEN_VTX 

Green-Gauss with vertex interpolation

CS_GRADIENT_GREEN_R 

Green-Gauss times a renormalization matrix (to be exact on affine functions)

Function Documentation

◆ cs_gradient_finalize()

void cs_gradient_finalize ( void )

Finalize gradient computation API.

◆ cs_gradient_free_quantities()

void cs_gradient_free_quantities ( void )

Free saved gradient quantities.

This is required when the mesh changes, so that the on-demand computation will be updated.

◆ cs_gradient_initialize()

void cs_gradient_initialize ( void )

Initialize gradient computation API.

◆ cs_gradient_porosity_balance()

void cs_gradient_porosity_balance ( int inc)

compute the steady balance due to porous modelling for the pressure gradient.

Parameters
[in]incif 0, solve on increment; 1 otherwise

◆ cs_gradient_scalar()

void cs_gradient_scalar ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int hyd_p_flag,
int w_stride,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
cs_real_3_t f_ext[],
const cs_field_bc_coeffs_t * bc_coeffs,
cs_real_t var[],
cs_real_t * c_weight,
const cs_internal_coupling_t * cpl,
cs_real_t grad[][3] )

◆ cs_gradient_scalar_cell()

void cs_gradient_scalar_cell ( const cs_mesh_t * m,
const cs_mesh_quantities_t * fvq,
cs_lnum_t c_id,
cs_halo_type_t halo_type,
const cs_field_bc_coeffs_t * bc_coeffs,
const cs_real_t var[],
const cs_real_t c_weight[],
cs_real_t grad[3] )

Compute the gradient of a scalar field at a given cell using least-squares reconstruction.

This assumes ghost cell values which might be used are already synchronized.

When boundary conditions are provided, both the bc_coeff_a and bc_coeff_b arrays must be given. If boundary values are known, bc_coeff_a can point to the boundary values array, and bc_coeff_b set to nullptr. If bc_coeff_a is nullptr, bc_coeff_b is ignored.

Parameters
[in]mpointer to associated mesh structure
[in]fvqpointer to associated finite volume quantities
[in]c_idcell id
[in]halo_typehalo type
[in]bc_coeffsboundary condition structure
[in]vargradient's base variable
[in]c_weightcell variable weight, or nullptr
[out]gradgradient

◆ cs_gradient_scalar_synced_input()

void cs_gradient_scalar_synced_input ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int hyd_p_flag,
int w_stride,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
cs_real_t f_ext[][3],
const cs_field_bc_coeffs_t * bc_coeffs,
const cs_real_t var[],
const cs_real_t * c_weight,
const cs_internal_coupling_t * cpl,
cs_real_t grad[][3] )

◆ cs_gradient_tensor()

void cs_gradient_tensor ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
const cs_field_bc_coeffs_t * bc_coeffs_ts,
cs_real_6_t * var,
cs_real_63_t * grad )

◆ cs_gradient_tensor_cell()

void cs_gradient_tensor_cell ( const cs_mesh_t * m,
const cs_mesh_quantities_t * fvq,
cs_lnum_t c_id,
cs_halo_type_t halo_type,
const cs_field_bc_coeffs_t * bc_coeffs_ts,
const cs_real_t var[][6],
const cs_real_t c_weight[],
cs_real_t grad[6][3] )

Compute the gradient of a tensor field at a given cell using least-squares reconstruction.

This assumes ghost cell values which might be used are already synchronized.

When boundary conditions are provided, both the bc_coeff_a and bc_coeff_b arrays must be given. If boundary values are known, bc_coeff_a can point to the boundary values array, and bc_coeff_b set to nullptr. If bc_coeff_a is nullptr, bc_coeff_b is ignored.

Parameters
[in]mpointer to associated mesh structure
[in]fvqpointer to associated finite volume quantities
[in]c_idcell id
[in]halo_typehalo type
[in]bc_coeffs_tsboundary condition structure
[in]vargradient's base variable
[in]c_weightcell variable weight, or nullptr
[out]gradgradient

◆ cs_gradient_tensor_synced_input()

void cs_gradient_tensor_synced_input ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
const cs_field_bc_coeffs_t * bc_coeffs_ts,
const cs_real_t var[][6],
const cs_real_t val_f[][6],
cs_real_63_t * grad )

Compute cell gradient of tensor.

This variant of the cs_gradient_tensor function assumes ghost cell values for input arrays (var and optionally c_weight) have already been synchronized.

Parameters
[in]var_namevariable name
[in]gradient_typegradient type
[in]halo_typehalo type
[in]incif 0, solve on increment; 1 otherwise
[in]n_r_sweepsif > 1, number of reconstruction sweeps (only used by CS_GRADIENT_GREEN_ITER)
[in]verbosityverbosity level
[in]clip_modeclipping mode
[in]epsilonprecision for iterative gradient calculation
[in]clip_coeffclipping coefficient
[in]bc_coeffs_tsboundary condition structure
[in,out]vargradient's base variable
[out]gradgradient ( $ \der{t_ij}{x_k} $ is grad[][ij][k])

◆ cs_gradient_type_by_imrgra()

void cs_gradient_type_by_imrgra ( int imrgra,
cs_gradient_type_t * gradient_type,
cs_halo_type_t * halo_type )

◆ cs_gradient_vector()

void cs_gradient_vector ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
const cs_field_bc_coeffs_t * bc_coeffs_v,
cs_real_t var[][3],
cs_real_t * c_weight,
const cs_internal_coupling_t * cpl,
cs_real_t gradv[][3][3] )

◆ cs_gradient_vector_cell()

void cs_gradient_vector_cell ( const cs_mesh_t * m,
const cs_mesh_quantities_t * fvq,
cs_lnum_t c_id,
cs_halo_type_t halo_type,
const cs_field_bc_coeffs_t * bc_coeffs_v,
const cs_real_t var[][3],
const cs_real_t c_weight[],
cs_real_t grad[3][3] )

Compute the gradient of a vector field at a given cell using least-squares reconstruction.

This assumes ghost cell values which might be used are already synchronized.

When boundary conditions are provided, both the bc_coeff_a and bc_coeff_b arrays must be given. If boundary values are known, bc_coeff_a can point to the boundary values array, and bc_coeff_b set to nullptr. If bc_coeff_a is nullptr, bc_coeff_b is ignored.

Parameters
[in]mpointer to associated mesh structure
[in]fvqpointer to associated finite volume quantities
[in]c_idcell id
[in]halo_typehalo type
[in]bc_coeffsboundary condition structure
[in]vargradient's base variable
[in]c_weightcell variable weight, or nullptr
[out]gradgradient

◆ cs_gradient_vector_synced_input()

void cs_gradient_vector_synced_input ( const char * var_name,
cs_gradient_type_t gradient_type,
cs_halo_type_t halo_type,
int inc,
int n_r_sweeps,
int verbosity,
cs_gradient_limit_t clip_mode,
double epsilon,
double clip_coeff,
const cs_field_bc_coeffs_t * bc_coeffs_v,
const cs_real_t var[][3],
const cs_real_t val_f[][3],
const cs_real_t c_weight[],
const cs_internal_coupling_t * cpl,
cs_real_t grad[][3][3] )

Compute cell gradient of vector field.

This variant of the cs_gradient_vector function assumes ghost cell values for input arrays (var and optionally c_weight) have already been synchronized.

Parameters
[in]var_namevariable name
[in]gradient_typegradient type
[in]halo_typehalo type
[in]incif 0, solve on increment; 1 otherwise
[in]n_r_sweepsif > 1, number of reconstruction sweeps (only used by CS_GRADIENT_GREEN_ITER)
[in]verbosityverbosity level
[in]clip_modeclipping mode
[in]epsilonprecision for iterative gradient calculation
[in]clip_coeffclipping coefficient
[in]bc_coeffs_vboundary condition structure
[in]vargradient's base variable
[in]c_weightcell variable weight, or nullptr
[in]cplassociated internal coupling, or nullptr
[out]gradgradient ( $ \der{u_i}{x_j} $ is gradv[][i][j])

Variable Documentation

◆ cs_gradient_type_name

const char* cs_gradient_type_name[]
extern