My Project
programmer's documentation
Macros | Typedefs | Functions
cs_source_term.h File Reference
#include "cs_base.h"
#include "cs_cdo_quantities.h"
#include "cs_cdo_local.h"
#include "cs_param.h"
#include "cs_quadrature.h"
#include "cs_xdef.h"
Include dependency graph for cs_source_term.h:

Go to the source code of this file.

Macros

#define CS_N_MAX_SOURCE_TERMS   8
 

Typedefs

typedef void() cs_source_term_cellwise_t(const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. More...
 

Functions

void cs_source_term_set_shared_pointers (const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect)
 Set shared pointers to main domain members. More...
 
cs_flag_t cs_source_term_set_default_flag (cs_param_space_scheme_t scheme)
 Set the default flag related to a source term according to the numerical scheme chosen for discretizing an equation. More...
 
void cs_source_term_set_reduction (cs_xdef_t *st, cs_flag_t flag)
 Set advanced parameters which are defined by default in a source term structure. More...
 
cs_flag_t cs_source_term_get_flag (const cs_xdef_t *st)
 Get metadata related to the given source term structure. More...
 
cs_flag_t cs_source_term_init (cs_param_space_scheme_t space_scheme, const int n_source_terms, cs_xdef_t *const *source_terms, cs_source_term_cellwise_t *compute_source[], cs_flag_t *sys_flag, cs_mask_t *source_mask[])
 Initialize data to build the source terms. More...
 
void cs_source_term_compute_cellwise (const int n_source_terms, cs_xdef_t *const *source_terms, const cs_cell_mesh_t *cm, const cs_mask_t *source_mask, cs_source_term_cellwise_t *compute_source[], cs_real_t time_eval, void *input, cs_cell_builder_t *cb, cs_real_t *result)
 Compute the local contributions of source terms in a cell. More...
 
void cs_source_term_compute_from_density (cs_flag_t loc, const cs_xdef_t *source, cs_real_t time_eval, double *p_values[])
 Compute the contribution related to a source term in the case of an input data which is a density. More...
 
void cs_source_term_compute_from_potential (cs_flag_t loc, const cs_xdef_t *source, cs_real_t time_eval, double *p_values[])
 Compute the contribution related to a source term in the case of an input data which is a potential. More...
 
void cs_source_term_pvsp_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices by a constant value. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure. More...
 
void cs_source_term_pvsp_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices by an analytical function. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure. More...
 
void cs_source_term_dcsd_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by a value. More...
 
void cs_source_term_dcvd_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density defined at dual cells by a value. More...
 
void cs_source_term_dcsd_by_array (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an array. More...
 
void cs_source_term_dcsd_bary_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. More...
 
void cs_source_term_dcsd_q1o1_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. More...
 
void cs_source_term_dcsd_q10o2_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a ten-point quadrature rule to evaluate the integral. Exact for quadratic function. More...
 
void cs_source_term_dcsd_q5o3_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a five-point quadrature rule to evaluate the integral. Exact for cubic function. This function may be expensive since many evaluations are needed. Please use it with care. More...
 
void cs_source_term_vcsp_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices and cells by a constant value. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure. More...
 
void cs_source_term_vcsp_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices and cells by an analytical function. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure. More...
 
void cs_source_term_pcsd_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density (sd) defined on primal cells by a value. Case of face-based schemes. More...
 
void cs_source_term_pcvd_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density (vd) defined on primal cells by a value. Case of face-based schemes. More...
 
void cs_source_term_pcsd_bary_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at primal cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. Case of face-based schemes. More...
 
void cs_source_term_pcsd_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by an analytic function. More...
 
void cs_source_term_pcsd_by_array (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at primal cells by an array. More...
 
void cs_source_term_pcvd_bary_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density defined at primal cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. Case of face-based schemes. More...
 
void cs_source_term_pcvd_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector average (va) defined on primal cells by an analytic function. More...
 
void cs_source_term_pcvd_by_array (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector density (vd) defined on primal cells by an array. More...
 
void cs_source_term_hhosd_by_value (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by a value. Case of HHO schemes. More...
 
void cs_source_term_hhosd_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by an analytic function. Case of HHO schemes. More...
 
void cs_source_term_hhovd_by_analytic (const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)
 Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector field (vd) defined on primal cells by an analytic function. Case of HHO schemes. More...
 

Macro Definition Documentation

◆ CS_N_MAX_SOURCE_TERMS

#define CS_N_MAX_SOURCE_TERMS   8

Typedef Documentation

◆ cs_source_term_cellwise_t

typedef void() cs_source_term_cellwise_t(const cs_xdef_t *source, const cs_cell_mesh_t *cm, cs_real_t time_eval, cs_cell_builder_t *cb, void *input, double *values)

Compute the contribution for a cell related to a source term and add it the given array of values.

Parameters
[in]sourcepointer to a cs_xdef_term_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

Function Documentation

◆ cs_source_term_compute_cellwise()

void cs_source_term_compute_cellwise ( const int  n_source_terms,
cs_xdef_t *const *  source_terms,
const cs_cell_mesh_t cm,
const cs_mask_t source_mask,
cs_source_term_cellwise_t compute_source[],
cs_real_t  time_eval,
void *  input,
cs_cell_builder_t cb,
cs_real_t result 
)

Compute the local contributions of source terms in a cell.

Parameters
[in]n_source_termsnumber of source terms
[in]source_termspointer to the definitions of source terms
[in]cmpointer to a cs_cell_mesh_t structure
[in]source_maskarray storing in a compact way which source term is defined in a given cell
[in]compute_sourcearray of function pointers
[in]time_evalphysical time at which one evaluates the term
[in,out]inputpointer to an element cast on-the-fly
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]resultarray storing the result of the evaluation

◆ cs_source_term_compute_from_density()

void cs_source_term_compute_from_density ( cs_flag_t  loc,
const cs_xdef_t source,
cs_real_t  time_eval,
double *  p_values[] 
)

Compute the contribution related to a source term in the case of an input data which is a density.

Parameters
[in]locdescribe where is located the associated DoF
[in]sourcepointer to a cs_xdef_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]p_valuespointer to the computed values (allocated if NULL)

◆ cs_source_term_compute_from_potential()

void cs_source_term_compute_from_potential ( cs_flag_t  loc,
const cs_xdef_t source,
cs_real_t  time_eval,
double *  p_values[] 
)

Compute the contribution related to a source term in the case of an input data which is a potential.

Parameters
[in]locdescribe where is located the associated DoF
[in]sourcepointer to a cs_xdef_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]p_valuespointer to the computed values (allocated if NULL)

◆ cs_source_term_dcsd_bary_by_analytic()

void cs_source_term_dcsd_bary_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcsd_by_array()

void cs_source_term_dcsd_by_array ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an array.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcsd_by_value()

void cs_source_term_dcsd_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by a value.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcsd_q10o2_by_analytic()

void cs_source_term_dcsd_q10o2_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a ten-point quadrature rule to evaluate the integral. Exact for quadratic function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcsd_q1o1_by_analytic()

void cs_source_term_dcsd_q1o1_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a the barycentric approximation as quadrature to evaluate the integral. Exact for linear function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcsd_q5o3_by_analytic()

void cs_source_term_dcsd_q5o3_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at dual cells by an analytical function. Use a five-point quadrature rule to evaluate the integral. Exact for cubic function. This function may be expensive since many evaluations are needed. Please use it with care.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_dcvd_by_value()

void cs_source_term_dcvd_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density defined at dual cells by a value.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_get_flag()

cs_flag_t cs_source_term_get_flag ( const cs_xdef_t st)

Get metadata related to the given source term structure.

Parameters
[in,out]stpointer to a cs_xdef_t structure
Returns
the value of the flag related to this source term

◆ cs_source_term_hhosd_by_analytic()

void cs_source_term_hhosd_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by an analytic function. Case of HHO schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_hhosd_by_value()

void cs_source_term_hhosd_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by a value. Case of HHO schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_hhovd_by_analytic()

void cs_source_term_hhovd_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector field (vd) defined on primal cells by an analytic function. Case of HHO schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_init()

cs_flag_t cs_source_term_init ( cs_param_space_scheme_t  space_scheme,
const int  n_source_terms,
cs_xdef_t *const *  source_terms,
cs_source_term_cellwise_t compute_source[],
cs_flag_t sys_flag,
cs_mask_t source_mask[] 
)

Initialize data to build the source terms.

Parameters
[in]space_schemescheme used to discretize in space
[in]n_source_termsnumber of source terms
[in]source_termspointer to the definitions of source terms
[in,out]compute_sourcearray of function pointers
[in,out]sys_flagmetadata about the algebraic system
[in,out]source_maskpointer to an array storing in a compact way which source term is defined in a given cell
Returns
a flag which indicates what to build in a cell mesh structure

◆ cs_source_term_pcsd_bary_by_analytic()

void cs_source_term_pcsd_bary_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at primal cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. Case of face-based schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in]inputNULL or pointer to a structure cast on-the-fly
[in,out]valuespointer to the computed values

◆ cs_source_term_pcsd_by_analytic()

void cs_source_term_pcsd_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by an analytic function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a scalar density (sd) defined on primal cells by an analytic function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_pcsd_by_array()

void cs_source_term_pcsd_by_array ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density defined at primal cells by an array.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_pcsd_by_value()

void cs_source_term_pcsd_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar density (sd) defined on primal cells by a value. Case of face-based schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_pcvd_bary_by_analytic()

void cs_source_term_pcvd_bary_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density defined at primal cells by an analytical function. Use the barycentric approximation as quadrature to evaluate the integral. Exact for linear function. Case of face-based schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in]inputNULL or pointer to a structure cast on-the-fly
[in,out]valuespointer to the computed values

