programmer's documentation
Typedefs | Functions
cs_cdo_diffusion.h File Reference
#include "cs_cdo.h"
#include "cs_cdo_bc.h"
#include "cs_cdo_connect.h"
#include "cs_cdo_local.h"
#include "cs_cdo_quantities.h"
#include "cs_hodge.h"
#include "cs_param.h"
Include dependency graph for cs_cdo_diffusion.h:

Go to the source code of this file.

Typedefs

typedef struct _cs_cdo_diff_t cs_cdo_diff_t
 

Functions

cs_cdo_diff_tcs_cdo_diffusion_builder_init (const cs_cdo_connect_t *connect, cs_space_scheme_t space_scheme, bool is_uniform, const cs_param_hodge_t h_info, const cs_param_bc_enforce_t bc_enforce)
 Initialize a builder structure used to build the stiffness matrix. More...
 
cs_cdo_diff_tcs_cdo_diffusion_builder_free (cs_cdo_diff_t *diff)
 Free a cs_cdo_diff_t structure. More...
 
cs_hodge_builder_tcs_cdo_diffusion_get_hodge_builder (cs_cdo_diff_t *diff)
 Get the related Hodge builder structure. More...
 
void cs_cdo_diffusion_get_tmp_buffers (const cs_cdo_diff_t *diff, cs_real_3_t **tmp_vec, cs_real_t **tmp_sca)
 Get temporary buffers attached to a cs_cdo_diff_t structure. More...
 
void cs_cdo_diffusion_build_c2bcbf (const cs_cdo_connect_t *connect, const cs_cdo_bc_list_t *dir_face, cs_lnum_t *p_c2bcbf_idx[], cs_lnum_t *p_c2bcbf_ids[])
 Define a cell –> Dirichlet boundary faces connectivity. More...
 
cs_locmat_tcs_cdo_diffusion_build_local (const cs_cdo_quantities_t *quant, const cs_cell_mesh_t *lm, const cs_real_3_t *tensor, cs_cdo_diff_t *diff)
 Define the local (cellwise) stiffness matrix. More...
 
void cs_cdo_diffusion_weak_bc (cs_lnum_t f_id, const cs_cell_mesh_t *cm, const cs_real_t matpty[3][3], cs_cdo_diff_t *diff, cs_cdo_locsys_t *csys)
 Define the local (cellwise) "normal trace gradient" matrix taking into account Dirichlet BCs by a weak enforcement using Nitsche technique (symmetrized or not) More...
 
void cs_cdo_diffusion_cellwise_flux (const cs_cell_mesh_t *cm, const cs_dface_t *dfaces, const cs_real_t pty_tens[3][3], const double *p_v, const double p_c, cs_cdo_diff_t *diff, double *c_flux)
 Compute the diffusive flux across dual faces for a given cell The computation takes into account a subdivision into tetrahedra of the current cell based on p_{ef,c}. More...
 
double cs_cdo_diffusion_face_flux (const cs_face_mesh_t *fm, const cs_real_t pty_tens[3][3], const double *p_v, const double p_f, const double p_c, cs_cdo_diff_t *diff)
 Compute the diffusive flux across a face (based on a subdivision into tetrahedra of the volume p_{f,c}) More...
 

Typedef Documentation

§ cs_cdo_diff_t

typedef struct _cs_cdo_diff_t cs_cdo_diff_t

Function Documentation

§ cs_cdo_diffusion_build_c2bcbf()

void cs_cdo_diffusion_build_c2bcbf ( const cs_cdo_connect_t connect,
const cs_cdo_bc_list_t dir_face,
cs_lnum_t p_c2bcbf_idx[],
cs_lnum_t p_c2bcbf_ids[] 
)

Define a cell –> Dirichlet boundary faces connectivity.

Parameters
[in]connectpointer to a cs_cdo_connect_t structure
[in]dir_facepointer to a cs_cdo_bc_list_t structure
[in,out]c2bcbf_idxpointer to the index to build
[in,out]c2bcbf_idspointer to the list of ids to build

§ cs_cdo_diffusion_build_local()

cs_locmat_t* cs_cdo_diffusion_build_local ( const cs_cdo_quantities_t quant,
const cs_cell_mesh_t cm,
const cs_real_3_t tensor,
cs_cdo_diff_t diff 
)

Define the local (cellwise) stiffness matrix.

Parameters
[in]quantpointer to a cs_cdo_quantities_t struct.
[in]lmcell-wise connectivity and quantitites
[in]tensor3x3 matrix attached to the diffusion property
[in,out]diffauxiliary structure used to build the diff. term
Returns
a pointer to a local stiffness matrix
Parameters
[in]quantpointer to a cs_cdo_quantities_t struct.
[in]cmcell-wise connectivity and quantitites
[in]tensor3x3 matrix attached to the diffusion property
[in,out]diffauxiliary structure used to build the diff. term
Returns
a pointer to a local stiffness matrix

