My Project
programmer's documentation
|
#include "cs_defs.h"
#include <assert.h>
#include "bft_error.h"
#include "bft_mem.h"
#include "cs_log.h"
#include "cs_time_step.h"
#include "cs_hho_builder.h"
Macros | |
#define | _dp3 cs_math_3_dot_product |
#define | _mv3 cs_math_33_3_product |
#define | CS_HHO_BUILDER_DBG 0 |
Functions | |
static void | _fill_vol_reco_op (cs_sdm_t *stiffness, cs_sdm_t *rhs_c_t, cs_hho_builder_t *hhob) |
Fill the volume-related part of the matrices when building the reconstruction operator. More... | |
static void | _add_tria_reduction (cs_real_t t_eval, const cs_xdef_analytic_input_t *anai, const cs_basis_func_t *fbf, const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const double surf, cs_cell_builder_t *cb, cs_real_t array[]) |
Compute the reduction onto the face polynomial space of a function defined by an analytical expression depending on the location and the current time. More... | |
static void | _add_tria_reduction_v (cs_real_t t_eval, const cs_xdef_analytic_input_t *anai, const cs_basis_func_t *fbf, const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const double surf, cs_cell_builder_t *cb, cs_real_t array[]) |
Compute the reduction onto the face polynomial space of a function defined by an analytical expression depending on the location and the current time. Vector case. More... | |
static void | _add_tetra_reduction (cs_real_t t_eval, const cs_xdef_analytic_input_t *anai, const cs_basis_func_t *cbf, const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const cs_real_3_t xv4, const double vol, cs_cell_builder_t *cb, cs_real_t array[]) |
Compute the reduction onto the cell polynomial space of a function defined by an analytical expression depending on the location and the current time. More... | |
static void | _add_tetra_reduction_v (cs_real_t t_eval, const cs_xdef_analytic_input_t *anai, const cs_basis_func_t *cbf, const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const cs_real_3_t xv4, const double vol, cs_cell_builder_t *cb, cs_real_t array[]) |
Compute the reduction onto the cell polynomial space of a function defined by an analytical expression depending on the location and the current time.Vector case. More... | |
static void | _add_tria_to_reco_op (const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const double surf, const cs_basis_func_t *fbf, const cs_real_t *kappa_nfc, cs_real_3_t *gpts, cs_sdm_t *rc, cs_sdm_t *rf, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the local stiffness matrix on the gradient basis and to the right-hand side. More... | |
static void | _add_tetra_to_reco_op (const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const cs_real_3_t xv4, const double vol, cs_sdm_t *stiffness, cs_real_3_t *gpts, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the local stiffness matrix on the gradient basis and to the right-hand side. More... | |
static void | _add_contrib_mcg (const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const cs_real_3_t xv4, const double vol, const cs_basis_func_t *cbf_kp1, cs_sdm_t *m_cg, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the M_cg matrix. More... | |
static void | _add_contrib_mf_cg (const cs_real_3_t xv1, const cs_real_3_t xv2, const cs_real_3_t xv3, const double surf, const cs_basis_func_t *fbf, const cs_basis_func_t *cbf_kp1, cs_sdm_t *mf_cg, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Routine for computing surfacic integrals over a triangle which are contributions to the Mf_cg matrix. More... | |
static cs_sdm_t * | _compute_mcg (const cs_cell_mesh_t *cm, cs_basis_func_t *cbf_kp1, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Compute the diffusion operator. The gradient reconstruction operator has to be built just before this call (cb->aux stores the rhs) More... | |
cs_hho_builder_t * | cs_hho_builder_create (int order, int n_fc) |
Allocate a cs_hho_builder_t structure. More... | |
void | cs_hho_builder_free (cs_hho_builder_t **p_builder) |
Free a cs_hho_builder_t structure. More... | |
void | cs_hho_builder_cellwise_setup (const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Set-up the basis functions related to a cell, its gradient and to the faces of this cell. Compute cell and face projection and its related modified Cholesky factorization. More... | |
void | cs_hho_builder_compute_grad_reco (const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Compute the gradient operator stemming from the relation stiffness * grad_op = rhs where stiffness is a square matrix of size grd_size rhs is matrix of size (n_fc*f_size + c_size) * grd_size Hence, grad_op a matrix grd_size * (n_fc*f_size + c_size) More... | |
void | cs_hho_builder_diffusion (const cs_cell_mesh_t *cm, cs_cell_builder_t *cb, cs_hho_builder_t *hhob) |
Compute the diffusion operator. The gradient reconstruction operator has to be built just before this call (cb->aux stores the rhs) More... | |
void | cs_hho_builder_reduction_from_analytic (const cs_xdef_t *def, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t red[]) |
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analytical expression depending on the location and the current time red array has to be allocated before calling this function. More... | |
void | cs_hho_builder_reduction_from_analytic_v (const cs_xdef_t *def, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t red[]) |
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analytical expression depending on the location and the current time This function handles the vector case. More... | |
void | cs_hho_builder_compute_dirichlet (const cs_xdef_t *def, short int f, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t res[]) |
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces. More... | |
void | cs_hho_builder_compute_dirichlet_v (const cs_xdef_t *def, short int f, const cs_cell_mesh_t *cm, cs_real_t t_eval, cs_cell_builder_t *cb, cs_hho_builder_t *hhob, cs_real_t res[]) |
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces. Vector case. More... | |
#define _dp3 cs_math_3_dot_product |
#define _mv3 cs_math_33_3_product |
#define CS_HHO_BUILDER_DBG 0 |
|
inlinestatic |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the M_cg matrix.
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | xv4 | fourth vertex |
[in] | vol | volume of the tetrahedron |
[in] | cbf_kp1 | pointer to a cell basis function order:k+1 |
[in,out] | m_cg | matrix to compute |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
|
inlinestatic |
Routine for computing surfacic integrals over a triangle which are contributions to the Mf_cg matrix.
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | surf | surface of the triangle to consider |
[in] | fbf | pointer to a set of face basis functions |
[in] | cbf_kp1 | pointer to a set of cell basis functions (k+1) |
[in,out] | mf_cg | matrix to compute |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
|
inlinestatic |
Compute the reduction onto the cell polynomial space of a function defined by an analytical expression depending on the location and the current time.
[in] | t_eval | time at which one performs the evaluation |
[in] | anai | pointer to an analytical definition |
[in] | cbf | pointer to a structure for face basis functions |
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | xv4 | third vertex |
[in] | vol | volume of the tetrahedron |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | array | array storing values to compute |
|
inlinestatic |
Compute the reduction onto the cell polynomial space of a function defined by an analytical expression depending on the location and the current time.Vector case.
[in] | t_eval | time at which one performs the evaluation |
[in] | anai | pointer to an analytical definition |
[in] | cbf | pointer to a structure for face basis functions |
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | xv4 | third vertex |
[in] | vol | volume of the tetrahedron |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | array | array storing values to compute |
|
static |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the local stiffness matrix on the gradient basis and to the right-hand side.
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | xv4 | fourth vertex |
[in] | vol | volume of the terahedron |
[in,out] | stiffness | stiffness matrix to compute |
[in,out] | gpts | coordinates of the Gauss points |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
|
inlinestatic |
Compute the reduction onto the face polynomial space of a function defined by an analytical expression depending on the location and the current time.
[in] | t_eval | time at which one performs the evaluation |
[in] | anai | pointer to an analytical definition |
[in] | fbf | pointer to a structure for face basis functions |
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | surf | area of the triangle |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | array | array storing values to compute |
|
inlinestatic |
Compute the reduction onto the face polynomial space of a function defined by an analytical expression depending on the location and the current time. Vector case.
[in] | t_eval | time at which one performs the evaluation |
[in] | anai | pointer to an analytical definition |
[in] | fbf | pointer to a structure for face basis functions |
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | surf | area of the triangle |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | array | array storing values to compute |
|
static |
Routine for computing volumetric integrals over a tetrahedron which are contributions to the local stiffness matrix on the gradient basis and to the right-hand side.
[in] | xv1 | first vertex |
[in] | xv2 | second vertex |
[in] | xv3 | third vertex |
[in] | surf | area of the triangle |
[in] | fbf | pointer to the related set of face basis functions |
[in] | kappa_nfc | permeability tensor times the related face normal |
[in,out] | gpts | coordinates of the Gauss points |
[in,out] | rc | right-hand side matrix to compute (cell part) |
[in,out] | rf | right-hand side matrix to compute (face part) |
[in,out] | gpts | coordinates of the Gauss points |
[in,out] | kappa_nfc | coordinates of the Gauss points |
[in,out] | cb | pointer to a cs_cell_builder_structure_t |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
|
static |
Compute the diffusion operator. The gradient reconstruction operator has to be built just before this call (cb->aux stores the rhs)
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | cbf_kp1 | pointer to a set of cell basis functions order=k+1 |
[in,out] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
|
static |
Fill the volume-related part of the matrices when building the reconstruction operator.
[in,out] | stiffness | pointer to the stiffness matrix |
[in,out] | rhs_c_t | pointer to the right-hand side (matrix) |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
void cs_hho_builder_cellwise_setup | ( | const cs_cell_mesh_t * | cm, |
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob | ||
) |
Set-up the basis functions related to a cell, its gradient and to the faces of this cell. Compute cell and face projection and its related modified Cholesky factorization.
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
void cs_hho_builder_compute_dirichlet | ( | const cs_xdef_t * | def, |
short int | f, | ||
const cs_cell_mesh_t * | cm, | ||
cs_real_t | t_eval, | ||
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob, | ||
cs_real_t | res[] | ||
) |
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces.
[in] | def | pointer to a cs_xdef_t structure |
[in] | f | local face id in the cellwise view of the mesh |
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in,out] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
[in,out] | res | vector containing the result |
void cs_hho_builder_compute_dirichlet_v | ( | const cs_xdef_t * | def, |
short int | f, | ||
const cs_cell_mesh_t * | cm, | ||
cs_real_t | t_eval, | ||
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob, | ||
cs_real_t | res[] | ||
) |
Compute the projection of the Dirichlet boundary conditions onto the polynomial spaces on faces. Vector case.
[in] | def | pointer to a cs_xdef_t structure |
[in] | f | local face id in the cellwise view of the mesh |
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in,out] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
[in,out] | res | vector containing the result |
void cs_hho_builder_compute_grad_reco | ( | const cs_cell_mesh_t * | cm, |
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob | ||
) |
Compute the gradient operator stemming from the relation stiffness * grad_op = rhs where stiffness is a square matrix of size grd_size rhs is matrix of size (n_fc*f_size + c_size) * grd_size Hence, grad_op a matrix grd_size * (n_fc*f_size + c_size)
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
cs_hho_builder_t* cs_hho_builder_create | ( | int | order, |
int | n_fc | ||
) |
Allocate a cs_hho_builder_t structure.
[in] | order | order of the polynomial basis function |
[in] | n_fc | max. number of faces in a cell |
void cs_hho_builder_diffusion | ( | const cs_cell_mesh_t * | cm, |
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob | ||
) |
Compute the diffusion operator. The gradient reconstruction operator has to be built just before this call (cb->aux stores the rhs)
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
void cs_hho_builder_free | ( | cs_hho_builder_t ** | p_builder | ) |
Free a cs_hho_builder_t structure.
[in,out] | p_builder | pointer of pointer on a cs_hho_builder_t struct. |
void cs_hho_builder_reduction_from_analytic | ( | const cs_xdef_t * | def, |
const cs_cell_mesh_t * | cm, | ||
cs_real_t | t_eval, | ||
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob, | ||
cs_real_t | red[] | ||
) |
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analytical expression depending on the location and the current time red array has to be allocated before calling this function.
[in] | def | pointer to a cs_xdef_t structure |
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in,out] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
[in,out] | red | vector containing the reduction |
void cs_hho_builder_reduction_from_analytic_v | ( | const cs_xdef_t * | def, |
const cs_cell_mesh_t * | cm, | ||
cs_real_t | t_eval, | ||
cs_cell_builder_t * | cb, | ||
cs_hho_builder_t * | hhob, | ||
cs_real_t | red[] | ||
) |
Compute the reduction onto the polynomial spaces (cell and faces) of a function defined by an analytical expression depending on the location and the current time This function handles the vector case.
red array has to be allocated before calling this function.
[in] | def | pointer to a cs_xdef_t structure |
[in] | cm | pointer to a cs_cell_mesh_t structure |
[in] | t_eval | time at which one performs the evaluation |
[in,out] | cb | pointer to a cell builder_t structure |
[in,out] | hhob | pointer to a cs_hho_builder_t structure |
[in,out] | red | vector containing the reduction |