◆ cs_source_term_pcvd_by_analytic()

void cs_source_term_pcvd_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector average (va) defined on primal cells by an analytic function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector average (va) defined on primal cells by an analytic function.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_pcvd_by_array()

void cs_source_term_pcvd_by_array ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution of a source term for a cell and add it to the given array of values. Case of a vector density (vd) defined on primal cells by an array.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_pcvd_by_value()

void cs_source_term_pcvd_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a vector-valued density (vd) defined on primal cells by a value. Case of face-based schemes.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed value

◆ cs_source_term_pvsp_by_analytic()

void cs_source_term_pvsp_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices by an analytical function. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_pvsp_by_value()

void cs_source_term_pvsp_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices by a constant value. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_set_default_flag()

cs_flag_t cs_source_term_set_default_flag ( cs_param_space_scheme_t  scheme)

Set the default flag related to a source term according to the numerical scheme chosen for discretizing an equation.

Parameters
[in]schemenumerical scheme used for the discretization
Returns
a default flag

◆ cs_source_term_set_reduction()

void cs_source_term_set_reduction ( cs_xdef_t st,
cs_flag_t  flag 
)

Set advanced parameters which are defined by default in a source term structure.

Parameters
[in,out]stpointer to a cs_xdef_t structure
[in]flagCS_FLAG_DUAL or CS_FLAG_PRIMAL

