13 #if !defined(KRATOS_QUADRILATERAL_GAUSS_LEGENDRE_INTEGRATION_POINTS_H_INCLUDED )
14 #define KRATOS_QUADRILATERAL_GAUSS_LEGENDRE_INTEGRATION_POINTS_H_INCLUDED
52 return s_integration_points;
57 std::stringstream buffer;
58 buffer <<
"Quadrilateral Gauss-Legendre quadrature 1 ";
92 return s_integration_points;
97 std::stringstream buffer;
98 buffer <<
"Quadrilateral Gauss-Legendre quadrature 2 ";
137 return s_integration_points;
142 std::stringstream buffer;
143 buffer <<
"Quadrilateral Gauss-Legendre quadrature 3 ";
169 IntegrationPointType( -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
170 IntegrationPointType( -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
171 IntegrationPointType( -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
172 IntegrationPointType( -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
173 IntegrationPointType( -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
174 IntegrationPointType( -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
175 IntegrationPointType( -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
176 IntegrationPointType( -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
177 IntegrationPointType( std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
178 IntegrationPointType( std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
179 IntegrationPointType( std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
180 IntegrationPointType( std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 + std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
181 IntegrationPointType( std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0)),
182 IntegrationPointType( std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), -std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
183 IntegrationPointType( std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 - 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 + std::sqrt(5.0/6.0)/6.0)),
184 IntegrationPointType( std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), std::sqrt( (3.0 + 2.0 * std::sqrt(6.0/5.0)) / 7.0 ), (0.5 - std::sqrt(5.0/6.0)/6.0)*(0.5 - std::sqrt(5.0/6.0)/6.0))
186 return s_integration_points;
191 std::stringstream buffer;
192 buffer <<
"Quadrilateral Gauss-Legendre quadrature 4 ";
217 const double a[] = {-0.906179845938664, -0.538469310105683, 0.000000000000000, 0.538469310105683, 0.906179845938664};
218 const double w[] = {0.236926885056189, 0.478628670499366, 0.568888888888889, 0.478628670499366, 0.236926885056189};
220 for(
unsigned int i = 0;
i < 5; ++
i) {
221 for(
unsigned int j = 0;
j < 5; ++
j) {
226 return s_integration_points;
231 std::stringstream buffer;
232 buffer <<
"Quadrilateral Gauss-Legendre quadrature 5 ";
Short class definition.
Definition: integration_point.h:52
Point class.
Definition: point.h:59
Definition: quadrilateral_gauss_legendre_integration_points.h:29
static const unsigned int Dimension
Definition: quadrilateral_gauss_legendre_integration_points.h:34
std::string Info() const
Definition: quadrilateral_gauss_legendre_integration_points.h:55
static const IntegrationPointsArrayType & IntegrationPoints()
Definition: quadrilateral_gauss_legendre_integration_points.h:47
std::size_t SizeType
Definition: quadrilateral_gauss_legendre_integration_points.h:32
KRATOS_CLASS_POINTER_DEFINITION(QuadrilateralGaussLegendreIntegrationPoints1)
static SizeType IntegrationPointsNumber()
Definition: quadrilateral_gauss_legendre_integration_points.h:42
IntegrationPointType::PointType PointType
Definition: quadrilateral_gauss_legendre_integration_points.h:40
IntegrationPoint< 2 > IntegrationPointType
Definition: quadrilateral_gauss_legendre_integration_points.h:36
std::array< IntegrationPointType, 1 > IntegrationPointsArrayType
Definition: quadrilateral_gauss_legendre_integration_points.h:38
Definition: quadrilateral_gauss_legendre_integration_points.h:66
static SizeType IntegrationPointsNumber()
Definition: quadrilateral_gauss_legendre_integration_points.h:79
std::array< IntegrationPointType, 4 > IntegrationPointsArrayType
Definition: quadrilateral_gauss_legendre_integration_points.h:75
static const IntegrationPointsArrayType & IntegrationPoints()
Definition: quadrilateral_gauss_legendre_integration_points.h:84
IntegrationPointType::PointType PointType
Definition: quadrilateral_gauss_legendre_integration_points.h:77
static const unsigned int Dimension
Definition: quadrilateral_gauss_legendre_integration_points.h:71
std::size_t SizeType
Definition: quadrilateral_gauss_legendre_integration_points.h:69
std::string Info() const
Definition: quadrilateral_gauss_legendre_integration_points.h:95
IntegrationPoint< 2 > IntegrationPointType
Definition: quadrilateral_gauss_legendre_integration_points.h:73
KRATOS_CLASS_POINTER_DEFINITION(QuadrilateralGaussLegendreIntegrationPoints2)
Definition: quadrilateral_gauss_legendre_integration_points.h:106
IntegrationPointType::PointType PointType
Definition: quadrilateral_gauss_legendre_integration_points.h:117
static const unsigned int Dimension
Definition: quadrilateral_gauss_legendre_integration_points.h:111
static SizeType IntegrationPointsNumber()
Definition: quadrilateral_gauss_legendre_integration_points.h:119
std::size_t SizeType
Definition: quadrilateral_gauss_legendre_integration_points.h:109
KRATOS_CLASS_POINTER_DEFINITION(QuadrilateralGaussLegendreIntegrationPoints3)
std::array< IntegrationPointType, 9 > IntegrationPointsArrayType
Definition: quadrilateral_gauss_legendre_integration_points.h:115
static const IntegrationPointsArrayType & IntegrationPoints()
Definition: quadrilateral_gauss_legendre_integration_points.h:124
IntegrationPoint< 2 > IntegrationPointType
Definition: quadrilateral_gauss_legendre_integration_points.h:113
std::string Info() const
Definition: quadrilateral_gauss_legendre_integration_points.h:140
Definition: quadrilateral_gauss_legendre_integration_points.h:151
static const IntegrationPointsArrayType & IntegrationPoints()
Definition: quadrilateral_gauss_legendre_integration_points.h:166
static SizeType IntegrationPointsNumber()
Definition: quadrilateral_gauss_legendre_integration_points.h:164
KRATOS_CLASS_POINTER_DEFINITION(QuadrilateralGaussLegendreIntegrationPoints4)
std::array< IntegrationPointType, 16 > IntegrationPointsArrayType
Definition: quadrilateral_gauss_legendre_integration_points.h:160
std::string Info() const
Definition: quadrilateral_gauss_legendre_integration_points.h:189
IntegrationPointType::PointType PointType
Definition: quadrilateral_gauss_legendre_integration_points.h:162
static const unsigned int Dimension
Definition: quadrilateral_gauss_legendre_integration_points.h:156
IntegrationPoint< 2 > IntegrationPointType
Definition: quadrilateral_gauss_legendre_integration_points.h:158
std::size_t SizeType
Definition: quadrilateral_gauss_legendre_integration_points.h:154
Definition: quadrilateral_gauss_legendre_integration_points.h:199
IntegrationPointType::PointType PointType
Definition: quadrilateral_gauss_legendre_integration_points.h:210
std::string Info() const
Definition: quadrilateral_gauss_legendre_integration_points.h:229
std::array< IntegrationPointType, 25 > IntegrationPointsArrayType
Definition: quadrilateral_gauss_legendre_integration_points.h:208
IntegrationPoint< 2 > IntegrationPointType
Definition: quadrilateral_gauss_legendre_integration_points.h:206
static SizeType IntegrationPointsNumber()
Definition: quadrilateral_gauss_legendre_integration_points.h:212
std::size_t SizeType
Definition: quadrilateral_gauss_legendre_integration_points.h:202
static const unsigned int Dimension
Definition: quadrilateral_gauss_legendre_integration_points.h:204
KRATOS_CLASS_POINTER_DEFINITION(QuadrilateralGaussLegendreIntegrationPoints5)
static IntegrationPointsArrayType IntegrationPoints()
Definition: quadrilateral_gauss_legendre_integration_points.h:214
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
w
Definition: generate_convection_diffusion_explicit_element.py:108
a
Definition: generate_stokes_twofluid_element.py:77
int j
Definition: quadrature.py:648
integer i
Definition: TensorModule.f:17