version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
15#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
16#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_PROJECTION_HH
17
18#include <vector>
19#include <algorithm>
20#include <string>
21
22#include <dune/common/timer.hh>
23#include <dune/common/exceptions.hh>
24#include <dune/common/fvector.hh>
25#include <dune/geometry/quadraturerules.hh>
26
27#include <dumux/common/tag.hh>
29#include <dumux/common/math.hh>
31
35
38
39namespace Dumux {
40
41// some implementation details of the projection coupling manager
42namespace Detail {
43
52template <class GridGeometry>
54{
55 using GridView = typename GridGeometry::GridView;
56 using GlobalPosition = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
57 using GridIndex = typename IndexTraits<GridView>::GridIndex;
58
59 struct Segment
60 {
61 GlobalPosition a, b;
62 double r;
63 GridIndex eIdx;
64 };
65
66public:
67 template <class RadiusFunction>
68 SegmentNetwork(const GridGeometry& gg, const RadiusFunction& radiusFunction)
69 {
70 const auto& gridView = gg.gridView();
71 segments_.resize(gridView.size(0));
72 for (const auto &element : elements(gridView))
73 {
74 const auto geometry = element.geometry();
75 const auto eIdx = gg.elementMapper().index(element);
76 const auto radius = radiusFunction(eIdx);
77
78 segments_[eIdx] = Segment{ geometry.corner(0), geometry.corner(1), radius, eIdx };
79 }
80
81 std::cout << "-- Extracted " << segments_.size() << " segments.\n";
82 }
83
85 auto radius(std::size_t eIdx) const
86 { return segments_[eIdx].r; }
87
89 const auto& segment(std::size_t eIdx) const
90 { return segments_[eIdx]; }
91
93 const auto& segments() const
94 { return segments_; }
95
98 auto findClosestSegmentSurface(const GlobalPosition& globalPos) const
99 {
100 GridIndex eIdx = 0;
101 double minSignedDist = std::numeric_limits<double>::max();
102
103 for (const auto &segment : segments_)
104 {
105 const auto signedDist = computeSignedDistance_(globalPos, segment);
106 if (signedDist < minSignedDist)
107 {
108 minSignedDist = signedDist;
109 eIdx = segment.eIdx;
110 }
111 }
112
113 using std::abs;
114 return std::make_tuple(eIdx, abs(minSignedDist));
115 }
116
118 template <class IndexRange>
119 auto findClosestSegmentSurface(const GlobalPosition& globalPos,
120 const IndexRange& segIndices) const
121 {
122 GridIndex eIdx = 0;
123 double minSignedDist = std::numeric_limits<double>::max();
124
125 for (const auto index : segIndices)
126 {
127 const auto &segment = segments_[index];
128 const auto signedDist = computeSignedDistance_(globalPos, segment);
129 if (signedDist < minSignedDist)
130 {
131 minSignedDist = signedDist;
132 eIdx = segment.eIdx;
133 }
134 }
135
136 using std::abs;
137 return std::make_tuple(eIdx, abs(minSignedDist));
138 }
139
140 // Compute the projection point of p onto the segment connecting a->b
141 GlobalPosition projectionPointOnSegment(const GlobalPosition& p, std::size_t idx) const
142 {
143 const auto& segment = segments_[idx];
144 return projectionPointOnSegment_(p, segment.a, segment.b);
145 }
146
147private:
148 // Compute the projection of p onto the segment
149 // returns the closest end point if the orthogonal projection onto the line implied by the segment
150 // is outside the segment
151 GlobalPosition projectionPointOnSegment_(const GlobalPosition& p, const GlobalPosition& a, const GlobalPosition& b) const
152 {
153 const auto v = b - a;
154 const auto w = p - a;
155
156 const auto proj1 = v*w;
157 if (proj1 <= 0.0)
158 return a;
159
160 const auto proj2 = v.two_norm2();
161 if (proj2 <= proj1)
162 return b;
163
164 const auto t = proj1 / proj2;
165 auto x = a;
166 x.axpy(t, v);
167 return x;
168 }
169
170 // Compute the signed (!) distance to a cylinder segment surface and the projection point onto the centerline
171 auto computeSignedDistance_(const GlobalPosition& p, const Segment& segment) const
172 {
173 const auto projPoint = projectionPointOnSegment_(p, segment.a, segment.b);
174 return (p-projPoint).two_norm()-segment.r;
175 }
176
177 std::vector<Segment> segments_;
178};
179
186template<class Network>
188{
189public:
190 NetworkIndicatorFunction(const Network& n)
191 : network_(n)
192 {}
193
194 template<class Pos>
195 auto operator()(const Pos& pos) const
196 {
197 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos);
198 return closestSegmentIdx;
199 }
200
201 template <class Pos, class HintContainer>
202 auto operator()(const Pos& pos, const HintContainer& hints) const
203 {
204 const auto& [closestSegmentIdx, _] = network_.findClosestSegmentSurface(pos, hints);
205 return closestSegmentIdx;
206 }
207
208private:
209 const Network& network_;
210};
211
217{
218public:
219 void push(const std::array<Dune::FieldVector<double, 3>, 3>& corners, double data)
220 {
221 vtkCells_.emplace_back(std::array<std::size_t, 3>{
222 { vtkPoints_.size(), vtkPoints_.size()+1, vtkPoints_.size()+2 }
223 });
224 for (const auto& c : corners)
225 vtkPoints_.push_back(c);
226 vtkCellData_.push_back(data);
227 }
228
229 void write(const std::string& name)
230 {
231 std::ofstream intersectionFile(name);
232 intersectionFile << "# vtk DataFile Version 2.0\n";
233 intersectionFile << "Really cool intersection data\n";
234 intersectionFile << "ASCII\n";
235 intersectionFile << "DATASET UNSTRUCTURED_GRID\n";
236 intersectionFile << "POINTS " << vtkPoints_.size() << " float\n";
237 for (const auto& p : vtkPoints_)
238 intersectionFile << p[0] << " " << p[1] << " " << p[2] << "\n";
239 intersectionFile << "\nCELLS " << vtkCells_.size() << " " << vtkCells_.size()*4 << "\n";
240 for (const auto& c : vtkCells_)
241 intersectionFile << "3 " << c[0] << " " << c[1] << " " << c[2] << "\n";
242 intersectionFile << "\nCELL_TYPES " << vtkCells_.size() << "\n";
243 for (int i = 0; i < vtkCells_.size(); ++i)
244 intersectionFile << "5\n";
245 intersectionFile << "\nCELL_DATA " << vtkCells_.size() << "\nSCALARS index float\nLOOKUP_TABLE default\n";
246 for (const auto& c : vtkCellData_)
247 intersectionFile << c << "\n";
248 }
249
250private:
251 std::vector<Dune::FieldVector<double, 3>> vtkPoints_;
252 std::vector<std::array<std::size_t, 3>> vtkCells_;
253 std::vector<double> vtkCellData_;
254};
255
256} // end namespace Detail
257
258namespace Embedded1d3dCouplingMode {
259struct Projection : public Utility::Tag<Projection> {
260 static std::string name() { return "projection"; }
261};
262
263inline constexpr Projection projection{};
264} // end namespace Embedded1d3dCouplingMode
265
266// forward declaration
267template<class MDTraits, class CouplingMode>
268class Embedded1d3dCouplingManager;
269
318template<class MDTraits>
319class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>
320: public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>>
321{
324
325 using Scalar = typename MDTraits::Scalar;
326 using SolutionVector = typename MDTraits::SolutionVector;
327 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
328
329 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
330 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
331
332 // the sub domain type aliases
333 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
334 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
335 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
336 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
337 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
338 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
339
340 static constexpr int lowDimDim = GridView<lowDimIdx>::dimension;
341 static constexpr int bulkDim = GridView<bulkIdx>::dimension;
342 static constexpr int dimWorld = GridView<bulkIdx>::dimensionworld;
343
344 template<std::size_t id>
345 static constexpr bool isBox()
346 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
347
348 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
349
351
352public:
353 static constexpr Embedded1d3dCouplingMode::Projection couplingMode{};
354
355 using ParentType::ParentType;
356
357 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
358 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
359 const SolutionVector& curSol)
360 {
361 ParentType::init(bulkProblem, lowDimProblem, curSol);
362 }
363
364 /* \brief Compute integration point point sources and associated data
365 * This method computes a source term on the explicitly given surface of the network grid
366 * by integrating over the elements of the surface, projecting the integration points onto
367 * the centerlines to evaluate 1d quantities and set point sources (minimum distance projection)
368 *
369 * \param order Unused in this algorithm! (passed from default 1D3D coupling manager)
370 * \param verbose If the source computation is verbose
371 */
372 void computePointSourceData(std::size_t order = 1, bool verbose = false)
373 {
374 // initialize the maps
375 // do some logging and profiling
376 Dune::Timer watch;
377 std::cout << "Initializing the coupling manager (projection)" << std::endl;
378
379 // debug VTK output
380 const bool enableIntersectionOutput = getParam<bool>("MixedDimension.Projection.EnableIntersectionOutput", false);
381 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkCoupledFaces =
382 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() : nullptr;
383 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkIntersections =
384 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() : nullptr;
385
386 // temporarily represent the network domain as a list of segments
387 // a more efficient data structure
388 const auto& lowDimProblem = this->problem(lowDimIdx);
389 const auto& lowDimFvGridGeometry = lowDimProblem.gridGeometry();
390 const auto& lowDimGridView = lowDimFvGridGeometry.gridView();
391 localAreaFactor_.resize(lowDimGridView.size(0), 0.0);
392 const auto radiusFunc = [&](const GridIndex<lowDimIdx> eIdx) {
393 return lowDimProblem.spatialParams().radius(eIdx);
394 };
395 const auto network = SegmentNetwork{
396 lowDimFvGridGeometry, radiusFunc
397 };
398
399 // precompute the vertex indices for efficiency (box method)
400 this->precomputeVertexIndices(bulkIdx);
401 this->precomputeVertexIndices(lowDimIdx);
402
403 const auto& bulkProblem = this->problem(bulkIdx);
404 const auto& bulkFvGridGeometry = bulkProblem.gridGeometry();
405 const auto& bulkGridView = bulkFvGridGeometry.gridView();
406 bulkElementMarker_.assign(bulkGridView.size(0), 0);
407 bulkVertexMarker_.assign(bulkGridView.size(bulkDim), 0);
408
409 // construct a bounding box that may be used to determine coupled faces
410 // the shrinking factor shrinks this bounding box
411 static const auto bBoxShrinking
412 = getParam<Scalar>("MixedDimension.Projection.CoupledBoundingBoxShrinkingFactor", 1e-2);
413 const GlobalPosition threshold(
414 bBoxShrinking*( bulkFvGridGeometry.bBoxMax()-bulkFvGridGeometry.bBoxMin() ).two_norm()
415 );
416 const auto bBoxMinSmall = bulkFvGridGeometry.bBoxMin() + threshold;
417 const auto bBoxMaxSmall = bulkFvGridGeometry.bBoxMax() - threshold;
418 auto insideBBox = [bBoxMin=bBoxMinSmall, bBoxMax=bBoxMaxSmall](const GlobalPosition& point) -> bool
419 {
420 static constexpr Scalar eps_ = 1.0e-7;
421 const Scalar eps0 = eps_*(bBoxMax[0] - bBoxMin[0]);
422 const Scalar eps1 = eps_*(bBoxMax[1] - bBoxMin[1]);
423 const Scalar eps2 = eps_*(bBoxMax[2] - bBoxMin[2]);
424 return (bBoxMin[0] - eps0 <= point[0] && point[0] <= bBoxMax[0] + eps0 &&
425 bBoxMin[1] - eps1 <= point[1] && point[1] <= bBoxMax[1] + eps1 &&
426 bBoxMin[2] - eps2 <= point[2] && point[2] <= bBoxMax[2] + eps2);
427 };
428
429 // setup simplex "quadrature" rule
430 static const auto quadSimplexRefineMaxLevel
431 = getParam<std::size_t>("MixedDimension.Projection.SimplexIntegrationRefine", 4);
432 const auto indicator = Detail::NetworkIndicatorFunction { network };
433
434 // estimate number of point sources for memory allocation
435 static const auto estimateNumberOfPointSources
436 = getParam<std::size_t>("MixedDimension.Projection.EstimateNumberOfPointSources", bulkFvGridGeometry.gridView().size(0));
437
438 this->pointSources(bulkIdx).reserve(estimateNumberOfPointSources);
439 this->pointSources(lowDimIdx).reserve(estimateNumberOfPointSources);
440 this->pointSourceData().reserve(estimateNumberOfPointSources);
441
442 // we go over the bulk surface intersection, check if we are coupled,
443 // and if this is the case, we add integration point sources using
444 // the local simplex refinement quadrature procedure (see paper Koch 2021)
445 for (const auto& element : elements(bulkGridView))
446 {
447 const auto bulkElementIdx = bulkFvGridGeometry.elementMapper().index(element);
448 const auto bulkGeometry = element.geometry();
449 const auto refElement = Dune::referenceElement(element);
450 for (const auto& intersection : intersections(bulkGridView, element))
451 {
452 // skip inner intersection
453 if (!intersection.boundary())
454 continue;
455
456 // check if we are close enough to the coupling interface
457 const auto [isCoupled, closestSegmentIdx, minDist] = isCoupled_(network, intersection, insideBBox);
458 if (!isCoupled)
459 continue;
460
461 // we found an intersection on the coupling interface: mark as coupled element
462 bulkElementMarker_[bulkElementIdx] = 1;
463
464 const auto interfaceGeometry = intersection.geometry();
465 for (int i = 0; i < interfaceGeometry.corners(); ++i)
466 {
467 const auto localVIdx = refElement.subEntity(intersection.indexInInside(), 1, i, bulkDim);
468 const auto vIdx = bulkFvGridGeometry.vertexMapper().subIndex(element, localVIdx, bulkDim);
469 bulkVertexMarker_[vIdx] = 1;
470 }
471
472 // debug graphical intersection output
473 if (enableIntersectionOutput)
474 vtkCoupledFaces->push({
475 { interfaceGeometry.corner(0), interfaceGeometry.corner(1), interfaceGeometry.corner(2)}
476 }, closestSegmentIdx);
477
478 const auto quadSimplex = LocalRefinementSimplexQuadrature{ interfaceGeometry, quadSimplexRefineMaxLevel, indicator };
479
480 // iterate over all quadrature points (children simplices)
481 for (const auto& qp : quadSimplex)
482 {
483 const auto surfacePoint = interfaceGeometry.global(qp.position());
484 const auto closestSegmentIdx = qp.indicator();
485 const auto projPoint = network.projectionPointOnSegment(surfacePoint, closestSegmentIdx);
486
487 addPointSourceAtIP_(bulkGeometry, interfaceGeometry, qp, bulkElementIdx, closestSegmentIdx, surfacePoint, projPoint);
488 addStencilEntries_(bulkElementIdx, closestSegmentIdx);
489 }
490 }
491 }
492
493 // make stencils unique
494 using namespace Dune::Hybrid;
495 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
496 {
497 for (auto&& stencil : this->couplingStencils(domainIdx))
498 {
499 std::sort(stencil.second.begin(), stencil.second.end());
500 stencil.second.erase(
501 std::unique(stencil.second.begin(), stencil.second.end()),
502 stencil.second.end()
503 );
504 }
505 });
506
507 // debug VTK output
508 if (enableIntersectionOutput)
509 vtkCoupledFaces->write("coupledfaces.vtk");
510
511 // compare theoretical cylinder surface to actual discrete surface
512 // of the coupled bulk interface facets
513 for (const auto& element : elements(lowDimGridView))
514 {
515 const auto length = element.geometry().volume();
516 const auto eIdx = lowDimFvGridGeometry.elementMapper().index(element);
517 const auto radius = this->problem(lowDimIdx).spatialParams().radius(eIdx);
518 const auto cylinderSurface = 2.0*M_PI*radius*length;
519 localAreaFactor_[eIdx] /= cylinderSurface;
520 localAreaFactor_[eIdx] -= 1.0;
521 }
522
523 this->pointSources(bulkIdx).shrink_to_fit();
524 this->pointSources(lowDimIdx).shrink_to_fit();
525 this->pointSourceData().shrink_to_fit();
526
527 std::cout << "-- Coupling at " << this->pointSourceData().size() << " integration points." << std::endl;
528 std::cout << "-- Finished initialization after " << watch.elapsed() << " seconds." << std::endl;
529 }
530
534 // \{
535
537 Scalar radius(std::size_t id) const
538 {
539 const auto& data = this->pointSourceData()[id];
540 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
541 }
542
543 // \}
544
546 const std::vector<int>& bulkElementMarker() const
547 { return bulkElementMarker_; }
548
550 const std::vector<int>& bulkVertexMarker() const
551 { return bulkVertexMarker_; }
552
556 const std::vector<Scalar>& localAreaFactor() const
557 { return localAreaFactor_; }
558
559private:
560 std::vector<int> bulkElementMarker_, bulkVertexMarker_;
561 std::vector<double> localAreaFactor_;
562
563 // add a point source containing all data needed to evaluate
564 // the coupling source term from both domains
565 template<class BulkGeometry, class InterfaceGeometry, class QP>
566 void addPointSourceAtIP_(const BulkGeometry& bulkGeometry,
567 const InterfaceGeometry& ifGeometry,
568 const QP& qp,
569 const GridIndex<bulkIdx> bulkElementIdx,
570 const GridIndex<lowDimIdx> lowDimElementIdx,
571 const GlobalPosition& surfacePoint,
572 const GlobalPosition& projPoint)
573 {
574 // add a new point source id (for the associated data)
575 const auto id = this->idCounter_++;
576
577 // add a bulk point source
578 const auto qpweight = qp.weight();
579 const auto integrationElement = ifGeometry.integrationElement(qp.position());
580 this->pointSources(bulkIdx).emplace_back(surfacePoint, id, qpweight, integrationElement, bulkElementIdx);
581
582 // find the corresponding projection point and 1d element
583 this->pointSources(lowDimIdx).emplace_back(projPoint, id, qpweight, integrationElement, lowDimElementIdx);
584 localAreaFactor_[lowDimElementIdx] += qpweight*integrationElement;
585
586 const auto& bulkFvGridGeometry = this->problem(bulkIdx).gridGeometry();
587 const auto& lowDimFvGridGeometry = this->problem(lowDimIdx).gridGeometry();
588
589 // pre compute additional data used for the evaluation of
590 // the actual solution dependent source term
591 PointSourceData psData;
592
593 if constexpr (isBox<lowDimIdx>())
594 {
595 std::vector<Dune::FieldVector<Scalar, 1> > shapeValues;
596 const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry();
597 this->getShapeValues(lowDimIdx, lowDimFvGridGeometry, lowDimGeometry, projPoint, shapeValues);
598 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
599 }
600 else
601 {
602 psData.addLowDimInterpolation(lowDimElementIdx);
603 }
604
605 if constexpr (isBox<bulkIdx>())
606 {
607 std::vector<Dune::FieldVector<Scalar, 1> > shapeValues;
608 this->getShapeValues(bulkIdx, bulkFvGridGeometry, bulkGeometry, surfacePoint, shapeValues);
609 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
610 }
611 else
612 {
613 psData.addBulkInterpolation(bulkElementIdx);
614 }
615
616 // publish point source data in the global vector
617 this->pointSourceData().emplace_back(std::move(psData));
618 }
619
620 void addStencilEntries_(const GridIndex<bulkIdx> bulkElementIdx, const GridIndex<lowDimIdx> lowDimElementIdx)
621 {
622 // export the bulk coupling stencil
623 if constexpr (isBox<lowDimIdx>())
624 {
625 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
626 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(
627 this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
628 vertices.begin(), vertices.end()
629 );
630
631 }
632 else
633 {
634 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
635 }
636
637 // export the lowdim coupling stencil
638 if constexpr (isBox<bulkIdx>())
639 {
640 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
641 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(
642 this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
643 vertices.begin(), vertices.end()
644 );
645 }
646 else
647 {
648 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
649 }
650 }
651
652 template<class BoundaryIntersection, class BBoxCheck>
653 auto isCoupled_(const SegmentNetwork& network, const BoundaryIntersection& intersection, const BBoxCheck& insideBBox) const
654 {
655 const auto center = intersection.geometry().center();
656 const auto [eIdx, minDist] = network.findClosestSegmentSurface(center);
657
658 // all elements in the bounding box are coupled if this is enabled
659 // this is a simplified check that can be enabled for box-shaped domains
660 static const bool considerBBox = getParam<bool>("MixedDimension.Projection.ConsiderFacesWithinBoundingBoxCoupled", false);
661 if (considerBBox && insideBBox(center))
662 return std::make_tuple(true, eIdx, minDist);
663
664 // general coupling face detection algorithm:
665 // if the distance is under a certain threshold we are coupled (and the angle is not very big)
666 static const auto rFactor = getParam<Scalar>("MixedDimension.Projection.CoupledRadiusFactor", 0.1);
667 static const auto spFactor = getParam<Scalar>("MixedDimension.Projection.CoupledAngleFactor", 0.3);
668 const auto projPoint = network.projectionPointOnSegment(center, eIdx);
669 const auto distVec = projPoint - center;
670 const bool isCoupled = minDist < rFactor * network.radius(eIdx) && distVec * intersection.centerUnitOuterNormal() > distVec.two_norm() * spFactor;
671 return std::make_tuple(isCoupled, eIdx, minDist);
672 }
673};
674
676template<class MDTraits>
677struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Projection>>
678: public std::true_type {};
679
680} // end namespace Dumux
681
682#endif
Simple legacy VTK writer for outputting debug data on the coupling interface.
Definition: couplingmanager1d3d_projection.hh:217
void push(const std::array< Dune::FieldVector< double, 3 >, 3 > &corners, double data)
Definition: couplingmanager1d3d_projection.hh:219
void write(const std::string &name)
Definition: couplingmanager1d3d_projection.hh:229
Get the closest segment for a given surface point.
Definition: couplingmanager1d3d_projection.hh:188
NetworkIndicatorFunction(const Network &n)
Definition: couplingmanager1d3d_projection.hh:190
auto operator()(const Pos &pos) const
Definition: couplingmanager1d3d_projection.hh:195
auto operator()(const Pos &pos, const HintContainer &hints) const
Definition: couplingmanager1d3d_projection.hh:202
Segment representation of a 1d network grid.
Definition: couplingmanager1d3d_projection.hh:54
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:119
const auto & segment(std::size_t eIdx) const
the segment with index eIdx
Definition: couplingmanager1d3d_projection.hh:89
auto radius(std::size_t eIdx) const
the segment radius
Definition: couplingmanager1d3d_projection.hh:85
GlobalPosition projectionPointOnSegment(const GlobalPosition &p, std::size_t idx) const
Definition: couplingmanager1d3d_projection.hh:141
const auto & segments() const
the segments
Definition: couplingmanager1d3d_projection.hh:93
SegmentNetwork(const GridGeometry &gg, const RadiusFunction &radiusFunction)
Definition: couplingmanager1d3d_projection.hh:68
auto findClosestSegmentSurface(const GlobalPosition &globalPos) const
Definition: couplingmanager1d3d_projection.hh:98
Manages the coupling between bulk elements and lower dimensional elements.
Definition: couplingmanager1d3d_projection.hh:321
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_projection.hh:372
const std::vector< int > & bulkElementMarker() const
Marker that is non-zero for bulk elements that are coupled.
Definition: couplingmanager1d3d_projection.hh:546
const std::vector< int > & bulkVertexMarker() const
Marker that is non-zero for bulk vertices that belong to coupled surface facets.
Definition: couplingmanager1d3d_projection.hh:550
const std::vector< Scalar > & localAreaFactor() const
Definition: couplingmanager1d3d_projection.hh:556
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_projection.hh:537
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_projection.hh:357
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:24
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:71
A quadrature rule using local refinement to approximate partitioned elements.
Definition: localrefinementquadrature.hh:41
A point source data class used for integration in multidimensional models.
Definition: pointsourcedata.hh:31
void addLowDimInterpolation(const ShapeValues &shapeValues, const std::vector< GridIndex< lowDimIdx > > &cornerIndices, GridIndex< lowDimIdx > eIdx)
Definition: pointsourcedata.hh:60
void addBulkInterpolation(const ShapeValues &shapeValues, const std::vector< GridIndex< bulkIdx > > &cornerIndices, GridIndex< bulkIdx > eIdx)
Definition: pointsourcedata.hh:50
Defines all properties used in Dumux.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
A function to compute the convex hull of a point cloud.
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Defines the index types used for grid and local indices.
A quadrature rule using local refinement to approximate partitioned elements.
Define some often used mathematical functions.
constexpr Box box
Definition: method.hh:147
constexpr Projection projection
Definition: couplingmanager1d3d_projection.hh:263
Definition: adapt.hh:17
A quadrature based on refinement.
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
Definition: couplingmanager1d3d_projection.hh:259
static std::string name()
Definition: couplingmanager1d3d_projection.hh:260
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
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.
Functionality to triangulate point clouds.