3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingmanager1d3d_projection.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
27#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
28#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
29
30#include <vector>
31#include <algorithm>
32#include <string>
33
34#include <dune/common/timer.hh>
35#include <dune/common/exceptions.hh>
36#include <dune/common/fvector.hh>
37#include <dune/geometry/quadraturerules.hh>
38
39#include <dumux/common/tag.hh>
41#include <dumux/common/math.hh>
43
47
50
51namespace Dumux {
52
53// some implementation details of the projection coupling manager
54namespace Detail {
55
64template <class GridGeometry>
66{
67 using GridView = typename GridGeometry::GridView;
68 using GlobalPosition = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
69 using GridIndex = typename IndexTraits<GridView>::GridIndex;
70
71 struct Segment
72 {
73 GlobalPosition a, b;
74 double r;
75 GridIndex eIdx;
76 };
77
78public:
79 template <class RadiusFunction>
80 SegmentNetwork(const GridGeometry& gg, const RadiusFunction& radiusFunction)
81 {
82 const auto& gridView = gg.gridView();
83 segments_.resize(gridView.size(0));
84 for (const auto &element : elements(gridView))
85 {
86 const auto geometry = element.geometry();
87 const auto eIdx = gg.elementMapper().index(element);
88 const auto radius = radiusFunction(eIdx);
89
90 segments_[eIdx] = Segment{ geometry.corner(0), geometry.corner(1), radius, eIdx };
91 }
92
93 std::cout << "-- Extracted " << segments_.size() << " segments.\n";
94 }
95
97 auto radius(std::size_t eIdx) const
98 { return segments_[eIdx].r; }
99
101 const auto& segment(std::size_t eIdx) const
102 { return segments_[eIdx]; }
103
105 const auto& segments() const
106 { return segments_; }
107
110 auto findClosestSegmentSurface(const GlobalPosition& globalPos) const
111 {
112 GridIndex eIdx = 0;
113 double minSignedDist = std::numeric_limits<double>::max();
114
115 for (const auto &segment : segments_)
116 {
117 const auto signedDist = computeSignedDistance_(globalPos, segment);
118 if (signedDist < minSignedDist)
119 {
120 minSignedDist = signedDist;
121 eIdx = segment.eIdx;
122 }
123 }
124
125 using std::abs;
126 return std::make_tuple(eIdx, abs(minSignedDist));
127 }
128
130 template <class IndexRange>
131 auto findClosestSegmentSurface(const GlobalPosition& globalPos,
132 const IndexRange& segIndices) const
133 {
134 GridIndex eIdx = 0;
135 double minSignedDist = std::numeric_limits<double>::max();
136
137 for (const auto index : segIndices)
138 {
139 const auto &segment = segments_[index];
140 const auto signedDist = computeSignedDistance_(globalPos, segment);
141 if (signedDist < minSignedDist)
142 {
143 minSignedDist = signedDist;
144 eIdx = segment.eIdx;
145 }
146 }
147
148 using std::abs;
149 return std::make_tuple(eIdx, abs(minSignedDist));
150 }
151
152 // Compute the projection point of p onto the segment connecting a->b
153 GlobalPosition projectionPointOnSegment(const GlobalPosition& p, std::size_t idx) const
154 {
155 const auto& segment = segments_[idx];
156 return projectionPointOnSegment_(p, segment.a, segment.b);
157 }
158
159private:
160 // Compute the projection of p onto the segment
161 // returns the closest end point if the orthogonal projection onto the line implied by the segment
162 // is outside the segment
163 GlobalPosition projectionPointOnSegment_(const GlobalPosition& p, const GlobalPosition& a, const GlobalPosition& b) const
164 {
165 const auto v = b - a;
166 const auto w = p - a;
167
168 const auto proj1 = v*w;
169 if (proj1 <= 0.0)
170 return a;
171
172 const auto proj2 = v.two_norm2();
173 if (proj2 <= proj1)
174 return b;
175
176 const auto t = proj1 / proj2;
177 auto x = a;
178 x.axpy(t, v);
179 return x;
180 }
181
182 // Compute the signed (!) distance to a cylinder segment surface and the projection point onto the centerline
183 auto computeSignedDistance_(const GlobalPosition& p, const Segment& segment) const
184 {
185 const auto projPoint = projectionPointOnSegment_(p, segment.a, segment.b);
186 return (p-projPoint).two_norm()-segment.r;
187 }
188
189 std::vector<Segment> segments_;
190};
191
198template<class Network>
200{
201public:
202 NetworkIndicatorFunction(const Network& n)
203 : network_(n)
204 {}
205
206 template<class Pos>
207 auto operator()(const Pos& pos) const
208 {
209 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos);
210 return closestSegmentIdx;
211 }
212
213 template <class Pos, class HintContainer>
214 auto operator()(const Pos& pos, const HintContainer& hints) const
215 {
216 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos, hints);
217 return closestSegmentIdx;
218 }
219
220private:
221 const Network& network_;
222};
223
229{
230public:
231 void push(const std::array<Dune::FieldVector<double, 3>, 3>& corners, double data)
232 {
233 vtkCells_.emplace_back(std::array<std::size_t, 3>{
234 { vtkPoints_.size(), vtkPoints_.size()+1, vtkPoints_.size()+2 }
235 });
236 for (const auto& c : corners)
237 vtkPoints_.push_back(c);
238 vtkCellData_.push_back(data);
239 }
240
241 void write(const std::string& name)
242 {
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";
260 }
261
262private:
263 std::vector<Dune::FieldVector<double, 3>> vtkPoints_;
264 std::vector<std::array<std::size_t, 3>> vtkCells_;
265 std::vector<double> vtkCellData_;
266};
267
268} // end namespace Detail
269
270namespace Embedded1d3dCouplingMode {
271struct Projection : public Utility::Tag<Projection> {
272 static std::string name() { return "projection"; }
273};
274
275inline constexpr Projection projection{};
276} // end namespace Embedded1d3dCouplingMode
277
278// forward declaration
279template<class MDTraits, class CouplingMode>
280class Embedded1d3dCouplingManager;
281
330template<class MDTraits>
331class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>
332: public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>>
333{
336
337 using Scalar = typename MDTraits::Scalar;
338 using SolutionVector = typename MDTraits::SolutionVector;
339 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
340
341 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
342 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
343
344 // the sub domain type aliases
345 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
346 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
347 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
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;
350 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
351
352 static constexpr int lowDimDim = GridView<lowDimIdx>::dimension;
353 static constexpr int bulkDim = GridView<bulkIdx>::dimension;
354 static constexpr int dimWorld = GridView<bulkIdx>::dimensionworld;
355
356 template<std::size_t id>
357 static constexpr bool isBox()
358 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
359
360 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
361
363
364public:
365 static constexpr Embedded1d3dCouplingMode::Projection couplingMode{};
366
367 using ParentType::ParentType;
368
369 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
370 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
371 const SolutionVector& curSol)
372 {
373 ParentType::init(bulkProblem, lowDimProblem, curSol);
374 }
375
376 /* \brief Compute integration point point sources and associated data
377 * This method computes a source term on the explicitly given surface of the network grid
378 * by integrating over the elements of the surface, projecting the integration points onto
379 * the centerlines to evaluate 1d quantities and set point sources (minimum distance projection)
380 *
381 * \param order Unused in this algorithm! (passed from default 1D3D coupling manager)
382 * \param verbose If the source computation is verbose
383 */
384 void computePointSourceData(std::size_t order = 1, bool verbose = false)
385 {
386 // initialize the maps
387 // do some logging and profiling
388 Dune::Timer watch;
389 std::cout << "Initializing the coupling manager (projection)" << std::endl;
390
391 // debug VTK output
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;
397
398 // temporarily represent the network domain as a list of segments
399 // a more efficient data structure
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);
406 };
407 const auto network = SegmentNetwork{
408 lowDimFvGridGeometry, radiusFunc
409 };
410
411 // precompute the vertex indices for efficiency (box method)
412 this->precomputeVertexIndices(bulkIdx);
413 this->precomputeVertexIndices(lowDimIdx);
414
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);
420
421 // construct a bounding box that may be used to determine coupled faces
422 // the shrinking factor shrinks this bounding box
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()
427 );
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
431 {
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);
439 };
440
441 // setup simplex "quadrature" rule
442 static const auto quadSimplexRefineMaxLevel
443 = getParam<std::size_t>("MixedDimension.Projection.SimplexIntegrationRefine", 4);
444 const auto indicator = Detail::NetworkIndicatorFunction { network };
445
446 // estimate number of point sources for memory allocation
447 static const auto estimateNumberOfPointSources
448 = getParam<std::size_t>("MixedDimension.Projection.EstimateNumberOfPointSources", bulkFvGridGeometry.gridView().size(0));
449
450 this->pointSources(bulkIdx).reserve(estimateNumberOfPointSources);
451 this->pointSources(lowDimIdx).reserve(estimateNumberOfPointSources);
452 this->pointSourceData().reserve(estimateNumberOfPointSources);
453
454 // we go over the bulk surface intersection, check if we are coupled,
455 // and if this is the case, we add integration point sources using
456 // the local simplex refinement quadrature procedure (see paper Koch 2021)
457 for (const auto& element : elements(bulkGridView))
458 {
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))
463 {
464 // skip inner intersection
465 if (!intersection.boundary())
466 continue;
467
468 // check if we are close enough to the coupling interface
469 const auto [isCoupled, closestSegmentIdx, minDist] = isCoupled_(network, intersection, insideBBox);
470 if (!isCoupled)
471 continue;
472
473 // we found an intersection on the coupling interface: mark as coupled element
474 bulkElementMarker_[bulkElementIdx] = 1;
475
476 const auto interfaceGeometry = intersection.geometry();
477 for (int i = 0; i < interfaceGeometry.corners(); ++i)
478 {
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;
482 }
483
484 // debug graphical intersection output
485 if (enableIntersectionOutput)
486 vtkCoupledFaces->push({
487 { interfaceGeometry.corner(0), interfaceGeometry.corner(1), interfaceGeometry.corner(2)}
488 }, closestSegmentIdx);
489
490 const auto quadSimplex = LocalRefinementSimplexQuadrature{ interfaceGeometry, quadSimplexRefineMaxLevel, indicator };
491
492 // iterate over all quadrature points (children simplices)
493 for (const auto& qp : quadSimplex)
494 {
495 const auto surfacePoint = interfaceGeometry.global(qp.position());
496 const auto closestSegmentIdx = qp.indicator();
497 const auto projPoint = network.projectionPointOnSegment(surfacePoint, closestSegmentIdx);
498
499 addPointSourceAtIP_(bulkGeometry, interfaceGeometry, qp, bulkElementIdx, closestSegmentIdx, surfacePoint, projPoint);
500 addStencilEntries_(bulkElementIdx, closestSegmentIdx);
501 }
502 }
503 }
504
505 // make stencils unique
506 using namespace Dune::Hybrid;
507 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
508 {
509 for (auto&& stencil : this->couplingStencils(domainIdx))
510 {
511 std::sort(stencil.second.begin(), stencil.second.end());
512 stencil.second.erase(
513 std::unique(stencil.second.begin(), stencil.second.end()),
514 stencil.second.end()
515 );
516 }
517 });
518
519 // debug VTK output
520 if (enableIntersectionOutput)
521 vtkCoupledFaces->write("coupledfaces.vtk");
522
523 // compare theoretical cylinder surface to actual discrete surface
524 // of the coupled bulk interface facets
525 for (const auto& element : elements(lowDimGridView))
526 {
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;
533 }
534
535 this->pointSources(bulkIdx).shrink_to_fit();
536 this->pointSources(lowDimIdx).shrink_to_fit();
537 this->pointSourceData().shrink_to_fit();
538
539 std::cout << "-- Coupling at " << this->pointSourceData().size() << " integration points." << std::endl;
540 std::cout << "-- Finished initialization after " << watch.elapsed() << " seconds." << std::endl;
541 }
542
546 // \{
547
549 Scalar radius(std::size_t id) const
550 {
551 const auto& data = this->pointSourceData()[id];
552 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
553 }
554
555 // \}
556
558 const std::vector<int>& bulkElementMarker() const
559 { return bulkElementMarker_; }
560
562 const std::vector<int>& bulkVertexMarker() const
563 { return bulkVertexMarker_; }
564
568 const std::vector<Scalar>& localAreaFactor() const
569 { return localAreaFactor_; }
570
571private:
572 std::vector<int> bulkElementMarker_, bulkVertexMarker_;
573 std::vector<double> localAreaFactor_;
574
575 // add a point source containing all data needed to evaluate
576 // the coupling source term from both domains
577 template<class BulkGeometry, class InterfaceGeometry, class QP>
578 void addPointSourceAtIP_(const BulkGeometry& bulkGeometry,
579 const InterfaceGeometry& ifGeometry,
580 const QP& qp,
581 const GridIndex<bulkIdx> bulkElementIdx,
582 const GridIndex<lowDimIdx> lowDimElementIdx,
583 const GlobalPosition& surfacePoint,
584 const GlobalPosition& projPoint)
585 {
586 // add a new point source id (for the associated data)
587 const auto id = this->idCounter_++;
588
589 // add a bulk point source
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);
593
594 // find the corresponding projection point and 1d element
595 this->pointSources(lowDimIdx).emplace_back(projPoint, id, qpweight, integrationElement, lowDimElementIdx);
596 localAreaFactor_[lowDimElementIdx] += qpweight*integrationElement;
597
598 const auto& bulkFvGridGeometry = this->problem(bulkIdx).gridGeometry();
599 const auto& lowDimFvGridGeometry = this->problem(lowDimIdx).gridGeometry();
600
601 // pre compute additional data used for the evaluation of
602 // the actual solution dependent source term
603 PointSourceData psData;
604
605 if constexpr (isBox<lowDimIdx>())
606 {
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);
611 }
612 else
613 {
614 psData.addLowDimInterpolation(lowDimElementIdx);
615 }
616
617 if constexpr (isBox<bulkIdx>())
618 {
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);
622 }
623 else
624 {
625 psData.addBulkInterpolation(bulkElementIdx);
626 }
627
628 // publish point source data in the global vector
629 this->pointSourceData().emplace_back(std::move(psData));
630 }
631
632 void addStencilEntries_(const GridIndex<bulkIdx> bulkElementIdx, const GridIndex<lowDimIdx> lowDimElementIdx)
633 {
634 // export the bulk coupling stencil
635 if constexpr (isBox<lowDimIdx>())
636 {
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()
641 );
642
643 }
644 else
645 {
646 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
647 }
648
649 // export the lowdim coupling stencil
650 if constexpr (isBox<bulkIdx>())
651 {
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()
656 );
657 }
658 else
659 {
660 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
661 }
662 }
663
664 template<class BoundaryIntersection, class BBoxCheck>
665 auto isCoupled_(const SegmentNetwork& network, const BoundaryIntersection& intersection, const BBoxCheck& insideBBox) const
666 {
667 const auto center = intersection.geometry().center();
668 const auto [eIdx, minDist] = network.findClosestSegmentSurface(center);
669
670 // all elements in the bounding box are coupled if this is enabled
671 // this is a simplified check that can be enabled for box-shaped domains
672 static const bool considerBBox = getParam<bool>("MixedDimension.Projection.ConsiderFacesWithinBoundingBoxCoupled", false);
673 if (considerBBox && insideBBox(center))
674 return std::make_tuple(true, eIdx, minDist);
675
676 // general coupling face detection algorithm:
677 // if the distance is under a certain threshold we are coupled (and the angle is not very big)
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);
684 }
685};
686
688template<class MDTraits>
689struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>>
690: public std::true_type {};
691
692} // end namespace Dumux
693
694#endif
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.