1#ifndef __CS_BASIS_FUNC_H__
2#define __CS_BASIS_FUNC_H__
49#define CS_BASIS_FUNC_MONOMIAL (1 << 0)
50#define CS_BASIS_FUNC_GRADIENT (1 << 1)
240static inline short int
void cs_basis_func_set_hho_flag(cs_flag_t face_flag, cs_flag_t cell_flag)
Set options for basis functions when using HHO schemes.
Definition: cs_basis_func.cpp:2703
void cs_basis_func_dump(const cs_basis_func_t *pbf)
Dump a cs_basis_func_t structure.
Definition: cs_basis_func.cpp:2736
static short int cs_basis_func_get_poly_order(const cs_basis_func_t *bf)
Get the polynomial order of the given basis.
Definition: cs_basis_func.h:241
void() cs_basis_func_project_t(const void *pbf, const cs_real_t *array, cs_real_t *dof)
Generic prototype for projecting an array on the polynomial basis function. This results from the app...
Definition: cs_basis_func.h:154
cs_basis_func_t * cs_basis_func_create(cs_flag_t flag, short int k, short int dim)
Allocate a cs_basis_func_t structure.
Definition: cs_basis_func.cpp:2367
cs_basis_func_t * cs_basis_func_grad_create(const cs_basis_func_t *ref)
Allocate a cs_basis_func_t structure which is associated to an existing set of basis functions....
Definition: cs_basis_func.cpp:2564
void() cs_basis_func_dump_proj_t(const void *pbf)
Generic prototype for printing the projector matrix related to the given basis function.
Definition: cs_basis_func.h:168
void cs_basis_func_copy_setup(const cs_basis_func_t *ref, cs_basis_func_t *rcv)
Copy the center and the different axis from the reference basis Up to now, only cell basis functions ...
Definition: cs_basis_func.cpp:2649
void cs_basis_func_get_hho_flag(cs_flag_t *face_flag, cs_flag_t *cell_flag)
Get options for basis functions when using HHO schemes.
Definition: cs_basis_func.cpp:2720
void() cs_basis_func_compute_proj_t(void *pbf, const cs_cell_mesh_t *cm, const short int id)
Generic prototype for computing the projector to the space spanned by the basis functions.
Definition: cs_basis_func.h:123
void() cs_basis_func_eval_all_at_point_t(const void *bf, const cs_real_t coords[3], cs_real_t *eval)
Generic prototype for evaluating all basis functions at a given point.
Definition: cs_basis_func.h:68
void() cs_basis_func_compute_facto_t(void *pbf)
Generic prototype for computing the Modified Choloesky factorization of the projection matrix (mass m...
Definition: cs_basis_func.h:137
cs_basis_func_t * cs_basis_func_free(cs_basis_func_t *pbf)
Free a cs_basis_func_t structure.
Definition: cs_basis_func.cpp:2675
void cs_basis_func_fprintf(FILE *fp, const char *fname, const cs_basis_func_t *pbf)
Print a cs_basis_func_t structure Print into the file f if given otherwise open a new file named fnam...
Definition: cs_basis_func.cpp:2777
void() cs_basis_func_eval_at_point_t(const void *bf, const cs_real_t coords[3], short int start, short int end, cs_real_t *eval)
Generic prototype for evaluating a set of basis functions at a given point.
Definition: cs_basis_func.h:86
void() cs_basis_func_setup_t(void *pbf, const cs_cell_mesh_t *cm, const short int id, const cs_real_t center[3], cs_cell_builder_t *cb)
Generic prototype for setting up a set of basis functions.
Definition: cs_basis_func.h:105
#define BEGIN_C_DECLS
Definition: cs_defs.h:542
double cs_real_t
Floating-point value.
Definition: cs_defs.h:342
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:359
#define END_C_DECLS
Definition: cs_defs.h:543
unsigned short int cs_flag_t
Definition: cs_defs.h:344
@ k
Definition: cs_field_pointer.h:72
void() cs_quadrature_tria_t(const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_real_3_t gpts[], double *weights)
Generic functoin pointer to compute the quadrature points for a triangle.
Definition: cs_quadrature.h:130
void() cs_quadrature_tet_t(const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, const cs_real_3_t v4, double vol, cs_real_3_t gpts[], double weights[])
Generic function to compute the quadrature points in a tetrehedra.
Definition: cs_quadrature.h:152
Definition: cs_basis_func.h:181
cs_nvec3_t * axis
Definition: cs_basis_func.h:189
cs_basis_func_compute_facto_t * compute_factorization
Definition: cs_basis_func.h:212
cs_basis_func_dump_proj_t * dump_projector
Definition: cs_basis_func.h:208
int size
Definition: cs_basis_func.h:186
short int * deg
Definition: cs_basis_func.h:195
int n_deg_elts
Definition: cs_basis_func.h:194
cs_basis_func_setup_t * setup
Definition: cs_basis_func.h:199
cs_basis_func_eval_all_at_point_t * eval_all_at_point
Definition: cs_basis_func.h:202
cs_basis_func_compute_proj_t * compute_projector
Definition: cs_basis_func.h:207
cs_basis_func_project_t * project
Definition: cs_basis_func.h:213
cs_real_t * facto
Definition: cs_basis_func.h:214
int n_gpts_tria
Definition: cs_basis_func.h:219
cs_flag_t flag
Definition: cs_basis_func.h:183
short int poly_order
Definition: cs_basis_func.h:184
cs_real_t phi0
Definition: cs_basis_func.h:188
short int dim
Definition: cs_basis_func.h:185
cs_basis_func_eval_at_point_t * eval_at_point
Definition: cs_basis_func.h:203
int n_gpts_tetra
Definition: cs_basis_func.h:221
cs_sdm_t * projector
Definition: cs_basis_func.h:206
cs_real_3_t center
Definition: cs_basis_func.h:190
cs_quadrature_tria_t * quadrature_tria
Definition: cs_basis_func.h:220
int facto_max_size
Definition: cs_basis_func.h:215
cs_quadrature_tet_t * quadrature_tetra
Definition: cs_basis_func.h:222
Set of local and temporary buffers.
Definition: cs_cdo_local.h:60
Set of local quantities and connectivities related to a mesh cell.
Definition: cs_cdo_local.h:203
Definition: cs_defs.h:400