25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_AVERAGE_HH
30#include <dune/common/timer.hh>
31#include <dune/geometry/quadraturerules.hh>
46namespace Embedded1d3dCouplingMode {
48 static std::string
name() {
return "average"; }
55template<
class MDTraits,
class CouplingMode>
56class Embedded1d3dCouplingManager;
64template<
class MDTraits>
67 CircleAveragePointSourceTraits<MDTraits>>
71 using Scalar =
typename MDTraits::Scalar;
72 using SolutionVector =
typename MDTraits::SolutionVector;
73 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
75 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
76 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
79 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
82 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
83 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
86 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
88 template<std::
size_t id>
89 static constexpr bool isBox()
95 bulkDim = GridView<bulkIdx>::dimension,
96 lowDimDim = GridView<lowDimIdx>::dimension,
97 dimWorld = GridView<bulkIdx>::dimensionworld
102 using ParentType::ParentType;
104 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
105 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
106 const SolutionVector& curSol)
108 ParentType::init(bulkProblem, lowDimProblem, curSol);
109 computeLowDimVolumeFractions();
118 template<std::
size_t id,
class JacobianPattern>
121 extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
131 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
133 const LocalAssemblerI& localAssemblerI,
134 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
135 JacobianMatrixDiagBlock& A,
136 GridVariables& gridVariables)
138 extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, A, gridVariables);
153 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
154 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
155 const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
160 std::cout <<
"Initializing the point sources..." << std::endl;
165 extendedSourceStencil_.clear();
168 this->precomputeVertexIndices(bulkIdx);
169 this->precomputeVertexIndices(lowDimIdx);
175 const auto& lowDimProblem = this->problem(lowDimIdx);
176 for (
const auto& is : intersections(this->glue()))
179 const auto& lowDimElement = is.targetEntity(0);
180 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
183 const auto intersectionGeometry = is.geometry();
185 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
194 for (
auto&& qp : quad)
197 const auto globalPos = intersectionGeometry.global(qp.position());
203 if (bulkElementIndices.empty())
210 static const auto numIp = getParam<int>(
"MixedDimension.NumCircleSegments");
211 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
212 const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
213 const auto circleAvgWeight = 2*M_PI*radius/numIp;
214 const auto integrationElement = intersectionGeometry.integrationElement(qp.position());
215 const auto qpweight = qp.weight();
218 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(
circlePoints.size());
219 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(
circlePoints.size());
222 std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
223 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
224 std::vector<ShapeValues> circleShapeValues;
227 int insideCirclePoints = 0;
231 if (circleBulkElementIndices.empty())
234 ++insideCirclePoints;
235 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices.size();
236 for (
const auto bulkElementIdx : circleBulkElementIndices)
238 circleStencil.push_back(bulkElementIdx);
239 circleIpWeight.push_back(localCircleAvgWeight);
242 if constexpr (isBox<bulkIdx>())
244 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
245 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
248 const auto bulkGeometry = bulkElement.geometry();
249 ShapeValues shapeValues;
250 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry,
circlePoints[k], shapeValues);
251 circleShapeValues.emplace_back(std::move(shapeValues));
258 if (circleStencil.empty())
262 if constexpr (isBox<bulkIdx>())
265 for (
const auto& vertices : circleCornerIndices)
267 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
268 vertices->begin(), vertices->end());
274 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
275 circleStencil.begin(), circleStencil.end());
279 const auto surfaceFraction = Scalar(insideCirclePoints)/Scalar(
circlePoints.size());
282 for (
auto bulkElementIdx : bulkElementIndices)
284 const auto id = this->idCounter_++;
286 this->pointSources(bulkIdx).emplace_back(globalPos,
id, qpweight, integrationElement*surfaceFraction, bulkElementIdx);
287 this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
288 this->pointSources(lowDimIdx).emplace_back(globalPos,
id, qpweight, integrationElement*surfaceFraction, lowDimElementIdx);
289 this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
293 PointSourceData psData;
295 if constexpr (isBox<lowDimIdx>())
297 ShapeValues shapeValues;
298 this->getShapeValues(lowDimIdx, lowDimGridGeometry, intersectionGeometry, globalPos, shapeValues);
299 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
303 psData.addLowDimInterpolation(lowDimElementIdx);
307 if constexpr (isBox<bulkIdx>())
309 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
311 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
312 ShapeValues shapeValues;
313 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, globalPos, shapeValues);
314 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
318 psData.addCircleInterpolation(circleIpWeight, circleStencil);
319 psData.addBulkInterpolation(bulkElementIdx);
323 this->pointSourceData().emplace_back(std::move(psData));
326 const auto outsideGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
330 if constexpr (isBox<lowDimIdx>())
332 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
333 this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
334 this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
339 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
343 if constexpr (isBox<bulkIdx>())
346 for (
const auto& vertices : circleCornerIndices)
348 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
349 vertices->begin(), vertices->end());
355 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
356 circleStencil.begin(), circleStencil.end());
363 for (
auto&& stencil : extendedSourceStencil_.stencil())
365 std::sort(stencil.second.begin(), stencil.second.end());
366 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
369 if constexpr (isBox<bulkIdx>())
371 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
372 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
373 [&](
auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
374 stencil.second.end());
379 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
380 [&](
auto i){ return i == stencil.first; }),
381 stencil.second.end());
386 using namespace Dune::Hybrid;
387 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
389 for (
auto&& stencil : this->couplingStencils(domainIdx))
391 std::sort(stencil.second.begin(), stencil.second.end());
392 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
396 std::cout <<
"took " << watch.elapsed() <<
" seconds." << std::endl;
403 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
405 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
406 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
409 for (
const auto& is : intersections(this->glue()))
412 const auto& inside = is.targetEntity(0);
413 const auto intersectionGeometry = is.geometry();
414 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
417 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
418 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
420 const auto& outside = is.domainEntity(outsideIdx);
421 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
422 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
435 const auto& data = this->pointSourceData()[id];
436 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
443 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
444 return lowDimVolumeInBulkElement_[eIdx];
451 const auto totalVolume = element.geometry().volume();
452 return lowDimVolume(element) / totalVolume;
462 std::vector<Scalar> lowDimVolumeInBulkElement_;
Helper class to create (named and comparable) tagged types.
Defines the index types used for grid and local indices.
Helper functions for distance queries.
Extended source stencil helper class for coupling managers.
Helper function to compute points on a circle.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Point source traits for average-based coupling modes.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
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:112
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:39
void circlePoints(std::vector< GlobalPosition > &points, const std::vector< Scalar > &sincos, const GlobalPosition ¢er, const GlobalPosition &normal, const Scalar radius)
Definition: circlepoints.hh:50
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
constexpr Box box
Definition: method.hh:139
constexpr Average average
Definition: couplingmanager1d3d_average.hh:51
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:57
Definition: common/properties.hh:102
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
Definition: couplingmanager1d3d_average.hh:47
static std::string name()
Definition: couplingmanager1d3d_average.hh:48
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_average.hh:68
void extendJacobianPattern(Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition: couplingmanager1d3d_average.hh:119
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_average.hh:104
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_average.hh:400
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_average.hh:150
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:449
void evalAdditionalDomainDerivatives(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: couplingmanager1d3d_average.hh:132
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_average.hh:433
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:441
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:79
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:46
Declares all properties used in Dumux.