13#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH 
   14#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH 
   18#include <dune/common/timer.hh> 
   19#include <dune/geometry/quadraturerules.hh> 
   36    static std::string 
name() { 
return "average"; }
 
 
 
   43template<
class MDTraits, 
class CouplingMode>
 
   44class Embedded1d3dCouplingManager;
 
   52template<
class MDTraits>
 
   53class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>
 
   54: 
public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>,
 
   55                                     CircleAveragePointSourceTraits<MDTraits>>
 
   57    using ThisType = Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Average>;
 
   58    using ParentType = EmbeddedCouplingManagerBase<MDTraits, ThisType, CircleAveragePointSourceTraits<MDTraits>>;
 
   59    using Scalar = 
typename MDTraits::Scalar;
 
   60    using SolutionVector = 
typename MDTraits::SolutionVector;
 
   61    using PointSourceData = 
typename ParentType::PointSourceTraits::PointSourceData;
 
   63    static constexpr auto bulkIdx = 
typename MDTraits::template SubDomain<0>::Index();
 
   64    static constexpr auto lowDimIdx = 
typename MDTraits::template SubDomain<1>::Index();
 
   67    template<std::
size_t id> 
using SubDomainTypeTag = 
typename MDTraits::template SubDomain<id>::TypeTag;
 
   68    template<std::
size_t id> 
using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
 
   69    template<std::
size_t id> 
using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
 
   70    template<std::
size_t id> 
using GridView = 
typename GridGeometry<id>::GridView;
 
   71    template<std::
size_t id> 
using Element = 
typename GridView<id>::template Codim<0>::Entity;
 
   72    template<std::
size_t id> 
using GridIndex = 
typename IndexTraits<GridView<id>>::GridIndex;
 
   74    using GlobalPosition = 
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
 
   76    template<std::
size_t id>
 
   77    static constexpr bool isBox()
 
   83        bulkDim = GridView<bulkIdx>::dimension,
 
   84        lowDimDim = GridView<lowDimIdx>::dimension,
 
   85        dimWorld = GridView<bulkIdx>::dimensionworld
 
   88    static constexpr Embedded1d3dCouplingMode::Average couplingMode{};
 
   90    using ParentType::ParentType;
 
   92    void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
 
   93              std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
 
   94              const SolutionVector& curSol)
 
   96        ParentType::init(bulkProblem, lowDimProblem, curSol);
 
   97        computeLowDimVolumeFractions();
 
  106    template<std::
size_t id, 
class JacobianPattern>
 
  107    void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern)
 const 
  109        extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
 
  119    template<std::
size_t i, 
class LocalAssemblerI, 
class JacobianMatrixDiagBlock, 
class Gr
idVariables>
 
  120    void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
 
  121                                         const LocalAssemblerI& localAssemblerI,
 
  122                                         const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
 
  123                                         JacobianMatrixDiagBlock& A,
 
  124                                         GridVariables& gridVariables)
 
  126        extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, A, gridVariables);
 
  138    void computePointSourceData(std::size_t order = 1, 
bool verbose = 
false)
 
  141        const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
 
  142        const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
 
  143        const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
 
  148        std::cout << 
"Initializing the point sources..." << std::endl;
 
  153        extendedSourceStencil_.clear();
 
  156        this->precomputeVertexIndices(bulkIdx);
 
  157        this->precomputeVertexIndices(lowDimIdx);
 
  163        const auto& lowDimProblem = this->problem(lowDimIdx);
 
  164        for (
const auto& is : intersections(this->glue()))
 
  167            const auto& lowDimElement = is.targetEntity(0);
 
  168            const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
 
  171            const auto intersectionGeometry = is.geometry();
 
  173            const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
 
  182            for (
auto&& qp : quad)
 
  185                const auto globalPos = intersectionGeometry.global(qp.position());
 
  191                if (bulkElementIndices.empty())
 
  198                static const auto numIp = 
