programmer's documentation
cs_quadrature.h
Go to the documentation of this file.
1 #ifndef __CS_QUADRATURE_H__
2 #define __CS_QUADRATURE_H__
3 
4 /*============================================================================
5  * Routines to handle quadrature rules
6  *============================================================================*/
7 
8 /*
9  This file is part of Code_Saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2018 EDF S.A.
12 
13  This program is free software; you can redistribute it and/or modify it under
14  the terms of the GNU General Public License as published by the Free Software
15  Foundation; either version 2 of the License, or (at your option) any later
16  version.
17 
18  This program is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21  details.
22 
23  You should have received a copy of the GNU General Public License along with
24  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25  Street, Fifth Floor, Boston, MA 02110-1301, USA.
26 */
27 
28 /*----------------------------------------------------------------------------
29  * Local headers
30  *----------------------------------------------------------------------------*/
31 
32 #include <bft_error.h>
33 
34 #include "cs_base.h"
35 #include "cs_defs.h"
36 #include "cs_math.h"
37 #include "cs_param.h"
38 
39 /*----------------------------------------------------------------------------*/
40 
42 
43 /*============================================================================
44  * Macro definitions
45  *============================================================================*/
46 
47 /*============================================================================
48  * Type definitions
49  *============================================================================*/
50 
51 typedef enum {
52 
54  CS_QUADRATURE_BARY, /* Value at the barycenter * meas */
55  CS_QUADRATURE_BARY_SUBDIV, /* Value at the barycenter * meas on a sub-mesh */
56  CS_QUADRATURE_HIGHER, /* Unique weight but several Gauss points */
57  CS_QUADRATURE_HIGHEST, /* Specific weight for each Gauss points */
59 
61 
62 /*----------------------------------------------------------------------------*/
73 /*----------------------------------------------------------------------------*/
74 
75 typedef void
77  const cs_real_3_t v2,
78  double len,
79  cs_real_3_t gpts[],
80  double *weights);
81 
82 /*----------------------------------------------------------------------------*/
94 /*----------------------------------------------------------------------------*/
95 
96 typedef void
98  const cs_real_3_t v2,
99  const cs_real_3_t v3,
100  double area,
101  cs_real_3_t gpts[],
102  double *weights);
103 
104 /*----------------------------------------------------------------------------*/
116 /*----------------------------------------------------------------------------*/
117 
118 typedef void
120  const cs_real_3_t v2,
121  const cs_real_3_t v3,
122  const cs_real_3_t v4,
123  double vol,
124  cs_real_3_t gpts[],
125  double weights[]);
126 
127 /*----------------------------------------------------------------------------*/
140 /*----------------------------------------------------------------------------*/
141 
142 typedef void
144  const cs_real_3_t v1,
145  const cs_real_3_t v2,
146  double len,
147  cs_analytic_func_t *ana,
148  void *input,
149  double results[]);
150 
151 /*----------------------------------------------------------------------------*/
165 /*----------------------------------------------------------------------------*/
166 
167 typedef void
169  const cs_real_3_t v1,
170  const cs_real_3_t v2,
171  const cs_real_3_t v3,
172  double area,
173  cs_analytic_func_t *ana,
174  void *input,
175  double results[]);
176 
177 /*----------------------------------------------------------------------------*/
192 /*----------------------------------------------------------------------------*/
193 
194 typedef void
196  const cs_real_3_t v1,
197  const cs_real_3_t v2,
198  const cs_real_3_t v3,
199  const cs_real_3_t v4,
200  double vol,
201  cs_analytic_func_t *ana,
202  void *input,
203  double results[]);
204 
205 
206 /*============================================================================
207  * Public function prototypes
208  *============================================================================*/
209 
210 /*----------------------------------------------------------------------------*/
214 /*----------------------------------------------------------------------------*/
215 
216 void
217 cs_quadrature_setup(void);
218 
219 /*----------------------------------------------------------------------------*/
227 /*----------------------------------------------------------------------------*/
228 
229 const char *
231 
232 /*----------------------------------------------------------------------------*/
243 /*----------------------------------------------------------------------------*/
244 
245 static inline void
247  const cs_real_3_t v2,
248  double len,
249  cs_real_3_t gpts[],
250  double *w)
251 {
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]);
255  w[0] = len;
256 }
257 
258 /*----------------------------------------------------------------------------*/
269 /*----------------------------------------------------------------------------*/
270 
271 void
273  const cs_real_3_t v2,
274  double len,
275  cs_real_3_t gpts[],
276  double *w);
277 
278 /*----------------------------------------------------------------------------*/
289 /*----------------------------------------------------------------------------*/
290 
291 void
293  const cs_real_3_t v2,
294  double len,
295  cs_real_3_t gpts[],
296  double w[]);
297 
298 /*----------------------------------------------------------------------------*/
310 /*----------------------------------------------------------------------------*/
311 
312 static inline void
314  const cs_real_3_t v2,
315  const cs_real_3_t v3,
316  double area,
317  cs_real_3_t gpts[],
318  double *w)
319 {
320  gpts[0][0] = cs_math_onethird * (v1[0] + v2[0] + v3[0]);
321  gpts[0][1] = cs_math_onethird * (v1[1] + v2[1] + v3[1]);
322  gpts[0][2] = cs_math_onethird * (v1[2] + v2[2] + v3[2]);
323  w[0] = area;
324 }
325 
326 /*----------------------------------------------------------------------------*/
338 /*----------------------------------------------------------------------------*/
339 
340 void
342  const cs_real_3_t v2,
343  const cs_real_3_t v3,
344  double area,
345  cs_real_3_t gpts[],
346  double *w);
347 
348 /*----------------------------------------------------------------------------*/
360 /*----------------------------------------------------------------------------*/
361 
362 void
364  const cs_real_3_t v2,
365  const cs_real_3_t v3,
366  double area,
367  cs_real_3_t gpts[],
368  double w[]);
369 
370 /*----------------------------------------------------------------------------*/
382 /*----------------------------------------------------------------------------*/
383 
384 void
386  const cs_real_3_t v2,
387  const cs_real_3_t v3,
388  double area,
389  cs_real_3_t gpts[],
390  double w[]);
391 
392 /*----------------------------------------------------------------------------*/
405 /*----------------------------------------------------------------------------*/
406 
407 static inline void
409  const cs_real_3_t v2,
410  const cs_real_3_t v3,
411  const cs_real_3_t v4,
412  double vol,
413  cs_real_3_t gpts[],
414  double weight[])
415 {
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]);
419  weight[0] = vol;
420 }
421 
422 /*----------------------------------------------------------------------------*/
435 /*----------------------------------------------------------------------------*/
436 
437 void
439  const cs_real_3_t v2,
440  const cs_real_3_t v3,
441  const cs_real_3_t v4,
442  double vol,
443  cs_real_3_t gpts[],
444  double weights[]);
445 
446 /*----------------------------------------------------------------------------*/
459 /*----------------------------------------------------------------------------*/
460 
461 void
463  const cs_real_3_t v2,
464  const cs_real_3_t v3,
465  const cs_real_3_t v4,
466  double vol,
467  cs_real_3_t gpts[],
468  double weights[]);
469 
470 /*----------------------------------------------------------------------------*/
483 /*----------------------------------------------------------------------------*/
484 
485 void
487  const cs_real_3_t v2,
488  const cs_real_3_t v3,
489  const cs_real_3_t v4,
490  double vol,
491  cs_real_3_t gpts[],
492  double weights[]);
493 
494 /*----------------------------------------------------------------------------*/
508 /*----------------------------------------------------------------------------*/
509 
510 static inline void
512  const cs_real_3_t v1,
513  const cs_real_3_t v2,
514  double len,
515  cs_analytic_func_t *ana,
516  void *input,
517  double results[])
518 {
519  cs_real_3_t xg;
520  double evaluation;
521 
522  // Copied from cs_quadrature_1pt
523  xg[0] = .5 * (v1[0] + v2[0]);
524  xg[1] = .5 * (v1[1] + v2[1]);
525  xg[2] = .5 * (v1[2] + v2[2]);
526 
527  ana(tcur, 1, NULL, xg, false, input, &evaluation);
528 
529  *results += len * evaluation;
530 }
531 
532 /*----------------------------------------------------------------------------*/
546 /*----------------------------------------------------------------------------*/
547 
548 static inline void
550  const cs_real_3_t v1,
551  const cs_real_3_t v2,
552  double len,
553  cs_analytic_func_t *ana,
554  void *input,
555  double results[])
556 {
557  cs_real_3_t gauss_pts[2];
558  double evaluation[2], weights[2];
559 
560  /* Compute Gauss points and its unique weight */
561  cs_quadrature_edge_2pts(v1, v2, len, gauss_pts, weights);
562 
563  ana(tcur, 2, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
564 
565  /* Return results */
566  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1];
567 }
568 
569 /*----------------------------------------------------------------------------*/
583 /*----------------------------------------------------------------------------*/
584 
585 static inline void
587  const cs_real_3_t v1,
588  const cs_real_3_t v2,
589  double len,
590  cs_analytic_func_t *ana,
591  void *input,
592  double results[])
593 {
594  cs_real_3_t gauss_pts[3];
595  double evaluation[3], weights[3];
596 
597  /* Compute Gauss points and its weights */
598  cs_quadrature_edge_3pts(v1, v2, len, gauss_pts, weights);
599 
600  ana(tcur, 3, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
601 
602  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
603  weights[2] * evaluation[2];
604 }
605 
606 
607 /*----------------------------------------------------------------------------*/
622 /*----------------------------------------------------------------------------*/
623 
624 static inline void
626  const cs_real_3_t v1,
627  const cs_real_3_t v2,
628  const cs_real_3_t v3,
629  double area,
630  cs_analytic_func_t *ana,
631  void *input,
632  double results[])
633 {
634  cs_real_3_t xg;
635  double evaluation;
636 
637  /* Copied from cs_quadrature_1pt */
638  xg[0] = cs_math_onethird * (v1[0] + v2[0] + v3[0]);
639  xg[1] = cs_math_onethird * (v1[1] + v2[1] + v3[1]);
640  xg[2] = cs_math_onethird * (v1[2] + v2[2] + v3[2]);
641 
642  ana(tcur, 1, NULL, xg, false, input, &evaluation);
643 
644  *results += area * evaluation;
645 }
646 
647 /*----------------------------------------------------------------------------*/
662 /*----------------------------------------------------------------------------*/
663 
664 static inline void
666  const cs_real_3_t v1,
667  const cs_real_3_t v2,
668  const cs_real_3_t v3,
669  double area,
670  cs_analytic_func_t *ana,
671  void *input,
672  double results[])
673 {
674  cs_real_3_t gauss_pts[3];
675  double evaluation[3], weights[3];
676 
677  /* Compute Gauss points and its unique weight */
678  cs_quadrature_tria_3pts(v1, v2, v3, area, gauss_pts, weights);
679 
680  ana(tcur, 3, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
681 
682  /* Return results */
683  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
684  weights[2] * evaluation[2];
685 }
686 
687 /*----------------------------------------------------------------------------*/
702 /*----------------------------------------------------------------------------*/
703 
704 static inline void
706  const cs_real_3_t v1,
707  const cs_real_3_t v2,
708  const cs_real_3_t v3,
709  double area,
710  cs_analytic_func_t *ana,
711  void *input,
712  double results[])
713 {
714  cs_real_3_t gauss_pts[4];
715  double evaluation[4], weights[4];
716 
717  /* Compute Gauss points and its weights */
718  cs_quadrature_tria_4pts(v1, v2, v3, area, gauss_pts, weights);
719 
720  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
721 
722  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
723  weights[2] * evaluation[2] + weights[3] * evaluation[3];
724 }
725 
726 /*----------------------------------------------------------------------------*/
740 /*----------------------------------------------------------------------------*/
741 
742 static inline void
744  const cs_real_3_t v1,
745  const cs_real_3_t v2,
746  const cs_real_3_t v3,
747  double area,
748  cs_analytic_func_t *ana,
749  void *input,
750  double results[])
751 {
752  cs_real_3_t xg;
753  double evaluation[3];
754 
755  /* Copied from cs_quadrature_1pt */
756  xg[0] = cs_math_onethird * (v1[0] + v2[0] + v3[0]);
757  xg[1] = cs_math_onethird * (v1[1] + v2[1] + v3[1]);
758  xg[2] = cs_math_onethird * (v1[2] + v2[2] + v3[2]);
759 
760  ana(tcur, 1, NULL, xg, false, input, evaluation);
761 
762  results[0] += area * evaluation[0];
763  results[1] += area * evaluation[1];
764  results[2] += area * evaluation[2];
765 }
766 
767 /*----------------------------------------------------------------------------*/
782 /*----------------------------------------------------------------------------*/
783 
784 static inline void
786  const cs_real_3_t v1,
787  const cs_real_3_t v2,
788  const cs_real_3_t v3,
789  double area,
790  cs_analytic_func_t *ana,
791  void *input,
792  double results[])
793 {
794  cs_real_3_t gauss_pts[3];
795  double evaluation[3*3], weights[3];
796 
797  /* Compute Gauss points and its unique weight */
798  cs_quadrature_tria_3pts(v1, v2, v3, area, gauss_pts, weights);
799 
800  ana(tcur, 3, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
801 
802  for (int p = 0; p < 3; p++) {
803  results[0] += weights[p] * evaluation[3*p ];
804  results[1] += weights[p] * evaluation[3*p+1];
805  results[2] += weights[p] * evaluation[3*p+2];
806  }
807 }
808 
809 /*----------------------------------------------------------------------------*/
824 /*----------------------------------------------------------------------------*/
825 
826 static inline void
828  const cs_real_3_t v1,
829  const cs_real_3_t v2,
830  const cs_real_3_t v3,
831  double area,
832  cs_analytic_func_t *ana,
833  void *input,
834  double results[])
835 {
836  cs_real_3_t gauss_pts[4];
837  double evaluation[3*4], weights[4];
838 
839  /* Compute Gauss points and its weights */
840  cs_quadrature_tria_4pts(v1, v2, v3, area, gauss_pts, weights);
841 
842  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
843 
844  for (int p = 0; p < 4; p++) {
845  results[0] += weights[p] * evaluation[3*p ];
846  results[1] += weights[p] * evaluation[3*p+1];
847  results[2] += weights[p] * evaluation[3*p+2];
848  }
849 }
850 
851 /*----------------------------------------------------------------------------*/
866 /*----------------------------------------------------------------------------*/
867 
868 static inline void
870  const cs_real_3_t v1,
871  const cs_real_3_t v2,
872  const cs_real_3_t v3,
873  double area,
874  cs_analytic_func_t *ana,
875  void *input,
876  double results[])
877 {
878  cs_real_3_t xg;
879  double evaluation[9];
880 
881  /* Copied from cs_quadrature_1pt */
882  xg[0] = cs_math_onethird * (v1[0] + v2[0] + v3[0]);
883  xg[1] = cs_math_onethird * (v1[1] + v2[1] + v3[1]);
884  xg[2] = cs_math_onethird * (v1[2] + v2[2] + v3[2]);
885 
886  ana(tcur, 1, NULL, xg, false, input, evaluation);
887 
888  for (short int ij = 0; ij < 9; ij++)
889  results[ij] += area * evaluation[ij];
890 }
891 
892 /*----------------------------------------------------------------------------*/
907 /*----------------------------------------------------------------------------*/
908 
909 static inline void
911  const cs_real_3_t v1,
912  const cs_real_3_t v2,
913  const cs_real_3_t v3,
914  double area,
915  cs_analytic_func_t *ana,
916  void *input,
917  double results[])
918 {
919  cs_real_3_t gauss_pts[3];
920  double evaluation[9*3], weights[3];
921 
922  /* Compute Gauss points and its unique weight */
923  cs_quadrature_tria_3pts(v1, v2, v3, area, gauss_pts, weights);
924 
925  ana(tcur, 3, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
926 
927  for (int p = 0; p < 3; p++) {
928  const double wp = weights[p];
929  double *eval_p = evaluation + 9*p;
930  for (short int ij = 0; ij < 9; ij++)
931  results[ij] += wp * eval_p[ij];
932  }
933 }
934 
935 /*----------------------------------------------------------------------------*/
950 /*----------------------------------------------------------------------------*/
951 
952 static inline void
954  const cs_real_3_t v1,
955  const cs_real_3_t v2,
956  const cs_real_3_t v3,
957  double area,
958  cs_analytic_func_t *ana,
959  void *input,
960  double results[])
961 {
962  cs_real_3_t gauss_pts[4];
963  double evaluation[9*4], weights[4];
964 
965  /* Compute Gauss points and its weights */
966  cs_quadrature_tria_4pts(v1, v2, v3, area, gauss_pts, weights);
967 
968  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
969 
970  for (int p = 0; p < 4; p++) {
971  const double wp = weights[p];
972  double *eval_p = evaluation + 9*p;
973  for (short int ij = 0; ij < 9; ij++)
974  results[ij] += wp * eval_p[ij];
975  }
976 }
977 
978 /*----------------------------------------------------------------------------*/
994 /*----------------------------------------------------------------------------*/
995 
996 static inline void
998  const cs_real_3_t v1,
999  const cs_real_3_t v2,
1000  const cs_real_3_t v3,
1001  const cs_real_3_t v4,
1002  double vol,
1003  cs_analytic_func_t *ana,
1004  void *input,
1005  double results[])
1006 {
1007  cs_real_3_t xg;
1008  double evaluation;
1009 
1010  /* Copied from cs_quadrature_tet_1pt */
1011  xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1012  xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1013  xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1014 
1015  ana(tcur, 1, NULL, xg, false, input, &evaluation);
1016 
1017  *results += vol * evaluation;
1018 }
1019 
1020 /*----------------------------------------------------------------------------*/
1036 /*----------------------------------------------------------------------------*/
1037 
1038 static inline void
1040  const cs_real_3_t v1,
1041  const cs_real_3_t v2,
1042  const cs_real_3_t v3,
1043  const cs_real_3_t v4,
1044  double vol,
1045  cs_analytic_func_t *ana,
1046  void *input,
1047  double results[])
1048 {
1049  cs_real_3_t gauss_pts[4];
1050  double evaluation[4], weights[4];
1051 
1052  /* Compute Gauss points and its unique weight */
1053  cs_quadrature_tet_4pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1054 
1055  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1056 
1057  /* Return results */
1058  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1059  weights[2] * evaluation[2] + weights[3] * evaluation[3];
1060 }
1061 
1062 /*----------------------------------------------------------------------------*/
1078 /*----------------------------------------------------------------------------*/
1079 
1080 static inline void
1082  const cs_real_3_t v1,
1083  const cs_real_3_t v2,
1084  const cs_real_3_t v3,
1085  const cs_real_3_t v4,
1086  double vol,
1087  cs_analytic_func_t *ana,
1088  void *input,
1089  double results[])
1090 {
1091  cs_real_3_t gauss_pts[5];
1092  double evaluation[5], weights[5];
1093 
1094  /* Compute Gauss points and its weights */
1095  cs_quadrature_tet_5pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1096 
1097  ana(tcur, 5, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1098 
1099  *results += weights[0] * evaluation[0] + weights[1] * evaluation[1] +
1100  weights[2] * evaluation[2] + weights[3] * evaluation[3] +
1101  weights[4] * evaluation[4];
1102 }
1103 
1104 /*----------------------------------------------------------------------------*/
1120 /*----------------------------------------------------------------------------*/
1121 
1122 static inline void
1124  const cs_real_3_t v1,
1125  const cs_real_3_t v2,
1126  const cs_real_3_t v3,
1127  const cs_real_3_t v4,
1128  double vol,
1129  cs_analytic_func_t *ana,
1130  void *input,
1131  double results[])
1132 {
1133  cs_real_3_t xg;
1134  double evaluation[3];
1135 
1136  /* Copied from cs_quadrature_tet_1pt */
1137  xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1138  xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1139  xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1140 
1141  ana(tcur, 1, NULL, xg, false, input, evaluation);
1142 
1143  results[0] += vol * evaluation[0];
1144  results[1] += vol * evaluation[1];
1145  results[2] += vol * evaluation[2];
1146 }
1147 
1148 /*----------------------------------------------------------------------------*/
1164 /*----------------------------------------------------------------------------*/
1165 
1166 static inline void
1168  const cs_real_3_t v1,
1169  const cs_real_3_t v2,
1170  const cs_real_3_t v3,
1171  const cs_real_3_t v4,
1172  double vol,
1173  cs_analytic_func_t *ana,
1174  void *input,
1175  double results[])
1176 {
1177  cs_real_3_t gauss_pts[4];
1178  double evaluation[3*4], weights[4];
1179 
1180  /* Compute Gauss points and its unique weight */
1181  cs_quadrature_tet_4pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1182 
1183  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1184 
1185  for (int p = 0; p < 4; p++) {
1186  results[0] += weights[p] * evaluation[3*p ];
1187  results[1] += weights[p] * evaluation[3*p+1];
1188  results[2] += weights[p] * evaluation[3*p+2];
1189  }
1190 }
1191 
1192 /*----------------------------------------------------------------------------*/
1208 /*----------------------------------------------------------------------------*/
1209 
1210 static inline void
1212  const cs_real_3_t v1,
1213  const cs_real_3_t v2,
1214  const cs_real_3_t v3,
1215  const cs_real_3_t v4,
1216  double vol,
1217  cs_analytic_func_t *ana,
1218  void *input,
1219  double results[])
1220 {
1221  cs_real_3_t gauss_pts[5];
1222  double evaluation[3*5], weights[5];
1223 
1224  /* Compute Gauss points and its weights */
1225  cs_quadrature_tet_5pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1226 
1227  ana(tcur, 5, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1228 
1229  for (int p = 0; p < 5; p++) {
1230  results[0] += weights[p] * evaluation[3*p ];
1231  results[1] += weights[p] * evaluation[3*p+1];
1232  results[2] += weights[p] * evaluation[3*p+2];
1233  }
1234 }
1235 
1236 /*----------------------------------------------------------------------------*/
1252 /*----------------------------------------------------------------------------*/
1253 
1254 static inline void
1256  const cs_real_3_t v1,
1257  const cs_real_3_t v2,
1258  const cs_real_3_t v3,
1259  const cs_real_3_t v4,
1260  double vol,
1261  cs_analytic_func_t *ana,
1262  void *input,
1263  double results[])
1264 {
1265  cs_real_3_t xg;
1266  double evaluation[9];
1267 
1268  /* Copied from cs_quadrature_tet_1pt */
1269  xg[0] = 0.25 * (v1[0] + v2[0] + v3[0] + v4[0]);
1270  xg[1] = 0.25 * (v1[1] + v2[1] + v3[1] + v4[1]);
1271  xg[2] = 0.25 * (v1[2] + v2[2] + v3[2] + v4[2]);
1272 
1273  ana(tcur, 1, NULL, xg, false, input, evaluation);
1274 
1275  for (short int ij = 0; ij < 9; ij++)
1276  results[ij] += vol * evaluation[ij];
1277 }
1278 
1279 /*----------------------------------------------------------------------------*/
1295 /*----------------------------------------------------------------------------*/
1296 
1297 static inline void
1299  const cs_real_3_t v1,
1300  const cs_real_3_t v2,
1301  const cs_real_3_t v3,
1302  const cs_real_3_t v4,
1303  double vol,
1304  cs_analytic_func_t *ana,
1305  void *input,
1306  double results[])
1307 {
1308  cs_real_3_t gauss_pts[4];
1309  double evaluation[9*4], weights[4];
1310 
1311  /* Compute Gauss points and its unique weight */
1312  cs_quadrature_tet_4pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1313 
1314  ana(tcur, 4, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1315 
1316  for (int p = 0; p < 4; p++) {
1317  const double wp = weights[p];
1318  double *eval_p = evaluation + 9*p;
1319  for (short int ij = 0; ij < 9; ij++)
1320  results[ij] += wp * eval_p[ij];
1321  }
1322 }
1323 
1324 /*----------------------------------------------------------------------------*/
1340 /*----------------------------------------------------------------------------*/
1341 
1342 static inline void
1344  const cs_real_3_t v1,
1345  const cs_real_3_t v2,
1346  const cs_real_3_t v3,
1347  const cs_real_3_t v4,
1348  double vol,
1349  cs_analytic_func_t *ana,
1350  void *input,
1351  double results[])
1352 {
1353  cs_real_3_t gauss_pts[5];
1354  double evaluation[9*5], weights[5];
1355 
1356  /* Compute Gauss points and its weights */
1357  cs_quadrature_tet_5pts(v1, v2, v3, v4, vol, gauss_pts, weights);
1358 
1359  ana(tcur, 5, NULL, (const cs_real_t *)gauss_pts, false, input, evaluation);
1360 
1361  for (int p = 0; p < 5; p++) {
1362  const double wp = weights[p];
1363  double *eval_p = evaluation + 9*p;
1364  for (short int ij = 0; ij < 9; ij++)
1365  results[ij] += wp * eval_p[ij];
1366  }
1367 }
1368 
1369 /*----------------------------------------------------------------------------*/
1379 /*----------------------------------------------------------------------------*/
1380 
1381 static inline cs_quadrature_tria_integral_t *
1383  cs_quadrature_type_t qtype)
1384 {
1385  switch (dim) {
1386 
1387  case 1: /* Scalar-valued integral */
1388 
1389  switch (qtype) {
1390 
1391  case CS_QUADRATURE_BARY:
1394  case CS_QUADRATURE_HIGHER:
1396  case CS_QUADRATURE_HIGHEST:
1398 
1399  default:
1400  bft_error(__FILE__, __LINE__, 0,
1401  " %s: Invalid quadrature type\n", __func__);
1402  }
1403  break;
1404 
1405  case 3: /* Vector-valued case */
1406 
1407  switch (qtype) {
1408 
1409  case CS_QUADRATURE_BARY:
1412  case CS_QUADRATURE_HIGHER:
1414  case CS_QUADRATURE_HIGHEST:
1416 
1417  default:
1418  bft_error(__FILE__, __LINE__, 0,
1419  " %s: Invalid quadrature type\n", __func__);
1420  }
1421  break;
1422 
1423  case 9: /* Tensor-valued case */
1424 
1425  switch (qtype) {
1426 
1427  case CS_QUADRATURE_BARY:
1430  case CS_QUADRATURE_HIGHER:
1432  case CS_QUADRATURE_HIGHEST:
1434 
1435  default:
1436  bft_error(__FILE__, __LINE__, 0,
1437  " %s: Invalid quadrature type\n", __func__);
1438  }
1439  break;
1440 
1441  default:
1442  bft_error(__FILE__, __LINE__, 0,
1443  " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
1444  __func__, dim);
1445 
1446  } /* switch on dim */
1447 
1448  return NULL; /* Should not go to this stage */
1449 }
1450 
1451 /*----------------------------------------------------------------------------*/
1461 /*----------------------------------------------------------------------------*/
1462 
1463 static inline cs_quadrature_tetra_integral_t*
1465  cs_quadrature_type_t qtype)
1466 {
1467  switch (dim) {
1468 
1469  case 1: /* Scalar-valued case */
1470 
1471  switch (qtype) {
1472 
1473  case CS_QUADRATURE_BARY:
1476  case CS_QUADRATURE_HIGHER:
1478  case CS_QUADRATURE_HIGHEST:
1480 
1481  default:
1482  bft_error(__FILE__, __LINE__, 0,
1483  " %s: Invalid quadrature type\n", __func__);
1484  }
1485  break;
1486 
1487  case 3: /* Vector-valued case */
1488 
1489  switch (qtype) {
1490 
1491  case CS_QUADRATURE_BARY:
1494  case CS_QUADRATURE_HIGHER:
1496  case CS_QUADRATURE_HIGHEST:
1498 
1499  default:
1500  bft_error(__FILE__, __LINE__, 0,
1501  " %s: Invalid quadrature type\n", __func__);
1502  }
1503  break;
1504 
1505  case 9: /* Tensor-valued case */
1506 
1507  switch (qtype) {
1508 
1509  case CS_QUADRATURE_BARY:
1512  case CS_QUADRATURE_HIGHER:
1514  case CS_QUADRATURE_HIGHEST:
1516 
1517  default:
1518  bft_error(__FILE__, __LINE__, 0,
1519  " %s: Invalid quadrature type\n", __func__);
1520  }
1521  break;
1522 
1523  default:
1524  bft_error(__FILE__, __LINE__, 0,
1525  " %s: Invalid dimension value %d. Only 1, 3 and 9 are valid.\n",
1526  __func__, dim);
1527 
1528  } /* Switch on dim */
1529 
1530  /* Avoid no return warning */
1531  return NULL;
1532 }
1533 
1534 /*----------------------------------------------------------------------------*/
1535 
1537 
1538 #endif /* __CS_QUADRATURE_H__ */
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
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:349
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:785
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:389
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, it enables to access in coords at the right location and the same thing to fill retval if compact is set to false.
Definition: cs_param.h:66
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
size_t len
Definition: mei_scanner.c:569
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:1081
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:953
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:1343
#define BEGIN_C_DECLS
Definition: cs_defs.h:461
Definition: cs_quadrature.h:55
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:1039
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
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_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:1298
Definition: cs_field_pointer.h:67
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:198
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
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...
Definition: cs_quadrature.h:58
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
double cs_real_t
Floating-point value.
Definition: cs_defs.h:297
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
Definition: cs_quadrature.h:57
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...
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:147
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
cs_quadrature_type_t
Definition: cs_quadrature.h:51
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:269
static cs_quadrature_tetra_integral_t * cs_quadrature_get_integral_function_tetra(short int dim, cs_quadrature_type_t qtype)
Retrieve the integral function according to the the quadrature type and the stride provided...
Definition: cs_quadrature.h:1464
void cs_quadrature_setup(void)
Compute constant weights for all quadratures used.
Definition: cs_quadrature.c:105
static int input(void)
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_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
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:1211
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_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:309
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_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:1123
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:910
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:306
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:997
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
const cs_real_t cs_math_onethird
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:1255
#define END_C_DECLS
Definition: cs_defs.h:462
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:433
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:869
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:827
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
Definition: cs_quadrature.h:56
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
Definition: cs_quadrature.h:53
static cs_quadrature_tria_integral_t * cs_quadrature_get_integral_function_tria(short int dim, cs_quadrature_type_t qtype)
Retrieve the integral function according to the the quadrature type and the stride provided...
Definition: cs_quadrature.h:1382
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:743
Definition: cs_quadrature.h:54
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:1167