57 using Scalar =
typename MDTraits::Scalar;
58 using SolutionVector =
typename MDTraits::SolutionVector;
59 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
61 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
62 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
65 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
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;
72 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
74 template<std::
size_t id>
75 static constexpr bool isBox()
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");
82 bulkDim = GridView<bulkIdx>::dimension,
83 lowDimDim = GridView<lowDimIdx>::dimension,
84 dimWorld = GridView<bulkIdx>::dimensionworld
88 template <
typename T,
typename ...Ts>
89 using VariableKernelWidthDetector =
decltype(std::declval<T>().kernelWidthFactor(std::declval<Ts>()...));
91 template<
class T,
typename ...Args>
92 static constexpr bool hasKernelWidthFactor()
93 {
return Dune::Std::is_detected<VariableKernelWidthDetector, T, Args...>::value; }
98 using ParentType::ParentType;
100 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
101 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
102 const SolutionVector& curSol)
104 ParentType::init(bulkProblem, lowDimProblem, curSol);
105 computeLowDimVolumeFractions();
107 const auto refinement = getParamFromGroup<int>(bulkProblem->paramGroup(),
"Grid.Refinement", 0);
109 DUNE_THROW(Dune::NotImplemented,
"The current intersection detection may likely fail for refined grids.");
118 template<std::
size_t id,
class JacobianPattern>
121 extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
131 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
133 const LocalAssemblerI& localAssemblerI,
134 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
135 JacobianMatrixDiagBlock& A,
136 GridVariables& gridVariables)
138 extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, A, gridVariables);
155 std::cout <<
"[coupling] Initializing the integration point source data structures..." << std::endl;
158 prepareDataStructures_();
159 std::cout <<
"[coupling] Resized data structures." << std::endl;
161 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
162 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
166 static const double characteristicRelativeLength = getParam<double>(
"MixedDimension.KernelIntegrationCRL", 0.1);
169 static const bool writeIntegrationPointsToFile = getParam<bool>(
"MixedDimension.WriteIntegrationPointsToFile",
false);
170 if (writeIntegrationPointsToFile)
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;
177 for (
const auto& is : intersections(this->glue()))
180 const auto& inside = is.targetEntity(0);
182 const auto intersectionGeometry = is.geometry();
183 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
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);
196 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
204 const auto id = this->idCounter_++;
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());
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());
218 PointSourceData psData;
221 if constexpr (isBox<lowDimIdx>())
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);
231 psData.addLowDimInterpolation(lowDimElementIdx);
235 if constexpr (isBox<bulkIdx>())
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);
245 psData.addBulkInterpolation(bulkElementIdx);
249 this->pointSourceData().emplace_back(std::move(psData));
252 this->averageDistanceToBulkCell().push_back(avgMinDist);
253 fluxScalingFactor_.push_back(this->problem(bulkIdx).fluxScalingFactor(avgMinDist, radius, kernelWidth));
257 if constexpr (isBox<bulkIdx>())
259 const auto& vertices = this->vertexIndices(bulkIdx, bulkElementIdx);
260 this->couplingStencils(lowDimIdx)[lowDimElementIdx].insert(this->couplingStencils(lowDimIdx)[lowDimElementIdx].end(),
261 vertices.begin(), vertices.end());
265 this->couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
271 makeUniqueStencil_();
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!");
276 std::cout <<
"[coupling] Finished preparing manager in " << watch.elapsed() <<
" seconds." << std::endl;
283 lowDimVolumeInBulkElement_.resize(this->gridView(bulkIdx).size(0));
285 const auto& lowDimGridGeometry = this->problem(lowDimIdx).gridGeometry();
286 const auto& bulkGridGeometry = this->problem(bulkIdx).gridGeometry();
289 for (
const auto& is : intersections(this->glue()))
292 const auto& inside = is.targetEntity(0);
293 const auto intersectionGeometry = is.geometry();
294 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
297 const auto radius = this->problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
298 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
300 const auto& outside = is.domainEntity(outsideIdx);
301 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
302 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*radius*radius;
315 const auto& data = this->pointSourceData()[id];
316 return this->problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
323 const auto eIdx = this->problem(bulkIdx).gridGeometry().elementMapper().index(element);
324 return lowDimVolumeInBulkElement_[eIdx];
331 const auto totalVolume = element.geometry().volume();
332 return lowDimVolume(element) / totalVolume;
336 const std::vector<std::size_t>&
bulkSourceIds(GridIndex<bulkIdx> eIdx,
int scvIdx = 0)
const
337 {
return bulkSourceIds_[eIdx][scvIdx]; }
341 {
return bulkSourceWeights_[eIdx][scvIdx]; }
345 {
return fluxScalingFactor_[id]; }
352 const typename ParentType::template CouplingStencils<bulkIdx>::mapped_type&
355 const auto& sourceStencils = extendedSourceStencil_.stencil();
356 if (
auto stencil = sourceStencils.find(eIdx); stencil != sourceStencils.end())
357 return stencil->second;
359 return this->emptyStencil(bulkIdx);
364 template<
class Line,
class CylIntegration>
365 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)
370 static const auto bulkParamGroup = this->problem(bulkIdx).paramGroup();
371 static const auto min = getParamFromGroup<GlobalPosition>(bulkParamGroup,
"Grid.LowerLeft");
372 static const auto max = getParamFromGroup<GlobalPosition>(bulkParamGroup,
"Grid.UpperRight");
373 static const auto cells = getParamFromGroup<std::array<int, bulkDim>>(bulkParamGroup,
"Grid.Cells");
374 const auto cylSamples = cylIntegration.size();
375 const auto& a = line.corner(0);
376 const auto& b = line.corner(1);
379 static const auto writeIntegrationPointsToFile = getParam<bool>(
"MixedDimension.WriteIntegrationPointsToFile",
false);
380 if (writeIntegrationPointsToFile)
382 std::ofstream ipPointFile(
"kernel_points.log", std::ios::app);
383 for (
int i = 0; i < cylSamples; ++i)
385 const auto& point = cylIntegration.integrationPoint(i);
387 ipPointFile << point[0] <<
"," << point[1] <<
"," << point[2] <<
'\n';
391 Scalar integral = 0.0;
392 for (
int i = 0; i < cylSamples; ++i)
394 const auto& point = cylIntegration.integrationPoint(i);
401 assert(bulkElementIdx < this->problem(bulkIdx).gridGeometry().gridView().size(0));
402 const auto localWeight = evalConstKernel_(a, b, point, radius, kernelWidth)*cylIntegration.integrationElement(i)/Scalar(embeddings);
403 integral += localWeight;
404 if (!bulkSourceIds_[bulkElementIdx][0].empty() &&
id == bulkSourceIds_[bulkElementIdx][0].back())
406 bulkSourceWeights_[bulkElementIdx][0].back() += localWeight;
410 bulkSourceIds_[bulkElementIdx][0].emplace_back(
id);
411 bulkSourceWeights_[bulkElementIdx][0].emplace_back(localWeight);
412 addBulkSourceStencils_(bulkElementIdx, lowDimElementIdx, coupledBulkElementIdx);
418 const auto length = (a-b).two_norm()/Scalar(embeddings);
419 return integral/length;
422 void prepareDataStructures_()
427 bulkSourceIds_.clear();
428 bulkSourceWeights_.clear();
429 extendedSourceStencil_.stencil().clear();
432 this->precomputeVertexIndices(bulkIdx);
433 this->precomputeVertexIndices(lowDimIdx);
435 bulkSourceIds_.resize(this->gridView(bulkIdx).size(0));
436 bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0));
442 this->pointSourceData().reserve(this->glue().size());
443 this->averageDistanceToBulkCell().reserve(this->glue().size());
444 fluxScalingFactor_.reserve(this->glue().size());
447 const auto numBulkElements = this->gridView(bulkIdx).size(0);
448 for (GridIndex<bulkIdx> bulkElementIdx = 0; bulkElementIdx < numBulkElements; ++bulkElementIdx)
450 this->couplingStencils(bulkIdx)[bulkElementIdx].reserve(10);
451 extendedSourceStencil_.stencil()[bulkElementIdx].reserve(10);
452 bulkSourceIds_[bulkElementIdx][0].reserve(10);
453 bulkSourceWeights_[bulkElementIdx][0].reserve(10);
458 void makeUniqueStencil_()
461 for (
auto&& stencil : extendedSourceStencil_.stencil())
463 std::sort(stencil.second.begin(), stencil.second.end());
464 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
467 if constexpr (isBox<bulkIdx>())
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());
477 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
478 [&](
auto i){ return i == stencil.first; }),
479 stencil.second.end());
484 using namespace Dune::Hybrid;
485 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
487 for (
auto&& stencil : this->couplingStencils(domainIdx))
489 std::sort(stencil.second.begin(), stencil.second.end());
490 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
496 void addBulkSourceStencils_(GridIndex<bulkIdx> bulkElementIdx, GridIndex<lowDimIdx> coupledLowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx)
499 if constexpr (isBox<lowDimIdx>())
501 const auto& vertices = this->vertexIndices(lowDimIdx, coupledLowDimElementIdx);
502 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
503 vertices.begin(), vertices.end());
508 auto& s = this->couplingStencils(bulkIdx)[bulkElementIdx];
509 s.push_back(coupledLowDimElementIdx);
514 if constexpr (isBox<bulkIdx>())
516 const auto& vertices = this->vertexIndices(bulkIdx, coupledBulkElementIdx);
517 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
518 vertices.begin(), vertices.end());
522 auto& s = extendedSourceStencil_.stencil()[bulkElementIdx];
523 s.push_back(coupledBulkElementIdx);
528 Scalar evalConstKernel_(
const GlobalPosition& a,
529 const GlobalPosition& b,
530 const GlobalPosition& point,
532 const Scalar rho)
const noexcept
535 const auto ab = b - a;
536 const auto t = (point - a)*ab/ab.two_norm2();
539 if (t < 0.0 || t > 1.0)
543 auto proj = a; proj.axpy(t, ab);
544 const auto radiusSquared = (proj - point).two_norm2();
546 if (radiusSquared > rho*rho)
549 return 1.0/(M_PI*rho*rho);
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); }
563 template<
class SpatialParams>
564 auto kernelWidthFactor_(
const SpatialParams& spatialParams,
unsigned int eIdx)
565 -> std::enable_if_t<!hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
567 static const Scalar kernelWidthFactor = getParam<Scalar>(
"MixedDimension.KernelWidthFactor");
568 return kernelWidthFactor;
572 EmbeddedCoupling::ExtendedSourceStencil<ThisType> extendedSourceStencil_;
574 std::vector<Scalar> lowDimVolumeInBulkElement_;
576 std::vector<std::array<std::vector<std::size_t>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceIds_;
578 std::vector<std::array<std::vector<Scalar>, isBox<bulkIdx>() ? 1<<bulkDim : 1>> bulkSourceWeights_;
581 std::vector<Scalar> fluxScalingFactor_;