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"; }
50 [[deprecated(
"Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
52 [[deprecated(
"Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
54 [[deprecated(
"Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
56 [[deprecated(
"Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
64template<
class MDTraits,
class CouplingMode>
65class Embedded1d3dCouplingManager;
73template<
class MDTraits>
76 CircleAveragePointSourceTraits<MDTraits>>
80 using Scalar =
typename MDTraits::Scalar;
81 using SolutionVector =
typename MDTraits::SolutionVector;
82 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
84 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
85 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
88 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
91 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
92 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
95 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
97 template<std::
size_t id>
98 static constexpr bool isBox()
104 bulkDim = GridView<bulkIdx>::dimension,
105 lowDimDim = GridView<lowDimIdx>::dimension,
106 dimWorld = GridView<bulkIdx>::dimensionworld
111 using ParentType::ParentType;
113 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
114 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
115 const SolutionVector& curSol)
117 ParentType::init(bulkProblem, lowDimProblem, curSol);
118 computeLowDimVolumeFractions();
127 template<std::
size_t id,
class JacobianPattern>
130 extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
140 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
142 const LocalAssemblerI& localAssemblerI,
143 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
144 JacobianMatrixDiagBlock& A,
145 GridVariables& gridVariables)
147 extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, this->curSol(), A, gridVariables);
162 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
163 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
164 const auto& bulkTree = bulkGridGeometry.boundingBoxTree();
169 std::cout <<
"Initializing the point sources..." << std::endl;
174 extendedSourceStencil_.clear();
177 this->precomputeVertexIndices(bulkIdx);
178 this->precomputeVertexIndices(lowDimIdx);
184 const auto& lowDimProblem = this->problem(lowDimIdx);
185 for (
const auto& is : intersections(this->glue()))
188 const auto& lowDimElement = is.targetEntity(0);
189 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement);
192 const auto intersectionGeometry = is.geometry();
194 const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order);
203 for (
auto&& qp : quad)
206 const auto globalPos = intersectionGeometry.global(qp.position());
212 if (bulkElementIndices.empty())
219 static const auto numIp = getParam<int>(
"MixedDimension.NumCircleSegments");
220 const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx);
221 const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0);
222 const auto circleAvgWeight = 2*M_PI*radius/numIp;
223 const auto integrationElement = intersectionGeometry.integrationElement(qp.position());
224 const auto qpweight = qp.weight();
227 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(
circlePoints.size());
228 std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(
circlePoints.size());
231 std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices;
232 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
233 std::vector<ShapeValues> circleShapeValues;
236 int insideCirclePoints = 0;
240 if (circleBulkElementIndices.empty())
243 ++insideCirclePoints;
244 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices.size();
245 for (
const auto bulkElementIdx : circleBulkElementIndices)
247 circleStencil.push_back(bulkElementIdx);
248 circleIpWeight.push_back(localCircleAvgWeight);
251 if constexpr (isBox<bulkIdx>())
253 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
254 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
257 const auto bulkGeometry = bulkElement.geometry();
258 ShapeValues shapeValues;
259 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry,
circlePoints[k], shapeValues);
260 circleShapeValues.emplace_back(std::move(shapeValues));
266 if constexpr (isBox<bulkIdx>())
269 for (
const auto& vertices : circleCornerIndices)
271 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
272 vertices->begin(), vertices->end());
278 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
279 circleStencil.begin(), circleStencil.end());
283 const auto surfaceFraction = Scalar(insideCirclePoints)/Scalar(
circlePoints.size());
286 for (
auto bulkElementIdx : bulkElementIndices)
288 const auto id = this->idCounter_++;
290 this->pointSources(bulkIdx).emplace_back(globalPos,
id, qpweight, integrationElement*surfaceFraction, bulkElementIdx);
291 this->pointSources(bulkIdx).back().setEmbeddings(bulkElementIndices.size());
292 this->pointSources(lowDimIdx).emplace_back(globalPos,
id, qpweight, integrationElement*surfaceFraction, lowDimElementIdx);
293 this->pointSources(lowDimIdx).back().setEmbeddings(bulkElementIndices.size());
297 PointSourceData psData;
299 if constexpr (isBox<lowDimIdx>())
301 ShapeValues shapeValues;
302 this->getShapeValues(lowDimIdx, lowDimGridGeometry, intersectionGeometry, globalPos, shapeValues);
303 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
307 psData.addLowDimInterpolation(lowDimElementIdx);
311 if constexpr (isBox<bulkIdx>())
313 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
315 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
316 ShapeValues shapeValues;
317 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, globalPos, shapeValues);
318 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
322 psData.addCircleInterpolation(circleIpWeight, circleStencil);
323 psData.addBulkInterpolation(bulkElementIdx);
327 this->pointSourceData().emplace_back(std::move(psData));
330 const auto outsideGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
334 if constexpr (isBox<lowDimIdx>())
336 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
337 this->vertexIndices(lowDimIdx, lowDimElementIdx).begin(),
338 this->vertexIndices(lowDimIdx, lowDimElementIdx).end());
343 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
347 if constexpr (isBox<bulkIdx>())
350 for (
const auto& vertices : circleCornerIndices)
352 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
353 vertices->begin(), vertices->end());
359 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
360 circleStencil.begin(), circleStencil.end());
367 for (
auto&& stencil : extendedSourceStencil_.stencil())
369 std::sort(stencil.second.begin(), stencil.second.end());
370 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
373 if constexpr (isBox<bulkIdx>())
375 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
376 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
377 [&](
auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
378 stencil.second.end());
383 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
384 [&](
auto i){ return i == stencil.first; }),
385 stencil.second.end());
390 using namespace Dune::Hybrid;
391 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
393 for (
auto&& stencil : this->couplingStencils(domainIdx))
395 std::sort(stencil.second.begin(), stencil.second.end());
396 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
400 std::cout <<
"took " << watch.elapsed() <<
" seconds." << std::endl;
407 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
409 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
410 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
413 for (
const auto& is : intersections(this->glue()))
416 const auto& inside = is.targetEntity(0);
417 const auto intersectionGeometry = is.geometry();
418 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
421 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
422 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
424 const auto& outside = is.domainEntity(outsideIdx);
425 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
426 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
439 const auto& data = this->pointSourceData()[id];
440 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
447 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
448 return lowDimVolumeInBulkElement_[eIdx];
455 const auto totalVolume = element.geometry().volume();
456 return lowDimVolume(element) / totalVolume;
466 std::vector<Scalar> lowDimVolumeInBulkElement_;
Helper function to compute points on a circle.
Point source traits for average-based coupling modes.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Extended source stencil helper class for coupling managers.
Helper class to create (named and comparable) tagged types.
Defines the index types used for grid and local indices.
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: geometry/distance.hh:36
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: geometry/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: geometry/intersectingentities.hh:100
EmbeddedCouplingMode
The coupling mode.
Definition: couplingmanager1d3d.hh:44
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
constexpr Average average
Definition: couplingmanager1d3d_average.hh:60
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:101
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
friend constexpr bool operator==(EmbeddedCouplingMode m, Average)
Definition: couplingmanager1d3d_average.hh:53
friend constexpr bool operator==(Average, EmbeddedCouplingMode m)
Definition: couplingmanager1d3d_average.hh:51
friend constexpr bool operator!=(Average, EmbeddedCouplingMode m)
Definition: couplingmanager1d3d_average.hh:55
static std::string name()
Definition: couplingmanager1d3d_average.hh:48
friend constexpr bool operator!=(EmbeddedCouplingMode m, Average)
Definition: couplingmanager1d3d_average.hh:57
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_average.hh:77
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:128
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_average.hh:113
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_average.hh:404
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_average.hh:159
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:453
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:141
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_average.hh:437
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_average.hh:445
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:78
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:46
Helper functions for distance queries.
Declares all properties used in Dumux.