§ cs_cdo_diffusion_builder_free()

cs_cdo_diff_t* cs_cdo_diffusion_builder_free ( cs_cdo_diff_t diff)

Free a cs_cdo_diff_t structure.

Parameters
[in,out]diffpointer to a cs_cdo_diff_t struc.
Returns
NULL

§ cs_cdo_diffusion_builder_init()

cs_cdo_diff_t* cs_cdo_diffusion_builder_init ( const cs_cdo_connect_t connect,
cs_space_scheme_t  space_scheme,
bool  is_uniform,
const cs_param_hodge_t  h_info,
const cs_param_bc_enforce_t  bc_enforce 
)

Initialize a builder structure used to build the stiffness matrix.

Parameters
[in]connectpointer to a cs_cdo_connect_t structure
[in]space_schemescheme used for discretizing in space
[in]is_uniformdiffusion tensor is uniform ? (true or false)
[in]h_infocs_param_hodge_t structure
[in]bc_enforcetype of boundary enforcement for Dirichlet values
Returns
a pointer to a new allocated cs_cdo_diff_t structure

§ cs_cdo_diffusion_cellwise_flux()

void cs_cdo_diffusion_cellwise_flux ( const cs_cell_mesh_t cm,
const cs_dface_t dface,
const cs_real_t  pty_tens[3][3],
const double *  p_v,
const double  p_c,
cs_cdo_diff_t diff,
double *  c_flux 
)

Compute the diffusive flux across dual faces for a given cell The computation takes into account a subdivision into tetrahedra of the current cell based on p_{ef,c}.

Parameters
[in]cmpointer to a cs_face_mesh_t structure
[in]dfacespointer to the dual faces related to cell edges
[in]pty_tens3x3 matrix related to the diffusion property
[in]p_varray of values attached to face vertices
[in]p_cvalue attached to the cell
[in,out]diffauxiliary structure dedicated to diffusion
[in,out]c_fluxflux across dual faces inside this cell
[in]cmpointer to a cs_face_mesh_t structure
[in]dfacepointer to the dual faces related to cell edges
[in]pty_tens3x3 matrix related to the diffusion property
[in]p_varray of values attached to face vertices
[in]p_cvalue attached to the cell
[in,out]diffauxiliary structure dedicated to diffusion
[in,out]c_fluxflux across dual faces inside this cell

§ cs_cdo_diffusion_face_flux()

double cs_cdo_diffusion_face_flux ( const cs_face_mesh_t fm,
const cs_real_t  pty_tens[3][3],
const double *  p_v,
const double  p_f,
const double  p_c,
cs_cdo_diff_t diff 
)

Compute the diffusive flux across a face (based on a subdivision into tetrahedra of the volume p_{f,c})

Parameters
[in]fmpointer to a cs_face_mesh_t structure
[in]pty_tens3x3 matrix related to the diffusion property
[in]p_varray of values attached to face vertices
[in]p_fvalue attached to the face
[in]p_cvalue attached to the cell
[in,out]diffauxiliary structure dedicated to diffusion
Returns
the value of the diffusive flux across the current face

§ cs_cdo_diffusion_get_hodge_builder()

cs_hodge_builder_t* cs_cdo_diffusion_get_hodge_builder ( cs_cdo_diff_t diff)

Get the related Hodge builder structure.

Parameters
[in]diffpointer to a cs_cdo_diff_t structure
Returns
a pointer to a cs_hodge_builder_t structure

§ cs_cdo_diffusion_get_tmp_buffers()

void cs_cdo_diffusion_get_tmp_buffers ( const cs_cdo_diff_t diff,
cs_real_3_t **  tmp_vec,
cs_real_t **  tmp_sca 
)

Get temporary buffers attached to a cs_cdo_diff_t structure.

Parameters
[in]diffpointer to a cs_cdo_diff_t structure
[in,out]tmp_vecpointer to a buffer of cs_real_3_t
[in,out]tmp_scapointer to a buffer of cs_real_t

§ cs_cdo_diffusion_weak_bc()

void cs_cdo_diffusion_weak_bc ( cs_lnum_t  f_id,
const cs_cell_mesh_t cm,
const cs_real_t  matpty[3][3],
cs_cdo_diff_t diff,
cs_cdo_locsys_t csys 
)

Define the local (cellwise) "normal trace gradient" matrix taking into account Dirichlet BCs by a weak enforcement using Nitsche technique (symmetrized or not)

Parameters
[in]f_idface id (a border face attached to a Dir. BC)
[in]cmpointer to a cs_cell_mesh_t struct.
[in]matpty3x3 matrix related to the diffusion property
[in,out]diffauxiliary structure used to build the diff. term
[in,out]csysstructure storing the cell-wise system