Branch data Line data Source code
1 : : // *****************************************************************************
2 : : /*!
3 : : \file src/PDE/Integrate/Quadrature.hpp
4 : : \copyright 2012-2015 J. Bakosi,
5 : : 2016-2018 Los Alamos National Security, LLC.,
6 : : 2019-2021 Triad National Security, LLC.
7 : : All rights reserved. See the LICENSE file for details.
8 : : \brief Quadrature coordinates and weights for numerical integration
9 : : \details This file contains functions that provide Gauss quadrature
10 : : coordinates and weights for numerical integration.
11 : : */
12 : : // *****************************************************************************
13 : : #ifndef Quadrature_h
14 : : #define Quadrature_h
15 : :
16 : : #include <array>
17 : : #include <vector>
18 : : #include <stdexcept>
19 : :
20 : : #include "Types.hpp"
21 : : #include "Exception.hpp"
22 : :
23 : : namespace tk {
24 : : //! Initialization of number of Gauss points for volume integration
25 : : //! \param[in] ndof Number of degrees of freedom
26 : 16410970 : constexpr std::size_t NGvol( const std::size_t ndof ) {
27 [ + + ][ + + ]: 16410970 : return ndof == 1 ? 1 :
[ - + ]
28 : : ndof == 4 ? 5 :
29 : : ndof == 10 ? 11 :
30 [ - - ]: 16410970 : throw std::logic_error("ndof must be one of 1,4,10");
31 : : }
32 : :
33 : :
34 : : //! Initialization of number of Gauss points for face integration
35 : : //! \param[in] ndof Number of degrees of freedom
36 : 52938101 : constexpr std::size_t NGfa( const std::size_t ndof ) {
37 [ + + ][ + + ]: 52938101 : return ndof == 1 ? 1 :
[ - + ]
38 : : ndof == 4 ? 3 :
39 : : ndof == 10 ? 6 :
40 [ - - ]: 52938101 : throw std::logic_error("ndof must be one of 1,4,10");
41 : : }
42 : :
43 : : //! \brief Initialization of number of Gauss points for volume integration
44 : : //! in error estimation
45 : : //! \param[in] ndof Number of degrees of freedom
46 : 2479338 : constexpr std::size_t NGdiag( const std::size_t ndof ) {
47 [ + + ][ + + ]: 2479338 : return ndof == 1 ? 1 :
[ - + ]
48 : : ndof == 4 ? 4 :
49 : : ndof == 10 ? 14 :
50 [ - - ]: 2479338 : throw std::logic_error("ndof must be one of 1,4,10");
51 : : }
52 : :
53 : : //! \brief Initialization of number of Gauss points for volume integration
54 : : //! in DG initialization
55 : : //! \param[in] ndof Number of degrees of freedom
56 : 1096 : constexpr std::size_t NGinit( const std::size_t ndof ) {
57 [ + + ][ + + ]: 1096 : return ndof == 1 ? 1 :
[ - + ]
58 : : ndof == 4 ? 14 :
59 : : ndof == 10 ? 14 :
60 [ - - ]: 1096 : throw std::logic_error("ndof must be one of 1,4,10");
61 : : }
62 : :
63 : : //! Diagonal mass matrix for Dubiner basis functions on a reference tetrahedron
64 : : //! \return Diagonal mass matrix
65 : : constexpr std::array< real, 10 > massMatrixDubiner() {
66 : : return {{
67 : : 1.0,
68 : : 1.0/10.0,
69 : : 3.0/10.0,
70 : : 3.0/5.0,
71 : : 1.0/35.0,
72 : : 1.0/21.0,
73 : : 1.0/14.0,
74 : : 1.0/7.0,
75 : : 3.0/14.0,
76 : : 3.0/7.0 }};
77 : : }
78 : :
79 : : //! Face centroid coords in standard-tet-reference frame
80 : : const std::array< std::array< tk::real, 3 > , 4 > fc_coord {{
81 : : { 1.0/3.0, 1.0/3.0, 1.0/3.0 },
82 : : { 0.0, 1.0/3.0, 1.0/3.0 },
83 : : { 1.0/3.0, 0.0, 1.0/3.0 },
84 : : { 1.0/3.0, 1.0/3.0, 0.0 } }};
85 : :
86 : : //! Initialize Gaussian quadrature points locations and weights for a tetrahedron
87 : : void
88 : : GaussQuadratureTet( std::size_t NG,
89 : : std::array< std::vector< real >, 3 >& coordgp,
90 : : std::vector< real >& wgp );
91 : :
92 : : //! Initialize Gaussian quadrature points locations and weights for a triangle
93 : : void
94 : : GaussQuadratureTri( std::size_t NG,
95 : : std::array< std::vector< real >, 2 >& coordgp,
96 : : std::vector< real >& wgp );
97 : : } // tk::
98 : :
99 : : #endif // Quadrature_h
|