25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_HH
30#include <dune/common/timer.hh>
31#include <dune/geometry/quadraturerules.hh>
44namespace Embedded1d3dCouplingMode {
46 static std::string
name() {
return "surface"; }
53template<
class MDTraits,
class CouplingMode>
54class Embedded1d3dCouplingManager;
63template<
class MDTraits>
66 CircleAveragePointSourceTraits<MDTraits>>
70 using Scalar =
typename MDTraits::Scalar;
71 using SolutionVector =
typename MDTraits::SolutionVector;
72 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
74 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
75 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
78 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
81 template<std::
size_t id>
using GridView =
typename GridGeometry<id>::GridView;
82 template<std::
size_t id>
using Element =
typename GridView<id>::template Codim<0>::Entity;
85 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
87 template<std::
size_t id>
88 static constexpr bool isBox()
92 bulkDim = GridView<bulkIdx>::dimension,
93 lowDimDim = GridView<lowDimIdx>::dimension,
94 dimWorld = GridView<bulkIdx>::dimensionworld
99 using ParentType::ParentType;
101 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
102 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
103 const SolutionVector& curSol)
105 ParentType::init(bulkProblem, lowDimProblem, curSol);
106 computeLowDimVolumeFractions();
115 template<std::
size_t id,
class JacobianPattern>
118 extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
128 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
130 const LocalAssemblerI& localAssemblerI,
131 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
132 JacobianMatrixDiagBlock& A,
133 GridVariables& gridVariables)
135 extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, A, gridVariables);
150 static const bool useCircleAverage = getParam<bool>(
"MixedDimension.UseCircleAverage",
true);
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_.stencil().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 integrationElement = intersectionGeometry.integrationElement(qp.position())*2*M_PI*radius/Scalar(numIp);
214 const auto qpweight = qp.weight()/(2*M_PI*radius);
215 const auto circleAvgWeight = 2*M_PI*radius/numIp;
218 std::vector<std::vector<std::size_t>> circleBulkElementIndices(
circlePoints.size());
219 std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(
circlePoints.size());
220 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;
230 if (circleBulkElementIndices[k].empty())
233 const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices[k].size();
234 for (
const auto bulkElementIdx : circleBulkElementIndices[k])
236 circleStencil.push_back(bulkElementIdx);
237 circleIpWeight.push_back(localCircleAvgWeight);
240 if constexpr (isBox<bulkIdx>())
242 const auto bulkElement = bulkGridGeometry.element(bulkElementIdx);
243 circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx)));
246 const auto bulkGeometry = bulkElement.geometry();
247 ShapeValues shapeValues;
248 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry,
circlePoints[k], shapeValues);
249 circleShapeValues.emplace_back(std::move(shapeValues));
255 if constexpr (isBox<bulkIdx>())
258 for (
const auto& vertices : circleCornerIndices)
260 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
261 vertices->begin(), vertices->end());
267 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
268 circleStencil.begin(), circleStencil.end());
276 if (circleBulkElementIndices[k].empty())
281 for (
const auto bulkElementIdx : circleBulkElementIndices[k])
283 const auto id = this->idCounter_++;
285 this->pointSources(bulkIdx).emplace_back(circlePos,
id, qpweight, integrationElement, bulkElementIdx);
286 this->pointSources(bulkIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
287 this->pointSources(lowDimIdx).emplace_back(globalPos,
id, qpweight, integrationElement, lowDimElementIdx);
288 this->pointSources(lowDimIdx).back().setEmbeddings(circleBulkElementIndices[k].size());
292 PointSourceData psData;
294 if constexpr (isBox<lowDimIdx>())
296 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
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 if (useCircleAverage)
310 psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil);
312 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
313 const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry();
314 ShapeValues shapeValues;
315 this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePos, shapeValues);
316 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
320 if (useCircleAverage)
321 psData.addCircleInterpolation(circleIpWeight, circleStencil);
323 psData.addBulkInterpolation(bulkElementIdx);
327 this->pointSourceData().emplace_back(std::move(psData));
331 if constexpr (isBox<lowDimIdx>())
333 const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx);
334 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
335 vertices.begin(), vertices.end());
340 this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx);
344 if (useCircleAverage)
346 if constexpr (isBox<bulkIdx>())
349 for (
const auto& vertices : circleCornerIndices)
351 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
352 vertices->begin(), vertices->end());
358 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
359 circleStencil.begin(), circleStencil.end());
368 for (
auto&& stencil : extendedSourceStencil_.stencil())
370 std::sort(stencil.second.begin(), stencil.second.end());
371 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
374 if constexpr (isBox<bulkIdx>())
376 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
377 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
378 [&](
auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
379 stencil.second.end());
384 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
385 [&](
auto i){ return i == stencil.first; }),
386 stencil.second.end());
391 using namespace Dune::Hybrid;
392 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
394 for (
auto&& stencil : this->couplingStencils(domainIdx))
396 std::sort(stencil.second.begin(), stencil.second.end());
397 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
401 std::cout <<
"took " << watch.elapsed() <<
" seconds." << std::endl;
408 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
410 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
411 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
414 for (
const auto& is : intersections(this->glue()))
417 const auto& inside = is.targetEntity(0);
418 const auto intersectionGeometry = is.geometry();
419 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
422 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
423 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
425 const auto& outside = is.domainEntity(outsideIdx);
426 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
427 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
440 const auto& data = this->pointSourceData()[id];
441 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
448 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
449 return lowDimVolumeInBulkElement_[eIdx];
456 const auto totalVolume = element.geometry().volume();
457 return lowDimVolume(element) / totalVolume;
467 std::vector<Scalar> lowDimVolumeInBulkElement_;
Defines the index types used for grid and local indices.
Helper class to create (named and comparable) tagged types.
Point source traits for average-based coupling modes.
Helper function to compute points on a circle.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Extended source stencil helper class for coupling managers.
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
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 Surface surface
Definition: couplingmanager1d3d_surface.hh:49
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_surface.hh:45
static std::string name()
Definition: couplingmanager1d3d_surface.hh:46
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_surface.hh:67
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_surface.hh:147
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_surface.hh:129
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_surface.hh:446
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_surface.hh:116
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_surface.hh:405
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_surface.hh:438
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_surface.hh:454
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_surface.hh:101
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.