KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
quadrilateral_gauss_legendre_integration_points.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Josep Maria Carbonell
11 //
12 
13 #if !defined(KRATOS_QUADRILATERAL_GAUSS_LEGENDRE_INTEGRATION_POINTS_H_INCLUDED )
14 #define KRATOS_QUADRILATERAL_GAUSS_LEGENDRE_INTEGRATION_POINTS_H_INCLUDED
15 
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "integration/quadrature.h"
23 
24 
25 namespace Kratos
26 {
27 
29 {
30 public:
32  typedef std::size_t SizeType;
33 
34  static const unsigned int Dimension = 2;
35 
37 
38  typedef std::array<IntegrationPointType, 1> IntegrationPointsArrayType;
39 
41 
43  {
44  return 1;
45  }
46 
48  {
49  static const IntegrationPointsArrayType s_integration_points{{
50  IntegrationPointType( 0.00 , 0.00 , 4.00 )
51  }};
52  return s_integration_points;
53  }
54 
55  std::string Info() const
56  {
57  std::stringstream buffer;
58  buffer << "Quadrilateral Gauss-Legendre quadrature 1 ";
59  return buffer.str();
60  }
61 
62 
63 }; // Class QuadrilateralGaussLegendreIntegrationPoints1
64 
66 {
67 public:
69  typedef std::size_t SizeType;
70 
71  static const unsigned int Dimension = 2;
72 
74 
75  typedef std::array<IntegrationPointType, 4> IntegrationPointsArrayType;
76 
78 
80  {
81  return 4;
82  }
83 
85  {
86  static const IntegrationPointsArrayType s_integration_points{{
87  IntegrationPointType( -1.00/std::sqrt(3.0) , -1.00/std::sqrt(3.0), 1.00 ),
88  IntegrationPointType( 1.00/std::sqrt(3.0) , -1.00/std::sqrt(3.0), 1.00 ),
89  IntegrationPointType( 1.00/std::sqrt(3.0) , 1.00/std::sqrt(3.0), 1.00 ),
90  IntegrationPointType( -1.00/std::sqrt(3.0) , 1.00/std::sqrt(3.0), 1.00 )
91  }};
92  return s_integration_points;
93  }
94 
95  std::string Info() const
96  {
97  std::stringstream buffer;
98  buffer << "Quadrilateral Gauss-Legendre quadrature 2 ";
99  return buffer.str();
100  }
101 
102 
103 }; // Class QuadrilateralGaussLegendreIntegrationPoints2
104 
106 {
107 public:
109  typedef std::size_t SizeType;
110 
111  static const unsigned int Dimension = 2;
112 
114 
115  typedef std::array<IntegrationPointType, 9> IntegrationPointsArrayType;
116 
118 
120  {
121  return 9;
122  }
123 
125  {
126  static const IntegrationPointsArrayType s_integration_points{{
127  IntegrationPointType( -std::sqrt(3.00/5.00) , -std::sqrt(3.00/5.00), 25.00/81.00 ),
128  IntegrationPointType( 0.00 , -std::sqrt(3.00/5.00), 40.00/81.00 ),
129  IntegrationPointType( std::sqrt(3.00/5.00) , -std::sqrt(3.00/5.00), 25.00/81.00 ),
130  IntegrationPointType( -std::sqrt(3.00/5.00), 0.00, 40.00/81.00 ),
131  IntegrationPointType( 0.00 , 0.00, 64.00/81.00 ),
132  IntegrationPointType( std::sqrt(3.00/5.00), 0.00, 40.00/81.00 ),
133  IntegrationPointType( -std::sqrt(3.00/5.00), std::sqrt(3.00/5.00), 25.00/81.00 ),
134  IntegrationPointType( 0.00, std::sqrt(3.00/5.00), 40.00/81.00 ),
135  IntegrationPointType( std::sqrt(3.00/5.00), std::sqrt(3.00/5.00), 25.00/81.00 )
136  }};
137  return s_integration_points;
138  }
139 
140  std::string Info() const
141  {
142  std::stringstream buffer;
143  buffer << "Quadrilateral Gauss-Legendre quadrature 3 ";
144  return buffer.str();
145  }
146 
147 
148 }; // Class QuadrilateralGaussLegendreIntegrationPoints3
149 
151 {
152 public:
154  typedef std::size_t SizeType;
155 
156  static const unsigned int Dimension = 2;
157 
159 
160  typedef std::array<IntegrationPointType, 16> IntegrationPointsArrayType;
161 
163 
164  static SizeType IntegrationPointsNumber() { return 16; }
165 
167  {
168  static const IntegrationPointsArrayType s_integration_points{{
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))
185  }};
186  return s_integration_points;
187  }
188 
189  std::string Info() const
190  {
191  std::stringstream buffer;
192  buffer << "Quadrilateral Gauss-Legendre quadrature 4 ";
193  return buffer.str();
194  }
195 
196 
197 }; // Class QuadrilateralGaussLegendreIntegrationPoints4
198 
200 public:
202  typedef std::size_t SizeType;
203 
204  static const unsigned int Dimension = 2;
205 
207 
208  typedef std::array<IntegrationPointType, 25> IntegrationPointsArrayType;
209 
211 
212  static SizeType IntegrationPointsNumber() {return 25;}
213 
215  {
216  static IntegrationPointsArrayType s_integration_points;
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};
219 
220  for(unsigned int i = 0; i < 5; ++i) {
221  for(unsigned int j = 0; j < 5; ++j) {
222  s_integration_points[5*i + j] = IntegrationPointType( a[i], a[j], w[i] * w[j]);
223  }
224  }
225 
226  return s_integration_points;
227  }
228 
229  std::string Info() const
230  {
231  std::stringstream buffer;
232  buffer << "Quadrilateral Gauss-Legendre quadrature 5 ";
233  return buffer.str();
234  }
235 
236 
237 }; // Class QuadrilateralGaussLegendreIntegrationPoints5
238 
239 
242 
243 
247 
248 
250 
251 
252 } // namespace Kratos.
253 
254 #endif // KRATOS_QUADRILATERAL_GAUSS_LEGENDRE_INTEGRATION_POINTS_H_INCLUDED defined
255 
256 
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