getParam<int>(
"MixedDimension.NumCircleSegments");
 
  199                const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
 
  200                const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
 
  201                const auto circleAvgWeight = 2*M_PI*radius/numIp;
 
  202                const auto integrationElement = intersectionGeometry.integrationElement(qp.position());
 
  203                const auto qpweight = qp.weight();
 
  206                std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(
circlePoints.size());
 
  207                std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(
circlePoints.size());
 
  210                std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
 
  211                using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
 
  212                std::vector<ShapeValues> circleShapeValues;
 
  215                int insideCirclePoints = 0;
 
  219                    if (circleBulkElementIndices.empty())
 
  222                    ++insideCirclePoints;
 
  223                    const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices.size();
 
  224                    for (
const auto bulkElementIdx : circleBulkElementIndices)
 
  226                        circleStencil.push_back(bulkElementIdx);
 
  227                        circleIpWeight.push_back(localCircleAvgWeight);
 
  230                        if constexpr (isBox<bulkIdx>())
 
  232                            const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
 
  233                            circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
 
  236                            const auto bulkGeometry = bulkElement.geometry();
 
  237                            ShapeValues shapeValues;
 
  238                            this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues);
 
  239                            circleShapeValues.emplace_back(std::move(shapeValues));
 
  246                if (circleStencil.empty())
 
  250                if constexpr (isBox<bulkIdx>())
 
  253                    for (
const auto& vertices : circleCornerIndices)
 
  255                        this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
 
  256                                                                                   vertices->begin(), vertices->end());
 
  262                    this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
 
  263                                                                               circleStencil.begin(), circleStencil.end());
 
  267                const auto surfaceFraction = Scalar(insideCirclePoints)/Scalar(
circlePoints.size());
 
  270                for (
auto bulkElementIdx : bulkElementIndices)
 
  272                    const auto id = this->idCounter_++;
 
  274                    this->pointSources(bulkIdx).emplace_back(globalPos, 
id, qpweight, integrationElement*surfaceFraction, bulkElementIdx);
 
  275                    this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
 
  276                    this->pointSources(lowDimIdx).emplace_back(globalPos, 
id, qpweight, integrationElement*surfaceFraction, lowDimElementIdx);
 
  277                    this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
 
  281                    PointSourceData psData;
 
  283                    if constexpr (isBox<lowDimIdx>())
 
  285                        ShapeValues shapeValues;
 
  286                        this->getShapeValues(lowDimIdx, lowDimGridGeometry, lowDimElement.geometry(), globalPos, shapeValues);
 
  287                        psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
 
  291                        psData.addLowDimInterpolation(lowDimElementIdx);
 
  295                    if constexpr (isBox<bulkIdx>())
 
  297                        psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
 
  299                        const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
 
  300                        ShapeValues shapeValues;
 
  301                        this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, globalPos, shapeValues);
 
  302                        psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
 
  306                        psData.addCircleInterpolation(circleIpWeight, circleStencil);
 
  307                        psData.addBulkInterpolation(bulkElementIdx);
 
  311                    this->pointSourceData().emplace_back(std::move(psData));
 
  314                    const auto outsideGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
 
  318                    if constexpr (isBox<lowDimIdx>())
 
  320                        this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
 
  321                                                                               this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
 
  322                                                                               this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
 
  327                        this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
 
  331                    if constexpr (isBox<bulkIdx>())
 
  334                        for (
const auto& vertices : circleCornerIndices)
 
  336                            extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
 
  337                                                                                    vertices->begin(), vertices->end());
 
  343                        extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
 
  344                                                                   circleStencil.begin(), circleStencil.end());
 
  351        for (
auto&& stencil : extendedSourceStencil_.stencil())
 
  353            std::sort(stencil.second.begin(), stencil.second.end());
 
  354            stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
 
  357            if constexpr (isBox<bulkIdx>())
 
  359                const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
 
  360                stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
 
  361                                                   [&](
auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
 
  362                                     stencil.second.end());
 
  367                stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
 
  368                                                   [&](
auto i){ return i == stencil.first; }),
 
  369                                     stencil.second.end());
 
  374        using namespace Dune::Hybrid;
 
  375        forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
 
  377            for (
auto&& stencil : this->couplingStencils(domainIdx))
 
  379                std::sort(stencil.second.begin(), stencil.second.end());
 
  380                stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
 
  384        std::cout << 
"took " << watch.elapsed() << 
" seconds." << std::endl;
 
  388    void computeLowDimVolumeFractions()
 
  391        lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
 
  393        const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
 
  394        const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
 
  397        for (
const auto& is : intersections(this->glue()))
 
  400            const auto& inside = is.targetEntity(0);
 
  401            const auto intersectionGeometry = is.geometry();
 
  402            const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
 
  405            const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
 
  406            for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
 
  408                const auto& outside = is.domainEntity(outsideIdx);
 
  409                const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
 
  410                lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
 
  421    Scalar radius(std::size_t 
id)
 const 
  423        const auto& data = this->pointSourceData()[id];
 
  424        return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
 
  429    Scalar lowDimVolume(
const Element<bulkIdx>& element)
 const 
  431        const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
 
  432        return lowDimVolumeInBulkElement_[eIdx];
 
  437    Scalar lowDimVolumeFraction(
const Element<bulkIdx>& element)
 const 
  439        const auto totalVolume = 
element.geometry().volume();
 
  440        return lowDimVolume(element) / totalVolume;
 
  448    const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type&
 
  449    extendedSourceStencil(std::size_t eIdx)
 const 
  451        const auto& sourceStencils = extendedSourceStencil_.stencil();
 
  452        if (
auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
 
  453            return stencil->second;
 
  455        return this->emptyStencil(bulkIdx);
 
  460    EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_;
 
  463    std::vector<Scalar> lowDimVolumeInBulkElement_;
 
  467template<
class MDTraits>
 
  469: 
public std::true_type {};
 
 
Point source traits for average-based coupling modes.
Helper function to compute points on a circle.
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition couplingmanager1d3d.hh:24
Defines all properties used in Dumux.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Helper functions for distance queries.
Extended source stencil helper class for coupling managers.
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
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition intersectingentities.hh:102
static Geometry::ctype averageDistancePointGeometry(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry, std::size_t integrationOrder=2)
Compute the average distance from a point to a geometry by integration.
Definition distance.hh:29
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
constexpr Box box
Definition method.hh:147
Definition couplingmanager1d3d_average.hh:34
constexpr Average average
Definition couplingmanager1d3d_average.hh:39
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition multistagemultidomainfvassembler.hh:78
Definition couplingmanager1d3d_average.hh:35
static std::string name()
Definition couplingmanager1d3d_average.hh:36
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition tag.hh:30
Helper class to create (named and comparable) tagged types.