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_KERNEL_HH | ||
14 | #define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_KERNEL_HH | ||
15 | |||
16 | #include <vector> | ||
17 | |||
18 | #include <dune/common/timer.hh> | ||
19 | #include <dune/geometry/quadraturerules.hh> | ||
20 | #include <dune/grid/common/capabilities.hh> | ||
21 | |||
22 | #include <dumux/common/tag.hh> | ||
23 | #include <dumux/common/properties.hh> | ||
24 | #include <dumux/common/indextraits.hh> | ||
25 | #include <dumux/geometry/distance.hh> | ||
26 | |||
27 | #include <dumux/multidomain/embedded/couplingmanagerbase.hh> | ||
28 | #include <dumux/multidomain/embedded/cylinderintegration.hh> | ||
29 | #include <dumux/multidomain/embedded/extendedsourcestencil.hh> | ||
30 | |||
31 | namespace Dumux { | ||
32 | |||
33 | namespace Embedded1d3dCouplingMode { | ||
34 | struct Kernel : public Utility::Tag<Kernel> { | ||
35 | static std::string name() { return "kernel"; } | ||
36 | }; | ||
37 | |||
38 | inline constexpr Kernel kernel{}; | ||
39 | } // end namespace Embedded1d3dCouplingMode | ||
40 | |||
41 | // forward declaration | ||
42 | template<class MDTraits, class CouplingMode> | ||
43 | class Embedded1d3dCouplingManager; | ||
44 | |||
45 | /*! | ||
46 | * \ingroup EmbeddedCoupling | ||
47 | * \brief Manages the coupling between bulk elements and lower dimensional elements | ||
48 | * Point sources on each integration point are computed by an AABB tree. | ||
49 | * \note Specialization for coupling method using a distributed kernel source with 3d quantities evaluated on the line | ||
50 | */ | ||
51 | template<class MDTraits> | ||
52 | class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel> | ||
53 | : public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>> | ||
54 | { | ||
55 | using ThisType = Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>; | ||
56 | using ParentType = EmbeddedCouplingManagerBase<MDTraits, ThisType>; | ||
57 | using Scalar = typename MDTraits::Scalar; | ||
58 | using SolutionVector = typename MDTraits::SolutionVector; | ||
59 | using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData; | ||
60 | |||
61 | static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index(); | ||
62 | static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index(); | ||
63 | |||
64 | // the sub domain type aliases | ||
65 | template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag; | ||
66 | template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>; | ||
67 | template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>; | ||
68 | template<std::size_t id> using GridView = typename GridGeometry<id>::GridView; | ||
69 | template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity; | ||
70 | template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex; | ||
71 | |||
72 | using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate; | ||
73 | |||
74 | template<std::size_t id> | ||
75 | static constexpr bool isBox() | ||
76 | { return GridGeometry<id>::discMethod == DiscretizationMethods::box; } | ||
77 | |||
78 | static_assert(!isBox<bulkIdx>() && !isBox<lowDimIdx>(), "The kernel coupling method is only implemented for the tpfa method"); | ||
79 | static_assert(Dune::Capabilities::isCartesian<typename GridView<bulkIdx>::Grid>::v, "The kernel coupling method is only implemented for structured grids"); | ||
80 | |||
81 | enum { | ||
82 | bulkDim = GridView<bulkIdx>::dimension, | ||
83 | lowDimDim = GridView<lowDimIdx>::dimension, | ||
84 | dimWorld = GridView<bulkIdx>::dimensionworld | ||
85 | }; | ||
86 | |||
87 | // detect if a class (the spatial params) has a kernelWidthFactor() function | ||
88 | template <typename T, typename ...Ts> | ||
89 | using VariableKernelWidthDetector = decltype(std::declval<T>().kernelWidthFactor(std::declval<Ts>()...)); | ||
90 | |||
91 | template<class T, typename ...Args> | ||
92 | static constexpr bool hasKernelWidthFactor() | ||
93 | { return Dune::Std::is_detected<VariableKernelWidthDetector, T, Args...>::value; } | ||
94 | |||
95 | public: | ||
96 | static constexpr Embedded1d3dCouplingMode::Kernel couplingMode{}; | ||
97 | |||
98 | using ParentType::ParentType; | ||
99 | |||
100 | 5 | void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem, | |
101 | std::shared_ptr<Problem<lowDimIdx>> lowDimProblem, | ||
102 | const SolutionVector& curSol) | ||
103 | { | ||
104 |
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); |
105 | 5 | computeLowDimVolumeFractions(); | |
106 | |||
107 | 15 | const auto refinement = getParamFromGroup<int>(bulkProblem->paramGroup(), "Grid.Refinement", 0); | |
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | if (refinement > 0) |
109 | ✗ | DUNE_THROW(Dune::NotImplemented, "The current intersection detection may likely fail for refined grids."); | |
110 | 5 | } | |
111 | |||
112 | /*! | ||
113 | * \brief extend the jacobian pattern of the diagonal block of domain i | ||
114 | * by those entries that are not already in the uncoupled pattern | ||
115 | * \note Such additional dependencies can arise from the coupling, e.g. if a coupling source | ||
116 | * term depends on a non-local average of a quantity of the same domain | ||
117 | */ | ||
118 | template<std::size_t id, class JacobianPattern> | ||
119 | ✗ | void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const | |
120 | { | ||
121 | ✗ | extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern); | |
122 | ✗ | } | |
123 | |||
124 | /*! | ||
125 | * \brief evaluate additional derivatives of the element residual of a domain with respect | ||
126 | * to dofs in the same domain that are not in the regular stencil (per default this is not the case) | ||
127 | * \note Such additional dependencies can arise from the coupling, e.g. if a coupling source | ||
128 | * term depends on a non-local average of a quantity of the same domain | ||
129 | * \note This is the same for box and cc | ||
130 | */ | ||
131 | template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables> | ||
132 | ✗ | void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI, | |
133 | const LocalAssemblerI& localAssemblerI, | ||
134 | const typename LocalAssemblerI::LocalResidual::ElementResidualVector&, | ||
135 | JacobianMatrixDiagBlock& A, | ||
136 | GridVariables& gridVariables) | ||
137 | { | ||
138 | 300595 | extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, A, gridVariables); | |
139 | ✗ | } | |
140 | |||
141 | /* \brief Compute integration point point sources and associated data | ||
142 | * | ||
143 | * This method uses grid glue to intersect the given grids. Over each intersection | ||
144 | * we later need to integrate a source term. This method places point sources | ||
145 | * at each quadrature point and provides the point source with the necessary | ||
146 | * information to compute integrals (quadrature weight and integration element) | ||
147 | * \param order The order of the quadrature rule for integration of sources over an intersection | ||
148 | * \param verbose If the point source computation is verbose | ||
149 | */ | ||
150 | 5 | void computePointSourceData(std::size_t order = 1, bool verbose = false) | |
151 | { | ||
152 | // initialize the maps | ||
153 | // do some logging and profiling | ||
154 | 5 | Dune::Timer watch; | |
155 | 10 | std::cout << "[coupling] Initializing the integration point source data structures..." << std::endl; | |
156 | |||
157 | // prepare the internal data structures | ||
158 | 5 | prepareDataStructures_(); | |
159 | 10 | std::cout << "[coupling] Resized data structures." << std::endl; | |
160 | |||
161 |
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(); |
162 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 5 times.
✗ Branch 3 not taken.
|
5 | const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry(); |
163 | |||
164 | // generate a bunch of random vectors and values for | ||
165 | // Monte-carlo integration on the cylinder defined by line and radius | ||
166 |
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 double characteristicRelativeLength = getParam<double>("MixedDimension.KernelIntegrationCRL", 0.1); |
167 | 5 | EmbeddedCoupling::CylinderIntegration<Scalar> cylIntegration(characteristicRelativeLength, 1); | |
168 | |||
169 |
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 writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false); |
170 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | if (writeIntegrationPointsToFile) |
171 | { | ||
172 | ✗ | std::ofstream ipPointFile("kernel_points.log", std::ios::trunc); | |
173 | ✗ | ipPointFile << "x,y,z\n"; | |
174 | ✗ | std::cout << "[coupling] Initialized kernel_points.log." << std::endl; | |
175 | } | ||
176 | |||
177 |
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())) |
178 | { | ||
179 | // all inside elements are identical... | ||
180 | 95 | const auto& inside = is.targetEntity(0); | |
181 | // get the intersection geometry for integrating over it | ||
182 |
1/2✓ Branch 1 taken 95 times.
✗ Branch 2 not taken.
|
95 | const auto intersectionGeometry = is.geometry(); |
183 | 190 | const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside); | |
184 | |||
185 | // for each intersection integrate kernel and add: | ||
186 | // * 1d: a new point source | ||
187 | // * 3d: a new kernel volume source | ||
188 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 95 times.
✓ Branch 3 taken 95 times.
✗ Branch 4 not taken.
|
95 | const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx); |
189 |
3/6✓ Branch 1 taken 95 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 95 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 95 times.
✗ Branch 8 not taken.
|
285 | const auto kernelWidthFactor = kernelWidthFactor_(this->problem(lowDimIdx).spatialParams(), lowDimElementIdx); |
190 | 95 | const auto kernelWidth = kernelWidthFactor*radius; | |
191 | 95 | const auto a = intersectionGeometry.corner(0); | |
192 | 95 | const auto b = intersectionGeometry.corner(1); | |
193 |
1/2✓ Branch 1 taken 95 times.
✗ Branch 2 not taken.
|
95 | cylIntegration.setGeometry(a, b, kernelWidth); |
194 | |||
195 | // we can have multiple 3d elements per 1d intersection | ||
196 |
4/4✓ Branch 0 taken 380 times.
✓ Branch 1 taken 95 times.
✓ Branch 2 taken 380 times.
✓ Branch 3 taken 95 times.
|
570 | for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx) |
197 | { | ||
198 | // compute source id | ||
199 | // for each id we have | ||
200 | // (1) a point source for the 1d domain | ||
201 | // (2) a bulk source weight for the each element in the 3d domain (the fraction of the total source/sink for the 1d element with the point source) | ||
202 | // TODO: i.e. one integration of the kernel should be enough (for each of the weights it's weight/embeddings) | ||
203 | // (3) the flux scaling factor for each outside element, i.e. each id | ||
204 | 380 | const auto id = this->idCounter_++; | |
205 | |||
206 | // compute the weights for the bulk volume sources | ||
207 |
1/2✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
|
380 | const auto& outside = is.domainEntity(outsideIdx); |
208 | 760 | const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside); | |
209 |
2/4✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
|
760 | const auto surfaceFactor = computeBulkSource(intersectionGeometry, radius, kernelWidth, id, lowDimElementIdx, bulkElementIdx, cylIntegration, is.numDomainNeighbors()); |
210 | |||
211 | // place a point source at the intersection center | ||
212 | 380 | const auto center = intersectionGeometry.center(); | |
213 |
3/6✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 380 times.
✗ Branch 8 not taken.
|
1140 | this->pointSources(lowDimIdx).emplace_back(center, id, surfaceFactor, intersectionGeometry.volume(), lowDimElementIdx); |
214 |
4/8✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 380 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 380 times.
✗ Branch 11 not taken.
|
1520 | this->pointSources(lowDimIdx).back().setEmbeddings(is.numDomainNeighbors()); |
215 | |||
216 | // pre compute additional data used for the evaluation of | ||
217 | // the actual solution dependent source term | ||
218 |
1/2✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
|
760 | PointSourceData psData; |
219 | |||
220 | // lowdim interpolation (evaluate at center) | ||
221 | if constexpr (isBox<lowDimIdx>()) | ||
222 | { | ||
223 | using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >; | ||
224 | const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry(); | ||
225 | ShapeValues shapeValues; | ||
226 | this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, center, shapeValues); | ||
227 | psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx); | ||
228 | } | ||
229 | else | ||
230 | { | ||
231 |
1/2✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
|
380 | psData.addLowDimInterpolation(lowDimElementIdx); |
232 | } | ||
233 | |||
234 | // bulk interpolation (evaluate at center) | ||
235 | if constexpr (isBox<bulkIdx>()) | ||
236 | { | ||
237 | using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >; | ||
238 | const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry(); | ||
239 | ShapeValues shapeValues; | ||
240 | this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, center, shapeValues); | ||
241 | psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx); | ||
242 | } | ||
243 | else | ||
244 | { | ||
245 |
1/2✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
|
380 | psData.addBulkInterpolation(bulkElementIdx); |
246 | } | ||
247 | |||
248 | // publish point source data in the global vector | ||
249 |
2/4✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
|
760 | this->pointSourceData().emplace_back(std::move(psData)); |
250 | |||
251 |
1/2✓ Branch 2 taken 380 times.
✗ Branch 3 not taken.
|
380 | const auto avgMinDist = averageDistanceSegmentGeometry(a, b, outside.geometry()); |
252 |
2/4✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
|
760 | this->averageDistanceToBulkCell().push_back(avgMinDist); |
253 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 380 times.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
|
380 | fluxScalingFactor_.push_back(this->problem(bulkIdx).fluxScalingFactor(avgMinDist, radius, kernelWidth)); |
254 | |||
255 | // export the lowdim coupling stencil | ||
256 | // we insert all vertices / elements and make it unique later | ||
257 | if constexpr (isBox<bulkIdx>()) | ||
258 | { | ||
259 | const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx); | ||
260 | this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(), | ||
261 | vertices.begin(), vertices.end()); | ||
262 | } | ||
263 | else | ||
264 | { | ||
265 |
3/6✓ Branch 1 taken 380 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 380 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 380 times.
✗ Branch 8 not taken.
|
760 | this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx); |
266 | } | ||
267 | } | ||
268 | } | ||
269 | |||
270 | // make the stencils unique | ||
271 | 5 | makeUniqueStencil_(); | |
272 | |||
273 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 5 times.
|
15 | if (!this->pointSources(bulkIdx).empty()) |
274 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Kernel method shouldn't have point sources in the bulk domain but only volume sources!"); | |
275 | |||
276 |
2/4✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 7 taken 5 times.
✗ Branch 8 not taken.
|
15 | std::cout << "[coupling] Finished preparing manager in " << watch.elapsed() << " seconds." << std::endl; |
277 | 5 | } | |
278 | |||
279 | //! Compute the low dim volume fraction in the bulk domain cells | ||
280 | 5 | void computeLowDimVolumeFractions() | |
281 | { | ||
282 | // resize the storage vector | ||
283 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0)); |
284 | // get references to the grid geometries | ||
285 |
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(); |
286 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry(); |
287 | |||
288 | // compute the low dim volume fractions | ||
289 |
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())) |
290 | { | ||
291 | // all inside elements are identical... | ||
292 | 95 | const auto& inside = is.targetEntity(0); | |
293 | 95 | const auto intersectionGeometry = is.geometry(); | |
294 | 190 | const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside); | |
295 | |||
296 | // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional | ||
297 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 95 times.
|
95 | const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx); |
298 |
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) |
299 | { | ||
300 | 380 | const auto& outside = is.domainEntity(outsideIdx); | |
301 | 760 | const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside); | |
302 | 1140 | lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius; | |
303 | } | ||
304 | } | ||
305 | 5 | } | |
306 | |||
307 | /*! | ||
308 | * \brief Methods to be accessed by the subproblems | ||
309 | */ | ||
310 | // \{ | ||
311 | |||
312 | //! Return a reference to the bulk problem | ||
313 | ✗ | Scalar radius(std::size_t id) const | |
314 | { | ||
315 |
4/12✗ Branch 0 not taken.
✓ Branch 1 taken 4180 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4180 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 138480 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 138480 times.
|
285320 | const auto& data = this->pointSourceData()[id]; |
316 |
6/10✗ Branch 0 not taken.
✓ Branch 1 taken 4180 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 138480 times.
✓ Branch 6 taken 5 times.
✓ Branch 7 taken 138475 times.
✓ Branch 8 taken 5 times.
✓ Branch 9 taken 138475 times.
|
142660 | return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx()); |
317 | } | ||
318 | |||
319 | //! The volume the lower dimensional domain occupies in the bulk domain element | ||
320 | // For one-dimensional low dim domain we assume radial tubes | ||
321 | Scalar lowDimVolume(const Element<bulkIdx>& element) const | ||
322 | { | ||
323 | const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element); | ||
324 | return lowDimVolumeInBulkElement_[eIdx]; | ||
325 | } | ||
326 | |||
327 | //! The volume fraction the lower dimensional domain occupies in the bulk domain element | ||
328 | // For one-dimensional low dim domain we assume radial tubes | ||
329 | Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const | ||
330 | { | ||
331 | const auto totalVolume = element.geometry().volume(); | ||
332 | return lowDimVolume(element) / totalVolume; | ||
333 | } | ||
334 | |||
335 | //! return all source ids for a bulk elements | ||
336 | const std::vector<std::size_t>& bulkSourceIds(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const | ||
337 | 2776860 | { return bulkSourceIds_[eIdx][scvIdx]; } | |
338 | |||
339 | //! return all source ids for a bulk elements | ||
340 | const std::vector<Scalar>& bulkSourceWeights(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const | ||
341 | 2776860 | { return bulkSourceWeights_[eIdx][scvIdx]; } | |
342 | |||
343 | //! The flux scaling factor for a source with id | ||
344 | Scalar fluxScalingFactor(std::size_t id) const | ||
345 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 138480 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 138480 times.
|
276960 | { return fluxScalingFactor_[id]; } |
346 | |||
347 | // \} | ||
348 | |||
349 | /*! | ||
350 | * \brief Extended source stencil (for the bulk domain) | ||
351 | */ | ||
352 | const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type& | ||
353 | 601000 | extendedSourceStencil(std::size_t eIdx) const | |
354 | { | ||
355 | 601000 | const auto& sourceStencils = extendedSourceStencil_.stencil(); | |
356 |
3/6✓ Branch 1 taken 601000 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 601000 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 601000 times.
✗ Branch 6 not taken.
|
2404000 | if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end()) |
357 | 1202000 | return stencil->second; | |
358 | |||
359 | ✗ | return this->emptyStencil(bulkIdx); | |
360 | } | ||
361 | |||
362 | private: | ||
363 | //! compute the kernel distributed sources and add stencils | ||
364 | template<class Line, class CylIntegration> | ||
365 | 380 | Scalar computeBulkSource(const Line& line, const Scalar radius, const Scalar kernelWidth, | |
366 | std::size_t id, GridIndex<lowDimIdx> lowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx, | ||
367 | const CylIntegration& cylIntegration, int embeddings) | ||
368 | { | ||
369 | // Monte-carlo integration on the cylinder defined by line and radius | ||
370 |
6/10✓ Branch 0 taken 5 times.
✓ Branch 1 taken 375 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 6 taken 5 times.
✓ Branch 8 taken 5 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 5 times.
✗ Branch 12 not taken.
|
380 | static const auto bulkParamGroup = this->problem(bulkIdx).paramGroup(); |
371 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 375 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
380 | static const auto min = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.LowerLeft"); |
372 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 375 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
380 | static const auto max = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.UpperRight"); |
373 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 375 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
380 | static const auto cells = getParamFromGroup<std::array<int, bulkDim>>(bulkParamGroup, "Grid.Cells"); |
374 | 380 | const auto cylSamples = cylIntegration.size(); | |
375 | 380 | const auto& a = line.corner(0); | |
376 | 380 | const auto& b = line.corner(1); | |
377 | |||
378 | // optionally write debugging / visual output of the integration points | ||
379 |
4/6✓ Branch 0 taken 5 times.
✓ Branch 1 taken 375 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
|
380 | static const auto writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false); |
380 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 380 times.
|
380 | if (writeIntegrationPointsToFile) |
381 | { | ||
382 | ✗ | std::ofstream ipPointFile("kernel_points.log", std::ios::app); | |
383 | ✗ | for (int i = 0; i < cylSamples; ++i) | |
384 | { | ||
385 | ✗ | const auto& point = cylIntegration.integrationPoint(i); | |
386 | ✗ | if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection) | |
387 | ✗ | ipPointFile << point[0] << "," << point[1] << "," << point[2] << '\n'; | |
388 | } | ||
389 | } | ||
390 | |||
391 | Scalar integral = 0.0; | ||
392 |
2/2✓ Branch 0 taken 11918700 times.
✓ Branch 1 taken 380 times.
|
11919080 | for (int i = 0; i < cylSamples; ++i) |
393 | { | ||
394 | 11918700 | const auto& point = cylIntegration.integrationPoint(i); | |
395 | // TODO: below only works for Cartesian grids with ijk numbering (e.g. level 0 YaspGrid (already fails when refined)) | ||
396 | // more general is the bounding box tree solution which always works, however it's much slower | ||
397 | //const auto bulkIndices = intersectingEntities(point, this->problem(bulkIdx).gridGeometry().boundingBoxTree(), true); | ||
398 |
1/2✓ Branch 1 taken 11918700 times.
✗ Branch 2 not taken.
|
11918700 | if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection) |
399 | { | ||
400 | 11918700 | const auto bulkElementIdx = intersectingEntityCartesianGrid(point, min, max, cells); | |
401 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 11918700 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 11918700 times.
|
11918700 | assert(bulkElementIdx < this->problem(bulkIdx).gridGeometry().gridView().size(0)); |
402 |
2/2✓ Branch 1 taken 11915200 times.
✓ Branch 2 taken 3500 times.
|
11918700 | const auto localWeight = evalConstKernel_(a, b, point, radius, kernelWidth)*cylIntegration.integrationElement(i)/Scalar(embeddings); |
403 | 11918700 | integral += localWeight; | |
404 |
14/14✓ Branch 0 taken 11915200 times.
✓ Branch 1 taken 3500 times.
✓ Branch 2 taken 11915200 times.
✓ Branch 3 taken 3500 times.
✓ Branch 4 taken 11915200 times.
✓ Branch 5 taken 3500 times.
✓ Branch 6 taken 11915200 times.
✓ Branch 7 taken 3500 times.
✓ Branch 8 taken 11904700 times.
✓ Branch 9 taken 10500 times.
✓ Branch 10 taken 11904700 times.
✓ Branch 11 taken 10500 times.
✓ Branch 12 taken 11904700 times.
✓ Branch 13 taken 10500 times.
|
47674800 | if (!bulkSourceIds_[bulkElementIdx][0].empty() && id == bulkSourceIds_[bulkElementIdx][0].back()) |
405 | { | ||
406 | 47618800 | bulkSourceWeights_[bulkElementIdx][0].back() += localWeight; | |
407 | } | ||
408 | else | ||
409 | { | ||
410 | 28000 | bulkSourceIds_[bulkElementIdx][0].emplace_back(id); | |
411 | 42000 | bulkSourceWeights_[bulkElementIdx][0].emplace_back(localWeight); | |
412 | 14000 | addBulkSourceStencils_(bulkElementIdx, lowDimElementIdx, coupledBulkElementIdx); | |
413 | } | ||
414 | } | ||
415 | } | ||
416 | |||
417 | // the surface factor (which fraction of the source is inside the domain and needs to be considered) | ||
418 | 380 | const auto length = (a-b).two_norm()/Scalar(embeddings); | |
419 | 380 | return integral/length; | |
420 | } | ||
421 | |||
422 | 5 | void prepareDataStructures_() | |
423 | { | ||
424 | // clear all internal members like pointsource vectors and stencils | ||
425 | // initializes the point source id counter | ||
426 | 5 | this->clear(); | |
427 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | bulkSourceIds_.clear(); |
428 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | bulkSourceWeights_.clear(); |
429 | 10 | extendedSourceStencil_.stencil().clear(); | |
430 | |||
431 | // precompute the vertex indices for efficiency for the box method | ||
432 | 5 | this->precomputeVertexIndices(bulkIdx); | |
433 | 5 | this->precomputeVertexIndices(lowDimIdx); | |
434 | |||
435 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | bulkSourceIds_.resize(this->gridView(bulkIdx).size(0)); |
436 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0)); |
437 | |||
438 | // intersect the bounding box trees | ||
439 | 5 | this->glueGrids(); | |
440 | |||
441 | // reserve memory for data | ||
442 | 20 | this->pointSourceData().reserve(this->glue().size()); | |
443 | 20 | this->averageDistanceToBulkCell().reserve(this->glue().size()); | |
444 | 15 | fluxScalingFactor_.reserve(this->glue().size()); | |
445 | |||
446 | // reserve memory for stencils | ||
447 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
5 | const auto numBulkElements = this->gridView(bulkIdx).size(0); |
448 |
2/2✓ Branch 0 taken 300500 times.
✓ Branch 1 taken 5 times.
|
300505 | for (GridIndex<bulkIdx> bulkElementIdx = 0; bulkElementIdx < numBulkElements; ++bulkElementIdx) |
449 | { | ||
450 | 601000 | this->couplingStencils(bulkIdx)[bulkElementIdx].reserve(10); | |
451 | 601000 | extendedSourceStencil_.stencil()[bulkElementIdx].reserve(10); | |
452 | 901500 | bulkSourceIds_[bulkElementIdx][0].reserve(10); | |
453 | 901500 | bulkSourceWeights_[bulkElementIdx][0].reserve(10); | |
454 | } | ||
455 | 5 | } | |
456 | |||
457 | //! Make the stencils unique | ||
458 | 5 | void makeUniqueStencil_() | |
459 | { | ||
460 | // make extra stencils unique | ||
461 |
2/2✓ Branch 0 taken 300500 times.
✓ Branch 1 taken 5 times.
|
300520 | for (auto&& stencil : extendedSourceStencil_.stencil()) |
462 | { | ||
463 | 901500 | std::sort(stencil.second.begin(), stencil.second.end()); | |
464 | 1502500 | stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end()); | |
465 | |||
466 | // remove the vertices element (box) | ||
467 | if constexpr (isBox<bulkIdx>()) | ||
468 | { | ||
469 | const auto& indices = this->vertexIndices(bulkIdx, stencil.first); | ||
470 | stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(), | ||
471 | [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }), | ||
472 | stencil.second.end()); | ||
473 | } | ||
474 | // remove the own element (cell-centered) | ||
475 | else | ||
476 | { | ||
477 | 1502500 | stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(), | |
478 |
3/10✓ Branch 0 taken 95 times.
✓ Branch 1 taken 3405 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 570 times.
✗ Branch 9 not taken.
|
4070 | [&](auto i){ return i == stencil.first; }), |
479 | stencil.second.end()); | ||
480 | } | ||
481 | } | ||
482 | |||
483 | // make stencils unique | ||
484 | using namespace Dune::Hybrid; | ||
485 | 45 | forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx) | |
486 | { | ||
487 |
4/4✓ Branch 0 taken 95 times.
✓ Branch 1 taken 5 times.
✓ Branch 3 taken 300500 times.
✓ Branch 4 taken 5 times.
|
300635 | for (auto&& stencil : this->couplingStencils(domainIdx)) |
488 | { | ||
489 | 901785 | std::sort(stencil.second.begin(), stencil.second.end()); | |
490 | 1502975 | stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end()); | |
491 | } | ||
492 | }); | ||
493 | 5 | } | |
494 | |||
495 | //! add additional stencil entries for the bulk element | ||
496 | 14000 | void addBulkSourceStencils_(GridIndex<bulkIdx> bulkElementIdx, GridIndex<lowDimIdx> coupledLowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx) | |
497 | { | ||
498 | // add the lowdim element to the coupling stencil of this bulk element | ||
499 | if constexpr (isBox<lowDimIdx>()) | ||
500 | { | ||
501 | const auto& vertices = this->vertexIndices(lowDimIdx, coupledLowDimElementIdx); | ||
502 | this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(), | ||
503 | vertices.begin(), vertices.end()); | ||
504 | |||
505 | } | ||
506 | else | ||
507 | { | ||
508 | 28000 | auto& s = this->couplingStencils(bulkIdx)[bulkElementIdx]; | |
509 | 14000 | s.push_back(coupledLowDimElementIdx); | |
510 | } | ||
511 | |||
512 | // the extended source stencil, every 3d element with a source is coupled to | ||
513 | // the element/dofs where the 3d quantities are measured | ||
514 | if constexpr (isBox<bulkIdx>()) | ||
515 | { | ||
516 | const auto& vertices = this->vertexIndices(bulkIdx, coupledBulkElementIdx); | ||
517 | extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(), | ||
518 | vertices.begin(), vertices.end()); | ||
519 | } | ||
520 | else | ||
521 | { | ||
522 | 28000 | auto& s = extendedSourceStencil_.stencil()[bulkElementIdx]; | |
523 | 14000 | s.push_back(coupledBulkElementIdx); | |
524 | } | ||
525 | 14000 | } | |
526 | |||
527 | //! a cylindrical kernel around the segment a->b | ||
528 | 11918700 | Scalar evalConstKernel_(const GlobalPosition& a, | |
529 | const GlobalPosition& b, | ||
530 | const GlobalPosition& point, | ||
531 | const Scalar R, | ||
532 | const Scalar rho) const noexcept | ||
533 | { | ||
534 | // projection of point onto line a + t*(b-a) | ||
535 | 11918700 | const auto ab = b - a; | |
536 | 11918700 | const auto t = (point - a)*ab/ab.two_norm2(); | |
537 | |||
538 | // return 0 if we are outside cylinder | ||
539 |
2/4✓ Branch 0 taken 11918700 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 11918700 times.
✗ Branch 3 not taken.
|
11918700 | if (t < 0.0 || t > 1.0) |
540 | return 0.0; | ||
541 | |||
542 | // compute distance | ||
543 | 11918700 | auto proj = a; proj.axpy(t, ab); | |
544 | 11918700 | const auto radiusSquared = (proj - point).two_norm2(); | |
545 | |||
546 |
1/2✓ Branch 0 taken 11918700 times.
✗ Branch 1 not taken.
|
11918700 | if (radiusSquared > rho*rho) |
547 | return 0.0; | ||
548 | |||
549 | 11918700 | return 1.0/(M_PI*rho*rho); | |
550 | } | ||
551 | |||
552 | /*! | ||
553 | * \brief Get the kernel width factor from the spatial params (if possible) | ||
554 | */ | ||
555 | template<class SpatialParams> | ||
556 | auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx) | ||
557 | -> std::enable_if_t<hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar> | ||
558 | { return spatialParams.kernelWidthFactor(eIdx); } | ||
559 | |||
560 | /*! | ||
561 | * \brief Get the kernel width factor (constant) from the input file (if possible) | ||
562 | */ | ||
563 | template<class SpatialParams> | ||
564 | ✗ | auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx) | |
565 | -> std::enable_if_t<!hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar> | ||
566 | { | ||
567 | ✗ | static const Scalar kernelWidthFactor = getParam<Scalar>("MixedDimension.KernelWidthFactor"); | |
568 | ✗ | return kernelWidthFactor; | |
569 | } | ||
570 | |||
571 | //! the extended source stencil object for kernel coupling | ||
572 | EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_; | ||
573 | //! vector for the volume fraction of the lowdim domain in the bulk domain cells | ||
574 | std::vector<Scalar> lowDimVolumeInBulkElement_; | ||
575 | //! kernel sources to integrate for each bulk element | ||
576 | std::vector<std::array<std::vector<std::size_t>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceIds_; | ||
577 | //! the integral of the kernel for each point source / integration point, i.e. weight for the source | ||
578 | std::vector<std::array<std::vector<Scalar>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceWeights_; | ||
579 | |||
580 | //! the flux scaling factor for the respective source id | ||
581 | std::vector<Scalar> fluxScalingFactor_; | ||
582 | }; | ||
583 | |||
584 | //! we support multithreaded assembly | ||
585 | template<class MDTraits> | ||
586 | struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>> | ||
587 | : public std::true_type {}; | ||
588 | |||
589 | } // end namespace Dumux | ||
590 | |||
591 | #endif | ||
592 |