Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup EmbeddedCoupling | ||
10 | * \brief Coupling manager for low-dimensional domains embedded in the bulk domain. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_HH | ||
14 | #define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_SURFACE_HH | ||
15 | |||
16 | #include <vector> | ||
17 | |||
18 | #include <dune/common/timer.hh> | ||
19 | #include <dune/geometry/quadraturerules.hh> | ||
20 | |||
21 | #include <dumux/common/tag.hh> | ||
22 | #include <dumux/common/properties.hh> | ||
23 | #include <dumux/common/indextraits.hh> | ||
24 | |||
25 | #include <dumux/multidomain/embedded/couplingmanagerbase.hh> | ||
26 | #include <dumux/multidomain/embedded/circlepoints.hh> | ||
27 | #include <dumux/multidomain/embedded/circleaveragepointsourcetraits.hh> | ||
28 | #include <dumux/multidomain/embedded/extendedsourcestencil.hh> | ||
29 | |||
30 | namespace Dumux { | ||
31 | |||
32 | namespace Embedded1d3dCouplingMode { | ||
33 | struct Surface : public Utility::Tag<Surface> { | ||
34 | static std::string name() { return "surface"; } | ||
35 | }; | ||
36 | |||
37 | inline constexpr Surface surface{}; | ||
38 | } // end namespace Embedded1d3dCouplingMode | ||
39 | |||
40 | // forward declaration | ||
41 | template<class MDTraits, class CouplingMode> | ||
42 | class Embedded1d3dCouplingManager; | ||
43 | |||
44 | /*! | ||
45 | * \ingroup EmbeddedCoupling | ||
46 | * \brief Manages the coupling between bulk elements and lower dimensional elements | ||
47 | * Point sources on each integration point are computed by an AABB tree. | ||
48 | * \note Specialization for coupling method using cylinder sources with 3d quantities evaluated on the cylinder surface | ||
49 | * \note the cylinder source is approximated by point sources on the cylinder surface | ||
50 | */ | ||
51 | template<class MDTraits> | ||
52 | class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface> | ||
53 | : public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>, | ||
54 | CircleAveragePointSourceTraits<MDTraits>> | ||
55 | { | ||
56 | using ThisType = Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>; | ||
57 | using ParentType = EmbeddedCouplingManagerBase<MDTraits, ThisType, CircleAveragePointSourceTraits<MDTraits>>; | ||
58 | using Scalar = typename MDTraits::Scalar; | ||
59 | using SolutionVector = typename MDTraits::SolutionVector; | ||
60 | using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData; | ||
61 | |||
62 | static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index(); | ||
63 | static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index(); | ||
64 | |||
65 | // the sub domain type aliases | ||
66 | template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag; | ||
67 | template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>; | ||
68 | template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>; | ||
69 | template<std::size_t id> using GridView = typename GridGeometry<id>::GridView; | ||
70 | template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity; | ||
71 | template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex; | ||
72 | |||
73 | using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate; | ||
74 | |||
75 | template<std::size_t id> | ||
76 | static constexpr bool isBox() | ||
77 | { return GridGeometry<id>::discMethod == DiscretizationMethods::box; } | ||
78 | |||
79 | enum { | ||
80 | bulkDim = GridView<bulkIdx>::dimension, | ||
81 | lowDimDim = GridView<lowDimIdx>::dimension, | ||
82 | dimWorld = GridView<bulkIdx>::dimensionworld | ||
83 | }; | ||
84 | public: | ||
85 | static constexpr Embedded1d3dCouplingMode::Surface couplingMode{}; | ||
86 | |||
87 |
4/8✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 5 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 5 times.
✗ Branch 8 not taken.
|
5 | using ParentType::ParentType; |
88 | |||
89 | 5 | void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem, | |
90 | std::shared_ptr<Problem<lowDimIdx>> lowDimProblem, | ||
91 | const SolutionVector& curSol) | ||
92 | { | ||
93 |
3/10✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 5 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 5 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
|
10 | ParentType::init(bulkProblem, lowDimProblem, curSol); |
94 | 5 | computeLowDimVolumeFractions(); | |
95 | 5 | } | |
96 | |||
97 | /*! | ||
98 | * \brief extend the jacobian pattern of the diagonal block of domain i | ||
99 | * by those entries that are not already in the uncoupled pattern | ||
100 | * \note Such additional dependencies can arise from the coupling, e.g. if a coupling source | ||
101 | * term depends on a non-local average of a quantity of the same domain | ||
102 | */ | ||
103 | template<std::size_t id, class JacobianPattern> | ||
104 | ✗ | void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const | |
105 | { | ||
106 | ✗ | extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern); | |
107 | ✗ | } | |
108 | |||
109 | /*! | ||
110 | * \brief evaluate additional derivatives of the element residual of a domain with respect | ||
111 | * to dofs in the same domain that are not in the regular stencil (per default this is not the case) | ||
112 | * \note Such additional dependencies can arise from the coupling, e.g. if a coupling source | ||
113 | * term depends on a non-local average of a quantity of the same domain | ||
114 | * \note This is the same for box and cc | ||
115 | */ | ||
116 | template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables> | ||
117 | ✗ | void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI, | |
118 | const LocalAssemblerI& localAssemblerI, | ||
119 | const typename LocalAssemblerI::LocalResidual::ElementResidualVector&, | ||
120 | JacobianMatrixDiagBlock& A, | ||
121 | GridVariables& gridVariables) | ||
122 | { | ||
123 | 300595 | extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables); | |
124 | ✗ | } | |
125 | |||
126 | /* \brief Compute integration point point sources and associated data | ||
127 | * | ||
128 | * This method uses grid glue to intersect the given grids. Over each intersection | ||
129 | * we later need to integrate a source term. This method places point sources | ||
130 | * at each quadrature point and provides the point source with the necessary | ||
131 | * information to compute integrals (quadrature weight and integration element) | ||
132 | * \param order The order of the quadrature rule for integration of sources over an intersection | ||
133 | * \param verbose If the point source computation is verbose | ||
134 | */ | ||
135 | 5 | void computePointSourceData(std::size_t order = 1, bool verbose = false) | |
136 | { | ||
137 | // if we use the circle average as the 3D values or a point evaluation | ||
138 |
3/6✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
5 | static const bool useCircleAverage = getParam<bool>("MixedDimension.UseCircleAverage", true); |
139 | |||
140 | // Initialize the bulk bounding box tree | ||
141 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
|
5 | const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry(); |
142 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry(); |
143 | 5 | const auto& bulkTree = bulkGridGeometry.boundingBoxTree(); | |
144 | |||
145 | // initialize the maps | ||
146 | // do some logging and profiling | ||
147 | 5 | Dune::Timer watch; | |
148 | 10 | std::cout << "Initializing the point sources..." << std::endl; | |
149 | |||
150 | // clear all internal members like pointsource vectors and stencils | ||
151 | // initializes the point source id counter | ||
152 | 5 | this->clear(); | |
153 | 10 | extendedSourceStencil_.stencil().clear(); | |
154 | |||
155 | // precompute the vertex indices for efficiency | ||
156 | 5 | this->precomputeVertexIndices(bulkIdx); | |
157 | 5 | this->precomputeVertexIndices(lowDimIdx); | |
158 | |||
159 | // intersect the bounding box trees | ||
160 | 5 | this->glueGrids(); | |
161 | |||
162 | // iterate over all intersection and add point sources | ||
163 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | const auto& lowDimProblem = this->problem(lowDimIdx); |
164 |
4/4✓ Branch 0 taken 95 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 95 times.
✓ Branch 3 taken 5 times.
|
110 | for (const auto& is : intersections(this->glue())) |
165 | { | ||
166 | // all inside elements are identical... | ||
167 | 95 | const auto& lowDimElement = is.targetEntity(0); | |
168 | 190 | const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(lowDimElement); | |
169 | |||
170 | // get the intersection geometry | ||
171 | 95 | const auto intersectionGeometry = is.geometry(); | |
172 | // get the Gaussian quadrature rule for the local intersection | ||
173 | 95 | const auto& quad = Dune::QuadratureRules<Scalar, lowDimDim>::rule(intersectionGeometry.type(), order); | |
174 | |||
175 | // apply the Gaussian quadrature rule and define point sources at each quadrature point | ||
176 | // note that the approximation is not optimal if | ||
177 | // (a) the one-dimensional elements are too large, | ||
178 | // (b) whenever a one-dimensional element is split between two or more elements, | ||
179 | // (c) when gradients of important quantities in the three-dimensional domain are large. | ||
180 | |||
181 | // iterate over all quadrature points | ||
182 |
4/4✓ Branch 0 taken 190 times.
✓ Branch 1 taken 95 times.
✓ Branch 2 taken 190 times.
✓ Branch 3 taken 95 times.
|
475 | for (auto&& qp : quad) |
183 | { | ||
184 | // global position of the quadrature point | ||
185 | 380 | const auto globalPos = intersectionGeometry.global(qp.position()); | |
186 | |||
187 |
1/4✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
380 | const auto bulkElementIndices = intersectingEntities(globalPos, bulkTree); |
188 | |||
189 | // do not add a point source if the qp is outside of the 3d grid | ||
190 | // this is equivalent to having a source of zero for that qp | ||
191 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 190 times.
|
380 | if (bulkElementIndices.empty()) |
192 | ✗ | continue; | |
193 | |||
194 | //////////////////////////////////////////////////////////////// | ||
195 | // get points on the cylinder surface at the integration point | ||
196 | //////////////////////////////////////////////////////////////// | ||
197 | |||
198 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 185 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
190 | static const auto numIp = getParam<int>("MixedDimension.NumCircleSegments"); |
199 | 190 | const auto radius = lowDimProblem.spatialParams().radius(lowDimElementIdx); | |
200 | 570 | const auto normal = intersectionGeometry.corner(1)-intersectionGeometry.corner(0); | |
201 |
1/2✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
|
190 | const auto integrationElement = intersectionGeometry.integrationElement(qp.position())*2*M_PI*radius/Scalar(numIp); |
202 | 190 | const auto qpweight = qp.weight()/(2*M_PI*radius); | |
203 | 190 | const auto circleAvgWeight = 2*M_PI*radius/numIp; | |
204 | |||
205 |
2/4✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 190 times.
✗ Branch 4 not taken.
|
380 | const auto circlePoints = EmbeddedCoupling::circlePoints(globalPos, normal, radius, numIp); |
206 |
4/10✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 190 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 190 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 190 times.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
760 | std::vector<std::vector<std::size_t>> circleBulkElementIndices(circlePoints.size()); |
207 |
3/8✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 190 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 190 times.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
760 | std::vector<Scalar> circleIpWeight; circleIpWeight.reserve(circlePoints.size()); |
208 |
4/10✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 190 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 190 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 190 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
|
760 | std::vector<GridIndex<bulkIdx>> circleStencil; circleStencil.reserve(circlePoints.size()); |
209 | // for box | ||
210 |
1/2✓ Branch 0 taken 190 times.
✗ Branch 1 not taken.
|
380 | std::vector<const std::vector<GridIndex<bulkIdx>>*> circleCornerIndices; |
211 | using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >; | ||
212 |
1/2✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
|
380 | std::vector<ShapeValues> circleShapeValues; |
213 | |||
214 | // go over all points of the average operator | ||
215 |
4/4✓ Branch 0 taken 48640 times.
✓ Branch 1 taken 190 times.
✓ Branch 2 taken 48640 times.
✓ Branch 3 taken 190 times.
|
48830 | for (int k = 0; k < circlePoints.size(); ++k) |
216 | { | ||
217 |
3/6✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 48640 times.
|
97280 | circleBulkElementIndices[k] = intersectingEntities(circlePoints[k], bulkTree); |
218 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 48640 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 48640 times.
|
145920 | if (circleBulkElementIndices[k].empty()) |
219 | ✗ | continue; | |
220 | |||
221 | 48640 | const auto localCircleAvgWeight = circleAvgWeight / circleBulkElementIndices[k].size(); | |
222 |
4/4✓ Branch 0 taken 48640 times.
✓ Branch 1 taken 48640 times.
✓ Branch 2 taken 48640 times.
✓ Branch 3 taken 48640 times.
|
194560 | for (const auto bulkElementIdx : circleBulkElementIndices[k]) |
223 | { | ||
224 |
1/2✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
|
48640 | circleStencil.push_back(bulkElementIdx); |
225 |
1/2✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
|
48640 | circleIpWeight.push_back(localCircleAvgWeight); |
226 | |||
227 | // precompute interpolation data for box scheme for each cut bulk element | ||
228 | if constexpr (isBox<bulkIdx>()) | ||
229 | { | ||
230 | const auto bulkElement = bulkGridGeometry.element(bulkElementIdx); | ||
231 | circleCornerIndices.push_back(&(this->vertexIndices(bulkIdx, bulkElementIdx))); | ||
232 | |||
233 | // evaluate shape functions at the integration point | ||
234 | const auto bulkGeometry = bulkElement.geometry(); | ||
235 | ShapeValues shapeValues; | ||
236 | this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePoints[k], shapeValues); | ||
237 | circleShapeValues.emplace_back(std::move(shapeValues)); | ||
238 | } | ||
239 | } | ||
240 | } | ||
241 | |||
242 | // export low dim circle stencil | ||
243 | if constexpr (isBox<bulkIdx>()) | ||
244 | { | ||
245 | // we insert all vertices and make it unique later | ||
246 | for (const auto& vertices : circleCornerIndices) | ||
247 | { | ||
248 | this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(), | ||
249 | vertices->begin(), vertices->end()); | ||
250 | |||
251 | } | ||
252 | } | ||
253 | else | ||
254 | { | ||
255 |
9/18✓ Branch 1 taken 190 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 190 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 190 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 190 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 190 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 190 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 190 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 190 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 190 times.
✗ Branch 26 not taken.
|
380 | this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(), |
256 | circleStencil.begin(), circleStencil.end()); | ||
257 | } | ||
258 | |||
259 |
4/4✓ Branch 0 taken 48640 times.
✓ Branch 1 taken 190 times.
✓ Branch 2 taken 48640 times.
✓ Branch 3 taken 190 times.
|
48830 | for (int k = 0; k < circlePoints.size(); ++k) |
260 | { | ||
261 |
1/2✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
|
48640 | const auto& circlePos = circlePoints[k]; |
262 | |||
263 | // find bulk elements intersection with the circle elements | ||
264 |
3/6✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 48640 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
|
145920 | if (circleBulkElementIndices[k].empty()) |
265 | continue; | ||
266 | |||
267 | // loop over the bulk elements at the integration points (usually one except when it is on a face or edge or vertex) | ||
268 | // and add a point source at every point on the circle | ||
269 |
4/4✓ Branch 0 taken 48640 times.
✓ Branch 1 taken 48640 times.
✓ Branch 2 taken 48640 times.
✓ Branch 3 taken 48640 times.
|
145920 | for (const auto bulkElementIdx : circleBulkElementIndices[k]) |
270 | { | ||
271 | 48640 | const auto id = this->idCounter_++; | |
272 | |||
273 |
2/4✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
|
97280 | this->pointSources(bulkIdx).emplace_back(circlePos, id, qpweight, integrationElement, bulkElementIdx); |
274 |
5/10✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 48640 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 48640 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 48640 times.
✗ Branch 14 not taken.
|
243200 | this->pointSources(bulkIdx).back().setEmbeddings(circleBulkElementIndices[k].size()); |
275 |
2/4✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
|
97280 | this->pointSources(lowDimIdx).emplace_back(globalPos, id, qpweight, integrationElement, lowDimElementIdx); |
276 |
5/10✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 48640 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 48640 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 48640 times.
✗ Branch 9 not taken.
|
243200 | this->pointSources(lowDimIdx).back().setEmbeddings(circleBulkElementIndices[k].size()); |
277 | |||
278 | // pre compute additional data used for the evaluation of | ||
279 | // the actual solution dependent source term | ||
280 |
1/2✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
|
97280 | PointSourceData psData; |
281 | |||
282 | if constexpr (isBox<lowDimIdx>()) | ||
283 | { | ||
284 | using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >; | ||
285 | ShapeValues shapeValues; | ||
286 | this->getShapeValues(lowDimIdx, lowDimGridGeometry, intersectionGeometry, globalPos, shapeValues); | ||
287 | psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx); | ||
288 | } | ||
289 | else | ||
290 | { | ||
291 |
1/2✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
|
48640 | psData.addLowDimInterpolation(lowDimElementIdx); |
292 | } | ||
293 | |||
294 | // add data needed to compute integral over the circle | ||
295 | if constexpr (isBox<bulkIdx>()) | ||
296 | { | ||
297 | if (useCircleAverage) | ||
298 | psData.addCircleInterpolation(circleCornerIndices, circleShapeValues, circleIpWeight, circleStencil); | ||
299 | |||
300 | using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >; | ||
301 | const auto bulkGeometry = bulkGridGeometry.element(bulkElementIdx).geometry(); | ||
302 | ShapeValues shapeValues; | ||
303 | this->getShapeValues(bulkIdx, bulkGridGeometry, bulkGeometry, circlePos, shapeValues); | ||
304 | psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx); | ||
305 | } | ||
306 | else | ||
307 | { | ||
308 |
1/2✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
|
48640 | if (useCircleAverage) |
309 |
1/2✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
|
48640 | psData.addCircleInterpolation(circleIpWeight, circleStencil); |
310 | |||
311 |
1/2✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
|
48640 | psData.addBulkInterpolation(bulkElementIdx); |
312 | } | ||
313 | |||
314 | // publish point source data in the global vector | ||
315 |
2/4✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
|
97280 | this->pointSourceData().emplace_back(std::move(psData)); |
316 | |||
317 | // export the bulk coupling stencil | ||
318 | // we insert all vertices / elements and make it unique later | ||
319 | if constexpr (isBox<lowDimIdx>()) | ||
320 | { | ||
321 | const auto& vertices = this->vertexIndices(lowDimIdx, lowDimElementIdx); | ||
322 | this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(), | ||
323 | vertices.begin(), vertices.end()); | ||
324 | |||
325 | } | ||
326 | else | ||
327 | { | ||
328 |
3/6✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 48640 times.
✗ Branch 8 not taken.
|
97280 | this->couplingStencils(bulkIdx)[bulkElementIdx].push_back(lowDimElementIdx); |
329 | } | ||
330 | |||
331 | // export bulk circle stencil (only needed for circle average) | ||
332 |
1/2✓ Branch 0 taken 48640 times.
✗ Branch 1 not taken.
|
48640 | if (useCircleAverage) |
333 | { | ||
334 | if constexpr (isBox<bulkIdx>()) | ||
335 | { | ||
336 | // we insert all vertices and make it unique later | ||
337 | for (const auto& vertices : circleCornerIndices) | ||
338 | { | ||
339 | extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(), | ||
340 | vertices->begin(), vertices->end()); | ||
341 | |||
342 | } | ||
343 | } | ||
344 | else | ||
345 | { | ||
346 |
9/18✓ Branch 1 taken 48640 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 48640 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 48640 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 48640 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 48640 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 48640 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 48640 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 48640 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 48640 times.
✗ Branch 26 not taken.
|
97280 | extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(), |
347 | circleStencil.begin(), circleStencil.end()); | ||
348 | } | ||
349 | } | ||
350 | } | ||
351 | } | ||
352 | } | ||
353 | } | ||
354 | |||
355 | // make the circle stencil unique (for source derivatives) | ||
356 |
2/2✓ Branch 0 taken 1820 times.
✓ Branch 1 taken 5 times.
|
1840 | for (auto&& stencil : extendedSourceStencil_.stencil()) |
357 | { | ||
358 | 5460 | std::sort(stencil.second.begin(), stencil.second.end()); | |
359 | 9100 | stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end()); | |
360 | |||
361 | // remove the vertices element (box) | ||
362 | if constexpr (isBox<bulkIdx>()) | ||
363 | { | ||
364 | const auto& indices = this->vertexIndices(bulkIdx, stencil.first); | ||
365 | stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(), | ||
366 | [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }), | ||
367 | stencil.second.end()); | ||
368 | } | ||
369 | // remove the own element (cell-centered) | ||
370 | else | ||
371 | { | ||
372 | 9100 | stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(), | |
373 |
3/10✓ Branch 0 taken 455 times.
✓ Branch 1 taken 5765 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 20330 times.
✗ Branch 9 not taken.
|
26550 | [&](auto i){ return i == stencil.first; }), |
374 | stencil.second.end()); | ||
375 | } | ||
376 | } | ||
377 | |||
378 | // make stencils unique | ||
379 | using namespace Dune::Hybrid; | ||
380 | 45 | forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx) | |
381 | { | ||
382 |
4/4✓ Branch 0 taken 95 times.
✓ Branch 1 taken 5 times.
✓ Branch 3 taken 1820 times.
✓ Branch 4 taken 5 times.
|
1955 | for (auto&& stencil : this->couplingStencils(domainIdx)) |
383 | { | ||
384 | 5745 | std::sort(stencil.second.begin(), stencil.second.end()); | |
385 | 9575 | stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end()); | |
386 | } | ||
387 | }); | ||
388 | |||
389 | 15 | std::cout << "took " << watch.elapsed() << " seconds." << std::endl; | |
390 | 5 | } | |
391 | |||
392 | //! Compute the low dim volume fraction in the bulk domain cells | ||
393 | 5 | void computeLowDimVolumeFractions() | |
394 | { | ||
395 | // resize the storage vector | ||
396 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0)); |
397 | // get references to the grid geometries | ||
398 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
|
5 | const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry(); |
399 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry(); |
400 | |||
401 | // compute the low dim volume fractions | ||
402 |
4/4✓ Branch 0 taken 95 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 95 times.
✓ Branch 3 taken 5 times.
|
110 | for (const auto& is : intersections(this->glue())) |
403 | { | ||
404 | // all inside elements are identical... | ||
405 | 95 | const auto& inside = is.targetEntity(0); | |
406 | 95 | const auto intersectionGeometry = is.geometry(); | |
407 | 190 | const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside); | |
408 | |||
409 | // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional | ||
410 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 95 times.
|
95 | const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx); |
411 |
4/4✓ Branch 0 taken 380 times.
✓ Branch 1 taken 95 times.
✓ Branch 2 taken 380 times.
✓ Branch 3 taken 95 times.
|
855 | for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx) |
412 | { | ||
413 | 380 | const auto& outside = is.domainEntity(outsideIdx); | |
414 | 760 | const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside); | |
415 | 1140 | lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius; | |
416 | } | ||
417 | } | ||
418 | 5 | } | |
419 | |||
420 | /*! | ||
421 | * \brief Methods to be accessed by the subproblems | ||
422 | */ | ||
423 | // \{ | ||
424 | |||
425 | //! Return a reference to the bulk problem | ||
426 | ✗ | Scalar radius(std::size_t id) const | |
427 | { | ||
428 |
4/12✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2009600 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2009600 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 1175040 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 1175040 times.
|
6369280 | const auto& data = this->pointSourceData()[id]; |
429 |
2/6✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2009600 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1175040 times.
|
3184640 | return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx()); |
430 | } | ||
431 | |||
432 | //! The volume the lower dimensional domain occupies in the bulk domain element | ||
433 | // For one-dimensional low dim domain we assume radial tubes | ||
434 | Scalar lowDimVolume(const Element<bulkIdx>& element) const | ||
435 | { | ||
436 | const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element); | ||
437 | return lowDimVolumeInBulkElement_[eIdx]; | ||
438 | } | ||
439 | |||
440 | //! The volume fraction the lower dimensional domain occupies in the bulk domain element | ||
441 | // For one-dimensional low dim domain we assume radial tubes | ||
442 | Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const | ||
443 | { | ||
444 | const auto totalVolume = element.geometry().volume(); | ||
445 | return lowDimVolume(element) / totalVolume; | ||
446 | } | ||
447 | |||
448 | // \} | ||
449 | |||
450 | /*! | ||
451 | * \brief Extended source stencil (for the bulk domain) | ||
452 | */ | ||
453 | const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type& | ||
454 | 601000 | extendedSourceStencil(std::size_t eIdx) const | |
455 | { | ||
456 | 601000 | const auto& sourceStencils = extendedSourceStencil_.stencil(); | |
457 |
4/4✓ Branch 1 taken 3640 times.
✓ Branch 2 taken 597360 times.
✓ Branch 3 taken 3640 times.
✓ Branch 4 taken 597360 times.
|
1803000 | if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end()) |
458 | 7280 | return stencil->second; | |
459 | |||
460 | 1194720 | return this->emptyStencil(bulkIdx); | |
461 | } | ||
462 | |||
463 | private: | ||
464 | //! the extended source stencil object | ||
465 | EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_; | ||
466 | |||
467 | //! vector for the volume fraction of the lowdim domain in the bulk domain cells | ||
468 | std::vector<Scalar> lowDimVolumeInBulkElement_; | ||
469 | }; | ||
470 | |||
471 | //! we support multithreaded assembly | ||
472 | template<class MDTraits> | ||
473 | struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Surface>> | ||
474 | : public std::true_type {}; | ||
475 | |||
476 | } // end namespace Dumux | ||
477 | |||
478 | #endif | ||
479 |