15#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_CYLINDERINTEGRATOR_HH 
   16#define DUMUX_MULTIDOMAIN_EMBEDDED_CYLINDERINTEGRATOR_HH 
   26#include <dune/common/fvector.hh> 
   49    using GlobalPosition = Dune::FieldVector<Scalar, 3>;
 
   59        rSamples_ = characteristicLengthToNumSamples_(rStep);
 
 
   72        using std::ceil; 
using std::clamp;
 
   73        rSamples_ = characteristicLengthToNumSamples_(rStep);
 
   74        zSamples_ = characteristicLengthToNumSamples_(zStep);
 
 
   86                     const GlobalPosition& topCenter,
 
   90        const auto zAxis = (topCenter-bottomCenter);
 
   91        const auto height = zAxis.two_norm();
 
   92        const auto zAxisUnitVector = zAxis/height;
 
   95        const auto rStep = radius/Scalar(rSamples_);
 
  100            zSamples_ = std::max<std::size_t>(1, ceil(height/rStep));
 
  101        const auto zStep = height/Scalar(zSamples_);
 
  104        auto kOffset = numRingSamples_;
 
  105        std::partial_sum(kOffset.begin(), kOffset.end(), kOffset.begin());
 
  108        numPoints_ = zSamples_*circleSamples_;
 
  110            std::cout << 
"CylinderIntegration -- number of integration points: " << numPoints_ << std::endl;
 
  113        integrationElement_.resize(numPoints_);
 
  114        points_.resize(numPoints_);
 
  117        std::vector<GlobalPosition> ringPoints;
 
  118        ringPoints.reserve(numRingSamples_.back());
 
  121        for (std::size_t i = 0; i < zSamples_; ++i)
 
  124            auto sliceCenter = bottomCenter;
 
  125            sliceCenter.axpy((Scalar(i)+0.5)*zStep, zAxisUnitVector);
 
  128            for (std::size_t j = 0; j < rSamples_; ++j)
 
  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];
 
  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)
 
  139                    const std::size_t idx = k + ringOffest + circleSamples_*i;
 
  140                    points_[idx] = ringPoints[k];
 
 
  149    { 
return integrationElement_[i]; }
 
 
  153    { 
return points_[i]; }
 
 
  157    { 
return numPoints_; }
 
 
  164        numRingSamples_.resize(rSamples_);
 
  165        for (
int j = 0; j < rSamples_; ++j)
 
  169            numRingSamples_[j] = floor(2*M_PI*(Scalar(j)+0.5));
 
  170            circleSamples_ += numRingSamples_[j];
 
  174        sincos_.resize(rSamples_);
 
  175        for (
int j = 0; j < rSamples_; ++j)
 
  177            const auto numPoints = numRingSamples_[j];
 
  178            sincos_[j].resize(2*numPoints);
 
  180            for (std::size_t i = 0; i < numPoints; ++i)
 
  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;
 
  193    std::size_t characteristicLengthToNumSamples_(
const Scalar cl)
 const 
  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));
 
  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_;
 
 
  212template<
class GlobalPosition>
 
  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)
 
  227    const auto da = (d*firstAxis);
 
  228    const auto db = (d*secondAxis);
 
  230    return (da*da/(a*a) + db*db/(b*b) < 1.0);
 
 
  234template<
class GlobalPosition>
 
  235inline std::pair<std::vector<GlobalPosition>, 
typename GlobalPosition::value_type>
 
  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)
 
  245    const auto aStep = characteristicLength;
 
  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;
 
  253    std::vector<GlobalPosition> points;
 
  254    points.reserve(aSamples*bSamples);
 
  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)
 
  263        posA.axpy(as*aStep, firstUnitAxis);
 
  264        for (std::size_t bs = 0; bs < bSamples; ++bs)
 
  267            pos.axpy(bs*bStep, secondUnitAxis);
 
  269                points.emplace_back(std::move(pos));
 
  273    return {points, aStep*bStep};
 
 
 
  288template<
