25#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_POROUSMEDIUM_COUPLINGMAPPER_HH
26#define DUMUX_MULTIDOMAIN_FREEFLOW_POROUSMEDIUM_COUPLINGMAPPER_HH
31#include <unordered_map>
35#include <dune/common/timer.hh>
36#include <dune/common/exceptions.hh>
37#include <dune/common/indices.hh>
38#include <dune/common/reservedvector.hh>
39#include <dune/geometry/axisalignedcubegeometry.hh>
55 using MapType = std::unordered_map<std::size_t, std::vector<std::size_t>>;
61 template<
class FreeFlowMomentumGr
idGeometry,
class FreeFlowMassGr
idGeometry,
class PoreNetworkGr
idGeometry>
62 void update(
const FreeFlowMomentumGridGeometry& ffMomentumGridGeometry,
63 const FreeFlowMassGridGeometry& ffMassGridGeometry,
64 const PoreNetworkGridGeometry& pnmGridGeometry)
67 resize_(ffMomentumGridGeometry, pnmGridGeometry);
69 std::cout <<
"Initializing the coupling map..." << std::endl;
71 auto ffFvGeometry =
localView(ffMomentumGridGeometry);
72 auto pnmFvGeometry =
localView(pnmGridGeometry);
74 using GlobalPosition =
typename FreeFlowMomentumGridGeometry::GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
76 for (
const auto& pnmElement : elements(pnmGridGeometry.gridView()))
78 const auto pnmElementIdx = pnmGridGeometry.elementMapper().index(pnmElement);
79 pnmFvGeometry.bindElement(pnmElement);
80 for (
const auto& pnmScv : scvs(pnmFvGeometry))
83 if (!pnmGridGeometry.dofOnBoundary(pnmScv.dofIndex()))
87 const auto pnmPos = pnmScv.dofPosition();
88 const auto pnmDofIdx = pnmScv.dofIndex();
89 const auto& otherPNMScv = pnmFvGeometry.scv(1 - pnmScv.indexInElement());
90 const auto otherPNMScvDofIdx = otherPNMScv.dofIndex();
93 const auto directlyCoupledFreeFlowElements =
intersectingEntities(pnmPos, ffMomentumGridGeometry.boundingBoxTree());
94 if (directlyCoupledFreeFlowElements.empty())
97 isCoupledPNMDof_[pnmDofIdx] =
true;
102 const std::size_t couplingNormalDirectionIndex = [&]
104 using Key = std::pair<std::size_t, bool>;
105 std::map<Key, std::size_t> result;
106 for (
const auto eIdx : directlyCoupledFreeFlowElements)
108 for (
const auto& intersection : intersections(ffMomentumGridGeometry.gridView(), ffMomentumGridGeometry.element(eIdx)))
112 const auto&
normal = intersection.centerUnitOuterNormal();
121 if (directlyCoupledFreeFlowElements.size() == 1 && result.size() > 1)
122 DUNE_THROW(Dune::InvalidStateException,
"Pore may not intersect with faces of different orientation when coupled to only one element");
124 return std::max_element(result.begin(), result.end(), [](
const auto& x,
const auto& y) { return x.second < y.second;})->first.first;
127 using Scalar =
typename FreeFlowMomentumGridGeometry::GridView::ctype;
128 const Scalar couplingPoreRadius = pnmGridGeometry.poreInscribedRadius(pnmDofIdx);
129 GlobalPosition lowerLeft = pnmPos - GlobalPosition(couplingPoreRadius);
130 lowerLeft[couplingNormalDirectionIndex] = pnmPos[couplingNormalDirectionIndex];
131 GlobalPosition upperRight = pnmPos + GlobalPosition(couplingPoreRadius);
132 upperRight[couplingNormalDirectionIndex] = pnmPos[couplingNormalDirectionIndex];
134 auto axes = std::move(std::bitset<FreeFlowMomentumGridGeometry::Grid::dimensionworld>{}.set());
135 axes.set(couplingNormalDirectionIndex,
false);
137 using PoreIntersectionGeometryType = Dune::AxisAlignedCubeGeometry<Scalar,
138 FreeFlowMomentumGridGeometry::GridView::dimension-1,
139 FreeFlowMomentumGridGeometry::GridView::dimensionworld>;
141 PoreIntersectionGeometryType poreIntersectionGeometry(lowerLeft, upperRight, axes);
142 const auto allCoupledFreeFlowElements =
intersectingEntities(std::move(poreIntersectionGeometry), ffMomentumGridGeometry.boundingBoxTree());
145 for (
const auto& ffElementInfo : allCoupledFreeFlowElements)
147 const auto freeFlowElementIndex = ffElementInfo.second();
148 pnmElementToFreeFlowElementsMap_[pnmElementIdx].push_back(freeFlowElementIndex);
149 freeFlowElementToPNMElementMap_[freeFlowElementIndex] = pnmElementIdx;
151 pnmToFreeFlowMassStencils_[pnmElementIdx].push_back(freeFlowElementIndex);
152 freeFlowMassToPNMStencils_[freeFlowElementIndex].push_back(pnmDofIdx);
154 ffFvGeometry.bindElement(ffMomentumGridGeometry.element(freeFlowElementIndex));
155 const auto coupledFreeFlowMomentumDofIndices = coupledFFMomentumDofs_(ffFvGeometry, ffMassGridGeometry, pnmPos, couplingPoreRadius, couplingNormalDirectionIndex);
157 pnmToFreeFlowMomentumStencils_[pnmElementIdx].push_back(coupledFreeFlowMomentumDofIndices.coupledFrontalDof);
158 freeFlowMomentumToPNMStencils_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof].push_back(pnmDofIdx);
159 freeFlowMomentumToPNMStencils_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof].push_back(otherPNMScvDofIdx);
161 isCoupledFreeFlowMomentumDof_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof] =
true;
162 isCoupledFreeFlowMomentumDofOnInterface_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof] =
true;
165 for (
const auto ffDofIdx : coupledFreeFlowMomentumDofIndices.coupledLateralDofs)
167 freeFlowMomentumToPNMStencils_[ffDofIdx].push_back(pnmDofIdx);
168 freeFlowMomentumToPNMStencils_[ffDofIdx].push_back(otherPNMScvDofIdx);
169 isCoupledFreeFlowMomentumDof_[ffDofIdx] =
true;
175 std::cout <<
"took " << watch.elapsed() <<
" seconds." << std::endl;
186 return pnmToFreeFlowMomentumStencils_.at(eIdxI);
188 return emptyStencil_;
199 return pnmToFreeFlowMassStencils_.at(eIdxI);
201 return emptyStencil_;
212 return freeFlowMassToPNMStencils_.at(eIdxI);
214 return emptyStencil_;
225 return freeFlowMomentumToPNMStencils_.at(dofIndex);
227 return emptyStencil_;
235 return static_cast<bool>(freeFlowElementToPNMElementMap_.count(eIdx));
243 return isCoupledFreeFlowMomentumDof_[dofIdx];
251 return static_cast<bool>(pnmElementToFreeFlowElementsMap_.count(eIdx));
259 return isCoupledPNMDof_[dofIdx];
264 return isCoupledFrontalFreeFlowMomentumScvf_.count(scvfIdx);
269 return isCoupledLateralFreeFlowMomentumScvf_.count(scvfIdx);
274 return isCoupledFreeFlowMassScvf_.count(scvfIdx);
278 {
return pnmElementToFreeFlowElementsMap_;}
281 {
return freeFlowElementToPNMElementMap_; }
287 pnmElementToFreeFlowElementsMap_.clear();
288 freeFlowElementToPNMElementMap_.clear();
289 isCoupledPNMDof_.clear();
290 isCoupledFreeFlowMomentumDof_.clear();
291 isCoupledFreeFlowMomentumDofOnInterface_.clear();
292 pnmToFreeFlowMassStencils_.clear();
293 pnmToFreeFlowMomentumStencils_.clear();
294 freeFlowMassToPNMStencils_.clear();
295 freeFlowMomentumToPNMStencils_.clear();
298 template<
class FreeFlowMomentumGr
idGeometry,
class PoreNetworkGr
idGeometry>
299 void resize_(
const FreeFlowMomentumGridGeometry& ffMomentumGridGeometry,
300 const PoreNetworkGridGeometry& pnmGridGeometry)
302 const auto numPNMDofs = pnmGridGeometry.numDofs();
303 const auto numFreeFlowMomentumDofs = ffMomentumGridGeometry.numDofs();
304 isCoupledPNMDof_.resize(numPNMDofs,
false);
305 isCoupledFreeFlowMomentumDof_.resize(numFreeFlowMomentumDofs,
false);
306 isCoupledFreeFlowMomentumDofOnInterface_.resize(numFreeFlowMomentumDofs,
false);
310 template<
class FVElementGeometry,
class FreeFlowMassGr
idGeometry,
class GlobalPosition,
class Scalar>
311 auto coupledFFMomentumDofs_(
const FVElementGeometry& fvGeometry,
312 const FreeFlowMassGridGeometry& ffMassGridGeometry,
313 const GlobalPosition& pnmPos,
314 const Scalar couplingPoreRadius,
315 const int couplingInterfaceDirectionIdx)
320 Dune::ReservedVector<std::size_t, FVElementGeometry::maxNumElementScvs> coupledLateralDofs;
321 std::size_t coupledFrontalDof;
326 for (
const auto& scv : scvs(fvGeometry))
328 const Scalar eps =
diameter(fvGeometry.geometry(scv))*1e-6;
329 assert(eps < couplingPoreRadius);
331 if (scv.dofAxis() == couplingInterfaceDirectionIdx)
333 if (abs(scv.dofPosition()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
335 result.coupledFrontalDof = scv.dofIndex();
338 for (
const auto& scvf : scvfs(fvGeometry, scv))
341 if (scvf.isLateral() && !scvf.boundary())
342 isCoupledLateralFreeFlowMomentumScvf_[scvf.index()] =
true;
343 else if (scvf.isFrontal() && scvf.boundary())
345 isCoupledFrontalFreeFlowMomentumScvf_[scvf.index()] =
true;
347 const auto&
element = ffMassGridGeometry.element(fvGeometry.elementIndex());
348 const auto ffMassFVGeometry =
localView(ffMassGridGeometry).bindElement(element);
349 for (
const auto& ffMassScvf : scvfs(ffMassFVGeometry))
351 if (abs(ffMassScvf.center()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
352 isCoupledFreeFlowMassScvf_[ffMassScvf.index()] =
true;
360 bool isCoupledDof =
false;
362 for (
int dimIdx = 0; dimIdx < GlobalPosition::dimension; ++dimIdx)
364 if (dimIdx == couplingInterfaceDirectionIdx || scv.boundary())
367 isCoupledDof = abs(scv.dofPosition()[dimIdx] - pnmPos[dimIdx]) < couplingPoreRadius + eps;
372 result.coupledLateralDofs.push_back(scv.dofIndex());
375 for (
const auto& scvf : scvfs(fvGeometry, scv))
377 if (scvf.isLateral() && scvf.boundary())
380 if (abs(scvf.ipGlobal()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
381 isCoupledLateralFreeFlowMomentumScvf_[scvf.index()] =
true;
391 std::vector<std::size_t> emptyStencil_;
393 std::unordered_map<std::size_t, std::vector<std::size_t>> pnmElementToFreeFlowElementsMap_;
394 std::unordered_map<std::size_t, std::size_t> freeFlowElementToPNMElementMap_;
396 std::vector<bool> isCoupledPNMDof_;
397 std::vector<bool> isCoupledFreeFlowMomentumDof_;
398 std::vector<bool> isCoupledFreeFlowMomentumDofOnInterface_;
400 std::unordered_map<std::size_t, bool> isCoupledLateralFreeFlowMomentumScvf_;
401 std::unordered_map<std::size_t, bool> isCoupledFrontalFreeFlowMomentumScvf_;
402 std::unordered_map<std::size_t, bool> isCoupledFreeFlowMassScvf_;
404 MapType pnmToFreeFlowMassStencils_;
405 MapType pnmToFreeFlowMomentumStencils_;
406 MapType freeFlowMassToPNMStencils_;
407 MapType freeFlowMomentumToPNMStencils_;
Algorithms that finds which geometric entities intersect.
A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
The available discretization methods in Dumux.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:40
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:38
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:114
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:38
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
static std::size_t normalAxis(const Vector &v)
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: normalaxis.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Coupling mapper for staggered free-flow and pore-network models.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:54
bool isCoupledFreeFlowMomentumLateralScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:267
const std::vector< std::size_t > & freeFlowMassToPoreNetworkCouplingStencil(const std::size_t eIdxI) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:209
void update(const FreeFlowMomentumGridGeometry &ffMomentumGridGeometry, const FreeFlowMassGridGeometry &ffMassGridGeometry, const PoreNetworkGridGeometry &pnmGridGeometry)
Main update routine.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:62
const std::vector< std::size_t > & freeFlowMomentumToPoreNetworkCouplingStencil(const std::size_t dofIndex) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:222
bool isCoupledFreeFlowMomentumScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:262
bool isCoupledPoreNetworkDof(std::size_t dofIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:257
bool isCoupledFreeFlowMassScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:272
const std::vector< std::size_t > & poreNetworkToFreeFlowMomentumCouplingStencil(const std::size_t eIdxI) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:183
const auto & freeFlowElementToPNMElementMap() const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:280
const auto & pnmElementToFreeFlowElementsMap() const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:277
bool isCoupledFreeFlowMomentumDof(std::size_t dofIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:241
const std::vector< std::size_t > & poreNetworkToFreeFlowMassCouplingStencil(const std::size_t eIdxI) const
returns an iterable container of all indices of degrees of freedom of domain j that couple with / inf...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:196
bool isCoupledFreeFlowElement(std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:233
bool isCoupledPoreNetworkElement(std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:249