24#ifndef DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
25#define DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/axisalignedcubegeometry.hh>
47template<
class Gr
idVariables,
class SolutionVector,
class LocalRes
idual>
50 using Scalar =
typename GridVariables::Scalar;
51 using GridGeometry =
typename GridVariables::GridGeometry;
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
54 using GridView =
typename GridGeometry::GridView;
55 using VolumeVariables =
typename GridVariables::VolumeVariables;
56 using Element =
typename GridView::template Codim<0>::Entity;
57 using NumEqVector =
typename LocalResidual::ElementResidualVector::value_type;
59 static constexpr auto dim = GridView::dimension;
60 static constexpr auto dimWorld = GridView::dimensionworld;
62 static_assert(dim > 1,
"Only implemented for dim > 1");
64 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
67 using SurfaceT = Dune::AxisAlignedCubeGeometry<Scalar, (dim == 2 ? 1 : 2), dimWorld>;
72 std::size_t normalDirectionIndex;
84 const SolutionVector& sol,
85 const LocalResidual& localResidual)
86 : gridVariables_(gridVariables)
88 , localResidual_(localResidual)
90 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(),
"FluxOverAxisAlignedSurface.Verbose",
false);
102 static_assert(std::is_same_v<std::decay_t<T>,
Surface>);
103 surfaces_.emplace(std::make_pair(
104 name, std::make_pair(
surface, NumEqVector(0.0))
116 const GlobalPosition& lowerLeft,
117 const GlobalPosition& upperRight)
120 const GlobalPosition v = upperRight - lowerLeft;
121 const auto it = std::find_if(v.begin(), v.end(), [](
const auto& x){ return abs(x) < 1e-20; });
123 DUNE_THROW(Dune::InvalidStateException,
"Surface is not axis-parallel!");
125 const std::size_t normalDirectionIndex =
std::distance(v.begin(), it);
126 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
127 inSurfaceAxes.set(normalDirectionIndex,
false);
130 surfaces_.emplace(std::make_pair(
133 std::move(
surface), normalDirectionIndex, NumEqVector(0.0)
146 const GlobalPosition& center,
147 const std::size_t normalDirectionIndex)
149 GlobalPosition lowerLeft = gridVariables_.gridGeometry().bBoxMin();
150 GlobalPosition upperRight = gridVariables_.gridGeometry().bBoxMax();
152 lowerLeft[normalDirectionIndex] = center[normalDirectionIndex];
153 upperRight[normalDirectionIndex] = center[normalDirectionIndex];
155 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
156 inSurfaceAxes.set(normalDirectionIndex,
false);
159 surfaces_.emplace(std::make_pair(
162 std::move(
surface), normalDirectionIndex, NumEqVector(0.0)
172 auto fluxType = [
this](
const auto& element,
173 const auto& fvGeometry,
174 const auto& elemVolVars,
176 const auto& elemFluxVarsCache)
178 return localResidual_.evalFlux(
179 problem_(), element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf
197 template<
class FluxType>
201 for (
auto&
surface : surfaces_)
204 snapSurfaceToClosestFace_();
205 calculateFluxes_(fluxType);
213 const auto&
flux(
const std::string& name)
const
215 return surfaces_.at(name).flux;
221 const std::map<std::string, SurfaceData>&
surfaces()
const
222 {
return surfaces_; }
229 for (
const auto& [name, data] : surfaces_)
230 std::cout <<
"Flux over surface " << name <<
": " << data.flux << std::endl;
235 template<
class FluxType>
236 void calculateFluxes_(
const FluxType& fluxType)
238 auto fvGeometry =
localView(problem_().gridGeometry());
239 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
240 auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache());
242 for (
const auto& element : elements(problem_().gridGeometry().gridView()))
244 fvGeometry.bindElement(element);
245 elemVolVars.bindElement(element, fvGeometry, sol_);
246 elemFluxVarsCache.bindElement(element, fvGeometry, elemVolVars);
248 for (
const auto& scvf : scvfs(fvGeometry))
252 for (
auto& [name, surfaceData] : surfaces_)
254 if (considerScvf_(scvf, surfaceData))
256 const auto result = fluxType(element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
257 surfaceData.flux += result;
260 std::cout <<
"At element " << problem_().gridGeometry().elementMapper().index(element)
261 <<
": Flux at face " << scvf.ipGlobal() <<
": " << result <<
" (" << name <<
")" << std::endl;
269 bool considerScvf_(
const SubControlVolumeFace& scvf,
const SurfaceData& SurfaceData)
const
274 if (scvf.boundary() || !std::signbit(scvf.unitOuterNormal()[SurfaceData.normalDirectionIndex]))
280 void snapSurfaceToClosestFace_()
283 const auto gridView = problem_().gridGeometry().gridView();
285 for (
auto& [name, surfaceData] : surfaces_)
287 GeometriesEntitySet entitySet({surfaceData.surface});
290 problem_().gridGeometry().boundingBoxTree(), geometriesTree
293 if (intersectingElements.empty())
295 std::cout <<
"surface boundaries: " << std::endl;
296 printSurfaceBoundaries_(surfaceData.surface);
298 DUNE_THROW(Dune::InvalidStateException,
"surface " << name <<
" does not intersect with any element");
301 std::vector<std::size_t> sortedResults;
302 sortedResults.reserve(gridView.size(0));
304 for (
const auto& i : intersectingElements)
305 sortedResults.push_back(i.first());
307 std::sort(sortedResults.begin(), sortedResults.end());
308 sortedResults.erase(std::unique(
309 sortedResults.begin(), sortedResults.end()
310 ), sortedResults.end());
313 GlobalPosition normalVector(0.0);
314 normalVector[surfaceData.normalDirectionIndex] = 1.0;
316 const auto& firstIntersectingElement = problem_().gridGeometry().element(sortedResults[0]);
317 Scalar
distance = std::numeric_limits<Scalar>::max();
318 bool snappingOcurred =
false;
320 GlobalPosition surfaceLowerLeft = surfaceData.surface.corner(0);
321 GlobalPosition surfaceUpperRight = surfaceData.surface.corner(3);
323 bool surfaceAlreadyOnFaces =
false;
324 for (
const auto& intersection : intersections(gridView, firstIntersectingElement))
326 if (surfaceAlreadyOnFaces)
330 if (abs(1.0 - abs(normalVector * intersection.centerUnitOuterNormal())) < 1e-8)
333 const auto getDistance = [](
const auto& p,
const auto& geo)
335 if constexpr (dim == 2)
341 const auto& geo = intersection.geometry();
342 if (
const Scalar d = getDistance(geo.center(), surfaceData.surface); d < 1e-8 *
diameter(geo))
345 surfaceAlreadyOnFaces =
true;
346 snappingOcurred =
false;
351 snappingOcurred =
true;
354 for (
int i = 0; i < surfaceData.surface.corners(); ++i)
356 const auto& faceCenter = geo.center();
357 surfaceLowerLeft[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
358 surfaceUpperRight[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
366 std::cout <<
"\n\nSurface '" << name <<
"' was automatically snapped to the closest faces" << std::endl;
367 std::cout <<
"Old surface boundaries: " << std::endl;
368 printSurfaceBoundaries_(surfaceData.surface);
371 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
372 inSurfaceAxes.set(surfaceData.normalDirectionIndex,
false);
373 surfaceData.surface =
Surface{surfaceLowerLeft, surfaceUpperRight, inSurfaceAxes};
375 std::cout <<
"New surface boundaries: " << std::endl;
376 printSurfaceBoundaries_(surfaceData.surface);
377 std::cout << std::endl;
384 for (
int i = 0; i <
surface.corners(); ++i)
385 std::cout <<
surface.corner(i) << std::endl;
388 const auto& problem_()
const {
return gridVariables_.curGridVolVars().problem(); }
390 std::map<std::string, SurfaceData> surfaces_;
391 const GridVariables& gridVariables_;
392 const SolutionVector& sol_;
393 const LocalResidual localResidual_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
Detect if a point intersects a geometry.
Helper functions for distance queries.
An interface for a set of geometric entities.
Algorithms that finds which geometric entites intersect.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:292
static Geometry::ctype distancePointPolygon(const typename Geometry::GlobalCoordinate &p, const Geometry &geometry)
Compute the shortest distance from a point to a given polygon geometry.
Definition: distance.hh:266
static Point::value_type distancePointSegment(const Point &p, const Point &a, const Point &b)
Compute the distance from a point to the segment connecting the points a and b.
Definition: distance.hh:145
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: intersectingentities.hh:112
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:36
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
constexpr Surface surface
Definition: couplingmanager1d3d_surface.hh:49
Class used to calculate fluxes over axis-aligned surfaces.
Definition: fluxoveraxisalignedsurface.hh:49
void addAxisAlignedPlane(const std::string &name, const GlobalPosition ¢er, const std::size_t normalDirectionIndex)
Add an axis-aligned plane (line in 2D) with a given name, specifying the planes's center and normal.
Definition: fluxoveraxisalignedsurface.hh:145
FluxOverAxisAlignedSurface(const GridVariables &gridVariables, const SolutionVector &sol, const LocalResidual &localResidual)
The constructor.
Definition: fluxoveraxisalignedsurface.hh:83
void printAllFluxes() const
Prints all fluxes.
Definition: fluxoveraxisalignedsurface.hh:227
void addAxisAlignedSurface(const std::string &name, T &&surface)
Add an axis-aligned surface with a given name.
Definition: fluxoveraxisalignedsurface.hh:100
const std::map< std::string, SurfaceData > & surfaces() const
Provides access to all surfaces.
Definition: fluxoveraxisalignedsurface.hh:221
const auto & flux(const std::string &name) const
Return the flux over given surface.
Definition: fluxoveraxisalignedsurface.hh:213
void calculateAllFluxes()
Calculate the fluxes over all surfaces.
Definition: fluxoveraxisalignedsurface.hh:170
SurfaceT Surface
Definition: fluxoveraxisalignedsurface.hh:78
void addAxisAlignedSurface(const std::string &name, const GlobalPosition &lowerLeft, const GlobalPosition &upperRight)
Add an axis-aligned surface (segment in 2D) with a given name, specifying the surface's corner points...
Definition: fluxoveraxisalignedsurface.hh:115
void calculateFluxes(const FluxType &fluxType)
Calculate the fluxes over all surfaces for a given flux type.
Definition: fluxoveraxisalignedsurface.hh:198
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:66
An interface for a set of geometric entities.
Definition: geometricentityset.hh:119