3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
couplingmanager1d3d_kernel.hh
Go to the documentation of this file.
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_KERNEL_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_COUPLINGMANAGER_1D3D_KERNEL_HH
27
28#include <vector>
29
30#include <dune/common/timer.hh>
31#include <dune/geometry/quadraturerules.hh>
32
33#include <dumux/common/tag.hh>
37
41
42namespace Dumux {
43
44namespace Embedded1d3dCouplingMode {
45struct Kernel : public Utility::Tag<Kernel> {
46 static std::string name() { return "kernel"; }
47
48 [[deprecated("Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
49 friend constexpr bool operator==(Kernel, EmbeddedCouplingMode m) { return m == EmbeddedCouplingMode::kernel; }
50 [[deprecated("Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
51 friend constexpr bool operator==(EmbeddedCouplingMode m, Kernel) { return m == EmbeddedCouplingMode::kernel; }
52 [[deprecated("Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
53 friend constexpr bool operator!=(Kernel, EmbeddedCouplingMode m) { return m != EmbeddedCouplingMode::kernel; }
54 [[deprecated("Comparison with enum is deprecated. Removed after 3.3. Use tags.")]]
55 friend constexpr bool operator!=(EmbeddedCouplingMode m, Kernel) { return m != EmbeddedCouplingMode::kernel; }
56};
57
58inline constexpr Kernel kernel{};
59} // end namespace Embedded1d3dCouplingMode
60
61// forward declaration
62template<class MDTraits, class CouplingMode>
63class Embedded1d3dCouplingManager;
64
71template<class MDTraits>
72class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>
73: public EmbeddedCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>>
74{
77 using Scalar = typename MDTraits::Scalar;
78 using SolutionVector = typename MDTraits::SolutionVector;
79 using PointSourceData = typename ParentType::PointSourceTraits::PointSourceData;
80
81 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
82 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
83
84 // the sub domain type aliases
85 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
86 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
87 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
88 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
89 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
90 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
91
92 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
93
94 template<std::size_t id>
95 static constexpr bool isBox()
96 { return GridGeometry<id>::discMethod == DiscretizationMethod::box; }
97
98 static_assert(!isBox<bulkIdx>() && !isBox<lowDimIdx>(), "The kernel coupling method is only implemented for the tpfa method");
99 static_assert(Dune::Capabilities::isCartesian<typename GridView<bulkIdx>::Grid>::v, "The kernel coupling method is only implemented for structured grids");
100
101 enum {
102 bulkDim = GridView<bulkIdx>::dimension,
103 lowDimDim = GridView<lowDimIdx>::dimension,
104 dimWorld = GridView<bulkIdx>::dimensionworld
105 };
106
107 // detect if a class (the spatial params) has a kernelWidthFactor() function
108 template <typename T, typename ...Ts>
109 using VariableKernelWidthDetector = decltype(std::declval<T>().kernelWidthFactor(std::declval<Ts>()...));
110
111 template<class T, typename ...Args>
112 static constexpr bool hasKernelWidthFactor()
113 { return Dune::Std::is_detected<VariableKernelWidthDetector, T, Args...>::value; }
114
115public:
116 static constexpr Embedded1d3dCouplingMode::Kernel couplingMode{};
117
118 using ParentType::ParentType;
119
120 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
121 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
122 const SolutionVector& curSol)
123 {
124 ParentType::init(bulkProblem, lowDimProblem, curSol);
125 computeLowDimVolumeFractions();
126
127 const auto refinement = getParamFromGroup<int>(bulkProblem->paramGroup(), "Grid.Refinement", 0);
128 if (refinement > 0)
129 DUNE_THROW(Dune::NotImplemented, "The current intersection detection may likely fail for refined grids.");
130 }
131
138 template<std::size_t id, class JacobianPattern>
139 void extendJacobianPattern(Dune::index_constant<id> domainI, JacobianPattern& pattern) const
140 {
141 extendedSourceStencil_.extendJacobianPattern(*this, domainI, pattern);
142 }
143
151 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
152 void evalAdditionalDomainDerivatives(Dune::index_constant<i> domainI,
153 const LocalAssemblerI& localAssemblerI,
154 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
155 JacobianMatrixDiagBlock& A,
156 GridVariables& gridVariables)
157 {
158 extendedSourceStencil_.evalAdditionalDomainDerivatives(*this, domainI, localAssemblerI, this->curSol(), A, gridVariables);
159 }
160
161 /* \brief Compute integration point point sources and associated data
162 *
163 * This method uses grid glue to intersect the given grids. Over each intersection
164 * we later need to integrate a source term. This method places point sources
165 * at each quadrature point and provides the point source with the necessary
166 * information to compute integrals (quadrature weight and integration element)
167 * \param order The order of the quadrature rule for integration of sources over an intersection
168 * \param verbose If the point source computation is verbose
169 */
170 void computePointSourceData(std::size_t order = 1, bool verbose = false)
171 {
172 // initialize the maps
173 // do some logging and profiling
174 Dune::Timer watch;
175 std::cout << "[coupling] Initializing the integration point source data structures..." << std::endl;
176
177 // prepare the internal data structures
178 prepareDataStructures_();
179 std::cout << "[coupling] Resized data structures." << std::endl;
180
181 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
182 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
183
184 // generate a bunch of random vectors and values for
185 // Monte-carlo integration on the cylinder defined by line and radius
186 static const double characteristicRelativeLength = getParam<double>("MixedDimension.KernelIntegrationCRL", 0.1);
187 EmbeddedCoupling::CylinderIntegration<Scalar> cylIntegration(characteristicRelativeLength, 1);
188
189 static const bool writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false);
190 if (writeIntegrationPointsToFile)
191 {
192 std::ofstream ipPointFile("kernel_points.log", std::ios::trunc);
193 ipPointFile << "x,y,z\n";
194 std::cout << "[coupling] Initialized kernel_points.log." << std::endl;
195 }
196
197 for (const auto& is : intersections(this->glue()))
198 {
199 // all inside elements are identical...
200 const auto& inside = is.targetEntity(0);
201 // get the intersection geometry for integrating over it
202 const auto intersectionGeometry = is.geometry();
203 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
204
205 // for each intersection integrate kernel and add:
206 // * 1d: a new point source
207 // * 3d: a new kernel volume source
208 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
209 const auto kernelWidthFactor = kernelWidthFactor_(this->problem(lowDimIdx).spatialParams(), lowDimElementIdx);
210 const auto kernelWidth = kernelWidthFactor*radius;
211 const auto a = intersectionGeometry.corner(0);
212 const auto b = intersectionGeometry.corner(1);
213 cylIntegration.setGeometry(a, b, kernelWidth);
214
215 // we can have multiple 3d elements per 1d intersection
216 for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
217 {
218 // compute source id
219 // for each id we have
220 // (1) a point source for the 1d domain
221 // (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)
222 // TODO: i.e. one integration of the kernel should be enough (for each of the weights it's weight/embeddings)
223 // (3) the flux scaling factor for each outside element, i.e. each id
224 const auto id = this->idCounter_++;
225
226 // compute the weights for the bulk volume sources
227 const auto& outside = is.domainEntity(outsideIdx);
228 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
229 const auto surfaceFactor = computeBulkSource(intersectionGeometry, radius, kernelWidth, id, lowDimElementIdx, bulkElementIdx, cylIntegration, is.numDomainNeighbors());
230
231 // place a point source at the intersection center
232 const auto center = intersectionGeometry.center();
233 this->pointSources(lowDimIdx).emplace_back(center, id, surfaceFactor, intersectionGeometry.volume(), lowDimElementIdx);
234 this->pointSources(lowDimIdx).back().setEmbeddings(is.numDomainNeighbors());
235
236 // pre compute additional data used for the evaluation of
237 // the actual solution dependent source term
238 PointSourceData psData;
239
240 // lowdim interpolation (evaluate at center)
241 if constexpr (isBox<lowDimIdx>())
242 {
243 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
244 const auto lowDimGeometry = this->problem(lowDimIdx).gridGeometry().element(lowDimElementIdx).geometry();
245 ShapeValues shapeValues;
246 this->getShapeValues(lowDimIdx, this->problem(lowDimIdx).gridGeometry(), lowDimGeometry, center, shapeValues);
247 psData.addLowDimInterpolation(shapeValues, this->vertexIndices(lowDimIdx, lowDimElementIdx), lowDimElementIdx);
248 }
249 else
250 {
251 psData.addLowDimInterpolation(lowDimElementIdx);
252 }
253
254 // bulk interpolation (evaluate at center)
255 if constexpr (isBox<bulkIdx>())
256 {
257 using ShapeValues = std::vector<Dune::FieldVector<Scalar, 1> >;
258 const auto bulkGeometry = this->problem(bulkIdx).gridGeometry().element(bulkElementIdx).geometry();
259 ShapeValues shapeValues;
260 this->getShapeValues(bulkIdx, this->problem(bulkIdx).gridGeometry(), bulkGeometry, center, shapeValues);
261 psData.addBulkInterpolation(shapeValues, this->vertexIndices(bulkIdx, bulkElementIdx), bulkElementIdx);
262 }
263 else
264 {
265 psData.addBulkInterpolation(bulkElementIdx);
266 }
267
268 // publish point source data in the global vector
269 this->pointSourceData().emplace_back(std::move(psData));
270
271 const auto avgMinDist = averageDistanceSegmentGeometry(a, b, outside.geometry());
272 this->averageDistanceToBulkCell().push_back(avgMinDist);
273 fluxScalingFactor_.push_back(this->problem(bulkIdx).fluxScalingFactor(avgMinDist, radius, kernelWidth));
274
275 // export the lowdim coupling stencil
276 // we insert all vertices / elements and make it unique later
277 if constexpr (isBox<bulkIdx>())
278 {
279 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
280 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
281 vertices.begin(), vertices.end());
282 }
283 else
284 {
285 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
286 }
287 }
288 }
289
290 // make the stencils unique
291 makeUniqueStencil_();
292
293 if (!this->pointSources(bulkIdx).empty())
294 DUNE_THROW(Dune::InvalidStateException, "Kernel method shouldn't have point sources in the bulk domain but only volume sources!");
295
296 std::cout << "[coupling] Finished preparing manager in " << watch.elapsed() << " seconds." << std::endl;
297 }
298
301 {
302 // resize the storage vector
303 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
304 // get references to the grid geometries
305 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
306 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
307
308 // compute the low dim volume fractions
309 for (const auto& is : intersections(this->glue()))
310 {
311 // all inside elements are identical...
312 const auto& inside = is.targetEntity(0);
313 const auto intersectionGeometry = is.geometry();
314 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
315
316 // compute the volume the low-dim domain occupies in the bulk domain if it were full-dimensional
317 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
318 for (int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
319 {
320 const auto& outside = is.domainEntity(outsideIdx);
321 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
322 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
323 }
324 }
325 }
326
330 // \{
331
333 Scalar radius(std::size_t id) const
334 {
335 const auto& data = this->pointSourceData()[id];
336 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
337 }
338
340 // For one-dimensional low dim domain we assume radial tubes
341 Scalar lowDimVolume(const Element<bulkIdx>& element) const
342 {
343 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
344 return lowDimVolumeInBulkElement_[eIdx];
345 }
346
348 // For one-dimensional low dim domain we assume radial tubes
349 Scalar lowDimVolumeFraction(const Element<bulkIdx>& element) const
350 {
351 const auto totalVolume = element.geometry().volume();
352 return lowDimVolume(element) / totalVolume;
353 }
354
356 const std::vector<std::size_t>& bulkSourceIds(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const
357 { return bulkSourceIds_[eIdx][scvIdx]; }
358
360 const std::vector<Scalar>& bulkSourceWeights(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const
361 { return bulkSourceWeights_[eIdx][scvIdx]; }
362
364 Scalar fluxScalingFactor(std::size_t id) const
365 { return fluxScalingFactor_[id]; }
366
367 // \}
368
369private:
371 template<class Line, class CylIntegration>
372 Scalar computeBulkSource(const Line& line, const Scalar radius, const Scalar kernelWidth,
373 std::size_t id, GridIndex<lowDimIdx> lowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx,
374 const CylIntegration& cylIntegration, int embeddings)
375 {
376 // Monte-carlo integration on the cylinder defined by line and radius
377 static const auto bulkParamGroup = this->problem(bulkIdx).paramGroup();
378 static const auto min = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.LowerLeft");
379 static const auto max = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.UpperRight");
380 static const auto cells = getParamFromGroup<std::array<int, bulkDim>>(bulkParamGroup, "Grid.Cells");
381 const auto cylSamples = cylIntegration.size();
382 const auto& a = line.corner(0);
383 const auto& b = line.corner(1);
384
385 // optionally write debugging / visual output of the integration points
386 static const auto writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false);
387 if (writeIntegrationPointsToFile)
388 {
389 std::ofstream ipPointFile("kernel_points.log", std::ios::app);
390 for (int i = 0; i < cylSamples; ++i)
391 {
392 const auto& point = cylIntegration.integrationPoint(i);
393 if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection)
394 ipPointFile << point[0] << "," << point[1] << "," << point[2] << '\n';
395 }
396 }
397
398 Scalar integral = 0.0;
399 for (int i = 0; i < cylSamples; ++i)
400 {
401 const auto& point = cylIntegration.integrationPoint(i);
402 // TODO: below only works for Cartesian grids with ijk numbering (e.g. level 0 YaspGrid (already fails when refined))
403 // more general is the bounding box tree solution which always works, however it's much slower
404 //const auto bulkIndices = intersectingEntities(point, this->problem(bulkIdx).gridGeometry().boundingBoxTree(), true);
405 if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection)
406 {
407 const auto bulkElementIdx = intersectingEntityCartesianGrid(point, min, max, cells);
408 assert(bulkElementIdx < this->problem(bulkIdx).gridGeometry().gridView().size(0));
409 const auto localWeight = evalConstKernel_(a, b, point, radius, kernelWidth)*cylIntegration.integrationElement(i)/Scalar(embeddings);
410 integral += localWeight;
411 if (!bulkSourceIds_[bulkElementIdx][0].empty() && id == bulkSourceIds_[bulkElementIdx][0].back())
412 {
413 bulkSourceWeights_[bulkElementIdx][0].back() += localWeight;
414 }
415 else
416 {
417 bulkSourceIds_[bulkElementIdx][0].emplace_back(id);
418 bulkSourceWeights_[bulkElementIdx][0].emplace_back(localWeight);
419 addBulkSourceStencils_(bulkElementIdx, lowDimElementIdx, coupledBulkElementIdx);
420 }
421 }
422 }
423
424 // the surface factor (which fraction of the source is inside the domain and needs to be considered)
425 const auto length = (a-b).two_norm()/Scalar(embeddings);
426 return integral/length;
427 }
428
429 void prepareDataStructures_()
430 {
431 // clear all internal members like pointsource vectors and stencils
432 // initializes the point source id counter
433 this->clear();
434 bulkSourceIds_.clear();
435 bulkSourceWeights_.clear();
436 extendedSourceStencil_.stencil().clear();
437
438 // precompute the vertex indices for efficiency for the box method
439 this->precomputeVertexIndices(bulkIdx);
440 this->precomputeVertexIndices(lowDimIdx);
441
442 bulkSourceIds_.resize(this->gridView(bulkIdx).size(0));
443 bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0));
444
445 // intersect the bounding box trees
446 this->glueGrids();
447
448 // reserve memory for data
449 this->pointSourceData().reserve(this->glue().size());
450 this->averageDistanceToBulkCell().reserve(this->glue().size());
451 fluxScalingFactor_.reserve(this->glue().size());
452
453 // reserve memory for stencils
454 const auto numBulkElements = this->gridView(bulkIdx).size(0);
455 for (GridIndex<bulkIdx> bulkElementIdx = 0; bulkElementIdx < numBulkElements; ++bulkElementIdx)
456 {
457 this->couplingStencils(bulkIdx)[bulkElementIdx].reserve(10);
458 extendedSourceStencil_.stencil()[bulkElementIdx].reserve(10);
459 bulkSourceIds_[bulkElementIdx][0].reserve(10);
460 bulkSourceWeights_[bulkElementIdx][0].reserve(10);
461 }
462 }
463
465 void makeUniqueStencil_()
466 {
467 // make extra stencils unique
468 for (auto&& stencil : extendedSourceStencil_.stencil())
469 {
470 std::sort(stencil.second.begin(), stencil.second.end());
471 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
472
473 // remove the vertices element (box)
474 if constexpr (isBox<bulkIdx>())
475 {
476 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
477 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
478 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
479 stencil.second.end());
480 }
481 // remove the own element (cell-centered)
482 else
483 {
484 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
485 [&](auto i){ return i == stencil.first; }),
486 stencil.second.end());
487 }
488 }
489
490 // make stencils unique
491 using namespace Dune::Hybrid;
492 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
493 {
494 for (auto&& stencil : this->couplingStencils(domainIdx))
495 {
496 std::sort(stencil.second.begin(), stencil.second.end());
497 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
498 }
499 });
500 }
501
503 void addBulkSourceStencils_(GridIndex<bulkIdx> bulkElementIdx, GridIndex<lowDimIdx> coupledLowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx)
504 {
505 // add the lowdim element to the coupling stencil of this bulk element
506 if constexpr (isBox<lowDimIdx>())
507 {
508 const auto& vertices = this->vertexIndices(lowDimIdx, coupledLowDimElementIdx);
509 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
510 vertices.begin(), vertices.end());
511
512 }
513 else
514 {
515 auto& s = this->couplingStencils(bulkIdx)[bulkElementIdx];
516 s.push_back(coupledLowDimElementIdx);
517 }
518
519 // the extended source stencil, every 3d element with a source is coupled to
520 // the element/dofs where the 3d quantities are measured
521 if constexpr (isBox<bulkIdx>())
522 {
523 const auto& vertices = this->vertexIndices(bulkIdx, coupledBulkElementIdx);
524 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
525 vertices.begin(), vertices.end());
526 }
527 else
528 {
529 auto& s = extendedSourceStencil_.stencil()[bulkElementIdx];
530 s.push_back(coupledBulkElementIdx);
531 }
532 }
533
535 Scalar evalConstKernel_(const GlobalPosition& a,
536 const GlobalPosition& b,
537 const GlobalPosition& point,
538 const Scalar R,
539 const Scalar rho) const noexcept
540 {
541 // projection of point onto line a + t*(b-a)
542 const auto ab = b - a;
543 const auto t = (point - a)*ab/ab.two_norm2();
544
545 // return 0 if we are outside cylinder
546 if (t < 0.0 || t > 1.0)
547 return 0.0;
548
549 // compute distance
550 auto proj = a; proj.axpy(t, ab);
551 const auto r = (proj - point).two_norm();
552
553 if (r > rho)
554 return 0.0;
555
556 return 1.0/(M_PI*rho*rho);
557 }
558
562 template<class SpatialParams>
563 auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx)
564 -> std::enable_if_t<hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
565 { return spatialParams.kernelWidthFactor(eIdx); }
566
570 template<class SpatialParams>
571 auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx)
572 -> std::enable_if_t<!hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
573 {
574 static const Scalar kernelWidthFactor = getParam<Scalar>("MixedDimension.KernelWidthFactor");
575 return kernelWidthFactor;
576 }
577
579 EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_;
581 std::vector<Scalar> lowDimVolumeInBulkElement_;
583 std::vector<std::array<std::vector<std::size_t>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceIds_;
585 std::vector<std::array<std::vector<Scalar>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceWeights_;
586
588 std::vector<Scalar> fluxScalingFactor_;
589};
590
591} // end namespace Dumux
592
593#endif
Defines the index types used for grid and local indices.
Helper class to create (named and comparable) tagged types.
Coupling manager for low-dimensional domains embedded in the bulk domain. Point sources on each integ...
Integration over cylindrical and elliptic cylindrical domains Lowest order integration formulas that ...
Extended source stencil helper class for coupling managers.
Geometry::ctype averageDistanceSegmentGeometry(const typename Geometry::GlobalCoordinate &a, const typename Geometry::GlobalCoordinate &b, const Geometry &geometry, std::size_t integrationOrder=2)
Compute the average distance from a segment to a geometry by integration.
Definition: geometry/distance.hh:121
std::size_t intersectingEntityCartesianGrid(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &min, const Dune::FieldVector< ctype, dimworld > &max, const std::array< int, std::size_t(dimworld)> &cells)
Compute the index of the intersecting element of a Cartesian grid with a point The grid is given by t...
Definition: geometry/intersectingentities.hh:376
EmbeddedCouplingMode
The coupling mode.
Definition: couplingmanager1d3d.hh:44
Definition: adapt.hh:29
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: geometry/boundingboxtree.hh:304
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 Kernel kernel
Definition: couplingmanager1d3d_kernel.hh:58
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_kernel.hh:45
friend constexpr bool operator!=(Kernel, EmbeddedCouplingMode m)
Definition: couplingmanager1d3d_kernel.hh:53
friend constexpr bool operator==(Kernel, EmbeddedCouplingMode m)
Definition: couplingmanager1d3d_kernel.hh:49
friend constexpr bool operator!=(EmbeddedCouplingMode m, Kernel)
Definition: couplingmanager1d3d_kernel.hh:55
static std::string name()
Definition: couplingmanager1d3d_kernel.hh:46
friend constexpr bool operator==(EmbeddedCouplingMode m, Kernel)
Definition: couplingmanager1d3d_kernel.hh:51
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_kernel.hh:74
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_kernel.hh:152
void computeLowDimVolumeFractions()
Compute the low dim volume fraction in the bulk domain cells.
Definition: couplingmanager1d3d_kernel.hh:300
Scalar radius(std::size_t id) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_kernel.hh:333
Scalar lowDimVolume(const Element< bulkIdx > &element) const
The volume the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_kernel.hh:341
const std::vector< Scalar > & bulkSourceWeights(GridIndex< bulkIdx > eIdx, int scvIdx=0) const
return all source ids for a bulk elements
Definition: couplingmanager1d3d_kernel.hh:360
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_kernel.hh:139
const std::vector< std::size_t > & bulkSourceIds(GridIndex< bulkIdx > eIdx, int scvIdx=0) const
return all source ids for a bulk elements
Definition: couplingmanager1d3d_kernel.hh:356
Scalar lowDimVolumeFraction(const Element< bulkIdx > &element) const
The volume fraction the lower dimensional domain occupies in the bulk domain element.
Definition: couplingmanager1d3d_kernel.hh:349
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_kernel.hh:170
Scalar fluxScalingFactor(std::size_t id) const
The flux scaling factor for a source with id.
Definition: couplingmanager1d3d_kernel.hh:364
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_kernel.hh:120
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:78
Helper class to integrate over a cylinder domain.
Definition: cylinderintegration.hh:60
void setGeometry(const GlobalPosition &bottomCenter, const GlobalPosition &topCenter, const Scalar radius, int verbosity=0)
Set the geometry of the cylinder.
Definition: cylinderintegration.hh:97
Helper functions for distance queries.
Declares all properties used in Dumux.