27#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
28#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
34#include <dune/common/timer.hh>
35#include <dune/common/exceptions.hh>
36#include <dune/common/fvector.hh>
37#include <dune/geometry/quadraturerules.hh>
64template <
class Gr
idGeometry>
67 using GridView =
typename GridGeometry::GridView;
68 using GlobalPosition =
typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
79 template <
class RadiusFunction>
82 const auto& gridView = gg.gridView();
83 segments_.resize(gridView.size(0));
84 for (
const auto &element : elements(gridView))
86 const auto geometry = element.geometry();
87 const auto eIdx = gg.elementMapper().index(element);
88 const auto radius = radiusFunction(eIdx);
90 segments_[eIdx] = Segment{ geometry.corner(0), geometry.corner(1),
radius, eIdx };
93 std::cout <<
"-- Extracted " << segments_.size() <<
" segments.\n";
98 {
return segments_[eIdx].r; }
102 {
return segments_[eIdx]; }
106 {
return segments_; }
113 double minSignedDist = std::numeric_limits<double>::max();
115 for (
const auto &
segment : segments_)
117 const auto signedDist = computeSignedDistance_(globalPos,
segment);
118 if (signedDist < minSignedDist)
120 minSignedDist = signedDist;
126 return std::make_tuple(eIdx, abs(minSignedDist));
130 template <
class IndexRange>
132 const IndexRange& segIndices)
const
135 double minSignedDist = std::numeric_limits<double>::max();
137 for (
const auto index : segIndices)
139 const auto &
segment = segments_[index];
140 const auto signedDist = computeSignedDistance_(globalPos,
segment);
141 if (signedDist < minSignedDist)
143 minSignedDist = signedDist;
149 return std::make_tuple(eIdx, abs(minSignedDist));
155 const auto&
segment = segments_[idx];
163 GlobalPosition projectionPointOnSegment_(
const GlobalPosition& p,
const GlobalPosition& a,
const GlobalPosition& b)
const
165 const auto v = b - a;
166 const auto w = p - a;
168 const auto proj1 = v*w;
172 const auto proj2 = v.two_norm2();
176 const auto t = proj1 / proj2;
183 auto computeSignedDistance_(
const GlobalPosition& p,
const Segment&
segment)
const
185 const auto projPoint = projectionPointOnSegment_(p,
segment.a,
segment.b);
186 return (p-projPoint).two_norm()-
segment.r;
189 std::vector<Segment> segments_;
198template<
class Network>
209 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos);
210 return closestSegmentIdx;
213 template <
class Pos,
class H
intContainer>
214 auto operator()(
const Pos& pos,
const HintContainer& hints)
const
216 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos, hints);
217 return closestSegmentIdx;
221 const Network& network_;
231 void push(
const std::array<Dune::FieldVector<double, 3>, 3>& corners,
double data)
233 vtkCells_.emplace_back(std::array<std::size_t, 3>{
234 { vtkPoints_.size(), vtkPoints_.size()+1, vtkPoints_.size()+2 }
236 for (
const auto& c : corners)
237 vtkPoints_.push_back(c);
238 vtkCellData_.push_back(data);
243 std::ofstream intersectionFile(name);
244 intersectionFile <<
"# vtk DataFile Version 2.0\n";
245 intersectionFile <<
"Really cool intersection data\n";
246 intersectionFile <<
"ASCII\n";
247 intersectionFile <<
"DATASET UNSTRUCTURED_GRID\n";
248 intersectionFile <<
"POINTS " << vtkPoints_.size() <<
" float\n";
249 for (
const auto& p : vtkPoints_)
250 intersectionFile << p[0] <<
" " << p[1] <<
" " << p[2] <<
"\n";
251 intersectionFile <<
"\nCELLS " << vtkCells_.size() <<
" " << vtkCells_.size()*4 <<
"\n";
252 for (
const auto& c : vtkCells_)
253 intersectionFile <<
"3 " << c[0] <<
" " << c[1] <<
" " << c[2] <<
"\n";
254 intersectionFile <<
"\nCELL_TYPES " << vtkCells_.size() <<
"\n";
255 for (
int i = 0; i < vtkCells_.size(); ++i)
256 intersectionFile <<
"5\n";
257 intersectionFile <<
"\nCELL_DATA " << vtkCells_.size() <<
"\nSCALARS index float\nLOOKUP_TABLE default\n";
258 for (
const auto& c : vtkCellData_)
259 intersectionFile << c <<
"\n";
263 std::vector<Dune::FieldVector<double, 3>> vtkPoints_;
264 std::vector<std::array<std::size_t, 3>> vtkCells_;
265 std::vector<double> vtkCellData_;
270namespace Embedded1d3dCouplingMode {
272 static std::string
name() {
return "projection"; }
279template<
class MDTraits,
class CouplingMode>
280class Embedded1d3dCouplingManager;
330template<
class MDTraits>
337 using Scalar =
typename MDTraits::Scalar;
338 using SolutionVector =
typename MDTraits::SolutionVector;
339 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
341 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
342 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
345 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
348 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
349 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
352 static constexpr int lowDimDim = GridView<lowDimIdx>::dimension;
353 static constexpr int bulkDim = GridView<bulkIdx>::dimension;
354 static constexpr int dimWorld = GridView<bulkIdx>::dimensionworld;
356 template<std::
size_t id>
357 static constexpr bool isBox()
360 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
367 using ParentType::ParentType;
369 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
370 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
371 const SolutionVector& curSol)
373 ParentType::init(bulkProblem, lowDimProblem, curSol);
389 std::cout <<
"Initializing the coupling manager (projection)" << std::endl;
392 const bool enableIntersectionOutput = getParam<bool>(
"MixedDimension.Projection.EnableIntersectionOutput",
false);
393 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkCoupledFaces =
394 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() :
nullptr;
395 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkIntersections =
396 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() :
nullptr;
400 const auto& lowDimProblem = this->problem(lowDimIdx);
401 const auto& lowDimFvGridGeometry = lowDimProblem.gridGeometry();
402 const auto& lowDimGridView = lowDimFvGridGeometry.gridView();
403 localAreaFactor_.resize(lowDimGridView.size(0), 0.0);
404 const auto radiusFunc = [&](
const GridIndex<lowDimIdx> eIdx) {
405 return lowDimProblem.spatialParams().radius(eIdx);
408 lowDimFvGridGeometry, radiusFunc
412 this->precomputeVertexIndices(bulkIdx);
413 this->precomputeVertexIndices(lowDimIdx);
415 const auto& bulkProblem = this->problem(bulkIdx);
416 const auto& bulkFvGridGeometry = bulkProblem.gridGeometry();
417 const auto& bulkGridView = bulkFvGridGeometry.gridView();
418 bulkElementMarker_.assign(bulkGridView.size(0), 0);
419 bulkVertexMarker_.assign(bulkGridView.size(bulkDim), 0);
423 static const auto bBoxShrinking
424 = getParam<Scalar>(
"MixedDimension.Projection.CoupledBoundingBoxShrinkingFactor", 1e-2);
425 const GlobalPosition threshold(
426 bBoxShrinking*( bulkFvGridGeometry.bBoxMax()-bulkFvGridGeometry.bBoxMin() ).two_norm()
428 const auto bBoxMinSmall = bulkFvGridGeometry.bBoxMin() + threshold;
429 const auto bBoxMaxSmall = bulkFvGridGeometry.bBoxMax() - threshold;
430 auto insideBBox = [bBoxMin=bBoxMinSmall, bBoxMax=bBoxMaxSmall](
const GlobalPosition& point) ->
bool
432 static constexpr Scalar eps_ = 1.0e-7;
433 const Scalar eps0 = eps_*(bBoxMax[0] - bBoxMin[0]);
434 const Scalar eps1 = eps_*(bBoxMax[1] - bBoxMin[1]);
435 const Scalar eps2 = eps_*(bBoxMax[2] - bBoxMin[2]);
436 return (bBoxMin[0] - eps0 <= point[0] && point[0] <= bBoxMax[0] + eps0 &&
437 bBoxMin[1] - eps1 <= point[1] && point[1] <= bBoxMax[1] + eps1 &&
438 bBoxMin[2] - eps2 <= point[2] && point[2] <= bBoxMax[2] + eps2);
442 static const auto quadSimplexRefineMaxLevel
443 = getParam<std::size_t>(
"MixedDimension.Projection.SimplexIntegrationRefine", 4);
447 static const auto estimateNumberOfPointSources
448 = getParam<std::size_t>(
"MixedDimension.Projection.EstimateNumberOfPointSources", bulkFvGridGeometry.gridView().size(0));
450 this->pointSources(bulkIdx).reserve(estimateNumberOfPointSources);
451 this->pointSources(lowDimIdx).reserve(estimateNumberOfPointSources);
452 this->pointSourceData().reserve(estimateNumberOfPointSources);
457 for (
const auto& element : elements(bulkGridView))
459 const auto bulkElementIdx = bulkFvGridGeometry.elementMapper().index(element);
460 const auto bulkGeometry = element.geometry();
461 const auto refElement = Dune::referenceElement(element);
462 for (
const auto& intersection : intersections(bulkGridView, element))
465 if (!intersection.boundary())
469 const auto [isCoupled, closestSegmentIdx, minDist] = isCoupled_(network, intersection, insideBBox);
474 bulkElementMarker_[bulkElementIdx] = 1;
476 const auto interfaceGeometry = intersection.geometry();
477 for (
int i = 0; i < interfaceGeometry.corners(); ++i)
479 const auto localVIdx = refElement.subEntity(intersection.indexInInside(), 1, i, bulkDim);
480 const auto vIdx = bulkFvGridGeometry.vertexMapper().subIndex(element, localVIdx, bulkDim);
481 bulkVertexMarker_[vIdx] = 1;
485 if (enableIntersectionOutput)
486 vtkCoupledFaces->push({
487 { interfaceGeometry.corner(0), interfaceGeometry.corner(1), interfaceGeometry.corner(2)}
488 }, closestSegmentIdx);
493 for (
const auto& qp : quadSimplex)
495 const auto surfacePoint = interfaceGeometry.global(qp.position());
496 const auto closestSegmentIdx = qp.indicator();
497 const auto projPoint = network.projectionPointOnSegment(surfacePoint, closestSegmentIdx);
499 addPointSourceAtIP_(bulkGeometry, interfaceGeometry, qp, bulkElementIdx, closestSegmentIdx, surfacePoint, projPoint);
500 addStencilEntries_(bulkElementIdx, closestSegmentIdx);
506 using namespace Dune::Hybrid;
507 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
509 for (
auto&& stencil : this->couplingStencils(domainIdx))
511 std::sort(stencil.second.begin(), stencil.second.end());
512 stencil.second.erase(
513 std::unique(stencil.second.begin(), stencil.second.end()),
520 if (enableIntersectionOutput)
521 vtkCoupledFaces->write(
"coupledfaces.vtk");
525 for (
const auto& element : elements(lowDimGridView))
527 const auto length = element.geometry().volume();
528 const auto eIdx = lowDimFvGridGeometry.elementMapper().index(element);
529 const auto radius = this->problem(lowDimIdx).spatialParams().radius(eIdx);
530 const auto cylinderSurface = 2.0*M_PI*radius*length;
531 localAreaFactor_[eIdx] /= cylinderSurface;
532 localAreaFactor_[eIdx] -= 1.0;
535 this->pointSources(bulkIdx).shrink_to_fit();
536 this->pointSources(lowDimIdx).shrink_to_fit();
537 this->pointSourceData().shrink_to_fit();
539 std::cout <<
"-- Coupling at " << this->pointSourceData().size() <<
" integration points." << std::endl;
540 std::cout <<
"-- Finished initialization after " << watch.elapsed() <<
" seconds." << std::endl;
551 const auto& data = this->pointSourceData()[id];
552 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
559 {
return bulkElementMarker_; }
563 {
return bulkVertexMarker_; }
569 {
return localAreaFactor_; }
572 std::vector<int> bulkElementMarker_, bulkVertexMarker_;
573 std::vector<double> localAreaFactor_;
577 template<
class BulkGeometry,
class InterfaceGeometry,
class QP>
578 void addPointSourceAtIP_(
const BulkGeometry& bulkGeometry,
579 const InterfaceGeometry& ifGeometry,
581 const GridIndex<bulkIdx> bulkElementIdx,
582 const GridIndex<lowDimIdx> lowDimElementIdx,
583 const GlobalPosition& surfacePoint,
584 const GlobalPosition& projPoint)
587 const auto id = this->idCounter_++;
590 const auto qpweight = qp.weight();
591 const auto integrationElement = ifGeometry.integrationElement(qp.position());
592 this->pointSources(bulkIdx).emplace_back(surfacePoint,
id, qpweight, integrationElement, bulkElementIdx);
595 this->pointSources(lowDimIdx).emplace_back(projPoint,
id, qpweight, integrationElement, lowDimElementIdx);
596 localAreaFactor_[lowDimElementIdx] += qpweight*integrationElement;
598 const auto& bulkFvGridGeometry = this->problem(bulkIdx).gridGeometry();
599 const auto& lowDimFvGridGeometry = this->problem(lowDimIdx).gridGeometry();
605 if constexpr (isBox<lowDimIdx>())
607 std::vector<Dune::FieldVector<Scalar, 1> > shapeValues;
608 const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry();
609 this->getShapeValues(lowDimIdx, lowDimFvGridGeometry, lowDimGeometry, projPoint, shapeValues);
610 psData.
addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
617 if constexpr (isBox<bulkIdx>())
619 std::vector<Dune::FieldVector<Scalar, 1> > shapeValues;
620 this->getShapeValues(bulkIdx, bulkFvGridGeometry, bulkGeometry, surfacePoint, shapeValues);
621 psData.
addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
629 this->pointSourceData().emplace_back(std::move(psData));
632 void addStencilEntries_(
const GridIndex<bulkIdx> bulkElementIdx,
const GridIndex<lowDimIdx> lowDimElementIdx)
635 if constexpr (isBox<lowDimIdx>())
637 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
638 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(
639 this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
640 vertices.begin(), vertices.end()
646 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
650 if constexpr (isBox<bulkIdx>())
652 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
653 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(
654 this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
655 vertices.begin(), vertices.end()
660 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
664 template<
class BoundaryIntersection,
class BBoxCheck>
665 auto isCoupled_(
const SegmentNetwork& network,
const BoundaryIntersection& intersection,
const BBoxCheck& insideBBox)
const
667 const auto center = intersection.geometry().center();
668 const auto [eIdx, minDist] = network.findClosestSegmentSurface(
center);
672 static const bool considerBBox = getParam<bool>(
"MixedDimension.Projection.ConsiderFacesWithinBoundingBoxCoupled",
false);
673 if (considerBBox && insideBBox(
center))
674 return std::make_tuple(
true, eIdx, minDist);
678 static const auto rFactor = getParam<Scalar>(
"MixedDimension.Projection.CoupledRadiusFactor", 0.1);
679 static const auto spFactor = getParam<Scalar>(
"MixedDimension.Projection.CoupledAngleFactor", 0.3);
680 const auto projPoint = network.projectionPointOnSegment(
center, eIdx);
681 const auto distVec = projPoint -
center;
682 const bool isCoupled = minDist < rFactor * network.radius(eIdx) && distVec * intersection.centerUnitOuterNormal() > distVec.two_norm() * spFactor;
683 return std::make_tuple(isCoupled, eIdx, minDist);
688template<
class MDTraits>
690:
public std::true_type {};
Helper class to create (named and comparable) tagged types.
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
A quadrature rule using local refinement to approximate partitioned elements.
A quadrature based on refinement.
A function to compute the convex hull of a point cloud.
Functionality to triangulate point clouds.
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:36
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
constexpr Projection projection
Definition: couplingmanager1d3d_projection.hh:275
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
Definition: common/properties.hh:100
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition: tag.hh:42
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:36
Segment representation of a 1d network grid.
Definition: couplingmanager1d3d_projection.hh:66
auto findClosestSegmentSurface(const GlobalPosition &globalPos, const IndexRange &segIndices) const
same as overload with taking a position but only look at the segments specified by the index vector
Definition: couplingmanager1d3d_projection.hh:131
const auto & segment(std::size_t eIdx) const
the segment with index eIdx
Definition: couplingmanager1d3d_projection.hh:101
auto radius(std::size_t eIdx) const
the segment radius
Definition: couplingmanager1d3d_projection.hh:97
GlobalPosition projectionPointOnSegment(const GlobalPosition &p, std::size_t idx) const
Definition: couplingmanager1d3d_projection.hh:153
const auto & segments() const
the segments
Definition: couplingmanager1d3d_projection.hh:105
SegmentNetwork(const GridGeometry &gg, const RadiusFunction &radiusFunction)
Definition: couplingmanager1d3d_projection.hh:80
auto findClosestSegmentSurface(const GlobalPosition &globalPos) const
Definition: couplingmanager1d3d_projection.hh:110
Get the closest segment for a given surface point.
Definition: couplingmanager1d3d_projection.hh:200
NetworkIndicatorFunction(const Network &n)
Definition: couplingmanager1d3d_projection.hh:202
auto operator()(const Pos &pos) const
Definition: couplingmanager1d3d_projection.hh:207
auto operator()(const Pos &pos, const HintContainer &hints) const
Definition: couplingmanager1d3d_projection.hh:214
Simple legacy VTK writer for outputting debug data on the coupling interface.
Definition: couplingmanager1d3d_projection.hh:229
void push(const std::array< Dune::FieldVector< double, 3 >, 3 > &corners, double data)
Definition: couplingmanager1d3d_projection.hh:231
void write(const std::string &name)
Definition: couplingmanager1d3d_projection.hh:241
Definition: couplingmanager1d3d_projection.hh:271
static std::string name()
Definition: couplingmanager1d3d_projection.hh:272
Manages the coupling between bulk elements and lower dimensional elements.
Definition: couplingmanager1d3d_projection.hh:333
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_projection.hh:384
const std::vector< int > & bulkElementMarker() const
Marker that is non-zero for bulk elements that are coupled.
Definition: couplingmanager1d3d_projection.hh:558
const std::vector< int > & bulkVertexMarker() const
Marker that is non-zero for bulk vertices that belong to coupled surface facets.
Definition: couplingmanager1d3d_projection.hh:562
const std::vector< Scalar > & localAreaFactor() const
Definition: couplingmanager1d3d_projection.hh:568
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_projection.hh:549
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_projection.hh:369
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:83
A quadrature rule using local refinement to approximate partitioned elements.
Definition: localrefinementquadrature.hh:53
A point source data class used for integration in multidimensional models.
Definition: pointsourcedata.hh:43
void addLowDimInterpolation(const ShapeValues &shapeValues, const std::vector< GridIndex< lowDimIdx > > &cornerIndices, GridIndex< lowDimIdx > eIdx)
Definition: pointsourcedata.hh:72
void addBulkInterpolation(const ShapeValues &shapeValues, const std::vector< GridIndex< bulkIdx > > &cornerIndices, GridIndex< bulkIdx > eIdx)
Definition: pointsourcedata.hh:62
trait that is specialized for coupling manager supporting multithreaded assembly
Definition: multidomain/fvassembler.hh:85
Declares all properties used in Dumux.