MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
intrules.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_INTRULES
13 #define MFEM_INTRULES
14 
15 #include "../config/config.hpp"
16 #include "../general/array.hpp"
17 
18 namespace mfem
19 {
20 
21 /* Classes for IntegrationPoint, IntegrationRule, and container class
22  IntegrationRules. Declares the global variable IntRules */
23 
24 /// Class for integration point with weight
26 {
27 public:
28  double x, y, z, weight;
29  int index;
30 
31  void Init(int const i)
32  {
33  x = y = z = weight = 0.0;
34  index = i;
35  }
36 
37  void Set(const double *p, const int dim)
38  {
39  MFEM_ASSERT(1 <= dim && dim <= 3, "invalid dim: " << dim);
40  x = p[0];
41  if (dim > 1)
42  {
43  y = p[1];
44  if (dim > 2)
45  {
46  z = p[2];
47  }
48  }
49  }
50 
51  void Get(double *p, const int dim) const
52  {
53  MFEM_ASSERT(1 <= dim && dim <= 3, "invalid dim: " << dim);
54  p[0] = x;
55  if (dim > 1)
56  {
57  p[1] = y;
58  if (dim > 2)
59  {
60  p[2] = z;
61  }
62  }
63  }
64 
65  void Set(const double x1, const double x2, const double x3, const double w)
66  { x = x1; y = x2; z = x3; weight = w; }
67 
68  void Set3w(const double *p) { x = p[0]; y = p[1]; z = p[2]; weight = p[3]; }
69 
70  void Set3(const double x1, const double x2, const double x3)
71  { x = x1; y = x2; z = x3; }
72 
73  void Set3(const double *p) { x = p[0]; y = p[1]; z = p[2]; }
74 
75  void Set2w(const double x1, const double x2, const double w)
76  { x = x1; y = x2; weight = w; }
77 
78  void Set2w(const double *p) { x = p[0]; y = p[1]; weight = p[2]; }
79 
80  void Set2(const double x1, const double x2) { x = x1; y = x2; }
81 
82  void Set2(const double *p) { x = p[0]; y = p[1]; }
83 
84  void Set1w(const double x1, const double w) { x = x1; weight = w; }
85 
86  void Set1w(const double *p) { x = p[0]; weight = p[1]; }
87 };
88 
89 /// Class for an integration rule - an Array of IntegrationPoint.
90 class IntegrationRule : public Array<IntegrationPoint>
91 {
92 private:
93  friend class IntegrationRules;
94  int Order;
95  /** @brief The quadrature weights gathered as a contiguous array. Created
96  by request with the method GetWeights(). */
97  mutable Array<double> weights;
98 
99  /// Sets the indices of each quadrature point on initialization.
100  void SetPointIndices();
101 
102  /// Define n-simplex rule (triangle/tetrahedron for n=2/3) of order (2s+1)
103  void GrundmannMollerSimplexRule(int s, int n = 3);
104 
105  void AddTriMidPoint(const int off, const double weight)
106  { IntPoint(off).Set2w(1./3., 1./3., weight); }
107 
108  void AddTriPoints3(const int off, const double a, const double b,
109  const double weight)
110  {
111  IntPoint(off + 0).Set2w(a, a, weight);
112  IntPoint(off + 1).Set2w(a, b, weight);
113  IntPoint(off + 2).Set2w(b, a, weight);
114  }
115 
116  void AddTriPoints3(const int off, const double a, const double weight)
117  { AddTriPoints3(off, a, 1. - 2.*a, weight); }
118 
119  void AddTriPoints3b(const int off, const double b, const double weight)
120  { AddTriPoints3(off, (1. - b)/2., b, weight); }
121 
122  void AddTriPoints3R(const int off, const double a, const double b,
123  const double c, const double weight)
124  {
125  IntPoint(off + 0).Set2w(a, b, weight);
126  IntPoint(off + 1).Set2w(c, a, weight);
127  IntPoint(off + 2).Set2w(b, c, weight);
128  }
129 
130  void AddTriPoints3R(const int off, const double a, const double b,
131  const double weight)
132  { AddTriPoints3R(off, a, b, 1. - a - b, weight); }
133 
134  void AddTriPoints6(const int off, const double a, const double b,
135  const double c, const double weight)
136  {
137  IntPoint(off + 0).Set2w(a, b, weight);
138  IntPoint(off + 1).Set2w(b, a, weight);
139  IntPoint(off + 2).Set2w(a, c, weight);
140  IntPoint(off + 3).Set2w(c, a, weight);
141  IntPoint(off + 4).Set2w(b, c, weight);
142  IntPoint(off + 5).Set2w(c, b, weight);
143  }
144 
145  void AddTriPoints6(const int off, const double a, const double b,
146  const double weight)
147  { AddTriPoints6(off, a, b, 1. - a - b, weight); }
148 
149  // add the permutations of (a,a,b)
150  void AddTetPoints3(const int off, const double a, const double b,
151  const double weight)
152  {
153  IntPoint(off + 0).Set(a, a, b, weight);
154  IntPoint(off + 1).Set(a, b, a, weight);
155  IntPoint(off + 2).Set(b, a, a, weight);
156  }
157 
158  // add the permutations of (a,b,c)
159  void AddTetPoints6(const int off, const double a, const double b,
160  const double c, const double weight)
161  {
162  IntPoint(off + 0).Set(a, b, c, weight);
163  IntPoint(off + 1).Set(a, c, b, weight);
164  IntPoint(off + 2).Set(b, c, a, weight);
165  IntPoint(off + 3).Set(b, a, c, weight);
166  IntPoint(off + 4).Set(c, a, b, weight);
167  IntPoint(off + 5).Set(c, b, a, weight);
168  }
169 
170  void AddTetMidPoint(const int off, const double weight)
171  { IntPoint(off).Set(0.25, 0.25, 0.25, weight); }
172 
173  // given a, add the permutations of (a,a,a,b), where 3*a + b = 1
174  void AddTetPoints4(const int off, const double a, const double weight)
175  {
176  IntPoint(off).Set(a, a, a, weight);
177  AddTetPoints3(off + 1, a, 1. - 3.*a, weight);
178  }
179 
180  // given b, add the permutations of (a,a,a,b), where 3*a + b = 1
181  void AddTetPoints4b(const int off, const double b, const double weight)
182  {
183  const double a = (1. - b)/3.;
184  IntPoint(off).Set(a, a, a, weight);
185  AddTetPoints3(off + 1, a, b, weight);
186  }
187 
188  // add the permutations of (a,a,b,b), 2*(a + b) = 1
189  void AddTetPoints6(const int off, const double a, const double weight)
190  {
191  const double b = 0.5 - a;
192  AddTetPoints3(off, a, b, weight);
193  AddTetPoints3(off + 3, b, a, weight);
194  }
195 
196  // given (a,b) or (a,c), add the permutations of (a,a,b,c), 2*a + b + c = 1
197  void AddTetPoints12(const int off, const double a, const double bc,
198  const double weight)
199  {
200  const double cb = 1. - 2*a - bc;
201  AddTetPoints3(off, a, bc, weight);
202  AddTetPoints3(off + 3, a, cb, weight);
203  AddTetPoints6(off + 6, a, bc, cb, weight);
204  }
205 
206  // given (b,c), add the permutations of (a,a,b,c), 2*a + b + c = 1
207  void AddTetPoints12bc(const int off, const double b, const double c,
208  const double weight)
209  {
210  const double a = (1. - b - c)/2.;
211  AddTetPoints3(off, a, b, weight);
212  AddTetPoints3(off + 3, a, c, weight);
213  AddTetPoints6(off + 6, a, b, c, weight);
214  }
215 
216 public:
218  Array<IntegrationPoint>(), Order(0) { }
219 
220  /// Construct an integration rule with given number of points
221  explicit IntegrationRule(int NP) :
222  Array<IntegrationPoint>(NP), Order(0)
223  {
224  for (int i = 0; i < this->Size(); i++)
225  {
226  (*this)[i].Init(i);
227  }
228  }
229 
230  /// Tensor product of two 1D integration rules
232 
233  /// Tensor product of three 1D integration rules
235  IntegrationRule &irz);
236 
237  /// Returns the order of the integration rule
238  int GetOrder() const { return Order; }
239 
240  /** @brief Sets the order of the integration rule. This is only for keeping
241  order information, it does not alter any data in the IntegrationRule. */
242  void SetOrder(const int order) { Order = order; }
243 
244  /// Returns the number of the points in the integration rule
245  int GetNPoints() const { return Size(); }
246 
247  /// Returns a reference to the i-th integration point
248  IntegrationPoint &IntPoint(int i) { return (*this)[i]; }
249 
250  /// Returns a const reference to the i-th integration point
251  const IntegrationPoint &IntPoint(int i) const { return (*this)[i]; }
252 
253  /// Return the quadrature weights in a contiguous array.
254  /** If a contiguous array is not required, the weights can be accessed with
255  a call like this: `IntPoint(i).weight`. */
256  const Array<double> &GetWeights() const;
257 
258  /// Destroys an IntegrationRule object
260 };
261 
262 /// A Class that defines 1-D numerical quadrature rules on [0,1].
264 {
265 public:
266  /** @name Methods for calculating quadrature rules.
267  These methods calculate the actual points and weights for the different
268  types of quadrature rules. */
269  ///@{
270  void GaussLegendre(const int np, IntegrationRule* ir);
271  void GaussLobatto(const int np, IntegrationRule *ir);
272  void OpenUniform(const int np, IntegrationRule *ir);
273  void ClosedUniform(const int np, IntegrationRule *ir);
274  void OpenHalfUniform(const int np, IntegrationRule *ir);
275  ///@}
276 
277  /// A helper function that will play nice with Poly_1D::OpenPoints and
278  /// Poly_1D::ClosedPoints
279  void GivePolyPoints(const int np, double *pts, const int type);
280 
281 private:
282  void CalculateUniformWeights(IntegrationRule *ir, const int type);
283 };
284 
285 /// A class container for 1D quadrature type constants.
287 {
288 public:
289  enum
290  {
291  Invalid = -1,
294  OpenUniform = 2, ///< aka open Newton-Cotes
295  ClosedUniform = 3, ///< aka closed Newton-Cotes
296  OpenHalfUniform = 4 ///< aka "open half" Newton-Cotes
297  };
298  /** @brief If the Quadrature1D type is not closed return Invalid; otherwise
299  return type. */
300  static int CheckClosed(int type);
301  /** @brief If the Quadrature1D type is not open return Invalid; otherwise
302  return type. */
303  static int CheckOpen(int type);
304 };
305 
306 /// Container class for integration rules
308 {
309 private:
310  /// Taken from the Quadrature1D class anonymous enum
311  /// Determines the type of numerical quadrature used for
312  /// segment, square, and cube geometries
313  const int quad_type;
314 
315  int own_rules, refined;
316 
317  /// Function that generates quadrature points and weights on [0,1]
318  QuadratureFunctions1D quad_func;
319 
320  Array<IntegrationRule *> PointIntRules;
321  Array<IntegrationRule *> SegmentIntRules;
322  Array<IntegrationRule *> TriangleIntRules;
323  Array<IntegrationRule *> SquareIntRules;
324  Array<IntegrationRule *> TetrahedronIntRules;
325  Array<IntegrationRule *> PrismIntRules;
326  Array<IntegrationRule *> CubeIntRules;
327 
328  void AllocIntRule(Array<IntegrationRule *> &ir_array, int Order)
329  {
330  if (ir_array.Size() <= Order)
331  {
332  ir_array.SetSize(Order + 1, NULL);
333  }
334  }
335  bool HaveIntRule(Array<IntegrationRule *> &ir_array, int Order)
336  {
337  return (ir_array.Size() > Order && ir_array[Order] != NULL);
338  }
339  int GetSegmentRealOrder(int Order) const
340  {
341  return Order | 1; // valid for all quad_type's
342  }
343 
344  /// The following methods allocate new IntegrationRule objects without
345  /// checking if they already exist. To avoid memory leaks use
346  /// IntegrationRules::Get(int GeomType, int Order) instead.
347  IntegrationRule *GenerateIntegrationRule(int GeomType, int Order);
348  IntegrationRule *PointIntegrationRule(int Order);
349  IntegrationRule *SegmentIntegrationRule(int Order);
350  IntegrationRule *TriangleIntegrationRule(int Order);
351  IntegrationRule *SquareIntegrationRule(int Order);
352  IntegrationRule *TetrahedronIntegrationRule(int Order);
353  IntegrationRule *PrismIntegrationRule(int Order);
354  IntegrationRule *CubeIntegrationRule(int Order);
355 
356  void DeleteIntRuleArray(Array<IntegrationRule *> &ir_array);
357 
358 public:
359  /// Sets initial sizes for the integration rule arrays, but rules
360  /// are defined the first time they are requested with the Get method.
361  explicit IntegrationRules(int Ref = 0,
362  int type = Quadrature1D::GaussLegendre);
363 
364  /// Returns an integration rule for given GeomType and Order.
365  const IntegrationRule &Get(int GeomType, int Order);
366 
367  void Set(int GeomType, int Order, IntegrationRule &IntRule);
368 
369  void SetOwnRules(int o) { own_rules = o; }
370 
371  /// Destroys an IntegrationRules object
373 };
374 
375 /// A global object with all integration rules (defined in intrules.cpp)
377 
378 /// A global object with all refined integration rules
380 
381 }
382 
383 #endif
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
void Set(const double *p, const int dim)
Definition: intrules.hpp:37
int Size() const
Logical size of the array.
Definition: array.hpp:124
void Get(double *p, const int dim) const
Definition: intrules.hpp:51
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:890
void ClosedUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:590
void Set1w(const double *p)
Definition: intrules.hpp:86
void GivePolyPoints(const int np, double *pts, const int type)
Definition: intrules.cpp:621
~IntegrationRules()
Destroys an IntegrationRules object.
Definition: intrules.cpp:981
Container class for integration rules.
Definition: intrules.hpp:307
void Set(const double x1, const double x2, const double x3, const double w)
Definition: intrules.hpp:65
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition: intrules.cpp:85
void Init(int const i)
Definition: intrules.hpp:31
void OpenUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:576
const IntegrationPoint & IntPoint(int i) const
Returns a const reference to the i-th integration point.
Definition: intrules.hpp:251
A class container for 1D quadrature type constants.
Definition: intrules.hpp:286
void GaussLobatto(const int np, IntegrationRule *ir)
Definition: intrules.cpp:454
void Set3(const double *p)
Definition: intrules.hpp:73
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
void SetOwnRules(int o)
Definition: intrules.hpp:369
double b
Definition: lissajous.cpp:42
~IntegrationRule()
Destroys an IntegrationRule object.
Definition: intrules.hpp:259
void Set3w(const double *p)
Definition: intrules.hpp:68
void Set2(const double x1, const double x2)
Definition: intrules.hpp:80
aka &quot;open half&quot; Newton-Cotes
Definition: intrules.hpp:296
void Set3(const double x1, const double x2, const double x3)
Definition: intrules.hpp:70
void Set2w(const double *p)
Definition: intrules.hpp:78
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
Definition: intrules.cpp:824
void OpenHalfUniform(const int np, IntegrationRule *ir)
Definition: intrules.cpp:608
int GetOrder() const
Returns the order of the integration rule.
Definition: intrules.hpp:238
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Definition: intrules.cpp:836
void SetOrder(const int order)
Sets the order of the integration rule. This is only for keeping order information, it does not alter any data in the IntegrationRule.
Definition: intrules.hpp:242
double a
Definition: lissajous.cpp:41
IntegrationRules(int Ref=0, int type=Quadrature1D::GaussLegendre)
Definition: intrules.cpp:856
Class for integration point with weight.
Definition: intrules.hpp:25
int dim
Definition: ex24.cpp:43
aka open Newton-Cotes
Definition: intrules.hpp:294
void Set2(const double *p)
Definition: intrules.hpp:82
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:379
void Set2w(const double x1, const double x2, const double w)
Definition: intrules.hpp:75
void pts(int iphi, int t, double x[])
IntegrationRule(int NP)
Construct an integration rule with given number of points.
Definition: intrules.hpp:221
void Set(int GeomType, int Order, IntegrationRule &IntRule)
Definition: intrules.cpp:936
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:376
A Class that defines 1-D numerical quadrature rules on [0,1].
Definition: intrules.hpp:263
void Set1w(const double x1, const double w)
Definition: intrules.hpp:84
void GaussLegendre(const int np, IntegrationRule *ir)
Definition: intrules.cpp:375
aka closed Newton-Cotes
Definition: intrules.hpp:295