MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
intrules.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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#if defined(MFEM_THREAD_SAFE) && defined(MFEM_USE_OPENMP)
18#include <omp.h>
19#endif
20
21#include <vector>
22#include <map>
23
24namespace mfem
25{
26
27class KnotVector;
28class Mesh;
29
30/* Classes for IntegrationPoint, IntegrationRule, and container class
31 IntegrationRules. Declares the global variable IntRules */
32
33/// Class for integration point with weight
35{
36public:
38 int index;
39
40 void Init(int const i)
41 {
42 x = y = z = weight = 0.0;
43 index = i;
44 }
45
46 void Set(const real_t *p, const int dim)
47 {
48 MFEM_ASSERT(1 <= dim && dim <= 3, "invalid dim: " << dim);
49 x = p[0];
50 if (dim > 1)
51 {
52 y = p[1];
53 if (dim > 2)
54 {
55 z = p[2];
56 }
57 }
58 }
59
60 void Get(real_t *p, const int dim) const
61 {
62 MFEM_ASSERT(1 <= dim && dim <= 3, "invalid dim: " << dim);
63 p[0] = x;
64 if (dim > 1)
65 {
66 p[1] = y;
67 if (dim > 2)
68 {
69 p[2] = z;
70 }
71 }
72 }
73
74 void Set(const real_t x1, const real_t x2, const real_t x3, const real_t w)
75 { x = x1; y = x2; z = x3; weight = w; }
76
77 void Set3w(const real_t *p) { x = p[0]; y = p[1]; z = p[2]; weight = p[3]; }
78
79 void Set3(const real_t x1, const real_t x2, const real_t x3)
80 { x = x1; y = x2; z = x3; }
81
82 void Set3(const real_t *p) { x = p[0]; y = p[1]; z = p[2]; }
83
84 void Set2w(const real_t x1, const real_t x2, const real_t w)
85 { x = x1; y = x2; weight = w; }
86
87 void Set2w(const real_t *p) { x = p[0]; y = p[1]; weight = p[2]; }
88
89 void Set2(const real_t x1, const real_t x2) { x = x1; y = x2; }
90
91 void Set2(const real_t *p) { x = p[0]; y = p[1]; }
92
93 void Set1w(const real_t x1, const real_t w) { x = x1; weight = w; }
94
95 void Set1w(const real_t *p) { x = p[0]; weight = p[1]; }
96};
97
98/// Class for an integration rule - an Array of IntegrationPoint.
99class IntegrationRule : public Array<IntegrationPoint>
100{
101private:
102 friend class IntegrationRules;
103 int Order = 0;
104 /** @brief The quadrature weights gathered as a contiguous array. Created
105 by request with the method GetWeights(). */
106 mutable Array<real_t> weights;
107
108 /// Define n-simplex rule (triangle/tetrahedron for n=2/3) of order (2s+1)
109 void GrundmannMollerSimplexRule(int s, int n = 3);
110
111 void AddTriMidPoint(const int off, const real_t weight)
112 { IntPoint(off).Set2w(1./3., 1./3., weight); }
113
114 void AddTriPoints3(const int off, const real_t a, const real_t b,
115 const real_t weight)
116 {
117 IntPoint(off + 0).Set2w(a, a, weight);
118 IntPoint(off + 1).Set2w(a, b, weight);
119 IntPoint(off + 2).Set2w(b, a, weight);
120 }
121
122 void AddTriPoints3(const int off, const real_t a, const real_t weight)
123 { AddTriPoints3(off, a, 1. - 2.*a, weight); }
124
125 void AddTriPoints3b(const int off, const real_t b, const real_t weight)
126 { AddTriPoints3(off, (1. - b)/2., b, weight); }
127
128 void AddTriPoints3R(const int off, const real_t a, const real_t b,
129 const real_t c, const real_t weight)
130 {
131 IntPoint(off + 0).Set2w(a, b, weight);
132 IntPoint(off + 1).Set2w(c, a, weight);
133 IntPoint(off + 2).Set2w(b, c, weight);
134 }
135
136 void AddTriPoints3R(const int off, const real_t a, const real_t b,
137 const real_t weight)
138 { AddTriPoints3R(off, a, b, 1. - a - b, weight); }
139
140 void AddTriPoints6(const int off, const real_t a, const real_t b,
141 const real_t c, const real_t weight)
142 {
143 IntPoint(off + 0).Set2w(a, b, weight);
144 IntPoint(off + 1).Set2w(b, a, weight);
145 IntPoint(off + 2).Set2w(a, c, weight);
146 IntPoint(off + 3).Set2w(c, a, weight);
147 IntPoint(off + 4).Set2w(b, c, weight);
148 IntPoint(off + 5).Set2w(c, b, weight);
149 }
150
151 void AddTriPoints6(const int off, const real_t a, const real_t b,
152 const real_t weight)
153 { AddTriPoints6(off, a, b, 1. - a - b, weight); }
154
155 // add the permutations of (a,a,b)
156 void AddTetPoints3(const int off, const real_t a, const real_t b,
157 const real_t weight)
158 {
159 IntPoint(off + 0).Set(a, a, b, weight);
160 IntPoint(off + 1).Set(a, b, a, weight);
161 IntPoint(off + 2).Set(b, a, a, weight);
162 }
163
164 // add the permutations of (a,b,c)
165 void AddTetPoints6(const int off, const real_t a, const real_t b,
166 const real_t c, const real_t weight)
167 {
168 IntPoint(off + 0).Set(a, b, c, weight);
169 IntPoint(off + 1).Set(a, c, b, weight);
170 IntPoint(off + 2).Set(b, c, a, weight);
171 IntPoint(off + 3).Set(b, a, c, weight);
172 IntPoint(off + 4).Set(c, a, b, weight);
173 IntPoint(off + 5).Set(c, b, a, weight);
174 }
175
176 void AddTetMidPoint(const int off, const real_t weight)
177 { IntPoint(off).Set(0.25, 0.25, 0.25, weight); }
178
179 // given a, add the permutations of (a,a,a,b), where 3*a + b = 1
180 void AddTetPoints4(const int off, const real_t a, const real_t weight)
181 {
182 IntPoint(off).Set(a, a, a, weight);
183 AddTetPoints3(off + 1, a, 1. - 3.*a, weight);
184 }
185
186 // given b, add the permutations of (a,a,a,b), where 3*a + b = 1
187 void AddTetPoints4b(const int off, const real_t b, const real_t weight)
188 {
189 const real_t a = (1. - b)/3.;
190 IntPoint(off).Set(a, a, a, weight);
191 AddTetPoints3(off + 1, a, b, weight);
192 }
193
194 // add the permutations of (a,a,b,b), 2*(a + b) = 1
195 void AddTetPoints6(const int off, const real_t a, const real_t weight)
196 {
197 const real_t b = 0.5 - a;
198 AddTetPoints3(off, a, b, weight);
199 AddTetPoints3(off + 3, b, a, weight);
200 }
201
202 // given (a,b) or (a,c), add the permutations of (a,a,b,c), 2*a + b + c = 1
203 void AddTetPoints12(const int off, const real_t a, const real_t bc,
204 const real_t weight)
205 {
206 const real_t cb = 1. - 2*a - bc;
207 AddTetPoints3(off, a, bc, weight);
208 AddTetPoints3(off + 3, a, cb, weight);
209 AddTetPoints6(off + 6, a, bc, cb, weight);
210 }
211
212 // given (b,c), add the permutations of (a,a,b,c), 2*a + b + c = 1
213 void AddTetPoints12bc(const int off, const real_t b, const real_t c,
214 const real_t weight)
215 {
216 const real_t a = (1. - b - c)/2.;
217 AddTetPoints3(off, a, b, weight);
218 AddTetPoints3(off + 3, a, c, weight);
219 AddTetPoints6(off + 6, a, b, c, weight);
220 }
221
222public:
225
226 /// Construct an integration rule with given number of points
227 explicit IntegrationRule(int NP) :
229 {
230 for (int i = 0; i < this->Size(); i++)
231 {
232 (*this)[i].Init(i);
233 }
234 }
235
236 /// Sets the indices of each quadrature point on initialization.
237 /** Note that most calls to IntegrationRule::SetSize should be paired with a
238 call to SetPointIndices in order for the indices to be set correctly. */
239 void SetPointIndices();
240
241 /// Tensor product of two 1D integration rules
243
244 /// Tensor product of three 1D integration rules
246 IntegrationRule &irz);
247
248 /// Returns the order of the integration rule
249 int GetOrder() const { return Order; }
250
251 /** @brief Sets the order of the integration rule. This is only for keeping
252 order information, it does not alter any data in the IntegrationRule. */
253 void SetOrder(const int order) { Order = order; }
254
255 /// Returns the number of the points in the integration rule
256 int GetNPoints() const { return Size(); }
257
258 /// Returns a reference to the i-th integration point
259 IntegrationPoint &IntPoint(int i) { return (*this)[i]; }
260
261 /// Returns a const reference to the i-th integration point
262 const IntegrationPoint &IntPoint(int i) const { return (*this)[i]; }
263
264 /// Return the quadrature weights in a contiguous array.
265 /** If a contiguous array is not required, the weights can be accessed with
266 a call like this: `IntPoint(i).weight`. */
267 const Array<real_t> &GetWeights() const;
268
269 /// @brief Return an integration rule for KnotVector @a kv, defined by
270 /// applying this rule on each knot interval.
272
273 /// Destroys an IntegrationRule object
275};
276
277/// Class for defining different integration rules on each NURBS patch.
279{
280public:
281 /// Construct a rule for each patch, using SetPatchRules1D.
282 NURBSMeshRules(const int numPatches, const int dim_) :
283 patchRules1D(numPatches, dim_),
284 npatches(numPatches), dim(dim_) { }
285
286 /// Returns a rule for the element.
287 IntegrationRule &GetElementRule(const int elem, const int patch,
288 const int *ijk,
289 Array<const KnotVector*> const& kv,
290 bool & deleteRule) const;
291
292 /// Add a rule to be used for individual elements. Returns the rule index.
293 std::size_t AddElementRule(IntegrationRule *ir_element)
294 {
295 elementRule.push_back(ir_element);
296 return elementRule.size() - 1;
297 }
298
299 /// @brief Set the integration rule for the element of the given index. This
300 /// rule is used instead of the rule for the patch containing the element.
301 void SetElementRule(const std::size_t element,
302 const std::size_t elementRuleIndex)
303 {
304 elementToRule[element] = elementRuleIndex;
305 }
306
307 /// @brief Set 1D integration rules to be used as a tensor product rule on
308 /// the patch with index @a patch. This class takes ownership of these rules.
309 void SetPatchRules1D(const int patch,
310 std::vector<const IntegrationRule*> & ir1D);
311
312 /// @brief For tensor product rules defined on each patch by
313 /// SetPatchRules1D(), return a pointer to the 1D rule in the specified
314 /// @a dimension.
315 const IntegrationRule* GetPatchRule1D(const int patch,
316 const int dimension) const
317 {
318 return patchRules1D(patch, dimension);
319 }
320
321 /// @brief For tensor product rules defined on each patch by
322 /// SetPatchRules1D(), return the integration point with index (i,j,k).
323 void GetIntegrationPointFrom1D(const int patch, int i, int j, int k,
324 IntegrationPoint & ip);
325
326 /// @brief Finalize() must be called before this class can be used for
327 /// assembly. In particular, it defines data used by GetPointElement().
328 void Finalize(Mesh const& mesh);
329
330 /// @brief For tensor product rules defined on each patch by
331 /// SetPatchRules1D(), returns the index of the element containing
332 /// integration point (i,j,k) for patch index @a patch. Finalize() must be
333 /// called first.
334 int GetPointElement(int patch, int i, int j, int k) const
335 {
336 return pointToElem[patch](i,j,k);
337 }
338
339 int GetDim() const { return dim; }
340
341 /// @brief For tensor product rules defined on each patch by
342 /// SetPatchRules1D(), returns an array of knot span indices for each
343 /// integration point in the specified @a dimension.
344 const Array<int>& GetPatchRule1D_KnotSpan(const int patch,
345 const int dimension) const
346 {
347 return patchRules1D_KnotSpan[patch][dimension];
348 }
349
351
352private:
353 /// Tensor-product rules defined on all patches independently.
355
356 /// Integration rules defined on elements.
357 std::vector<IntegrationRule*> elementRule;
358
359 std::map<std::size_t, std::size_t> elementToRule;
360
361 std::vector<Array3D<int>> pointToElem;
362 std::vector<std::vector<Array<int>>> patchRules1D_KnotSpan;
363
364 const int npatches;
365 const int dim;
366};
367
368/// A Class that defines 1-D numerical quadrature rules on [0,1].
370{
371public:
372 /** @name Methods for calculating quadrature rules.
373 These methods calculate the actual points and weights for the different
374 types of quadrature rules. */
375 ///@{
376 static void GaussLegendre(const int np, IntegrationRule* ir);
377 static void GaussLobatto(const int np, IntegrationRule *ir);
378 static void OpenUniform(const int np, IntegrationRule *ir);
379 static void ClosedUniform(const int np, IntegrationRule *ir);
380 static void OpenHalfUniform(const int np, IntegrationRule *ir);
381 static void ClosedGL(const int np, IntegrationRule *ir);
382 ///@}
383
384 /// A helper function that will play nice with Poly_1D::OpenPoints and
385 /// Poly_1D::ClosedPoints
386 static void GivePolyPoints(const int np, real_t *pts, const int type);
387
388private:
389 static void CalculateUniformWeights(IntegrationRule *ir, const int type);
390};
391
392/// A class container for 1D quadrature type constants.
394{
395public:
396 enum
397 {
401 OpenUniform = 2, ///< aka open Newton-Cotes
402 ClosedUniform = 3, ///< aka closed Newton-Cotes
403 OpenHalfUniform = 4, ///< aka "open half" Newton-Cotes
404 ClosedGL = 5 ///< aka closed Gauss Legendre
405 };
406 /** @brief If the Quadrature1D type is not closed return Invalid; otherwise
407 return type. */
408 static int CheckClosed(int type);
409 /** @brief If the Quadrature1D type is not open return Invalid; otherwise
410 return type. */
411 static int CheckOpen(int type);
412};
413
414/// Container class for integration rules
416{
417private:
418 /// Taken from the Quadrature1D class anonymous enum
419 /// Determines the type of numerical quadrature used for
420 /// segment, square, and cube geometries
421 const int quad_type;
422
423 int own_rules, refined;
424
425 Array<IntegrationRule *> PointIntRules;
426 Array<IntegrationRule *> SegmentIntRules;
427 Array<IntegrationRule *> TriangleIntRules;
428 Array<IntegrationRule *> SquareIntRules;
429 Array<IntegrationRule *> TetrahedronIntRules;
430 Array<IntegrationRule *> PyramidIntRules;
431 Array<IntegrationRule *> PrismIntRules;
432 Array<IntegrationRule *> CubeIntRules;
433
434#if defined(MFEM_THREAD_SAFE) && defined(MFEM_USE_OPENMP)
435 Array<omp_lock_t> IntRuleLocks;
436#endif
437
438 void AllocIntRule(Array<IntegrationRule *> &ir_array, int Order) const
439 {
440 if (ir_array.Size() <= Order)
441 {
442 ir_array.SetSize(Order + 1, NULL);
443 }
444 }
445 bool HaveIntRule(Array<IntegrationRule *> &ir_array, int Order) const
446 {
447 return (ir_array.Size() > Order && ir_array[Order] != NULL);
448 }
449 int GetSegmentRealOrder(int Order) const
450 {
451 return Order | 1; // valid for all quad_type's
452 }
453 void DeleteIntRuleArray(Array<IntegrationRule *> &ir_array) const;
454
455 /// The following methods allocate new IntegrationRule objects without
456 /// checking if they already exist. To avoid memory leaks use
457 /// IntegrationRules::Get(int GeomType, int Order) instead.
458 IntegrationRule *GenerateIntegrationRule(int GeomType, int Order);
459 IntegrationRule *PointIntegrationRule(int Order);
460 IntegrationRule *SegmentIntegrationRule(int Order);
461 IntegrationRule *TriangleIntegrationRule(int Order);
462 IntegrationRule *SquareIntegrationRule(int Order);
463 IntegrationRule *TetrahedronIntegrationRule(int Order);
464 IntegrationRule *PyramidIntegrationRule(int Order);
465 IntegrationRule *PrismIntegrationRule(int Order);
466 IntegrationRule *CubeIntegrationRule(int Order);
467
468public:
469 /// Sets initial sizes for the integration rule arrays, but rules
470 /// are defined the first time they are requested with the Get method.
471 explicit IntegrationRules(int ref = 0,
472 int type = Quadrature1D::GaussLegendre);
473
474 /// Returns an integration rule for given GeomType and Order.
475 const IntegrationRule &Get(int GeomType, int Order);
476
477 void Set(int GeomType, int Order, IntegrationRule &IntRule);
478
479 void SetOwnRules(int o) { own_rules = o; }
480
481 /// Destroys an IntegrationRules object
483};
484
485/// A global object with all integration rules (defined in intrules.cpp)
486extern MFEM_EXPORT IntegrationRules IntRules;
487
488/// A global object with all refined integration rules
490
491}
492
493#endif
Dynamic 2D array using row-major layout.
Definition array.hpp:372
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
Class for integration point with weight.
Definition intrules.hpp:35
void Get(real_t *p, const int dim) const
Definition intrules.hpp:60
void Set(const real_t *p, const int dim)
Definition intrules.hpp:46
void Set2w(const real_t *p)
Definition intrules.hpp:87
void Init(int const i)
Definition intrules.hpp:40
void Set2w(const real_t x1, const real_t x2, const real_t w)
Definition intrules.hpp:84
void Set2(const real_t *p)
Definition intrules.hpp:91
void Set3(const real_t *p)
Definition intrules.hpp:82
void Set1w(const real_t x1, const real_t w)
Definition intrules.hpp:93
void Set3w(const real_t *p)
Definition intrules.hpp:77
void Set2(const real_t x1, const real_t x2)
Definition intrules.hpp:89
void Set1w(const real_t *p)
Definition intrules.hpp:95
void Set(const real_t x1, const real_t x2, const real_t x3, const real_t w)
Definition intrules.hpp:74
void Set3(const real_t x1, const real_t x2, const real_t x3)
Definition intrules.hpp:79
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
~IntegrationRule()
Destroys an IntegrationRule object.
Definition intrules.hpp:274
IntegrationRule(int NP)
Construct an integration rule with given number of points.
Definition intrules.hpp:227
int GetOrder() const
Returns the order of the integration rule.
Definition intrules.hpp:249
const IntegrationPoint & IntPoint(int i) const
Returns a const reference to the i-th integration point.
Definition intrules.hpp:262
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationRule * ApplyToKnotIntervals(KnotVector const &kv) const
Return an integration rule for KnotVector kv, defined by applying this rule on each knot interval.
Definition intrules.cpp:181
const Array< real_t > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition intrules.cpp:86
void SetOrder(const int order)
Sets the order of the integration rule. This is only for keeping order information,...
Definition intrules.hpp:253
void SetPointIndices()
Sets the indices of each quadrature point on initialization.
Definition intrules.cpp:99
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
Container class for integration rules.
Definition intrules.hpp:416
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
IntegrationRules(int ref=0, int type=Quadrature1D::GaussLegendre)
Definition intrules.cpp:961
void SetOwnRules(int o)
Definition intrules.hpp:479
void Set(int GeomType, int Order, IntegrationRule &IntRule)
~IntegrationRules()
Destroys an IntegrationRules object.
Mesh data type.
Definition mesh.hpp:56
Class for defining different integration rules on each NURBS patch.
Definition intrules.hpp:279
const IntegrationRule * GetPatchRule1D(const int patch, const int dimension) const
For tensor product rules defined on each patch by SetPatchRules1D(), return a pointer to the 1D rule ...
Definition intrules.hpp:315
void Finalize(Mesh const &mesh)
Finalize() must be called before this class can be used for assembly. In particular,...
NURBSMeshRules(const int numPatches, const int dim_)
Construct a rule for each patch, using SetPatchRules1D.
Definition intrules.hpp:282
int GetPointElement(int patch, int i, int j, int k) const
For tensor product rules defined on each patch by SetPatchRules1D(), returns the index of the element...
Definition intrules.hpp:334
IntegrationRule & GetElementRule(const int elem, const int patch, const int *ijk, Array< const KnotVector * > const &kv, bool &deleteRule) const
Returns a rule for the element.
void SetPatchRules1D(const int patch, std::vector< const IntegrationRule * > &ir1D)
Set 1D integration rules to be used as a tensor product rule on the patch with index patch....
const Array< int > & GetPatchRule1D_KnotSpan(const int patch, const int dimension) const
For tensor product rules defined on each patch by SetPatchRules1D(), returns an array of knot span in...
Definition intrules.hpp:344
void GetIntegrationPointFrom1D(const int patch, int i, int j, int k, IntegrationPoint &ip)
For tensor product rules defined on each patch by SetPatchRules1D(), return the integration point wit...
void SetElementRule(const std::size_t element, const std::size_t elementRuleIndex)
Set the integration rule for the element of the given index. This rule is used instead of the rule fo...
Definition intrules.hpp:301
std::size_t AddElementRule(IntegrationRule *ir_element)
Add a rule to be used for individual elements. Returns the rule index.
Definition intrules.hpp:293
A class container for 1D quadrature type constants.
Definition intrules.hpp:394
@ ClosedUniform
aka closed Newton-Cotes
Definition intrules.hpp:402
@ ClosedGL
aka closed Gauss Legendre
Definition intrules.hpp:404
@ OpenHalfUniform
aka "open half" Newton-Cotes
Definition intrules.hpp:403
@ OpenUniform
aka open Newton-Cotes
Definition intrules.hpp:401
static int CheckOpen(int type)
If the Quadrature1D type is not open return Invalid; otherwise return type.
Definition intrules.cpp:941
static int CheckClosed(int type)
If the Quadrature1D type is not closed return Invalid; otherwise return type.
Definition intrules.cpp:928
A Class that defines 1-D numerical quadrature rules on [0,1].
Definition intrules.hpp:370
static void GaussLegendre(const int np, IntegrationRule *ir)
Definition intrules.cpp:424
static void ClosedUniform(const int np, IntegrationRule *ir)
Definition intrules.cpp:660
static void OpenUniform(const int np, IntegrationRule *ir)
Definition intrules.cpp:644
static void ClosedGL(const int np, IntegrationRule *ir)
Definition intrules.cpp:695
static void GivePolyPoints(const int np, real_t *pts, const int type)
Definition intrules.cpp:717
static void OpenHalfUniform(const int np, IntegrationRule *ir)
Definition intrules.cpp:680
static void GaussLobatto(const int np, IntegrationRule *ir)
Definition intrules.cpp:512
int dim
Definition ex24.cpp:53
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition hooke.cpp:45
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
float real_t
Definition config.hpp:43
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition intrules.hpp:489
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
real_t p(const Vector &x, real_t t)
RefCoord s[3]
void pts(int iphi, int t, real_t x[])