377 std::cout <<
"Initializing the coupling manager (projection)" << std::endl;
380 const bool enableIntersectionOutput = getParam<bool>(
"MixedDimension.Projection.EnableIntersectionOutput",
false);
381 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkCoupledFaces =
382 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() :
nullptr;
383 std::unique_ptr<Detail::DebugIntersectionVTKOutput> vtkIntersections =
384 enableIntersectionOutput ? std::make_unique<Detail::DebugIntersectionVTKOutput>() :
nullptr;
388 const auto& lowDimProblem = this->problem(lowDimIdx);
389 const auto& lowDimFvGridGeometry = lowDimProblem.gridGeometry();
390 const auto& lowDimGridView = lowDimFvGridGeometry.gridView();
391 localAreaFactor_.resize(lowDimGridView.size(0), 0.0);
392 const auto radiusFunc = [&](
const GridIndex<lowDimIdx> eIdx) {
393 return lowDimProblem.spatialParams().radius(eIdx);
396 lowDimFvGridGeometry, radiusFunc
400 this->precomputeVertexIndices(bulkIdx);
401 this->precomputeVertexIndices(lowDimIdx);
403 const auto& bulkProblem = this->problem(bulkIdx);
404 const auto& bulkFvGridGeometry = bulkProblem.gridGeometry();
405 const auto& bulkGridView = bulkFvGridGeometry.gridView();
406 bulkElementMarker_.assign(bulkGridView.size(0), 0);
407 bulkVertexMarker_.assign(bulkGridView.size(bulkDim), 0);
411 static const auto bBoxShrinking
412 = getParam<Scalar>(
"MixedDimension.Projection.CoupledBoundingBoxShrinkingFactor", 1e-2);
413 const GlobalPosition threshold(
414 bBoxShrinking*( bulkFvGridGeometry.bBoxMax()-bulkFvGridGeometry.bBoxMin() ).two_norm()
416 const auto bBoxMinSmall = bulkFvGridGeometry.bBoxMin() + threshold;
417 const auto bBoxMaxSmall = bulkFvGridGeometry.bBoxMax() - threshold;
418 auto insideBBox = [bBoxMin=bBoxMinSmall, bBoxMax=bBoxMaxSmall](
const GlobalPosition& point) ->
bool
420 static constexpr Scalar eps_ = 1.0e-7;
421 const Scalar eps0 = eps_*(bBoxMax[0] - bBoxMin[0]);
422 const Scalar eps1 = eps_*(bBoxMax[1] - bBoxMin[1]);
423 const Scalar eps2 = eps_*(bBoxMax[2] - bBoxMin[2]);
424 return (bBoxMin[0] - eps0 <= point[0] && point[0] <= bBoxMax[0] + eps0 &&
425 bBoxMin[1] - eps1 <= point[1] && point[1] <= bBoxMax[1] + eps1 &&
426 bBoxMin[2] - eps2 <= point[2] && point[2] <= bBoxMax[2] + eps2);
430 static const auto quadSimplexRefineMaxLevel
431 = getParam<std::size_t>(
"MixedDimension.Projection.SimplexIntegrationRefine", 4);
435 static const auto estimateNumberOfPointSources
436 = getParam<std::size_t>(
"MixedDimension.Projection.EstimateNumberOfPointSources", bulkFvGridGeometry.gridView().size(0));
438 this->pointSources(bulkIdx).reserve(estimateNumberOfPointSources);
439 this->pointSources(lowDimIdx).reserve(estimateNumberOfPointSources);
440 this->pointSourceData().reserve(estimateNumberOfPointSources);
445 for (
const auto& element : elements(bulkGridView))
447 const auto bulkElementIdx = bulkFvGridGeometry.elementMapper().index(element);
448 const auto bulkGeometry = element.geometry();
449 const auto refElement = Dune::referenceElement(element);
450 for (
const auto& intersection : intersections(bulkGridView, element))
453 if (!intersection.boundary())
457 const auto [isCoupled, closestSegmentIdx, minDist] = isCoupled_(network, intersection, insideBBox);
462 bulkElementMarker_[bulkElementIdx] = 1;
464 const auto interfaceGeometry = intersection.geometry();
465 for (
int i = 0; i < interfaceGeometry.corners(); ++i)
467 const auto localVIdx = refElement.subEntity(intersection.indexInInside(), 1, i, bulkDim);
468 const auto vIdx = bulkFvGridGeometry.vertexMapper().subIndex(element, localVIdx, bulkDim);
469 bulkVertexMarker_[vIdx] = 1;
473 if (enableIntersectionOutput)
474 vtkCoupledFaces->push({
475 { interfaceGeometry.corner(0), interfaceGeometry.corner(1), interfaceGeometry.corner(2)}
476 }, closestSegmentIdx);
481 for (
const auto& qp : quadSimplex)
483 const auto surfacePoint = interfaceGeometry.global(qp.position());
484 const auto closestSegmentIdx = qp.indicator();
485 const auto projPoint = network.projectionPointOnSegment(surfacePoint, closestSegmentIdx);
487 addPointSourceAtIP_(bulkGeometry, interfaceGeometry, qp, bulkElementIdx, closestSegmentIdx, surfacePoint, projPoint);
488 addStencilEntries_(bulkElementIdx, closestSegmentIdx);
494 using namespace Dune::Hybrid;
495 forEach(integralRange(Dune::index_constant<2>{}), [&](
const auto domainIdx)
497 for (
auto&& stencil : this->couplingStencils(domainIdx))
499 std::sort(stencil.second.begin(), stencil.second.end());
500 stencil.second.erase(
501 std::unique(stencil.second.begin(), stencil.second.end()),
508 if (enableIntersectionOutput)
509 vtkCoupledFaces->write(
"coupledfaces.vtk");
513 for (
const auto& element : elements(lowDimGridView))
515 const auto length = element.geometry().volume();
516 const auto eIdx = lowDimFvGridGeometry.elementMapper().index(element);
517 const auto radius = this->problem(lowDimIdx).spatialParams().radius(eIdx);
518 const auto cylinderSurface = 2.0*M_PI*radius*length;
519 localAreaFactor_[eIdx] /= cylinderSurface;
520 localAreaFactor_[eIdx] -= 1.0;
523 this->pointSources(bulkIdx).shrink_to_fit();
524 this->pointSources(lowDimIdx).shrink_to_fit();
525 this->pointSourceData().shrink_to_fit();
527 std::cout <<
"-- Coupling at " << this->pointSourceData().size() <<
" integration points." << std::endl;
528 std::cout <<
"-- Finished initialization after " << watch.elapsed() <<
" seconds." << std::endl;