|
My Project
programmer's documentation
|
Go to the documentation of this file. 1 #ifndef __CS_QUADRATURE_H__
2 #define __CS_QUADRATURE_H__
252 gpts[0][0] = 0.5*(v1[0] + v2[0]);
253 gpts[0][1] = 0.5*(v1[1] + v2[1]);
254 gpts[0][2] = 0.5*(v1[2] + v2[2]);
416 gpts[0][0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
417 gpts[0][1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
418 gpts[0][2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
523 xg[0] = .5 * (v1[0] + v2[0]);
524 xg[1] = .5 * (v1[1] + v2[1]);
525 xg[2] = .5 * (v1[2] + v2[2]);
527 ana(tcur, 1, NULL, xg,
false,
input, &evaluation);
529 *results +=
len * evaluation;
558 double evaluation[2], weights[2];
563 ana(tcur, 2, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
566 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1];
595 double evaluation[3], weights[3];
600 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
602 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
603 weights[2] * evaluation[2];
642 ana(tcur, 1, NULL, xg,
false,
input, &evaluation);
644 *results += area * evaluation;
675 double evaluation[3], weights[3];
680 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
683 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
684 weights[2] * evaluation[2];
715 double evaluation[4], weights[4];
720 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
722 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
723 weights[2] * evaluation[2] + weights[3] * evaluation[3];
754 double evaluation[3];
761 ana(tcur, 1, NULL, xg,
false,
input, evaluation);
763 results[0] += area * evaluation[0];
764 results[1] += area * evaluation[1];
765 results[2] += area * evaluation[2];
796 double evaluation[3*3], weights[3];
801 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
803 for (
int p = 0;
p < 3;
p++) {
804 results[0] += weights[
p] * evaluation[3*
p ];
805 results[1] += weights[
p] * evaluation[3*
p+1];
806 results[2] += weights[
p] * evaluation[3*
p+2];
838 double evaluation[3*4], weights[4];
843 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
845 for (
int p = 0;
p < 4;
p++) {
846 results[0] += weights[
p] * evaluation[3*
p ];
847 results[1] += weights[
p] * evaluation[3*
p+1];
848 results[2] += weights[
p] * evaluation[3*
p+2];
880 double evaluation[9];
887 ana(tcur, 1, NULL, xg,
false,
input, evaluation);
889 for (
short int ij = 0; ij < 9; ij++)
890 results[ij] += area * evaluation[ij];
921 double evaluation[9*3], weights[3];
926 ana(tcur, 3, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
928 for (
int p = 0;
p < 3;
p++) {
929 const double wp = weights[
p];
930 double *eval_p = evaluation + 9*
p;
931 for (
short int ij = 0; ij < 9; ij++)
932 results[ij] += wp * eval_p[ij];
964 double evaluation[9*4], weights[4];
969 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
971 for (
int p = 0;
p < 4;
p++) {
972 const double wp = weights[
p];
973 double *eval_p = evaluation + 9*
p;
974 for (
short int ij = 0; ij < 9; ij++)
975 results[ij] += wp * eval_p[ij];
1012 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1013 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1014 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1016 ana(tcur, 1, NULL, xg,
false,
input, &evaluation);
1018 *results += vol * evaluation;
1051 double evaluation[4], weights[4];
1056 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1059 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1060 weights[2] * evaluation[2] + weights[3] * evaluation[3];
1093 double evaluation[5], weights[5];
1098 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1100 *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1101 weights[2] * evaluation[2] + weights[3] * evaluation[3] +
1102 weights[4] * evaluation[4];
1135 double evaluation[3];
1138 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1139 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1140 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1142 ana(tcur, 1, NULL, xg,
false,
input, evaluation);
1144 results[0] += vol * evaluation[0];
1145 results[1] += vol * evaluation[1];
1146 results[2] += vol * evaluation[2];
1179 double evaluation[3*4], weights[4];
1184 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1186 for (
int p = 0;
p < 4;
p++) {
1187 results[0] += weights[
p] * evaluation[3*
p ];
1188 results[1] += weights[
p] * evaluation[3*
p+1];
1189 results[2] += weights[
p] * evaluation[3*
p+2];
1223 double evaluation[3*5], weights[5];
1228 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1230 for (
int p = 0;
p < 5;
p++) {
1231 results[0] += weights[
p] * evaluation[3*
p ];
1232 results[1] += weights[
p] * evaluation[3*
p+1];
1233 results[2] += weights[
p] * evaluation[3*
p+2];
1267 double evaluation[9];
1270 xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1271 xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1272 xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1274 ana(tcur, 1, NULL, xg,
false,
input, evaluation);
1276 for (
short int ij = 0; ij < 9; ij++)
1277 results[ij] += vol * evaluation[ij];
1310 double evaluation[9*4], weights[4];
1315 ana(tcur, 4, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1317 for (
int p = 0;
p < 4;
p++) {
1318 const double wp = weights[
p];
1319 double *eval_p = evaluation + 9*
p;
1320 for (
short int ij = 0; ij < 9; ij++)
1321 results[ij] += wp * eval_p[ij];
1355 double evaluation[9*5], weights[5];
1360 ana(tcur, 5, NULL, (
const cs_real_t *)gauss_pts,
false,
input, evaluation);
1362 for (
int p = 0;
p < 5;
p++) {
1363 const double wp = weights[
p];
1364 double *eval_p = evaluation + 9*
p;
1365 for (
short int ij = 0; ij < 9; ij++)
1366 results[ij] += wp * eval_p[ij];
1402 " %s: Invalid quadrature type\n", __func__);
1420 " %s: Invalid quadrature type\n", __func__);
1438 " %s: Invalid quadrature type\n", __func__);
1444 " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
1484 " %s: Invalid quadrature type\n", __func__);
1502 " %s: Invalid quadrature type\n", __func__);
1520 " %s: Invalid quadrature type\n", __func__);
1526 " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
Definition: cs_quadrature.h:56
Definition: cs_quadrature.h:58
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:1212
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:549
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:911
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.c:200
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:246
size_t len
Definition: mei_scanner.c:560
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:1383
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:1465
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:315
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:744
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:1256
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.c:308
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:870
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:1124
#define END_C_DECLS
Definition: cs_defs.h:468
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
#define BEGIN_C_DECLS
Definition: cs_defs.h:467
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:954
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:828
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:76
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
Definition: cs_quadrature.h:53
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:1299
Definition: cs_quadrature.h:54
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:625
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 sho...
Definition: cs_quadrature.c:498
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.c:351
Definition: cs_quadrature.h:57
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:408
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:705
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:168
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:1344
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:195
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.
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:1082
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.c:435
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:511
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:1040
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:97
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 compact, void *input, cs_real_t *retval)
Generic function pointer for an analytic function elt_ids is optional. If not NULL,...
Definition: cs_param.h:66
const cs_real_t cs_math_1ov3
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:586
void cs_quadrature_setup(void)
Compute constant weights for all quadratures used.
Definition: cs_quadrature.c:107
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:143
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_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:665
unsigned short int cs_flag_t
Definition: cs_defs.h:304
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:998
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:1168
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:119
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:786
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.c:391
cs_quadrature_type_t
Definition: cs_quadrature.h:51
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:313
Definition: cs_quadrature.h:55
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.c:149
void const cs_int_t * type
Definition: cs_measures_util.h:425
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.c:271
Definition: cs_field_pointer.h:67