version 3.11-dev
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
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>
26
30
31namespace Dumux {
32
33namespace Embedded1d3dCouplingMode {
34struct Kernel : public Utility::Tag<Kernel> {
35 static std::string name() { return "kernel"; }
36};
37
38inline constexpr Kernel kernel{};
39} // end namespace Embedded1d3dCouplingMode
40
41// forward declaration
42template<class MDTraits, class CouplingMode>
43class Embedded1d3dCouplingManager;
44
51template<class MDTraits>
52class Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>
53: public Embedded1d3dCouplingManagerBase<MDTraits, Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>>
54{
58 using Scalar = typename MDTraits::Scalar;
59 using SolutionVector = typename MDTraits::SolutionVector;
60 using PointSourceData = typename BaseType::PointSourceTraits::PointSourceData;
61
62 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
63 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
64
65 // the sub domain type aliases
66 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
67 template<std::size_t id> using Problem = GetPropType<SubDomainTypeTag<id>, Properties::Problem>;
68 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
69 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
70 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
71 template<std::size_t id> using GridIndex = typename IndexTraits<GridView<id>>::GridIndex;
72
73 using GlobalPosition = typename Element<bulkIdx>::Geometry::GlobalCoordinate;
74
75 template<std::size_t id>
76 static constexpr bool isBox()
77 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
78
79 static_assert(!isBox<bulkIdx>() && !isBox<lowDimIdx>(), "The kernel coupling method is only implemented for the tpfa method");
80 static_assert(Dune::Capabilities::isCartesian<typename GridView<bulkIdx>::Grid>::v, "The kernel coupling method is only implemented for structured grids");
81
82 enum {
83 bulkDim = GridView<bulkIdx>::dimension,
84 lowDimDim = GridView<lowDimIdx>::dimension,
85 dimWorld = GridView<bulkIdx>::dimensionworld
86 };
87
88 // detect if a class (the spatial params) has a kernelWidthFactor() function
89 template <typename T, typename ...Ts>
90 using VariableKernelWidthDetector = decltype(std::declval<T>().kernelWidthFactor(std::declval<Ts>()...));
91
92 template<class T, typename ...Args>
93 static constexpr bool hasKernelWidthFactor()
94 { return Dune::Std::is_detected<VariableKernelWidthDetector, T, Args...>::value; }
95
96public:
97 static constexpr Embedded1d3dCouplingMode::Kernel couplingMode{};
98
99 using ParentType::ParentType;
100
101 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
102 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
103 const SolutionVector& curSol)
104 {
105 ParentType::init(bulkProblem, lowDimProblem, curSol);
106
107 const auto refinement = getParamFromGroup<int>(bulkProblem->paramGroup(), "Grid.Refinement", 0);
108 if (refinement > 0)
109 DUNE_THROW(Dune::NotImplemented, "The current intersection detection may likely fail for refined grids.");
110 }
111
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
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 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 void computePointSourceData(std::size_t order = 1, bool verbose = false)
151 {
152 // initialize the maps
153 // do some logging and profiling
154 Dune::Timer watch;
155 std::cout << "[coupling] Initializing the integration point source data structures..." << std::endl;
156
157 // prepare the internal data structures
158 prepareDataStructures_();
159 std::cout << "[coupling] Resized data structures." << std::endl;
160
161 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
162 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 static const double characteristicRelativeLength = getParam<double>("MixedDimension.KernelIntegrationCRL", 0.1);
167 EmbeddedCoupling::CylinderIntegration<Scalar> cylIntegration(characteristicRelativeLength, 1);
168
169 static const bool writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false);
170 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 for (const auto& is : intersections(this->glue()))
178 {
179 // all inside elements are identical...
180 const auto& inside = is.targetEntity(0);
181 // get the intersection geometry for integrating over it
182 const auto intersectionGeometry = is.geometry();
183 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 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
189 const auto kernelWidthFactor = kernelWidthFactor_(this->problem(lowDimIdx).spatialParams(), lowDimElementIdx);
190 const auto kernelWidth = kernelWidthFactor*radius;
191 const auto a = intersectionGeometry.corner(0);
192 const auto b = intersectionGeometry.corner(1);
193 cylIntegration.setGeometry(a, b, kernelWidth);
194
195 // we can have multiple 3d elements per 1d intersection
196 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 const auto id = this->idCounter_++;
205
206 // compute the weights for the bulk volume sources
207 const auto& outside = is.domainEntity(outsideIdx);
208 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
209 const auto surfaceFactor = computeBulkSource(intersectionGeometry, radius, kernelWidth, id, lowDimElementIdx, bulkElementIdx, cylIntegration, is.numDomainNeighbors());
210
211 // place a point source at the intersection center
212 const auto center = intersectionGeometry.center();
213 this->pointSources(lowDimIdx).emplace_back(center, id, surfaceFactor, intersectionGeometry.volume(), lowDimElementIdx);
214 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 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 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 psData.addBulkInterpolation(bulkElementIdx);
246 }
247
248 // publish point source data in the global vector
249 this->pointSourceData().emplace_back(std::move(psData));
250
251 const auto avgMinDist = averageDistanceSegmentGeometry(a, b, outside.geometry());
252 this->averageDistanceToBulkCell().push_back(avgMinDist);
253 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 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
266 }
267 }
268 }
269
270 // make the stencils unique
271 makeUniqueStencil_();
272
273 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 std::cout << "[coupling] Finished preparing manager in " << watch.elapsed() << " seconds." << std::endl;
277 }
278
282 // \{
283
284
286 const std::vector<std::size_t>& bulkSourceIds(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const
287 { return bulkSourceIds_[eIdx][scvIdx]; }
288
290 const std::vector<Scalar>& bulkSourceWeights(GridIndex<bulkIdx> eIdx, int scvIdx = 0) const
291 { return bulkSourceWeights_[eIdx][scvIdx]; }
292
294 Scalar fluxScalingFactor(std::size_t id) const
295 { return fluxScalingFactor_[id]; }
296
297 // \}
298
302 const typename BaseType::template CouplingStencils<bulkIdx>::mapped_type&
303 extendedSourceStencil(std::size_t eIdx) const
304 {
305 const auto& sourceStencils = extendedSourceStencil_.stencil();
306 if (auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
307 return stencil->second;
308
309 return this->emptyStencil(bulkIdx);
310 }
311
312private:
314 template<class Line, class CylIntegration>
315 Scalar computeBulkSource(const Line& line, const Scalar radius, const Scalar kernelWidth,
316 std::size_t id, GridIndex<lowDimIdx> lowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx,
317 const CylIntegration& cylIntegration, int embeddings)
318 {
319 // Monte-carlo integration on the cylinder defined by line and radius
320 static const auto bulkParamGroup = this->problem(bulkIdx).paramGroup();
321 static const auto min = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.LowerLeft");
322 static const auto max = getParamFromGroup<GlobalPosition>(bulkParamGroup, "Grid.UpperRight");
323 static const auto cells = getParamFromGroup<std::array<int, bulkDim>>(bulkParamGroup, "Grid.Cells");
324 const auto cylSamples = cylIntegration.size();
325 const auto& a = line.corner(0);
326 const auto& b = line.corner(1);
327
328 // optionally write debugging / visual output of the integration points
329 static const auto writeIntegrationPointsToFile = getParam<bool>("MixedDimension.WriteIntegrationPointsToFile", false);
330 if (writeIntegrationPointsToFile)
331 {
332 std::ofstream ipPointFile("kernel_points.log", std::ios::app);
333 for (int i = 0; i < cylSamples; ++i)
334 {
335 const auto& point = cylIntegration.integrationPoint(i);
336 if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection)
337 ipPointFile << point[0] << "," << point[1] << "," << point[2] << '\n';
338 }
339 }
340
341 Scalar integral = 0.0;
342 for (int i = 0; i < cylSamples; ++i)
343 {
344 const auto& point = cylIntegration.integrationPoint(i);
345 // TODO: below only works for Cartesian grids with ijk numbering (e.g. level 0 YaspGrid (already fails when refined))
346 // more general is the bounding box tree solution which always works, however it's much slower
347 //const auto bulkIndices = intersectingEntities(point, this->problem(bulkIdx).gridGeometry().boundingBoxTree(), true);
348 if (const bool hasIntersection = intersectsPointBoundingBox(point, min, max); hasIntersection)
349 {
350 const auto bulkElementIdx = intersectingEntityCartesianGrid(point, min, max, cells);
351 assert(bulkElementIdx < this->problem(bulkIdx).gridGeometry().gridView().size(0));
352 const auto localWeight = evalConstKernel_(a, b, point, radius, kernelWidth)*cylIntegration.integrationElement(i)/Scalar(embeddings);
353 integral += localWeight;
354 if (!bulkSourceIds_[bulkElementIdx][0].empty() && id == bulkSourceIds_[bulkElementIdx][0].back())
355 {
356 bulkSourceWeights_[bulkElementIdx][0].back() += localWeight;
357 }
358 else
359 {
360 bulkSourceIds_[bulkElementIdx][0].emplace_back(id);
361 bulkSourceWeights_[bulkElementIdx][0].emplace_back(localWeight);
362 addBulkSourceStencils_(bulkElementIdx, lowDimElementIdx, coupledBulkElementIdx);
363 }
364 }
365 }
366
367 // the surface factor (which fraction of the source is inside the domain and needs to be considered)
368 const auto length = (a-b).two_norm()/Scalar(embeddings);
369 return integral/length;
370 }
371
372 void prepareDataStructures_()
373 {
374 // clear all internal members like pointsource vectors and stencils
375 // initializes the point source id counter
376 this->clear();
377 bulkSourceIds_.clear();
378 bulkSourceWeights_.clear();
379 extendedSourceStencil_.stencil().clear();
380
381 // precompute the vertex indices for efficiency for the box method
382 this->precomputeVertexIndices(bulkIdx);
383 this->precomputeVertexIndices(lowDimIdx);
384
385 bulkSourceIds_.resize(this->gridView(bulkIdx).size(0));
386 bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0));
387
388 // intersect the bounding box trees
389 this->glueGrids();
390
391 // reserve memory for data
392 this->pointSourceData().reserve(this->glue().size());
393 this->averageDistanceToBulkCell().reserve(this->glue().size());
394 fluxScalingFactor_.reserve(this->glue().size());
395
396 // reserve memory for stencils
397 const auto numBulkElements = this->gridView(bulkIdx).size(0);
398 for (GridIndex<bulkIdx> bulkElementIdx = 0; bulkElementIdx < numBulkElements; ++bulkElementIdx)
399 {
400 this->couplingStencils(bulkIdx)[bulkElementIdx].reserve(10);
401 extendedSourceStencil_.stencil()[bulkElementIdx].reserve(10);
402 bulkSourceIds_[bulkElementIdx][0].reserve(10);
403 bulkSourceWeights_[bulkElementIdx][0].reserve(10);
404 }
405 }
406
408 void makeUniqueStencil_()
409 {
410 // make extra stencils unique
411 for (auto&& stencil : extendedSourceStencil_.stencil())
412 {
413 std::sort(stencil.second.begin(), stencil.second.end());
414 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
415
416 // remove the vertices element (box)
417 if constexpr (isBox<bulkIdx>())
418 {
419 const auto& indices = this->vertexIndices(bulkIdx, stencil.first);
420 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
421 [&](auto i){ return std::find(indices.begin(), indices.end(), i) != indices.end(); }),
422 stencil.second.end());
423 }
424 // remove the own element (cell-centered)
425 else
426 {
427 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
428 [&](auto i){ return i == stencil.first; }),
429 stencil.second.end());
430 }
431 }
432
433 // make stencils unique
434 using namespace Dune::Hybrid;
435 forEach(integralRange(Dune::index_constant<2>{}), [&](const auto domainIdx)
436 {
437 for (auto&& stencil : this->couplingStencils(domainIdx))
438 {
439 std::sort(stencil.second.begin(), stencil.second.end());
440 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
441 }
442 });
443 }
444
446 void addBulkSourceStencils_(GridIndex<bulkIdx> bulkElementIdx, GridIndex<lowDimIdx> coupledLowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx)
447 {
448 // add the lowdim element to the coupling stencil of this bulk element
449 if constexpr (isBox<lowDimIdx>())
450 {
451 const auto& vertices = this->vertexIndices(lowDimIdx, coupledLowDimElementIdx);
452 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
453 vertices.begin(), vertices.end());
454
455 }
456 else
457 {
458 auto& s = this->couplingStencils(bulkIdx)[bulkElementIdx];
459 s.push_back(coupledLowDimElementIdx);
460 }
461
462 // the extended source stencil, every 3d element with a source is coupled to
463 // the element/dofs where the 3d quantities are measured
464 if constexpr (isBox<bulkIdx>())
465 {
466 const auto& vertices = this->vertexIndices(bulkIdx, coupledBulkElementIdx);
467 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
468 vertices.begin(), vertices.end());
469 }
470 else
471 {
472 auto& s = extendedSourceStencil_.stencil()[bulkElementIdx];
473 s.push_back(coupledBulkElementIdx);
474 }
475 }
476
478 Scalar evalConstKernel_(const GlobalPosition& a,
479 const GlobalPosition& b,
480 const GlobalPosition& point,
481 const Scalar R,
482 const Scalar rho) const noexcept
483 {
484 // projection of point onto line a + t*(b-a)
485 const auto ab = b - a;
486 const auto t = (point - a)*ab/ab.two_norm2();
487
488 // return 0 if we are outside cylinder
489 if (t < 0.0 || t > 1.0)
490 return 0.0;
491
492 // compute distance
493 auto proj = a; proj.axpy(t, ab);
494 const auto radiusSquared = (proj - point).two_norm2();
495
496 if (radiusSquared > rho*rho)
497 return 0.0;
498
499 return 1.0/(M_PI*rho*rho);
500 }
501
505 template<class SpatialParams>
506 auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx)
507 -> std::enable_if_t<hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
508 { return spatialParams.kernelWidthFactor(eIdx); }
509
513 template<class SpatialParams>
514 auto kernelWidthFactor_(const SpatialParams& spatialParams, unsigned int eIdx)
515 -> std::enable_if_t<!hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
516 {
517 static const Scalar kernelWidthFactor = getParam<Scalar>("MixedDimension.KernelWidthFactor");
518 return kernelWidthFactor;
519 }
520
522 EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_;
523
525 std::vector<std::array<std::vector<std::size_t>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceIds_;
527 std::vector<std::array<std::vector<Scalar>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceWeights_;
528
530 std::vector<Scalar> fluxScalingFactor_;
531};
532
534template<class MDTraits>
535struct CouplingManagerSupportsMultithreadedAssembly<Embedded1d3dCouplingManager<MDTraits, Embedded1d3dCouplingMode::Kernel>>
536: public std::true_type {};
537
538} // end namespace Dumux
539
540#endif
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d_kernel.hh:54
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:132
const BaseType::template CouplingStencils< bulkIdx >::mapped_type & extendedSourceStencil(std::size_t eIdx) const
Extended source stencil (for the bulk domain)
Definition: couplingmanager1d3d_kernel.hh:303
const std::vector< Scalar > & bulkSourceWeights(GridIndex< bulkIdx > eIdx, int scvIdx=0) const
return all source ids for a bulk elements
Definition: couplingmanager1d3d_kernel.hh:290
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:119
const std::vector< std::size_t > & bulkSourceIds(GridIndex< bulkIdx > eIdx, int scvIdx=0) const
Methods to be accessed by the subproblems.
Definition: couplingmanager1d3d_kernel.hh:286
void computePointSourceData(std::size_t order=1, bool verbose=false)
Definition: couplingmanager1d3d_kernel.hh:150
Scalar fluxScalingFactor(std::size_t id) const
The flux scaling factor for a source with id.
Definition: couplingmanager1d3d_kernel.hh:294
void init(std::shared_ptr< Problem< bulkIdx > > bulkProblem, std::shared_ptr< Problem< lowDimIdx > > lowDimProblem, const SolutionVector &curSol)
Definition: couplingmanager1d3d_kernel.hh:101
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3dbase.hh:35
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanager1d3d.hh:24
Helper class to integrate over a cylinder domain.
Definition: cylinderintegration.hh:48
void setGeometry(const GlobalPosition &bottomCenter, const GlobalPosition &topCenter, const Scalar radius, int verbosity=0)
Set the geometry of the cylinder.
Definition: cylinderintegration.hh:85
Manages the coupling between bulk elements and lower dimensional elements Point sources on each integ...
Definition: couplingmanagerbase.hh:71
Defines all properties used in Dumux.
Coupling manager for low-dimensional domains embedded in the bulk domain.
Integration over cylindrical and elliptic cylindrical domains Lowest order integration formulas that ...
Helper functions for distance queries.
Extended source stencil helper class for coupling managers.
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:439
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
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:265
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Defines the index types used for grid and local indices.
constexpr Box box
Definition: method.hh:147
constexpr Kernel kernel
Definition: couplingmanager1d3d_kernel.hh:38
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: adapt.hh:17
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:305
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
Definition: couplingmanager1d3d_kernel.hh:34
static std::string name()
Definition: couplingmanager1d3d_kernel.hh:35
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:26
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition: tag.hh:30
Helper class to create (named and comparable) tagged types.