#include "cs_cdo_connect.h"
#include "cs_cdo_local.h"
#include "cs_cdo_quantities.h"
#include "cs_param_cdo.h"
#include "cs_property.h"
#include "cs_sdm.h"
Go to the source code of this file.
Data Structures | |
struct | cs_hodge_param_t |
Structure storing all metadata/parameters related to the usage of a discrete Hodge operator. More... | |
struct | cs_hodge_t |
Structure associated to a discrete Hodge operator. More... | |
Typedefs | |
typedef void() | cs_hodge_compute_t(const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a discrete Hodge operator or a related operator (such as the stiffmess matrix) for a given cell. The discrete Hodge operator is stored in hodge->matrix whereas the associated operator is stored in cb->loc. More... | |
Enumerations | |
enum | cs_hodge_type_t { CS_HODGE_TYPE_VPCD, CS_HODGE_TYPE_EPFD, CS_HODGE_TYPE_FPED, CS_HODGE_TYPE_EDFP, CS_HODGE_TYPE_CPVD, CS_HODGE_TYPE_FB, CS_HODGE_TYPE_VC, CS_HODGE_N_TYPES } |
Type of discrete Hodge operators. More... | |
enum | cs_hodge_algo_t { CS_HODGE_ALGO_VORONOI, CS_HODGE_ALGO_WBS, CS_HODGE_ALGO_COST, CS_HODGE_ALGO_OCS2, CS_HODGE_ALGO_BUBBLE, CS_HODGE_ALGO_AUTO, CS_HODGE_N_ALGOS } |
Type of algorithm to build a discrete Hodge operator. More... | |
Functions | |
cs_hodge_t * | cs_hodge_create (const cs_cdo_connect_t *connect, const cs_property_t *property, const cs_hodge_param_t *hp, bool need_tensor, bool need_eigen) |
Create and initialize a pointer to a cs_hodge_t structure. More... | |
cs_hodge_t ** | cs_hodge_init_context (const cs_cdo_connect_t *connect, const cs_property_t *property, const cs_hodge_param_t *hp, bool need_tensor, bool need_eigen) |
Create and initialize an array of pointers to a cs_hodge_t structures. This array is of size the number of OpenMP threads. Only the one associated to the current thread is set. More... | |
void | cs_hodge_free (cs_hodge_t **p_hodge) |
Free a cs_hodge_t structure. More... | |
void | cs_hodge_free_context (cs_hodge_t ***p_hodges) |
Free a set of cs_hodge_t structures. More... | |
cs_hodge_compute_t * | cs_hodge_get_func (const char *calling_func, const cs_hodge_param_t hp) |
Retrieve a function pointer to compute a discrete Hodge operator. More... | |
void | cs_hodge_param_log (const char *prefix, const cs_property_t *property, const cs_hodge_param_t hp) |
Output the settings related to a cs_hodge_param_t structure. More... | |
void | cs_hodge_copy_parameters (const cs_hodge_param_t *h_ref, cs_hodge_param_t *h_cpy) |
Copy the set of parameters associated to a discrete Hodge operator to another one. More... | |
void | cs_hodge_set_property_value (const cs_lnum_t c_id, const cs_real_t t_eval, const cs_flag_t c_flag, cs_hodge_t *hodge) |
Set the property value (scalar- or tensor-valued) related to a discrete Hodge operator inside a cell and if needed other related quantities. More... | |
void | cs_hodge_set_property_value_cw (const cs_cell_mesh_t *cm, const cs_real_t t_eval, const cs_flag_t c_flag, cs_hodge_t *hodge) |
Set the property value (scalar- or tensor-valued) related to a discrete Hodge operator inside a cell and if needed ohter related quantities. Cell-wise variant (usage of cs_cell_mesh_t structure) More... | |
void | cs_hodge_fb_cost_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes. More... | |
void | cs_hodge_fb_bubble_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic COST algo. with the usage of bubble stabilization. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes. More... | |
void | cs_hodge_fb_voro_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes. More... | |
void | cs_hodge_vb_cost_get_iso_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and an isotropic property. More... | |
void | cs_hodge_vb_cost_get_aniso_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and an anistropic property. More... | |
void | cs_hodge_vb_bubble_get_iso_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic Bubble algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and isotropic material property. More... | |
void | cs_hodge_vb_bubble_get_aniso_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic Bubble algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and anisotropic material property. More... | |
void | cs_hodge_vb_ocs2_get_aniso_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the Orthogonal Consistent/Sub-Stabilization decomposition (OCS2) with a subdivision of pvol_{e,c}. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case Vb schemes and an anisotropic material property. More... | |
void | cs_hodge_vb_cost_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes. More... | |
void | cs_hodge_vb_voro_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes. More... | |
void | cs_hodge_vb_wbs_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdivision (WBS) algo. The computed matrix is stored in cb->loc. More... | |
void | cs_hodge_vcb_get_stiffness (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdivision (WBS) algo. The computed matrix is stored in cb->loc. More... | |
void | cs_hodge_fb_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator on a given cell which is equivalent of a mass matrix. It relies on a CO+ST algo. and is specific to CDO-Fb schemes. The discrete Hodge operator is stored in hodge->matrix. More... | |
void | cs_hodge_vcb_voro_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the Voronoi algo. This leads to a diagonal operator. This function is specific for vertex+cell-based schemes The discrete Hodge operator is stored in hodge->matrix. More... | |
void | cs_hodge_vcb_wbs_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the WBS algo. This function is specific for vertex+cell-based schemes The discrete Hodge operator is stored in hodge->matrix. More... | |
void | cs_hodge_vpcd_wbs_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using WBS algo. Hodge op. from primal vertices to dual cells. This function is specific for vertex-based schemes The discrete Hodge operator is stored in hodge->matrix. More... | |
void | cs_hodge_vpcd_voro_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal vertices to dual cells. This function is specific for vertex-based schemes The discrete Hodge operator is stored in hodge->matrix. More... | |
void | cs_hodge_epfd_voro_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal edges to dual faces. The discrete Hodge operator is stored in hodge->matrix This function is specific for vertex-based schemes. More... | |
void | cs_hodge_epfd_cost_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal edges to dual faces. The discrete Hodge operator is stored in hodge->matrix This function is specific for vertex-based schemes. More... | |
void | cs_hodge_epfd_bubble_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the COST algo. with a bubble stabilization. The discrete Hodge operator is stored in hodge->matrix Hodge op. from primal edges to dual faces. This function is specific for vertex-based schemes. More... | |
void | cs_hodge_epfd_ocs2_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the Orthogonal Consistent/Sub-Stabilization decomposition (OCS2) with a subdivision of pvol_{e,c}. The discrete Hodge operator is stored in hodge->matrix Hodge op. from primal edges to dual faces. This function is specific for vertex-based schemes. More... | |
void | cs_hodge_fped_voro_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes. More... | |
void | cs_hodge_fped_cost_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes. More... | |
void | cs_hodge_fped_bubble_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the Bubble algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes. More... | |
void | cs_hodge_edfp_voro_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from dual edges to primal faces. The discrete Hodge operator is stored in hodge->matrix This function is related to face-based schemes. More... | |
void | cs_hodge_edfp_cost_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to primal faces. This function is related to face-based schemes. More... | |
void | cs_hodge_edfp_bubble_get (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the Bubble algo. Hodge op. from dual edges to primal faces. The discrete Hodge operator is stored in hodge->matrix This function is related to face-based schemes. More... | |
void | cs_hodge_edfp_cost_get_opt (const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to primal faces. This function is related to face-based schemes. More... | |
void | cs_hodge_matvec (const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_hodge_param_t hodgep, const cs_property_t *pty, const cs_real_t in_vals[], cs_real_t t_eval, cs_real_t result[]) |
Compute cellwise a discrete hodge operator and multiple it with a vector. More... | |
void | cs_hodge_circulation_from_flux (const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, cs_real_t t_eval, const cs_hodge_param_t hodgep, const cs_property_t *pty, const cs_real_t flux[], cs_real_t circul[]) |
Compute cellwise a discrete hodge operator in order to define a circulation array from a flux array. More... | |
void | cs_hodge_compute_wbs_surfacic (const cs_face_mesh_t *fm, cs_sdm_t *hf) |
Compute the hodge operator related to a face (i.e. a mass matrix with unity property) using a Whitney Barycentric Subdivision (WBS) algorithm. More... | |
typedef void() cs_hodge_compute_t(const cs_cell_mesh_t *cm, cs_hodge_t *hodge, cs_cell_builder_t *cb) |
Build a discrete Hodge operator or a related operator (such as the stiffmess matrix) for a given cell. The discrete Hodge operator is stored in hodge->matrix whereas the associated operator is stored in cb->loc.
[in] | cm | pointer to a cs_cell_mesh_t struct. |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
enum cs_hodge_algo_t |
Type of algorithm to build a discrete Hodge operator.
Enumerator | |
---|---|
CS_HODGE_ALGO_VORONOI | Assume that an orthogonality condition holds between geometrical entities This leads to a diagonal discrete Hodge operator (and is related to a two-point flux approximation) |
CS_HODGE_ALGO_WBS | WBS: Whitney Barycentric Subdivision Rely on a subdivision into tetrahedra of a polyhedral cell and then apply algorithms close to what exists in the FE litterature |
CS_HODGE_ALGO_COST | Orthogonal decomposition between the consistency (CO) and the stabilization (ST) parts. Several discretizations share this way to build discrete Hodge operators like SUSHI, DGA and GCR (these discretizations only differ from the scaling coefficient in front of the stabilization term) |
CS_HODGE_ALGO_OCS2 | This algorithm is close to CS_HODGE_ALGO_COST but rely on subdivision of each polyhedral cell which is a refinement of what is considered in CS_HODGE_ALGO_COST for the building the stabilization |
CS_HODGE_ALGO_BUBBLE | This algorithm also relies on an orthogonal decomposition between the consistency (CO) and the stabilization (ST) parts but the stabilization part relies on a bubble (similar to what can be found in the FE literature) |
CS_HODGE_ALGO_AUTO | Automatic switch between the above-mentioned algorithms. |
CS_HODGE_N_ALGOS |
enum cs_hodge_type_t |
Type of discrete Hodge operators.
void cs_hodge_circulation_from_flux | ( | const cs_cdo_connect_t * | connect, |
const cs_cdo_quantities_t * | quant, | ||
cs_real_t | t_eval, | ||
const cs_hodge_param_t | hodgep, | ||
const cs_property_t * | pty, | ||
const cs_real_t | flux[], | ||
cs_real_t | circul[] | ||
) |
Compute cellwise a discrete hodge operator in order to define a circulation array from a flux array.
[in] | connect | pointer to a cs_cdo_connect_t structure |
[in] | quant | pointer to a cs_cdo_quantities_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in] | hodgep | cs_hodge_param_t structure |
[in] | pty | pointer to a cs_property_t structure or NULL |
[in] | flux | vector to multiply with the discrete Hodge op. |
[in,out] | circul | array storing the resulting matrix-vector product |
void cs_hodge_compute_wbs_surfacic | ( | const cs_face_mesh_t * | fm, |
cs_sdm_t * | hf | ||
) |
Compute the hodge operator related to a face (i.e. a mass matrix with unity property) using a Whitney Barycentric Subdivision (WBS) algorithm.
[in] | fm | pointer to a cs_face_mesh_t structure |
[in,out] | hf | pointer to a cs_sdm_t structure to define |
void cs_hodge_copy_parameters | ( | const cs_hodge_param_t * | h_ref, |
cs_hodge_param_t * | h_cpy | ||
) |
Copy the set of parameters associated to a discrete Hodge operator to another one.
[in] | h_ref | reference set of parameters |
[in,out] | h_cpy | set of parameters to update |
cs_hodge_t* cs_hodge_create | ( | const cs_cdo_connect_t * | connect, |
const cs_property_t * | property, | ||
const cs_hodge_param_t * | hp, | ||
bool | need_tensor, | ||
bool | need_eigen | ||
) |
Create and initialize a pointer to a cs_hodge_t structure.
[in] | connect | pointer to cs_cdo_connect_t structure |
[in] | property | pointer to a property structure |
[in] | hp | pointer to a cs_hodge_param_t structure |
[in] | need_tensor | true if one needs a tensor otherwise false |
[in] | need_eigen | true if one needs to compute eigen valuese |
void cs_hodge_edfp_bubble_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the Bubble algo. Hodge op. from dual edges to primal faces. The discrete Hodge operator is stored in hodge->matrix This function is related to face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_edfp_cost_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to primal faces. This function is related to face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t struct. |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to primal faces. This function is related to face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_edfp_cost_get_opt | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from dual edges to primal faces. This function is related to face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t struct. |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_edfp_voro_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from dual edges to primal faces. The discrete Hodge operator is stored in hodge->matrix This function is related to face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_epfd_bubble_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the COST algo. with a bubble stabilization. The discrete Hodge operator is stored in hodge->matrix Hodge op. from primal edges to dual faces. This function is specific for vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_epfd_cost_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal edges to dual faces. The discrete Hodge operator is stored in hodge->matrix This function is specific for vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_epfd_ocs2_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the Orthogonal Consistent/Sub-Stabilization decomposition (OCS2) with a subdivision of pvol_{e,c}. The discrete Hodge operator is stored in hodge->matrix Hodge op. from primal edges to dual faces. This function is specific for vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_epfd_voro_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal edges to dual faces. The discrete Hodge operator is stored in hodge->matrix This function is specific for vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fb_bubble_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic COST algo. with the usage of bubble stabilization. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fb_cost_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fb_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator on a given cell which is equivalent of a mass matrix. It relies on a CO+ST algo. and is specific to CDO-Fb schemes. The discrete Hodge operator is stored in hodge->matrix.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fb_voro_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO face-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fped_bubble_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the Bubble algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fped_cost_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the COST algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_fped_voro_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal faces to dual edges. The discrete Hodge operator is stored in hodge->matrix This function is related to cell-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_free | ( | cs_hodge_t ** | p_hodge | ) |
Free a cs_hodge_t structure.
[in,out] | p_hodge | double pointer to a cs_hodge_t structure |
void cs_hodge_free_context | ( | cs_hodge_t *** | p_hodges | ) |
Free a set of cs_hodge_t structures.
[in,out] | p_hodges | triple pointer to cs_hodge_t structures |
cs_hodge_compute_t* cs_hodge_get_func | ( | const char * | calling_func, |
const cs_hodge_param_t | hp | ||
) |
Retrieve a function pointer to compute a discrete Hodge operator.
[in] | calling_func | name of the calling function |
[in] | hp | a cs_hodge_param_t structure |
cs_hodge_t** cs_hodge_init_context | ( | const cs_cdo_connect_t * | connect, |
const cs_property_t * | property, | ||
const cs_hodge_param_t * | hp, | ||
bool | need_tensor, | ||
bool | need_eigen | ||
) |
Create and initialize an array of pointers to a cs_hodge_t structures. This array is of size the number of OpenMP threads. Only the one associated to the current thread is set.
[in] | connect | pointer to cs_cdo_connect_t structure |
[in] | property | pointer to a property structure |
[in] | hp | pointer to a cs_hodge_param_t structure |
[in] | need_tensor | true if one needs a tensor otherwise false |
[in] | need_eigen | true if one needs to compute eigen valuese |
void cs_hodge_matvec | ( | const cs_cdo_connect_t * | connect, |
const cs_cdo_quantities_t * | quant, | ||
const cs_hodge_param_t | hodgep, | ||
const cs_property_t * | pty, | ||
const cs_real_t | in_vals[], | ||
cs_real_t | t_eval, | ||
cs_real_t | result[] | ||
) |
Compute cellwise a discrete hodge operator and multiple it with a vector.
[in] | connect | pointer to a cs_cdo_connect_t structure |
[in] | quant | pointer to a cs_cdo_quantities_t structure |
[in] | hodgep | cs_hodge_param_t structure |
[in] | pty | pointer to a cs_property_t structure or NULL |
[in] | in_vals | vector to multiply with the discrete Hodge op. |
[in] | t_eval | time at which one performs the evaluation |
[in,out] | result | array storing the resulting matrix-vector product |
void cs_hodge_param_log | ( | const char * | prefix, |
const cs_property_t * | property, | ||
const cs_hodge_param_t | hp | ||
) |
Output the settings related to a cs_hodge_param_t structure.
[in] | prefix | optional string |
[in] | property | optional pointer to a property structure |
[in] | hp | a cs_hodge_param_t structure |
void cs_hodge_set_property_value | ( | const cs_lnum_t | c_id, |
const cs_real_t | t_eval, | ||
const cs_flag_t | c_flag, | ||
cs_hodge_t * | hodge | ||
) |
Set the property value (scalar- or tensor-valued) related to a discrete Hodge operator inside a cell and if needed other related quantities.
[in] | c_id | id of the cell to deal with |
[in] | t_eval | time at which one performs the evaluation |
[in] | c_flag | flag related to this cell |
[in,out] | hodge | pointer to a cs_hodge_t structure |
Set the property value (scalar- or tensor-valued) related to a discrete Hodge operator inside a cell and if needed other related quantities.
[in] | c_id | id of the cell to deal with |
[in] | t_eval | time at which one performs the evaluation |
[in] | c_flag | flag related to this cell |
[in,out] | hodge | pointer to a cs_hodge_t structure |
void cs_hodge_set_property_value_cw | ( | const cs_cell_mesh_t * | cm, |
const cs_real_t | t_eval, | ||
const cs_flag_t | c_flag, | ||
cs_hodge_t * | hodge | ||
) |
Set the property value (scalar- or tensor-valued) related to a discrete Hodge operator inside a cell and if needed ohter related quantities. Cell-wise variant (usage of cs_cell_mesh_t structure)
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in] | c_flag | flag related to this cell |
[in,out] | hodge | pointer to a cs_hodge_t structure |
void cs_hodge_vb_bubble_get_aniso_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic Bubble algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and anisotropic material property.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_bubble_get_iso_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic Bubble algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and isotropic material property.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_cost_get_aniso_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and an anistropic property.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_cost_get_iso_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes and an isotropic property.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_cost_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic COST algo. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_ocs2_get_aniso_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the Orthogonal Consistent/Sub-Stabilization decomposition (OCS2) with a subdivision of pvol_{e,c}. The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case Vb schemes and an anisotropic material property.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_voro_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the Voronoi algorithm The computed matrix is stored in cb->loc and the related discrete hodge operator in hodge->matrix Case of CDO vertex-based schemes.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vb_wbs_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdivision (WBS) algo. The computed matrix is stored in cb->loc.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vcb_get_stiffness | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local stiffness matrix using the generic WBS algo. WBS standing for Whitney Barycentric Subdivision (WBS) algo. The computed matrix is stored in cb->loc.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vcb_voro_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the Voronoi algo. This leads to a diagonal operator. This function is specific for vertex+cell-based schemes The discrete Hodge operator is stored in hodge->matrix.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vcb_wbs_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using the WBS algo. This function is specific for vertex+cell-based schemes The discrete Hodge operator is stored in hodge->matrix.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vpcd_voro_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using VORONOI algo. Hodge op. from primal vertices to dual cells. This function is specific for vertex-based schemes The discrete Hodge operator is stored in hodge->matrix.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |
void cs_hodge_vpcd_wbs_get | ( | const cs_cell_mesh_t * | cm, |
cs_hodge_t * | hodge, | ||
cs_cell_builder_t * | cb | ||
) |
Build a local Hodge operator for a given cell using WBS algo. Hodge op. from primal vertices to dual cells. This function is specific for vertex-based schemes The discrete Hodge operator is stored in hodge->matrix.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in,out] | hodge | pointer to a cs_hodge_t structure |
[in,out] | cb | pointer to a cs_cell_builder_t structure |