1 #ifndef __CS_QUADRATURE_H__
2 #define __CS_QUADRATURE_H__
284 gpts[0][0] = 0.5*(v1[0] + v2[0]);
285 gpts[0][1] = 0.5*(v1[1] + v2[1]);
286 gpts[0][2] = 0.5*(v1[2] + v2[2]);
354 gpts[0][0] = c_1ov3 * (v1[0] + v2[0] + v3[0]);
355 gpts[0][1] = c_1ov3 * (v1[1] + v2[1] + v3[1]);
356 gpts[0][2] = c_1ov3 * (v1[2] + v2[2] + v3[2]);
456 gpts[0][0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
457 gpts[0][1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
458 gpts[0][2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
563 xg[0] = .5 * (v1[0] + v2[0]);
564 xg[1] = .5 * (v1[1] + v2[1]);
565 xg[2] = .5 * (v1[2] + v2[2]);
568 ana(tcur, 1, NULL, xg,
false, input, &feval);
571 *results += len * feval;
600 double feval[2], weights[2];
606 ana(tcur, 2, NULL, (
const cs_real_t *)gauss_pts,
false, input, feval);
609 *results += weights[0] * feval[0] + weights[1] * feval[1];
638 double feval[3], weights[3];
644 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false, input, feval);
647 *results += weights[0]*feval[0] + weights[1]*feval[1] + weights[2]*feval[2];
679 xg[0] = .5 * (v1[0] + v2[0]);
680 xg[1] = .5 * (v1[1] + v2[1]);
681 xg[2] = .5 * (v1[2] + v2[2]);
684 ana(tcur, 1, NULL, xg,
false, input, feval);
687 results[0] += len * feval[0];
688 results[1] += len * feval[1];
689 results[2] += len * feval[2];
718 double feval[6], weights[2];
724 ana(tcur, 2, NULL, (
const cs_real_t *)gauss_pts,
false, input, feval);
727 results[0] += weights[0] * feval[0] + weights[1] * feval[3];
728 results[1] += weights[0] * feval[1] + weights[1] * feval[4];
729 results[2] += weights[0] * feval[2] + weights[1] * feval[5];
758 double feval[9], weights[3];
764 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false, input, feval);
767 for (
int p = 0;
p < 3;
p++) {
768 results[0] += weights[
p] * feval[3*
p ];
769 results[1] += weights[
p] * feval[3*
p+1];
770 results[2] += weights[
p] * feval[3*
p+2];
807 xg[0] = c_1ov3 * (v1[0] + v2[0] + v3[0]);
808 xg[1] = c_1ov3 * (v1[1] + v2[1] + v3[1]);
809 xg[2] = c_1ov3 * (v1[2] + v2[2] + v3[2]);
816 ana(tcur, 1, NULL, xg,
false, input, &evaluation);
818 *results += area * evaluation;
849 double evaluation[3], weights[3];
854 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
857 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
858 weights[2] * evaluation[2];
889 double evaluation[4], weights[4];
894 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
896 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
897 weights[2] * evaluation[2] + weights[3] * evaluation[3];
928 double evaluation[7], weights[7];
933 ana(tcur, 7, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
935 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
936 weights[2] * evaluation[2] + weights[3] * evaluation[3] +
937 weights[4] * evaluation[4] + weights[5] * evaluation[5] +
938 weights[6] * evaluation[6] ;
969 double evaluation[3];
974 xg[0] = c_1ov3 * (v1[0] + v2[0] + v3[0]);
975 xg[1] = c_1ov3 * (v1[1] + v2[1] + v3[1]);
976 xg[2] = c_1ov3 * (v1[2] + v2[2] + v3[2]);
983 ana(tcur, 1, NULL, xg,
false, input, evaluation);
985 results[0] += area * evaluation[0];
986 results[1] += area * evaluation[1];
987 results[2] += area * evaluation[2];
1018 double evaluation[3*3], weights[3];
1023 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1025 for (
int p = 0;
p < 3;
p++) {
1026 results[0] += weights[
p] * evaluation[3*
p ];
1027 results[1] += weights[
p] * evaluation[3*
p+1];
1028 results[2] += weights[
p] * evaluation[3*
p+2];
1060 double evaluation[3*4], weights[4];
1065 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1067 for (
int p = 0;
p < 4;
p++) {
1068 results[0] += weights[
p] * evaluation[3*
p ];
1069 results[1] += weights[
p] * evaluation[3*
p+1];
1070 results[2] += weights[
p] * evaluation[3*
p+2];
1102 double evaluation[3*7], weights[7];
1107 ana(tcur, 7, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1109 for (
int p = 0;
p < 7;
p++) {
1110 results[0] += weights[
p] * evaluation[3*
p ];
1111 results[1] += weights[
p] * evaluation[3*
p+1];
1112 results[2] += weights[
p] * evaluation[3*
p+2];
1144 double evaluation[9];
1149 xg[0] = c_1ov3 * (v1[0] + v2[0] + v3[0]);
1150 xg[1] = c_1ov3 * (v1[1] + v2[1] + v3[1]);
1151 xg[2] = c_1ov3 * (v1[2] + v2[2] + v3[2]);
1158 ana(tcur, 1, NULL, xg,
false, input, evaluation);
1160 for (
short int ij = 0; ij < 9; ij++)
1161 results[ij] += area * evaluation[ij];
1192 double evaluation[9*3], weights[3];
1197 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1199 for (
int p = 0;
p < 3;
p++) {
1200 const double wp = weights[
p];
1201 double *eval_p = evaluation + 9*
p;
1202 for (
short int ij = 0; ij < 9; ij++)
1203 results[ij] += wp * eval_p[ij];
1235 double evaluation[9*4], weights[4];
1240 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1242 for (
int p = 0;
p < 4;
p++) {
1243 const double wp = weights[
p];
1244 double *eval_p = evaluation + 9*
p;
1245 for (
short int ij = 0; ij < 9; ij++)
1246 results[ij] += wp * eval_p[ij];
1278 double evaluation[9*7], weights[7];
1283 ana(tcur, 7, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1285 for (
int p = 0;
p < 7;
p++) {
1286 const double wp = weights[
p];
1287 double *eval_p = evaluation + 9*
p;
1288 for (
short int ij = 0; ij < 9; ij++)
1289 results[ij] += wp * eval_p[ij];
1326 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1327 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1328 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1330 ana(tcur, 1, NULL, xg,
false, input, &evaluation);
1332 *results += vol * evaluation;
1365 double evaluation[4], weights[4];
1370 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1373 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1374 weights[2] * evaluation[2] + weights[3] * evaluation[3];
1407 double evaluation[5], weights[5];
1412 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1414 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1415 weights[2] * evaluation[2] + weights[3] * evaluation[3] +
1416 weights[4] * evaluation[4];
1449 double evaluation[3];
1452 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1453 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1454 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1456 ana(tcur, 1, NULL, xg,
false, input, evaluation);
1458 results[0] += vol * evaluation[0];
1459 results[1] += vol * evaluation[1];
1460 results[2] += vol * evaluation[2];
1493 double evaluation[3*4], weights[4];
1498 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1500 for (
int p = 0;
p < 4;
p++) {
1501 results[0] += weights[
p] * evaluation[3*
p ];
1502 results[1] += weights[
p] * evaluation[3*
p+1];
1503 results[2] += weights[
p] * evaluation[3*
p+2];
1537 double evaluation[3*5], weights[5];
1542 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1544 for (
int p = 0;
p < 5;
p++) {
1545 results[0] += weights[
p] * evaluation[3*
p ];
1546 results[1] += weights[
p] * evaluation[3*
p+1];
1547 results[2] += weights[
p] * evaluation[3*
p+2];
1581 double evaluation[9];
1584 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1585 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1586 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1588 ana(tcur, 1, NULL, xg,
false, input, evaluation);
1590 for (
short int ij = 0; ij < 9; ij++)
1591 results[ij] += vol * evaluation[ij];
1624 double evaluation[9*4], weights[4];
1629 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1631 for (
int p = 0;
p < 4;
p++) {
1632 const double wp = weights[
p];
1633 double *eval_p = evaluation + 9*
p;
1634 for (
short int ij = 0; ij < 9; ij++)
1635 results[ij] += wp * eval_p[ij];
1669 double evaluation[9*5], weights[5];
1674 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false, input, evaluation);
1676 for (
int p = 0;
p < 5;
p++) {
1677 const double wp = weights[
p];
1678 double *eval_p = evaluation + 9*
p;
1679 for (
short int ij = 0; ij < 9; ij++)
1680 results[ij] += wp * eval_p[ij];
1717 " %s: Invalid quadrature type\n", __func__);
1735 " %s: Invalid quadrature type\n", __func__);
1741 " %s: Invalid dimension value %d. Only 1 and 3 are valid.\n",
1781 " %s: Invalid quadrature type\n", __func__);
1799 " %s: Invalid quadrature type\n", __func__);
1817 " %s: Invalid quadrature type\n", __func__);
1823 " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
1863 " %s: Invalid quadrature type\n", __func__);
1881 " %s: Invalid quadrature type\n", __func__);
1899 " %s: Invalid quadrature type\n", __func__);
1905 " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
void bft_error(const char *const file_name, const int line_num, const int sys_error_code, const char *const format,...)
Calls the error handler (set by bft_error_handler_set() or default).
Definition: bft_error.c:193
#define BEGIN_C_DECLS
Definition: cs_defs.h:528
double cs_real_t
Floating-point value.
Definition: cs_defs.h:332
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:347
#define END_C_DECLS
Definition: cs_defs.h:529
unsigned short int cs_flag_t
Definition: cs_defs.h:334
@ p
Definition: cs_field_pointer.h:67
unsigned int cs_eflag_t
Definition: cs_flag.h:190
const cs_real_t cs_math_1ov3
void() cs_analytic_func_t(cs_real_t time, cs_lnum_t n_elts, const cs_lnum_t *elt_ids, const cs_real_t *coords, bool dense_output, void *input, cs_real_t *retval)
Generic function pointer for an evaluation relying on an analytic function.
Definition: cs_param_types.h:127
static cs_quadrature_edge_integral_t * cs_quadrature_get_edge_integral(int dim, cs_quadrature_type_t qtype)
Retrieve the integral function according to the quadrature type and the stride provided Case of integ...
Definition: cs_quadrature.h:1698
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).
Definition: cs_quadrature.cxx:398
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-valu...
Definition: cs_quadrature.h:551
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 weig...
Definition: cs_quadrature.h:1482
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).
Definition: cs_quadrature.h:448
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...
Definition: cs_quadrature.h:1396
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 ...
Definition: cs_quadrature.h:629
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...
Definition: cs_quadrature.h:1658
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 ad...
Definition: cs_quadrature.h:879
static void cs_quadrature_tria_7pts_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 7 Gauss points and 7 weights and ad...
Definition: cs_quadrature.h:918
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.
Definition: cs_quadrature.cxx:314
cs_quadrature_type_t
Type of quadrature to use when computing an integral quantity. This rationale is used for integrals a...
Definition: cs_quadrature.h:84
@ CS_QUADRATURE_HIGHEST
Definition: cs_quadrature.h:90
@ CS_QUADRATURE_NONE
Definition: cs_quadrature.h:86
@ CS_QUADRATURE_BARY_SUBDIV
Definition: cs_quadrature.h:88
@ CS_QUADRATURE_N_TYPES
Definition: cs_quadrature.h:91
@ CS_QUADRATURE_HIGHER
Definition: cs_quadrature.h:89
@ CS_QUADRATURE_BARY
Definition: cs_quadrature.h:87
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.
Definition: cs_quadrature.h:201
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....
Definition: cs_quadrature.h:1570
static void cs_quadrature_edge_2pts_vect(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...
Definition: cs_quadrature.h:709
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 ad...
Definition: cs_quadrature.h:1225
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
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 or...
Definition: cs_quadrature.h:278
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 ad...
Definition: cs_quadrature.h:1050
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....
Definition: cs_quadrature.h:1438
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...
Definition: cs_quadrature.h:1526
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 or...
Definition: cs_quadrature.cxx:206
void cs_quadrature_setup(void)
Compute constant weights for all quadratures used.
Definition: cs_quadrature.cxx:108
static void cs_quadrature_tria_7pts_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 7 Gauss points and 7 weights and ad...
Definition: cs_quadrature.h:1092
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 or...
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...
Definition: cs_quadrature.h:591
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....
Definition: cs_quadrature.h:1312
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 ...
Definition: cs_quadrature.h:1008
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 o...
Definition: cs_quadrature.h:792
const char * cs_quadrature_get_type_name(const cs_quadrature_type_t type)
Return th name associated to a type of quadrature.
Definition: cs_quadrature.cxx:155
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 weig...
Definition: cs_quadrature.h:1354
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 (baryc...
Definition: cs_quadrature.h:345
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.
Definition: cs_quadrature.cxx:277
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.
Definition: cs_quadrature.h:228
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.
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).
Definition: cs_quadrature.cxx:444
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 o...
Definition: cs_quadrature.h:959
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 ...
Definition: cs_quadrature.h:839
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.
Definition: cs_quadrature.h:109
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.
Definition: cs_quadrature.h:1762
static void cs_quadrature_edge_1pt_vect(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 using a barycentric quadrature rule and add it to results Case of a...
Definition: cs_quadrature.h:667
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
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.
Definition: cs_quadrature.h:176
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.
Definition: cs_quadrature.h:1844
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 ...
Definition: cs_quadrature.h:1182
static void cs_quadrature_tria_7pts_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 7 Gauss points and 7 weights and ad...
Definition: cs_quadrature.h:1268
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).
Definition: cs_quadrature.cxx:356
cs_eflag_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 sho...
Definition: cs_quadrature.cxx:507
static void cs_quadrature_edge_3pts_vect(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 ...
Definition: cs_quadrature.h:749
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....
Definition: cs_quadrature.h:1134
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 weig...
Definition: cs_quadrature.h:1613