version 3.8
cylinderintegration.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
15#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_CYLINDERINTEGRATOR_HH
16#define DUMUX_MULTIDOMAIN_EMBEDDED_CYLINDERINTEGRATOR_HH
17
18#include <cmath>
19#include <algorithm>
20#include <numeric>
21#include <iterator>
22#include <utility>
23#include <tuple>
24#include <limits>
25
26#include <dune/common/fvector.hh>
27
28#include <dumux/common/math.hh>
30
32
46template<class Scalar>
48{
49 using GlobalPosition = Dune::FieldVector<Scalar, 3>;
50
51public:
57 CylinderIntegration(const Scalar rStep)
58 {
59 rSamples_ = characteristicLengthToNumSamples_(rStep);
60 initialize_();
61 }
62
69 CylinderIntegration(const Scalar rStep, const Scalar zStep)
70 : zStepFixed_(true)
71 {
72 using std::ceil; using std::clamp;
73 rSamples_ = characteristicLengthToNumSamples_(rStep);
74 zSamples_ = characteristicLengthToNumSamples_(zStep);
75 initialize_();
76 }
77
85 void setGeometry(const GlobalPosition& bottomCenter,
86 const GlobalPosition& topCenter,
87 const Scalar radius,
88 int verbosity = 0)
89 {
90 const auto zAxis = (topCenter-bottomCenter);
91 const auto height = zAxis.two_norm();
92 const auto zAxisUnitVector = zAxis/height;
93
94 // compute step size in r-direction
95 const auto rStep = radius/Scalar(rSamples_);
96
97 // compute zSamples from r samples if not specified by the user
98 using std::ceil;
99 if (!zStepFixed_)
100 zSamples_ = std::max<std::size_t>(1, ceil(height/rStep));
101 const auto zStep = height/Scalar(zSamples_);
102
103 // compute offsets for index calculation
104 auto kOffset = numRingSamples_;
105 std::partial_sum(kOffset.begin(), kOffset.end(), kOffset.begin());
106
107 // compute total number of samples
108 numPoints_ = zSamples_*circleSamples_;
109 if (verbosity > 0)
110 std::cout << "CylinderIntegration -- number of integration points: " << numPoints_ << std::endl;
111
112 // resize integration points and integration elements
113 integrationElement_.resize(numPoints_);
114 points_.resize(numPoints_);
115
116 // reserve enough memory for the points on each ring
117 std::vector<GlobalPosition> ringPoints;
118 ringPoints.reserve(numRingSamples_.back());
119
120 // compute integration points and integration elements
121 for (std::size_t i = 0; i < zSamples_; ++i)
122 {
123 // for each cylinder slice i
124 auto sliceCenter = bottomCenter;
125 sliceCenter.axpy((Scalar(i)+0.5)*zStep, zAxisUnitVector);
126
127 // generate circle points for each ring j
128 for (std::size_t j = 0; j < rSamples_; ++j)
129 {
130 const auto r = (Scalar(j)+0.5)*rStep;
131 circlePoints(ringPoints, sincos_[j], sliceCenter, zAxisUnitVector, r);
132 const auto ringOffest = j > 0 ? kOffset[j-1] : 0;
133 const auto ringSamples = numRingSamples_[j];
134 // volume of each element in ring
135 // total ring volume is given by M_PI*rStep*rStep*zStep*((j+1)^2 - j^2)
136 const auto integrationElement = M_PI*rStep*rStep*zStep*(1.0 + 2.0*Scalar(j))/ringSamples;
137 for (std::size_t k = 0; k < ringSamples; ++k)
138 {
139 const std::size_t idx = k + ringOffest + circleSamples_*i;
140 points_[idx] = ringPoints[k];
141 integrationElement_[idx] = integrationElement;
142 }
143 }
144 }
145 }
146
148 Scalar integrationElement(std::size_t i) const
149 { return integrationElement_[i]; }
150
152 const GlobalPosition& integrationPoint(std::size_t i) const
153 { return points_[i]; }
154
156 std::size_t size() const
157 { return numPoints_; }
158
159private:
160 void initialize_()
161 {
162 // precompute number of cells in ring and total circle samples
163 circleSamples_ = 0;
164 numRingSamples_.resize(rSamples_);
165 for (int j = 0; j < rSamples_; ++j)
166 {
167 // number of cells in ring j
168 using std::floor;
169 numRingSamples_[j] = floor(2*M_PI*(Scalar(j)+0.5));
170 circleSamples_ += numRingSamples_[j];
171 }
172
173 // optimization because calling too many sin/cos can be really expensive
174 sincos_.resize(rSamples_);
175 for (int j = 0; j < rSamples_; ++j)
176 {
177 const auto numPoints = numRingSamples_[j];
178 sincos_[j].resize(2*numPoints);
179 Scalar t = 0 + 0.1;
180 for (std::size_t i = 0; i < numPoints; ++i)
181 {
182 using std::sin; using std::cos;
183 sincos_[j][2*i] = sin(t);
184 sincos_[j][2*i + 1] = cos(t);
185 t += 2*M_PI/numPoints;
186 if (t > 2*M_PI) t -= 2*M_PI;
187 }
188 }
189 }
190
191 // invert characteristic length (between 0 and 1)
192 // and make sure there is no integer wrap around for small ratios
193 std::size_t characteristicLengthToNumSamples_(const Scalar cl) const
194 {
195 using std::clamp; using std::min; using std::max; using std::ceil; using std::nexttoward;
196 const Scalar clampedCl = clamp(cl, std::numeric_limits<Scalar>::min(), 1.0);
197 const Scalar floatNumSamples = ceil(1.0/clampedCl);
198 const Scalar largestValidFloat = nexttoward(std::numeric_limits<std::size_t>::max(), 0.0);
199 return static_cast<std::size_t>(min(floatNumSamples, largestValidFloat));
200 }
201
202 bool zStepFixed_ = false;
203 std::size_t zSamples_, rSamples_, numPoints_, circleSamples_;
204 std::vector<Scalar> integrationElement_;
205 std::vector<GlobalPosition> points_;
206 std::vector<std::size_t> numRingSamples_;
207 std::vector<std::vector<Scalar>> sincos_;
208};
209
210namespace Detail {
212template<class GlobalPosition>
213inline bool pointInEllipse(const GlobalPosition& p,
214 const GlobalPosition& center,
215 const GlobalPosition& firstAxis,
216 const GlobalPosition& secondAxis,
217 const GlobalPosition& normal,
218 const typename GlobalPosition::value_type a,
219 const typename GlobalPosition::value_type b)
220{
221 const auto d = p-center;
222 // early return if point is not on the ellipse plane
223 if (d*normal > 1e-7*a)
224 return false;
225
226 // check if it's actually within the ellipse
227 const auto da = (d*firstAxis);
228 const auto db = (d*secondAxis);
229
230 return (da*da/(a*a) + db*db/(b*b) < 1.0);
231}
232
234template<class GlobalPosition>
235inline std::pair<std::vector<GlobalPosition>, typename GlobalPosition::value_type>
236ellipseIntegrationPoints(const GlobalPosition& center,
237 const GlobalPosition& firstUnitAxis,
238 const GlobalPosition& secondUnitAxis,
239 typename GlobalPosition::value_type a,
240 typename GlobalPosition::value_type b,
241 const GlobalPosition& normal,
242 typename GlobalPosition::value_type characteristicLength)
243{
244 // choose step size in a-direction
245 const auto aStep = characteristicLength;
246
247 using std::floor;
248 const std::size_t aSamples = std::max<std::size_t>(floor(2*a/aStep), 1);
249 const std::size_t bSamples = std::max<std::size_t>(floor(2*b/aStep), 1);
250 const auto bStep = 2*b / bSamples;
251
252 // reserve upper limit estimate for memory needed
253 std::vector<GlobalPosition> points;
254 points.reserve(aSamples*bSamples);
255
256 // walk over lattice grid and reject points outside the ellipse
257 auto startAB = center;
258 startAB.axpy(-a + 0.5*aStep, firstUnitAxis);
259 startAB.axpy(-b + 0.5*bStep, secondUnitAxis);
260 for (std::size_t as = 0; as < aSamples; ++as)
261 {
262 auto posA = startAB;
263 posA.axpy(as*aStep, firstUnitAxis);
264 for (std::size_t bs = 0; bs < bSamples; ++bs)
265 {
266 auto pos = posA;
267 pos.axpy(bs*bStep, secondUnitAxis);
268 if (pointInEllipse(pos, center, firstUnitAxis, secondUnitAxis, normal, a, b))
269 points.emplace_back(std::move(pos));
270 }
271 }
272 // return points and integration element
273 return {points, aStep*bStep};
274}
275} // end namespace Detail
276
288template<class Scalar>
290{
291 using GlobalPosition = Dune::FieldVector<Scalar, 3>;
292
293public:
299 explicit EllipticCylinderIntegration(const Scalar relCharLength)
300 : relCharLength_(relCharLength)
301 {}
302
311 void setGeometry(const GlobalPosition& bottomCenter,
312 const GlobalPosition& topCenter,
313 const GlobalPosition& firstAxis,
314 const GlobalPosition& secondAxis,
315 int verbosity = 0)
316 {
317 const auto a = firstAxis.two_norm();
318 const auto b = secondAxis.two_norm();
319 const auto characteristicLength = relCharLength_*a;
320
321 auto firstUnitAxis = firstAxis; firstUnitAxis /= a;
322 auto secondUnitAxis = secondAxis; secondUnitAxis /= b;
323 const auto normal = crossProduct(firstUnitAxis, secondUnitAxis);
324
325 auto zAxis = (topCenter-bottomCenter);
326 const auto length = zAxis.two_norm();
327 zAxis /= length;
328 const auto height = (zAxis*normal)*length;
329 using std::floor;
330 const std::size_t zSamples = std::max<std::size_t>(floor(height/characteristicLength), 1);
331 const auto zStep = length / Scalar(zSamples);
332 const auto zStepFactor = length/height;
333
334 // walk over lattice grid and reject points outside the ellipse
335 auto startZ = bottomCenter;
336 startZ.axpy(0.5*zStep, zAxis);
337 auto [ellipseIPs, ellipseIntegrationElement]
338 = Detail::ellipseIntegrationPoints(startZ, firstUnitAxis, secondUnitAxis, a, b, normal, characteristicLength);
339
340 // compute total number of points
341 const auto abPoints = ellipseIPs.size();
342 numPoints_ = abPoints*zSamples;
343
344 // first move in the first layer of ellipse integration points
345 points_.clear();
346 points_.reserve(numPoints_);
347 std::move(ellipseIPs.begin(), ellipseIPs.end(), std::back_inserter(points_));
348
349 // then extrude along the z axis
350 auto zStepVector = zAxis; zStepVector *= zStep;
351 for (std::size_t zs = 1; zs < zSamples; ++zs)
352 std::transform(points_.end() - abPoints, points_.end(), std::back_inserter(points_),
353 [&](const auto& p){ return p + zStepVector; });
354
355 // computation and error correction for integral element to obtain exact integral of constant functions
356 integrationElement_ = ellipseIntegrationElement*zStep/zStepFactor;
357 const auto meanLocalError = (height*a*b*M_PI - integrationElement_*numPoints_)/Scalar(numPoints_);
358 integrationElement_ += meanLocalError;
359
360 if (verbosity > 0)
361 {
362 std::cout << "EllipticCylinderIntegration -- number of integration points: " << numPoints_ <<"\n";
363 if (verbosity > 1)
364 std::cout << " -- a: " << a << ", b: " << b << ", h: " << height << "\n"
365 << " -- volume: " << integrationElement_*numPoints_ << " (expected " << height*a*b*M_PI << ")\n"
366 << " -- corrected a volume error of: " << meanLocalError/integrationElement_*100 << "%\n";
367 }
368 }
369
371 Scalar integrationElement(std::size_t i) const
372 { return integrationElement_; }
373
375 const GlobalPosition& integrationPoint(std::size_t i) const
376 { return points_[i]; }
377
379 std::size_t size() const
380 { return numPoints_; }
381
382private:
383 const Scalar relCharLength_;
384 std::size_t numPoints_;
385 Scalar integrationElement_;
386 std::vector<GlobalPosition> points_;
387};
388
399template<class Scalar>
401{
402 using GlobalPosition = Dune::FieldVector<Scalar, 3>;
403
404public:
410 explicit EllipseIntegration(const Scalar relCharLength)
411 : relCharLength_(relCharLength)
412 {}
413
421 void setGeometry(const GlobalPosition& center,
422 const GlobalPosition& firstAxis,
423 const GlobalPosition& secondAxis,
424 int verbosity = 0)
425 {
426 const auto a = firstAxis.two_norm();
427 const auto b = secondAxis.two_norm();
428 const auto characteristicLength = relCharLength_*a;
429
430 auto firstUnitAxis = firstAxis; firstUnitAxis /= a;
431 auto secondUnitAxis = secondAxis; secondUnitAxis /= b;
432 const auto normal = crossProduct(firstUnitAxis, secondUnitAxis);
433
434 // generate integration points
435 std::tie(points_, integrationElement_)
436 = Detail::ellipseIntegrationPoints(center, firstUnitAxis, secondUnitAxis, a, b, normal, characteristicLength);
437
438 // store number of sample points
439 const auto abPoints = points_.size();
440 numPoints_ = abPoints;
441
442 // computation and error correction for integral element to obtain exact integral of constant functions
443 const auto meanLocalError = (a*b*M_PI - integrationElement_*numPoints_)/Scalar(numPoints_);
444 integrationElement_ += meanLocalError;
445
446 if (verbosity > 0)
447 {
448 std::cout << "EllipseIntegration -- number of integration points: " << numPoints_ << "\n";
449 if (verbosity > 1)
450 std::cout << " -- a: " << a << ", b: " << b << "\n"
451 << " -- area: " << integrationElement_*numPoints_ << " (expected " << a*b*M_PI << ")\n"
452 << " -- corrected an area error of: " << meanLocalError/integrationElement_*100 << "%\n";
453 }
454 }
455
457 Scalar integrationElement(std::size_t i) const
458 { return integrationElement_; }
459
461 const GlobalPosition& integrationPoint(std::size_t i) const
462 { return points_[i]; }
463
465 std::size_t size() const
466 { return numPoints_; }
467
468private:
469 const Scalar relCharLength_;
470 std::size_t numPoints_;
471 Scalar integrationElement_;
472 std::vector<GlobalPosition> points_;
473};
474
475} // end namespace Dumux::EmbeddedCoupling
476
477#endif
Helper function to compute points on a circle.
Helper class to integrate over a cylinder domain.
Definition: cylinderintegration.hh:48
CylinderIntegration(const Scalar rStep, const Scalar zStep)
Constructor.
Definition: cylinderintegration.hh:69
void setGeometry(const GlobalPosition &bottomCenter, const GlobalPosition &topCenter, const Scalar radius, int verbosity=0)
Set the geometry of the cylinder.
Definition: cylinderintegration.hh:85
std::size_t size() const
The number of integration points.
Definition: cylinderintegration.hh:156
Scalar integrationElement(std::size_t i) const
The integration element of the ith integration point.
Definition: cylinderintegration.hh:148
const GlobalPosition & integrationPoint(std::size_t i) const
The ith integration point.
Definition: cylinderintegration.hh:152
CylinderIntegration(const Scalar rStep)
Constructor.
Definition: cylinderintegration.hh:57
Helper class to integrate over an elliptic domain.
Definition: cylinderintegration.hh:401
const GlobalPosition & integrationPoint(std::size_t i) const
obtain ith integration point
Definition: cylinderintegration.hh:461
EllipseIntegration(const Scalar relCharLength)
Constructor.
Definition: cylinderintegration.hh:410
Scalar integrationElement(std::size_t i) const
obtain ith integration element
Definition: cylinderintegration.hh:457
void setGeometry(const GlobalPosition &center, const GlobalPosition &firstAxis, const GlobalPosition &secondAxis, int verbosity=0)
set geometry of an ellipse
Definition: cylinderintegration.hh:421
std::size_t size() const
number of integration points
Definition: cylinderintegration.hh:465
Helper class to integrate over an elliptic cylinder domain.
Definition: cylinderintegration.hh:290
void setGeometry(const GlobalPosition &bottomCenter, const GlobalPosition &topCenter, const GlobalPosition &firstAxis, const GlobalPosition &secondAxis, int verbosity=0)
Set the geometry of elliptic cylinder.
Definition: cylinderintegration.hh:311
const GlobalPosition & integrationPoint(std::size_t i) const
obtain ith integration point
Definition: cylinderintegration.hh:375
Scalar integrationElement(std::size_t i) const
obtain ith integration element
Definition: cylinderintegration.hh:371
EllipticCylinderIntegration(const Scalar relCharLength)
Constructor.
Definition: cylinderintegration.hh:299
std::size_t size() const
number of samples points
Definition: cylinderintegration.hh:379
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:642
void circlePoints(std::vector< GlobalPosition > &points, const std::vector< Scalar > &sincos, const GlobalPosition &center, const GlobalPosition &normal, const Scalar radius)
Definition: circlepoints.hh:38
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Define some often used mathematical functions.
std::pair< std::vector< GlobalPosition >, typename GlobalPosition::value_type > ellipseIntegrationPoints(const GlobalPosition &center, const GlobalPosition &firstUnitAxis, const GlobalPosition &secondUnitAxis, typename GlobalPosition::value_type a, typename GlobalPosition::value_type b, const GlobalPosition &normal, typename GlobalPosition::value_type characteristicLength)
construct evenly distributed integration points on an ellipse
Definition: cylinderintegration.hh:236
bool pointInEllipse(const GlobalPosition &p, const GlobalPosition &center, const GlobalPosition &firstAxis, const GlobalPosition &secondAxis, const GlobalPosition &normal, const typename GlobalPosition::value_type a, const typename GlobalPosition::value_type b)
check if a point is in an ellipse
Definition: cylinderintegration.hh:213
Definition: circlepoints.hh:24