My Project
programmer's documentation
Data Structures | Macros | Typedefs | Functions
cs_basis_func.h File Reference
#include "cs_cdo_local.h"
#include "cs_quadrature.h"
Include dependency graph for cs_basis_func.h:

Go to the source code of this file.

Data Structures

struct  cs_basis_func_t
 

Macros

#define CS_BASIS_FUNC_MONOMIAL   (1 << 0) /* 1: Use mononial basis functions */
 
#define CS_BASIS_FUNC_GRADIENT   (1 << 1) /* 2: Basis functions for gradient */
 

Typedefs

typedef 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. More...
 
typedef 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. More...
 
typedef 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. More...
 
typedef 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. More...
 
typedef void() cs_basis_func_compute_facto_t(void *pbf)
 Generic prototype for computing the Modified Choloesky factorization of the projection matrix (mass matrix) related to the basis function. More...
 
typedef 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 application of a Modified Choloesky factorization which should be performed before calling this function. The input array is defined as follows: int_elem v.phi_i for all i in the basis. v is a function to estimate. More...
 
typedef void() cs_basis_func_dump_proj_t(const void *pbf)
 Generic prototype for printing the projector matrix related to the given basis function. More...
 

Functions

static short int cs_basis_func_get_poly_order (const cs_basis_func_t *bf)
 Get the polynomial order of the given basis. More...
 
cs_basis_func_tcs_basis_func_create (cs_flag_t flag, short int k, short int dim)
 Allocate a cs_basis_func_t structure. More...
 
cs_basis_func_tcs_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. Up to now, only cell basis functions are handled. Building a projection matrix is not possible in this case. More...
 
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 are handled. More...
 
cs_basis_func_tcs_basis_func_free (cs_basis_func_t *pbf)
 Free a cs_basis_func_t structure. More...
 
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. More...
 
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. More...
 
void cs_basis_func_dump (const cs_basis_func_t *pbf)
 Dump a cs_basis_func_t structure. More...
 
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 fname if given otherwise print into the standard output. More...
 

Macro Definition Documentation

◆ CS_BASIS_FUNC_GRADIENT

#define CS_BASIS_FUNC_GRADIENT   (1 << 1) /* 2: Basis functions for gradient */

◆ CS_BASIS_FUNC_MONOMIAL

#define CS_BASIS_FUNC_MONOMIAL   (1 << 0) /* 1: Use mononial basis functions */

Typedef Documentation

◆ cs_basis_func_compute_facto_t

typedef void() cs_basis_func_compute_facto_t(void *pbf)

Generic prototype for computing the Modified Choloesky factorization of the projection matrix (mass matrix) related to the basis function.

Parameters
[in,out]pbfpointer to be cast into a cs_basis_func_t structure

◆ cs_basis_func_compute_proj_t

typedef 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.

Parameters
[in,out]pbfpointer to be cast into a cs_basis_func_t structure
[in]cmpointer to a cs_cell_mesh_t
[in]idid of the element to consider

◆ cs_basis_func_dump_proj_t

typedef void() cs_basis_func_dump_proj_t(const void *pbf)

Generic prototype for printing the projector matrix related to the given basis function.

Parameters
[in]pbfpointer to be cast into a cs_basis_func_t structure

◆ cs_basis_func_eval_all_at_point_t

typedef 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.

Parameters
[in]bfpointer to be cast into a cs_basis_func_t structure
[in]coordspoint coordinates
[in,out]evalvector containing the evaluations of the functions

◆ cs_basis_func_eval_at_point_t

typedef 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.

Parameters
[in]bfpointer to be cast into a cs_basis_func_t structure
[in]coordspoint coordinates
[in]startstarts evaluating basis function from this id-1
[in]endends evaluating basis function at this id
[in,out]evalvector containing the evaluations of the functions

◆ cs_basis_func_project_t

typedef 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 application of a Modified Choloesky factorization which should be performed before calling this function. The input array is defined as follows: int_elem v.phi_i for all i in the basis. v is a function to estimate.

Parameters
[in]pbfpointer to be cast into a cs_basis_func_t structure
[in]arrayarray to project
[in,out]dofprojection of the array (= DoF values in this basis)

◆ cs_basis_func_setup_t

typedef 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.

Parameters
[in,out]pbfpointer to be cast into a cs_basis_func_t structure
[in]cmpointer to a cs_cell_mesh_t
[in]idid of the element to consider
[in]centerpoint used for centering the set of basis functions
[in,out]cbpointer to a cs_cell_builder_t structure

Function Documentation

◆ cs_basis_func_copy_setup()

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 are handled.

Parameters
[in]refset of basis function used as a reference
[in,out]rcvset of basis function where members are set

◆ cs_basis_func_create()

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.

Parameters
[in]flagmetadata related to the way of building basis functions
[in]korder
[in]dim2 or 3 geometrical dimension
Returns
a pointer to the new cs_basis_func_t
Parameters
[in]flagmetadata related to the way of building basis functions
[in]kpolynomial order
[in]dim2 or 3 w.r.t. the geometrical dimension
Returns
a pointer to the new cs_basis_func_t

◆ cs_basis_func_dump()

void cs_basis_func_dump ( const cs_basis_func_t pbf)

Dump a cs_basis_func_t structure.

Parameters
[in]pbfpointer to the cs_basis_func_t structure to dump

◆ cs_basis_func_fprintf()

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 fname if given otherwise print into the standard output.

Parameters
[in]fppointer to a file structure or NULL
[in]fnamefilename or NULL
[in]pbfpointer to the cs_basis_func_t structure to dump

◆ cs_basis_func_free()

cs_basis_func_t* cs_basis_func_free ( cs_basis_func_t pbf)

Free a cs_basis_func_t structure.

Parameters
[in,out]pbfpointer to the cs_basis_func_t structure to free
Returns
a NULL pointer

◆ cs_basis_func_get_hho_flag()

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.

Parameters
[out]face_flagpointer to options related to face basis functinos
[out]cell_flagpointer to options related to cell basis functinos

◆ cs_basis_func_get_poly_order()

static short int cs_basis_func_get_poly_order ( const cs_basis_func_t bf)
inlinestatic

Get the polynomial order of the given basis.

Parameters
[in]bfset of basis functions
Returns
the polynomial order

◆ cs_basis_func_grad_create()

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. Up to now, only cell basis functions are handled. Building a projection matrix is not possible in this case.

Parameters
[in]refset of basis function used as a reference
Returns
a pointer to the new cs_basis_func_t for gradient of the current basis functions

◆ cs_basis_func_set_hho_flag()

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.

Parameters
[in]face_flagoptions related to face basis functinos
[in]cell_flagoptions related to cell basis functinos