69 using Scalar =
typename MDTraits::Scalar;
70 using SolutionVector =
typename MDTraits::SolutionVector;
71 using PointSourceData =
typename ParentType::PointSourceTraits::PointSourceData;
73 static constexpr auto bulkIdx =
typename MDTraits::template SubDomain<0>::Index();
74 static constexpr auto lowDimIdx =
typename MDTraits::template SubDomain<1>::Index();
77 template<std::
size_t id>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
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;
84 using GlobalPosition =
typename Element<bulkIdx>::Geometry::GlobalCoordinate;
86 template<std::
size_t id>
87 static constexpr bool isBox()
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");
94 bulkDim = GridView<bulkIdx>::dimension,
95 lowDimDim = GridView<lowDimIdx>::dimension,
96 dimWorld = GridView<bulkIdx>::dimensionworld
100 template <
typename T,
typename ...Ts>
101 using VariableKernelWidthDetector =
decltype(std::declval<T>().kernelWidthFactor(std::declval<Ts>()...));
103 template<
class T,
typename ...Args>
104 static constexpr bool hasKernelWidthFactor()
105 {
return Dune::Std::is_detected<VariableKernelWidthDetector, T, Args...>::value; }
110 using ParentType::ParentType;
112 void init(std::shared_ptr<Problem<bulkIdx>> bulkProblem,
113 std::shared_ptr<Problem<lowDimIdx>> lowDimProblem,
114 const SolutionVector&
curSol)
121 DUNE_THROW(Dune::NotImplemented,
"The current intersection detection may likely fail for refined grids.");
130 template<std::
size_t id,
class JacobianPattern>
133 extendedSourceStencil_.extendJacobianPattern(*
this, domainI, pattern);
143 template<std::
size_t i,
class LocalAssemblerI,
class JacobianMatrixDiagBlock,
class Gr
idVariables>
145 const LocalAssemblerI& localAssemblerI,
146 const typename LocalAssemblerI::LocalResidual::ElementResidualVector&,
147 JacobianMatrixDiagBlock& A,
148 GridVariables& gridVariables)
150 extendedSourceStencil_.evalAdditionalDomainDerivatives(*
this, domainI, localAssemblerI, this->
curSol(), A, gridVariables);
167 std::cout <<
"[coupling] Initializing the integration point source data structures..." << std::endl;
170 prepareDataStructures_();
171 std::cout <<
"[coupling] Resized data structures." << std::endl;
173 const auto& bulkGridGeometry = this->
problem(bulkIdx).gridGeometry();
174 const auto& lowDimGridGeometry = this->
problem(lowDimIdx).gridGeometry();
178 static const double characteristicRelativeLength =
getParam<double>(
"MixedDimension.KernelIntegrationCRL", 0.1);
181 static const bool writeIntegrationPointsToFile =
getParam<bool>(
"MixedDimension.WriteIntegrationPointsToFile",
false);
182 if (writeIntegrationPointsToFile)
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;
189 for (
const auto& is : intersections(this->
glue()))
192 const auto& inside = is.targetEntity(0);
194 const auto intersectionGeometry = is.geometry();
195 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
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);
208 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
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());
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());
230 PointSourceData psData;
233 if constexpr (isBox<lowDimIdx>())
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);
243 psData.addLowDimInterpolation(lowDimElementIdx);
247 if constexpr (isBox<bulkIdx>())
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);
257 psData.addBulkInterpolation(bulkElementIdx);
269 if constexpr (isBox<bulkIdx>())
271 const auto& vertices = this->
vertexIndices(bulkIdx, bulkElementIdx);
273 vertices.begin(), vertices.end());
277 this->
couplingStencils(lowDimIdx)[lowDimElementIdx].push_back(bulkElementIdx);
283 makeUniqueStencil_();
286 DUNE_THROW(Dune::InvalidStateException,
"Kernel method shouldn't have point sources in the bulk domain but only volume sources!");
288 std::cout <<
"[coupling] Finished preparing manager in " << watch.elapsed() <<
" seconds." << std::endl;
295 lowDimVolumeInBulkElement_.resize(this->
gridView(bulkIdx).size(0));
297 const auto& lowDimGridGeometry = this->
problem(lowDimIdx).gridGeometry();
298 const auto& bulkGridGeometry = this->
problem(bulkIdx).gridGeometry();
301 for (
const auto& is : intersections(this->
glue()))
304 const auto& inside = is.targetEntity(0);
305 const auto intersectionGeometry = is.geometry();
306 const auto lowDimElementIdx = lowDimGridGeometry.elementMapper().index(inside);
309 const auto radius = this->
problem(lowDimIdx).spatialParams().radius(lowDimElementIdx);
310 for (
int outsideIdx = 0; outsideIdx < is.numDomainNeighbors(); ++outsideIdx)
312 const auto& outside = is.domainEntity(outsideIdx);
313 const auto bulkElementIdx = bulkGridGeometry.elementMapper().index(outside);
314 lowDimVolumeInBulkElement_[bulkElementIdx] += intersectionGeometry.volume()*M_PI*
radius*
radius;
328 return this->
problem(lowDimIdx).spatialParams().radius(data.lowDimElementIdx());
335 const auto eIdx = this->
problem(bulkIdx).gridGeometry().elementMapper().index(element);
336 return lowDimVolumeInBulkElement_[eIdx];
343 const auto totalVolume = element.geometry().volume();
348 const std::vector<std::size_t>&
bulkSourceIds(GridIndex<bulkIdx> eIdx,
int scvIdx = 0)
const
349 {
return bulkSourceIds_[eIdx][scvIdx]; }
353 {
return bulkSourceWeights_[eIdx][scvIdx]; }
357 {
return fluxScalingFactor_[id]; }
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)
369 static const auto bulkParamGroup = this->problem(bulkIdx).paramGroup();
373 const auto cylSamples = cylIntegration.size();
374 const auto& a = line.corner(0);
375 const auto& b = line.corner(1);
378 static const auto writeIntegrationPointsToFile =
getParam<bool>(
"MixedDimension.WriteIntegrationPointsToFile",
false);
379 if (writeIntegrationPointsToFile)
381 std::ofstream ipPointFile(
"kernel_points.log", std::ios::app);
382 for (
int i = 0; i < cylSamples; ++i)
384 const auto& point = cylIntegration.integrationPoint(i);
386 ipPointFile << point[0] <<
"," << point[1] <<
"," << point[2] <<
'\n';
390 Scalar integral = 0.0;
391 for (
int i = 0; i < cylSamples; ++i)
393 const auto& point = cylIntegration.integrationPoint(i);
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())
405 bulkSourceWeights_[bulkElementIdx][0].back() += localWeight;
409 bulkSourceIds_[bulkElementIdx][0].emplace_back(
id);
410 bulkSourceWeights_[bulkElementIdx][0].emplace_back(localWeight);
411 addBulkSourceStencils_(bulkElementIdx, lowDimElementIdx, coupledBulkElementIdx);
417 const auto length = (a-b).two_norm()/Scalar(embeddings);
418 return integral/length;
421 void prepareDataStructures_()
426 bulkSourceIds_.clear();
427 bulkSourceWeights_.clear();
428 extendedSourceStencil_.stencil().clear();
431 this->precomputeVertexIndices(bulkIdx);
432 this->precomputeVertexIndices(lowDimIdx);
434 bulkSourceIds_.resize(this->gridView(bulkIdx).size(0));
435 bulkSourceWeights_.resize(this->gridView(bulkIdx).size(0));
441 this->pointSourceData().reserve(this->glue().size());
442 this->averageDistanceToBulkCell().reserve(this->glue().size());
443 fluxScalingFactor_.reserve(this->glue().size());
446 const auto numBulkElements = this->gridView(bulkIdx).size(0);
447 for (GridIndex<bulkIdx> bulkElementIdx = 0; bulkElementIdx < numBulkElements; ++bulkElementIdx)
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);
457 void makeUniqueStencil_()
460 for (
auto&& stencil : extendedSourceStencil_.stencil())
462 std::sort(stencil.second.begin(), stencil.second.end());
463 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
466 if constexpr (isBox<bulkIdx>())
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());
476 stencil.second.erase(std::remove_if(stencil.second.begin(), stencil.second.end(),
477 [&](
auto i){ return i == stencil.first; }),
478 stencil.second.end());
483 using namespace Dune::Hybrid;
484 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
486 for (
auto&& stencil : this->couplingStencils(domainIdx))
488 std::sort(stencil.second.begin(), stencil.second.end());
489 stencil.second.erase(std::unique(stencil.second.begin(), stencil.second.end()), stencil.second.end());
495 void addBulkSourceStencils_(GridIndex<bulkIdx> bulkElementIdx, GridIndex<lowDimIdx> coupledLowDimElementIdx, GridIndex<bulkIdx> coupledBulkElementIdx)
498 if constexpr (isBox<lowDimIdx>())
500 const auto& vertices = this->vertexIndices(lowDimIdx, coupledLowDimElementIdx);
501 this->couplingStencils(bulkIdx)[bulkElementIdx].insert(this->couplingStencils(bulkIdx)[bulkElementIdx].end(),
502 vertices.begin(), vertices.end());
507 auto& s = this->couplingStencils(bulkIdx)[bulkElementIdx];
508 s.push_back(coupledLowDimElementIdx);
513 if constexpr (isBox<bulkIdx>())
515 const auto& vertices = this->vertexIndices(bulkIdx, coupledBulkElementIdx);
516 extendedSourceStencil_.stencil()[bulkElementIdx].insert(extendedSourceStencil_.stencil()[bulkElementIdx].end(),
517 vertices.begin(), vertices.end());
521 auto& s = extendedSourceStencil_.stencil()[bulkElementIdx];
522 s.push_back(coupledBulkElementIdx);
527 Scalar evalConstKernel_(
const GlobalPosition& a,
528 const GlobalPosition& b,
529 const GlobalPosition& point,
531 const Scalar rho)
const noexcept
534 const auto ab = b - a;
535 const auto t = (point - a)*ab/ab.two_norm2();
538 if (t < 0.0 || t > 1.0)
542 auto proj = a; proj.axpy(t, ab);
543 const auto radiusSquared = (proj - point).two_norm2();
545 if (radiusSquared > rho*rho)
548 return 1.0/(M_PI*rho*rho);
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); }
562 template<
class SpatialParams>
563 auto kernelWidthFactor_(
const SpatialParams& spatialParams,
unsigned int eIdx)
564 -> std::enable_if_t<!hasKernelWidthFactor<SpatialParams, unsigned int>(), Scalar>
566 static const Scalar kernelWidthFactor =
getParam<Scalar>(
"MixedDimension.KernelWidthFactor");
567 return kernelWidthFactor;
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_;
580 std::vector<Scalar> fluxScalingFactor_;