class Scalar>
 
  291    using GlobalPosition = Dune::FieldVector<Scalar, 3>;
 
  300    : relCharLength_(relCharLength)
 
 
  312                     const GlobalPosition& topCenter,
 
  313                     const GlobalPosition& firstAxis,
 
  314                     const GlobalPosition& secondAxis,
 
  317        const auto a = firstAxis.two_norm();
 
  318        const auto b = secondAxis.two_norm();
 
  319        const auto characteristicLength = relCharLength_*a;
 
  321        auto firstUnitAxis = firstAxis; firstUnitAxis /= a;
 
  322        auto secondUnitAxis = secondAxis; secondUnitAxis /= b;
 
  325        auto zAxis = (topCenter-bottomCenter);
 
  326        const auto length = zAxis.two_norm();
 
  328        const auto height = (zAxis*
normal)*length;
 
  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;
 
  335        auto startZ = bottomCenter;
 
  336        startZ.axpy(0.5*zStep, zAxis);
 
  337        auto [ellipseIPs, ellipseIntegrationElement]
 
  341        const auto abPoints = ellipseIPs.size();
 
  342        numPoints_ = abPoints*zSamples;
 
  346        points_.reserve(numPoints_);
 
  347        std::move(ellipseIPs.begin(), ellipseIPs.end(), std::back_inserter(points_));
 
  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; });
 
  356        integrationElement_ = ellipseIntegrationElement*zStep/zStepFactor;
 
  357        const auto meanLocalError = (height*a*b*M_PI - integrationElement_*numPoints_)/Scalar(numPoints_);
 
  358        integrationElement_ += meanLocalError;
 
  362            std::cout << 
"EllipticCylinderIntegration -- number of integration points: " << numPoints_ <<
"\n";
 
  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";
 
 
  372    { 
return integrationElement_; }
 
 
  376    { 
return points_[i]; }
 
 
  380    { 
return numPoints_; }
 
 
  383    const Scalar relCharLength_;
 
  384    std::size_t numPoints_;
 
  385    Scalar integrationElement_;
 
  386    std::vector<GlobalPosition> points_;
 
 
  399template<
class Scalar>
 
  402    using GlobalPosition = Dune::FieldVector<Scalar, 3>;
 
  411    : relCharLength_(relCharLength)
 
 
  422                     const GlobalPosition& firstAxis,
 
  423                     const GlobalPosition& secondAxis,
 
  426        const auto a = firstAxis.two_norm();
 
  427        const auto b = secondAxis.two_norm();
 
  428        const auto characteristicLength = relCharLength_*a;
 
  430        auto firstUnitAxis = firstAxis; firstUnitAxis /= a;
 
  431        auto secondUnitAxis = secondAxis; secondUnitAxis /= b;
 
  435        std::tie(points_, integrationElement_)
 
  439        const auto abPoints = points_.size();
 
  440        numPoints_ = abPoints;
 
  443        const auto meanLocalError = (a*b*M_PI - integrationElement_*numPoints_)/Scalar(numPoints_);
 
  444        integrationElement_ += meanLocalError;
 
  448            std::cout << 
"EllipseIntegration -- number of integration points: " << numPoints_ << 
"\n";
 
  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";
 
 
  458    { 
return integrationElement_; }
 
 
  462    { 
return points_[i]; }
 
 
  466    { 
return numPoints_; }
 
 
  469    const Scalar relCharLength_;
 
  470    std::size_t numPoints_;
 
  471    Scalar integrationElement_;
 
  472    std::vector<GlobalPosition> points_;
 
 
Helper function to compute points on a circle.
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
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 ¢er, 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
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:671
void circlePoints(std::vector< GlobalPosition > &points, const std::vector< Scalar > &sincos, const GlobalPosition ¢er, 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.
Definition cylinderintegration.hh:210
std::pair< std::vector< GlobalPosition >, typename GlobalPosition::value_type > ellipseIntegrationPoints(const GlobalPosition ¢er, 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 ¢er, 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