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