|
programmer's documentation
|
#include <bft_error.h>#include "cs_base.h"#include "cs_defs.h"#include "cs_math.h"#include "cs_param.h"
Go to the source code of this file.
Typedefs | |
| typedef void() | cs_quadrature_edge_t(const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_real_3_t gpts[], double *weights) |
| Generic function pointer to compute the quadrature points for an edge from v1 -> v2. More... | |
| typedef 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. More... | |
| typedef 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. More... | |
| typedef void() | cs_quadrature_edge_integral_t(double tcur, const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over an edge based on a specified quadrature rule and add it to results. More... | |
| typedef void() | cs_quadrature_tria_integral_t(double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle based on a specified quadrature rule and add it to results. More... | |
| typedef void() | cs_quadrature_tetra_integral_t(double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron based on a specified quadrature rule and add it to results. More... | |
Enumerations | |
| enum | cs_quadrature_type_t { CS_QUADRATURE_NONE, CS_QUADRATURE_BARY, CS_QUADRATURE_BARY_SUBDIV, CS_QUADRATURE_HIGHER, CS_QUADRATURE_HIGHEST, CS_QUADRATURE_N_TYPES } |
Functions | |
| void | cs_quadrature_setup (void) |
| Compute constant weights for all quadratures used. More... | |
| const char * | cs_quadrature_get_type_name (const cs_quadrature_type_t type) |
| Return th name associated to a type of quadrature. More... | |
| static void | cs_quadrature_edge_1pt (const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_real_3_t gpts[], double *w) |
| Compute quadrature points for an edge from v1 -> v2 (2 points) Exact for polynomial function up to order 1. More... | |
| void | cs_quadrature_edge_2pts (const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_real_3_t gpts[], double *w) |
| Compute quadrature points for an edge from v1 -> v2 (2 points) Exact for polynomial function up to order 3. More... | |
| void | cs_quadrature_edge_3pts (const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_real_3_t gpts[], double w[]) |
| Compute quadrature points for an edge from v1 -> v2 (3 points) Exact for polynomial function up to order 5. More... | |
| static void | cs_quadrature_tria_1pt (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 *w) |
| Compute quadrature points for a triangle (1 point) Exact for polynomial function up to order 1 (barycentric approx.) More... | |
| void | cs_quadrature_tria_3pts (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 *w) |
| Compute quadrature points for a triangle (3 points) Exact for polynomial function up to order 2. More... | |
| void | cs_quadrature_tria_4pts (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 w[]) |
| Compute quadrature points for a triangle (4 points) Exact for polynomial function up to order 3. More... | |
| void | cs_quadrature_tria_7pts (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 w[]) |
| Compute quadrature points for a triangle (7 points) Exact for polynomial function up to order 5. More... | |
| static void | cs_quadrature_tet_1pt (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 weight[]) |
| Compute the quadrature in a tetrehedra. Exact for 1st order polynomials (order 2). More... | |
| void | cs_quadrature_tet_4pts (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[]) |
| Compute the quadrature in a tetrehedra. Exact for 2nd order polynomials (order 3). More... | |
| void | cs_quadrature_tet_5pts (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[]) |
| Compute the quadrature in a tetrehedra. Exact for 3rd order polynomials (order 4). More... | |
| void | cs_quadrature_tet_15pts (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[]) |
| Compute the quadrature in a tetrehedra. Exact for 5th order polynomials (order 6). More... | |
| static void | cs_quadrature_edge_1pt_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over an edge with the mid-point rule and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_edge_2pts_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over an edge with a quadrature rule using 2 Gauss points and a unique weight and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_edge_3pts_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over an edge with a quadrature rule using 3 Gauss points and weights and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tria_1pt_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tria_3pts_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tria_4pts_scal (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tria_1pt_vect (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results Case of a vector-valued function. More... | |
| static void | cs_quadrature_tria_3pts_vect (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results Case of a vector-valued function. More... | |
| static void | cs_quadrature_tria_4pts_vect (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results. Case of a vector-valued function. More... | |
| static void | cs_quadrature_tria_1pt_tens (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results. Case of a tensor-valued function. More... | |
| static void | cs_quadrature_tria_3pts_tens (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results. Case of a tensor-valued function. More... | |
| static void | cs_quadrature_tria_4pts_tens (double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results. Case of a tensor-valued function. More... | |
| static void | cs_quadrature_tet_1pt_scal (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tet_4pts_scal (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tet_5pts_scal (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results. Case of a scalar-valued function. More... | |
| static void | cs_quadrature_tet_1pt_vect (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a vector-valued function. More... | |
| static void | cs_quadrature_tet_4pts_vect (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a vector-valued function. More... | |
| static void | cs_quadrature_tet_5pts_vect (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results. Case of a vector-valued function. More... | |
| static void | cs_quadrature_tet_1pt_tens (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a tensor-valued function. More... | |
| static void | cs_quadrature_tet_4pts_tens (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a tensor-valued function. More... | |
| static void | cs_quadrature_tet_5pts_tens (double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results Case of a tensor-valued function. More... | |
| static cs_quadrature_tria_integral_t * | cs_quadrature_get_tria_integral (int dim, cs_quadrature_type_t qtype) |
| Retrieve the integral function according to the quadrature type and the stride provided. More... | |
| static cs_quadrature_tetra_integral_t * | cs_quadrature_get_tetra_integral (int dim, cs_quadrature_type_t qtype) |
| Retrieve the integral function according to the quadrature type and the stride provided. More... | |
| cs_flag_t | cs_quadrature_get_flag (const cs_quadrature_type_t qtype, const cs_flag_t loc) |
Get the flags adapted to the given quadrature type qtype and the location on which the quadrature should be performed. More... | |
| typedef void() cs_quadrature_edge_integral_t(double tcur, const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over an edge based on a specified quadrature rule and add it to results.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the edge |
| [in] | v2 | second point of the edge |
| [in] | len | length of the edge |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| typedef void() cs_quadrature_edge_t(const cs_real_3_t v1, const cs_real_3_t v2, double len, cs_real_3_t gpts[], double *weights) |
Generic function pointer to compute the quadrature points for an edge from v1 -> v2.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | len | length of edge [v1, v2] |
| [in,out] | gpts | gauss points |
| [in,out] | w | weight (same weight for the two points) |
| typedef 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.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | v4 | fourth vertex |
| [in] | vol | volume of tetrahedron {v1, v2, v3, v4} |
| [in,out] | gpts | Gauss points |
| [in,out] | weights | weigths related to each Gauss point |
| typedef void() cs_quadrature_tetra_integral_t(double tcur, 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_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a tetrahedron based on a specified quadrature rule and add it to results.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| typedef void() cs_quadrature_tria_integral_t(double tcur, const cs_real_3_t v1, const cs_real_3_t v2, const cs_real_3_t v3, double area, cs_analytic_func_t *ana, void *input, double results[]) |
Compute the integral over a triangle based on a specified quadrature rule and add it to results.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| typedef 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.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | area | area of triangle {v1, v2, v3} |
| [in,out] | gpts | Gauss points |
| [in,out] | weights | weights values |
| enum cs_quadrature_type_t |
|
inlinestatic |
Compute quadrature points for an edge from v1 -> v2 (2 points) Exact for polynomial function up to order 1.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | len | length of edge [v1, v2] |
| [in,out] | gpts | gauss points |
| [in,out] | w | weight (same weight for the two points) |
|
inlinestatic |
Compute the integral over an edge with the mid-point rule and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the edge |
| [in] | v2 | second point of the edge |
| [in] | len | length of the edge |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_edge_2pts | ( | const cs_real_3_t | v1, |
| const cs_real_3_t | v2, | ||
| double | len, | ||
| cs_real_3_t | gpts[], | ||
| double * | w | ||
| ) |
Compute quadrature points for an edge from v1 -> v2 (2 points) Exact for polynomial function up to order 3.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | len | length of edge [v1, v2] |
| [in,out] | gpts | gauss points |
| [in,out] | w | weight (same weight for the two points) |
|
inlinestatic |
Compute the integral over an edge with a quadrature rule using 2 Gauss points and a unique weight and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the edge |
| [in] | v2 | second point of the edge |
| [in] | len | length of the edge |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_edge_3pts | ( | const cs_real_3_t | v1, |
| const cs_real_3_t | v2, | ||
| double | len, | ||
| cs_real_3_t | gpts[], | ||
| double | w[] | ||
| ) |
Compute quadrature points for an edge from v1 -> v2 (3 points) Exact for polynomial function up to order 5.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | len | length of edge [v1, v2] |
| [in,out] | gpts | gauss points |
| [in,out] | w | weights |
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | len | length of edge [v1, v2] |
| [in,out] | gpts | gauss points |
| [in,out] | w | weights |
|
inlinestatic |
Compute the integral over an edge with a quadrature rule using 3 Gauss points and weights and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the edge |
| [in] | v2 | second point of the edge |
| [in] | len | length of the edge |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| cs_flag_t cs_quadrature_get_flag | ( | const cs_quadrature_type_t | qtype, |
| const cs_flag_t | loc | ||
| ) |
Get the flags adapted to the given quadrature type qtype and the location on which the quadrature should be performed.
| [in] | qtype | cs_quadrature_type_t |
| [in] | loc | It could be CS_FLAG_CELL, CS_FLAG_FACE or CS_FLAG_EDGE plus CS_FLAG_PRIMAL or CS_FLAG_DUAL |
|
inlinestatic |
Retrieve the integral function according to the quadrature type and the stride provided.
| [in] | dim | dimension of the function to integrate |
| [in] | qtype | quadrature type |
|
inlinestatic |
Retrieve the integral function according to the quadrature type and the stride provided.
| [in] | dim | dimension of the function to integrate |
| [in] | qtype | quadrature type |
| const char* cs_quadrature_get_type_name | ( | const cs_quadrature_type_t | type | ) |
Return th name associated to a type of quadrature.
| [in] | type | cs_quadrature_type_t |
| void cs_quadrature_setup | ( | void | ) |
Compute constant weights for all quadratures used.
| void cs_quadrature_tet_15pts | ( | 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[] | ||
| ) |
Compute the quadrature in a tetrehedra. Exact for 5th order polynomials (order 6).
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | v4 | fourth vertex |
| [in] | vol | volume of tetrahedron {v1, v2, v3, v4} |
| [in,out] | gpts | 15 Gauss points (size = 3*15) |
| [in,out] | weights | 15 weigths related to each Gauss point |
|
inlinestatic |
Compute the quadrature in a tetrehedra. Exact for 1st order polynomials (order 2).
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | v4 | fourth vertex |
| [in] | vol | volume of tetrahedron |
| [in,out] | gpts | 1 Gauss point |
| [in,out] | weight | weight related to this point (= volume) |
|
inlinestatic |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron using a barycentric quadrature rule and add it to results. Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_tet_4pts | ( | 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[] | ||
| ) |
Compute the quadrature in a tetrehedra. Exact for 2nd order polynomials (order 3).
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | v4 | fourth vertex |
| [in] | vol | volume of tetrahedron {v1, v2, v3, v4} |
| [in,out] | gpts | 4 Gauss points (size = 3*4) |
| [in,out] | weights | weight (same value for all points) |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 4 Gauss points and a unique weight and add it to results. Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_tet_5pts | ( | 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[] | ||
| ) |
Compute the quadrature in a tetrehedra. Exact for 3rd order polynomials (order 4).
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | v4 | fourth vertex |
| [in] | vol | volume of tetrahedron {v1, v2, v3, v4} |
| [in,out] | gpts | 5 Gauss points (size = 3*5) |
| [in,out] | weights | 5 weigths related to each Gauss point |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results. Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a tetrahedron with a quadrature rule using 5 Gauss points and 5 weights and add it to results. Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the tetrahedron |
| [in] | v2 | second point of the tetrahedron |
| [in] | v3 | third point of the tetrahedron |
| [in] | v4 | fourth point of the tetrahedron |
| [in] | vol | volume of the tetrahedron |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute quadrature points for a triangle (1 point) Exact for polynomial function up to order 1 (barycentric approx.)
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | area | area of triangle {v1, v2, v3} |
| [in,out] | gpts | gauss points |
| [in,out] | w | weight |
|
inlinestatic |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results. Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle using a barycentric quadrature rule and add it to results Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | 1st point of the triangle |
| [in] | v2 | 2nd point of the triangle |
| [in] | v3 | 3rd point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_tria_3pts | ( | 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 * | w | ||
| ) |
Compute quadrature points for a triangle (3 points) Exact for polynomial function up to order 2.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | area | area of triangle {v1, v2, v3} |
| [in,out] | gpts | gauss points |
| [in,out] | w | weight (same weight for the three points) |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results. Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 3 Gauss points and a unique weight and add it to results Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_tria_4pts | ( | 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 | w[] | ||
| ) |
Compute quadrature points for a triangle (4 points) Exact for polynomial function up to order 3.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | area | area of triangle {v1, v2, v3} |
| [in,out] | gpts | gauss points |
| [in,out] | w | weights |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results Case of a scalar-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results. Case of a tensor-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
|
inlinestatic |
Compute the integral over a triangle with a quadrature rule using 4 Gauss points and 4 weights and add it to results. Case of a vector-valued function.
| [in] | tcur | current physical time of the simulation |
| [in] | v1 | first point of the triangle |
| [in] | v2 | second point of the triangle |
| [in] | v3 | third point of the triangle |
| [in] | area | area of the triangle |
| [in] | ana | pointer to the analytic function |
| [in] | input | NULL or pointer to a structure cast on-the-fly |
| [in,out] | results | array of double |
| void cs_quadrature_tria_7pts | ( | 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 | w[] | ||
| ) |
Compute quadrature points for a triangle (7 points) Exact for polynomial function up to order 5.
| [in] | v1 | first vertex |
| [in] | v2 | second vertex |
| [in] | v3 | third vertex |
| [in] | area | area of triangle {v1, v2, v3} |
| [in,out] | gpts | gauss points |
| [in,out] | w | weights |
1.8.13