#include <ple_locator.h>
#include "cs_defs.h"
#include "cs_base.h"
#include "cs_field.h"
#include "cs_matrix_assembler.h"
#include "cs_mesh.h"
#include "cs_mesh_quantities.h"
#include "cs_parameters.h"
Go to the source code of this file.
Data Structures | |
struct | cs_internal_coupling_t |
Functions | |
int | cs_internal_coupling_n_couplings (void) |
Return number of defined internal couplings. More... | |
void | cs_internal_coupling_add (const char criteria_cells[], const char criteria_faces[]) |
Define coupling volume using given selection criteria. More... | |
void | cs_internal_coupling_add_volume (const char criteria_cells[]) |
Define coupling volume using given criteria. Then, this volume will be separated from the rest of the domain with thin walls. More... | |
void | cs_internal_coupling_add_volume_zone (const cs_zone_t *z) |
Define coupling volume using a cs_zone_t. Then, this volume will be separated from the rest of the domain with thin walls. More... | |
void | cs_internal_coupling_add_volume_zones (int n_zones, const int zone_ids[]) |
Define coupling volume using given cs_zone_t. Then, this volume will be separated from the rest of the domain with thin walls. More... | |
void | cs_internal_coupling_add_boundary_groups (cs_internal_coupling_t *cpl, const char *interior_name, const char *exterior_name) |
Define internal coupling volume boundary group names. More... | |
void | cs_internal_coupling_bcs (int bc_type[]) |
Impose wall BCs to internal coupled faces if not yet defined. More... | |
void | cs_internal_coupling_finalize (void) |
Destruction of all internal coupling related structures. More... | |
cs_internal_coupling_t * | cs_internal_coupling_by_id (int coupling_id) |
Return the coupling associated with a given coupling_id. More... | |
void | cs_internal_coupling_exchange_var (const cs_internal_coupling_t *cpl, int stride, cs_real_t distant[], cs_real_t local[]) |
Exchange quantities from distant to local (update local using distant). More... | |
void | cs_internal_coupling_exchange_by_cell_id (const cs_internal_coupling_t *cpl, int stride, const cs_real_t tab[], cs_real_t local[]) |
Exchange variable between groups using cell id. More... | |
void | cs_internal_coupling_exchange_by_face_id (const cs_internal_coupling_t *cpl, int stride, const cs_real_t tab[], cs_real_t local[]) |
Exchange variable between groups using face id. More... | |
void | cs_internal_coupling_it_cocg_contribution (const cs_internal_coupling_t *cpl, cs_real_33_t cocg[]) |
void | cs_internal_coupling_setup (void) |
Setup internal coupling related parameters. More... | |
void | cs_internal_coupling_initialize (void) |
Initialize internal coupling related structures. More... | |
void | cs_internal_coupling_iterative_scalar_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], cs_real_3_t *restrict grad, const cs_real_t pvar[], cs_real_3_t rhs[]) |
Add internal coupling rhs contribution for iterative gradient calculation. More... | |
void | cs_internal_coupling_iterative_vector_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], cs_real_33_t *restrict grad, const cs_real_3_t pvar[], cs_real_33_t rhs[]) |
Add internal coupling rhs contribution for iterative vector gradient calculation. More... | |
void | cs_internal_coupling_iterative_tensor_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], cs_real_63_t *restrict grad, const cs_real_6_t pvar[], cs_real_63_t rhs[]) |
Add internal coupling rhs contribution for iterative tensor gradient calculation. More... | |
void | cs_internal_coupling_reconstruct_scalar_gradient (const cs_internal_coupling_t *cpl, cs_real_3_t *restrict r_grad, cs_real_3_t grad[]) |
Add internal coupling contribution for reconstruction of the gradient of a scalar. More... | |
void | cs_internal_coupling_reconstruct_vector_gradient (const cs_internal_coupling_t *cpl, cs_real_33_t *restrict r_grad, cs_real_33_t grad[]) |
Add internal coupling contribution for reconstruction of the gradient of a vector. More... | |
void | cs_internal_coupling_reconstruct_tensor_gradient (const cs_internal_coupling_t *cpl, cs_real_63_t *restrict r_grad, cs_real_63_t grad[]) |
Add internal coupling contribution for reconstruction of the gradient of a symmetric tensor. More... | |
void | cs_internal_coupling_update_bc_coeff_s (cs_field_bc_coeffs_t *bc_coeffs, const cs_internal_coupling_t *cpl, cs_halo_type_t halo_type, int w_stride, double clip_coeff, cs_real_t *bc_coeff_a, cs_real_t *bc_coeff_b, const cs_real_t *var, const cs_real_t *c_weight) |
Compute scalar boundary condition coefficients for internal coupling. More... | |
void | cs_internal_coupling_update_bc_coeff_v (cs_field_bc_coeffs_t *bc_coeffs, const cs_internal_coupling_t *cpl, cs_halo_type_t halo_type, double clip_coeff, cs_real_t bc_coeff_a[][3], cs_real_t bc_coeff_b[][3][3], const cs_real_3_t *var, const cs_real_t *c_weight) |
Update vector boundary condition coefficients for internal coupling. More... | |
void | cs_internal_coupling_spmv_contribution (bool exclude_diag, const cs_field_t *f, const cs_real_t *restrict x, cs_real_t *restrict y) |
void | cs_internal_coupling_matrix_add_ids (int coupling_id, const cs_gnum_t *r_g_id, cs_matrix_assembler_t *ma) |
void | cs_internal_coupling_matrix_add_values (const cs_field_t *f, cs_lnum_t db_size, cs_lnum_t eb_size, const cs_gnum_t r_g_id[], cs_matrix_assembler_values_t *mav) |
void | cs_internal_coupling_coupled_faces (const cs_internal_coupling_t *cpl, cs_lnum_t *n_local, const cs_lnum_t *faces_local[], cs_lnum_t *n_distant, const cs_lnum_t *faces_distant[]) |
void | cs_internal_coupling_log (const cs_internal_coupling_t *cpl) |
void | cs_internal_coupling_dump (void) |
void | cs_internal_coupling_preprocess (cs_mesh_t *mesh) |
void | cs_internal_coupling_map (cs_mesh_t *mesh) |
void | cs_internal_coupling_add_entity (int f_id) |
void | cs_internal_coupling_initialize_scalar_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], const cs_real_t pvar[], cs_real_3_t *restrict grad) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative scalar gradient calculation. More... | |
void | cs_internal_coupling_initialize_vector_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], const cs_real_3_t pvar[], cs_real_33_t *restrict grad) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative vector gradient calculation. More... | |
void | cs_internal_coupling_initialize_tensor_gradient (const cs_internal_coupling_t *cpl, const cs_real_t c_weight[], const cs_real_6_t pvar[], cs_real_63_t *restrict grad) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative symmetric tensor gradient calculation. More... | |
void | cs_ic_field_set_exchcoeff (const cs_field_t *f, const cs_real_t *hbnd) |
Update internal coupling coefficients of the field of the given id using given boundary exchange coefficients passed by face id. More... | |
void | cs_ic_field_dist_data_by_face_id (const int field_id, int stride, const cs_real_t tab_distant[], cs_real_t tab_local[]) |
Get distant data using face id at all coupling faces for a given field id. More... | |
void cs_ic_field_dist_data_by_face_id | ( | const int | field_id, |
int | stride, | ||
const cs_real_t | tab_distant[], | ||
cs_real_t | tab_local[] | ||
) |
Get distant data using face id at all coupling faces for a given field id.
[in] | field_id | field id |
[in] | stride | number of values (interlaced) by entity |
[in] | tab_distant | exchanged data by face id |
[out] | tab_local | local data by face id |
void cs_ic_field_set_exchcoeff | ( | const cs_field_t * | f, |
const cs_real_t * | hbnd | ||
) |
Update internal coupling coefficients of the field of the given id using given boundary exchange coefficients passed by face id.
[in] | f | pointer to field |
[in] | hbnd | boundary exchange coefficients passed by face id |
void cs_internal_coupling_add | ( | const char | criteria_cells[], |
const char | criteria_faces[] | ||
) |
Define coupling volume using given selection criteria.
Then, this volume must be separated from the rest of the domain with a wall.
[in] | criteria_cells | criteria for the first group of cells |
[in] | criteria_faces | criteria for faces to be joined |
void cs_internal_coupling_add_boundary_groups | ( | cs_internal_coupling_t * | cpl, |
const char * | interior_name, | ||
const char * | exterior_name | ||
) |
Define internal coupling volume boundary group names.
This is used only for internal couplings based on a separation of volumes (cs_internal_coupling_add_volume, cs_internal_coupling_add_volume_zone, cs_internal_coupling_add_volume_zones).
The interior name is used for faces adjacent to the main volume, and the exterio name for faces adjacent to the selected (exterior) volume.
This allows filtering faces on each side of the boundary in a simpler manner.
[in,out] | cpl | pointer to mesh structure to modify |
[in] | criteria_cells | criteria for the first group of cells |
void cs_internal_coupling_add_entity | ( | int | f_id | ) |
void cs_internal_coupling_add_volume | ( | const char | criteria_cells[] | ) |
Define coupling volume using given criteria. Then, this volume will be separated from the rest of the domain with thin walls.
[in] | criteria_cells | criteria for the first group of cells |
void cs_internal_coupling_add_volume_zone | ( | const cs_zone_t * | z | ) |
Define coupling volume using a cs_zone_t. Then, this volume will be separated from the rest of the domain with thin walls.
[in] | z | pointer to cs_volume_zone_t |
Define coupling volume using a cs_zone_t. Then, this volume will be separated from the rest of the domain with thin walls.
[in] | z | pointer to cs_volume_zone_t |
void cs_internal_coupling_add_volume_zones | ( | int | n_zones, |
const int | zone_ids[] | ||
) |
Define coupling volume using given cs_zone_t. Then, this volume will be separated from the rest of the domain with thin walls.
[in] | n_zones | number of associated volume zones |
[in] | zone_ids | ids of associated volume zones |
void cs_internal_coupling_bcs | ( | int | bc_type[] | ) |
Impose wall BCs to internal coupled faces if not yet defined.
[in,out] | bc_type | face boundary condition type |
cs_internal_coupling_t* cs_internal_coupling_by_id | ( | int | coupling_id | ) |
Return the coupling associated with a given coupling_id.
[in] | coupling_id | associated with a coupling entity |
void cs_internal_coupling_coupled_faces | ( | const cs_internal_coupling_t * | cpl, |
cs_lnum_t * | n_local, | ||
const cs_lnum_t * | faces_local[], | ||
cs_lnum_t * | n_distant, | ||
const cs_lnum_t * | faces_distant[] | ||
) |
void cs_internal_coupling_dump | ( | void | ) |
void cs_internal_coupling_exchange_by_cell_id | ( | const cs_internal_coupling_t * | cpl, |
int | stride, | ||
const cs_real_t | tab[], | ||
cs_real_t | local[] | ||
) |
Exchange variable between groups using cell id.
[in] | cpl | pointer to coupling entity |
[in] | stride | number of values (non interlaced) by entity |
[in] | tab | variable exchanged |
[out] | local | local data |
void cs_internal_coupling_exchange_by_face_id | ( | const cs_internal_coupling_t * | cpl, |
int | stride, | ||
const cs_real_t | tab[], | ||
cs_real_t | local[] | ||
) |
Exchange variable between groups using face id.
[in] | cpl | pointer to coupling entity |
[in] | stride | number of values (non interlaced) by entity |
[in] | tab | variable exchanged |
[out] | local | local data |
void cs_internal_coupling_exchange_var | ( | const cs_internal_coupling_t * | cpl, |
int | stride, | ||
cs_real_t | distant[], | ||
cs_real_t | local[] | ||
) |
Exchange quantities from distant to local (update local using distant).
[in] | cpl | pointer to coupling entity |
[in] | stride | stride (e.g. 1 for double, 3 for interleaved coordinates) |
[in] | distant | distant values, size coupling->n_distant |
[out] | local | local values, size coupling->n_local |
void cs_internal_coupling_finalize | ( | void | ) |
Destruction of all internal coupling related structures.
void cs_internal_coupling_initialize | ( | void | ) |
Initialize internal coupling related structures.
void cs_internal_coupling_initialize_scalar_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
const cs_real_t | pvar[], | ||
cs_real_3_t *restrict | grad | ||
) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative scalar gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in] | pvar | variable |
[in,out] | grad | gradient |
void cs_internal_coupling_initialize_tensor_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
const cs_real_6_t | pvar[], | ||
cs_real_63_t *restrict | grad | ||
) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative symmetric tensor gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in,out] | pvar | variable |
[in,out] | grad | gradient |
void cs_internal_coupling_initialize_vector_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
const cs_real_3_t | pvar[], | ||
cs_real_33_t *restrict | grad | ||
) |
Add contribution from coupled faces (internal coupling) to initialisation for iterative vector gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in] | pvar | variable |
[in,out] | grad | gradient |
void cs_internal_coupling_it_cocg_contribution | ( | const cs_internal_coupling_t * | cpl, |
cs_real_33_t | cocg[] | ||
) |
void cs_internal_coupling_iterative_scalar_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
cs_real_3_t *restrict | grad, | ||
const cs_real_t | pvar[], | ||
cs_real_3_t | rhs[] | ||
) |
Add internal coupling rhs contribution for iterative gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in] | grad | pointer to gradient |
[in] | pvar | pointer to variable |
[in,out] | rhs | pointer to rhs contribution |
void cs_internal_coupling_iterative_tensor_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
cs_real_63_t *restrict | grad, | ||
const cs_real_6_t | pvar[], | ||
cs_real_63_t | rhs[] | ||
) |
Add internal coupling rhs contribution for iterative tensor gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in] | grad | pointer to gradient |
[in] | pvar | pointer to variable |
[in,out] | rhs | pointer to rhs contribution |
void cs_internal_coupling_iterative_vector_gradient | ( | const cs_internal_coupling_t * | cpl, |
const cs_real_t | c_weight[], | ||
cs_real_33_t *restrict | grad, | ||
const cs_real_3_t | pvar[], | ||
cs_real_33_t | rhs[] | ||
) |
Add internal coupling rhs contribution for iterative vector gradient calculation.
[in] | cpl | pointer to coupling entity |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
[in] | grad | pointer to gradient |
[in] | pvar | pointer to variable |
[in,out] | rhs | pointer to rhs contribution |
void cs_internal_coupling_log | ( | const cs_internal_coupling_t * | cpl | ) |
void cs_internal_coupling_map | ( | cs_mesh_t * | mesh | ) |
void cs_internal_coupling_matrix_add_ids | ( | int | coupling_id, |
const cs_gnum_t * | r_g_id, | ||
cs_matrix_assembler_t * | ma | ||
) |
void cs_internal_coupling_matrix_add_values | ( | const cs_field_t * | f, |
cs_lnum_t | db_size, | ||
cs_lnum_t | eb_size, | ||
const cs_gnum_t | r_g_id[], | ||
cs_matrix_assembler_values_t * | mav | ||
) |
int cs_internal_coupling_n_couplings | ( | void | ) |
Return number of defined internal couplings.
void cs_internal_coupling_preprocess | ( | cs_mesh_t * | mesh | ) |
void cs_internal_coupling_reconstruct_scalar_gradient | ( | const cs_internal_coupling_t * | cpl, |
cs_real_3_t *restrict | r_grad, | ||
cs_real_3_t | grad[] | ||
) |
Add internal coupling contribution for reconstruction of the gradient of a scalar.
[in] | cpl | pointer to coupling entity |
[in] | r_grad | pointer to reconstruction gradient |
[in,out] | grad | pointer to gradient to be reconstructed var |
void cs_internal_coupling_reconstruct_tensor_gradient | ( | const cs_internal_coupling_t * | cpl, |
cs_real_63_t *restrict | r_grad, | ||
cs_real_63_t | grad[] | ||
) |
Add internal coupling contribution for reconstruction of the gradient of a symmetric tensor.
[in] | cpl | pointer to coupling entity |
[in] | r_grad | pointer to reconstruction gradient |
[in,out] | grad | pointer to gradient to be reconstructed var |
void cs_internal_coupling_reconstruct_vector_gradient | ( | const cs_internal_coupling_t * | cpl, |
cs_real_33_t *restrict | r_grad, | ||
cs_real_33_t | grad[] | ||
) |
Add internal coupling contribution for reconstruction of the gradient of a vector.
[in] | cpl | pointer to coupling entity |
[in] | r_grad | pointer to reconstruction gradient |
[in,out] | grad | pointer to gradient to be reconstructed var |
void cs_internal_coupling_setup | ( | void | ) |
Setup internal coupling related parameters.
void cs_internal_coupling_spmv_contribution | ( | bool | exclude_diag, |
const cs_field_t * | f, | ||
const cs_real_t *restrict | x, | ||
cs_real_t *restrict | y | ||
) |
void cs_internal_coupling_update_bc_coeff_s | ( | cs_field_bc_coeffs_t * | bc_coeffs, |
const cs_internal_coupling_t * | cpl, | ||
cs_halo_type_t | halo_type, | ||
int | w_stride, | ||
double | clip_coeff, | ||
cs_real_t * | bc_coeff_a, | ||
cs_real_t * | bc_coeff_b, | ||
const cs_real_t * | var, | ||
const cs_real_t * | c_weight | ||
) |
Compute scalar boundary condition coefficients for internal coupling.
[in] | bc_coeffs | associated BC coefficients structure |
[in] | cpl | structure associated with internal coupling |
[in] | halo_type | halo type |
[in] | w_stride | stride for weighting coefficient |
[in] | clip_coeff | clipping coefficient |
[out] | bc_coeff_a | boundary condition term a |
[out] | bc_coeff_b | boundary condition term b |
[in] | var | gradient's base variable |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
Compute scalar boundary condition coefficients for internal coupling.
[in] | bc_coeffs | associated BC coefficients structure |
[in] | cpl | structure associated with internal coupling |
[in] | halo_type | halo type |
[in] | w_stride | stride for weighting coefficient |
[in] | clip_coeff | clipping coefficient |
[out] | bc_coeff_a | boundary condition term a |
[out] | bc_coeff_b | boundary condition term b |
[in] | var | gradient's base variable |
[in] | c_weight | weighted gradient coefficient variable, or NULL |
void cs_internal_coupling_update_bc_coeff_v | ( | cs_field_bc_coeffs_t * | bc_coeffs, |
const cs_internal_coupling_t * | cpl, | ||
cs_halo_type_t | halo_type, | ||
double | clip_coeff, | ||
cs_real_t | bc_coeff_a[][3], | ||
cs_real_t | bc_coeff_b[][3][3], | ||
const cs_real_3_t * | var, | ||
const cs_real_t * | c_weight | ||
) |
Update vector boundary condition coefficients for internal coupling.
[in] | bc_coeffs | associated BC coefficients structure |
[in] | cpl | structure associated with internal coupling |
[in] | halo_type | halo type |
[in] | clip_coeff | clipping coefficient |
[out] | bc_coeff_a | boundary condition term a |
[out] | bc_coeff_b | boundary condition term b |
[in] | var | gradient's base variable |
[in] | c_weight | weighted gradient coefficient variable, or NULL |