◆ cs_source_term_set_shared_pointers()

void cs_source_term_set_shared_pointers ( const cs_cdo_quantities_t quant,
const cs_cdo_connect_t connect 
)

Set shared pointers to main domain members.

Parameters
[in]quantadditional mesh quantities struct.
[in]connectpointer to a cs_cdo_connect_t struct.

◆ cs_source_term_vcsp_by_analytic()

void cs_source_term_vcsp_by_analytic ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices and cells by an analytical function. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values

◆ cs_source_term_vcsp_by_value()

void cs_source_term_vcsp_by_value ( const cs_xdef_t source,
const cs_cell_mesh_t cm,
cs_real_t  time_eval,
cs_cell_builder_t cb,
void *  input,
double *  values 
)

Compute the contribution for a cell related to a source term and add it the given array of values. Case of a scalar potential defined at primal vertices and cells by a constant value. A discrete Hodge operator has to be computed before this call and stored inside a cs_cell_builder_t structure.

Parameters
[in]sourcepointer to a cs_xdef_t structure
[in]cmpointer to a cs_cell_mesh_t structure
[in]time_evalphysical time at which one evaluates the term
[in,out]cbpointer to a cs_cell_builder_t structure
[in,out]inputpointer to an element cast on-the-fly (or NULL)
[in,out]valuespointer to the computed values