12#ifndef DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
13#define DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/axisalignedcubegeometry.hh>
35template<
class Gr
idVariables,
class SolutionVector,
class LocalRes
idual>
38 using Scalar =
typename GridVariables::Scalar;
39 using GridGeometry =
typename GridVariables::GridGeometry;
40 using FVElementGeometry =
typename GridGeometry::LocalView;
41 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
42 using GridView =
typename GridGeometry::GridView;
43 using VolumeVariables =
typename GridVariables::VolumeVariables;
44 using Element =
typename GridView::template Codim<0>::Entity;
45 using NumEqVector =
typename LocalResidual::ElementResidualVector::value_type;
47 static constexpr auto dim = GridView::dimension;
48 static constexpr auto dimWorld = GridView::dimensionworld;
50 static_assert(dim > 1,
"Only implemented for dim > 1");
52 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
55 using SurfaceT = Dune::AxisAlignedCubeGeometry<Scalar, (dim == 2 ? 1 : 2), dimWorld>;
60 std::size_t normalDirectionIndex;
72 const SolutionVector& sol,
73 const LocalResidual& localResidual,
74 bool nonIntersectingSurfaceIsError =
false)
75 : gridVariables_(gridVariables)
77 , localResidual_(localResidual)
78 , nonIntersectingSurfaceIsError_(nonIntersectingSurfaceIsError)
80 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(),
"FluxOverAxisAlignedSurface.Verbose",
false);
92 static_assert(std::is_same_v<std::decay_t<T>,
Surface>);
93 surfaces_.emplace(std::make_pair(
94 name, std::make_pair(
surface, NumEqVector(0.0))
106 const GlobalPosition& lowerLeft,
107 const GlobalPosition& upperRight)
110 const GlobalPosition v = upperRight - lowerLeft;
111 const auto it = std::find_if(v.begin(), v.end(), [](
const auto& x){ return abs(x) < 1e-20; });
113 DUNE_THROW(Dune::InvalidStateException,
"Surface is not axis-parallel!");
115 const std::size_t normalDirectionIndex =
std::distance(v.begin(), it);
116 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
117 inSurfaceAxes.set(normalDirectionIndex,
false);
120 surfaces_.emplace(std::make_pair(
123 std::move(
surface), normalDirectionIndex, NumEqVector(0.0)
136 const GlobalPosition&
center,
137 const std::size_t normalDirectionIndex)
139 GlobalPosition lowerLeft = gridVariables_.gridGeometry().bBoxMin();
140 GlobalPosition upperRight = gridVariables_.gridGeometry().bBoxMax();
142 lowerLeft[normalDirectionIndex] =
center[normalDirectionIndex];
143 upperRight[normalDirectionIndex] =
center[normalDirectionIndex];
145 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
146 inSurfaceAxes.set(normalDirectionIndex,
false);
149 surfaces_.emplace(std::make_pair(
152 std::move(
surface), normalDirectionIndex, NumEqVector(0.0)
162 auto fluxType = [
this](
const auto& element,
163 const auto& fvGeometry,
164 const auto& elemVolVars,
166 const auto& elemFluxVarsCache)
168 return localResidual_.evalFlux(
169 problem_(), element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf
187 template<
class FluxType>
191 for (
auto&
surface : surfaces_)
194 snapSurfaceToClosestFace_();
195 calculateFluxes_(fluxType);
203 const auto&
flux(
const std::string& name)
const
205 return surfaces_.at(name).flux;
211 const std::map<std::string, SurfaceData>&
surfaces()
const
212 {
return surfaces_; }
219 for (
const auto& [name, data] : surfaces_)
220 std::cout <<
"Flux over surface " << name <<
": " << data.flux << std::endl;
227 { nonIntersectingSurfaceIsError_ = isError; }
231 template<
class FluxType>
232 void calculateFluxes_(
const FluxType& fluxType)
234 auto fvGeometry =
localView(problem_().gridGeometry());
235 auto elemVolVars =
localView(gridVariables_.curGridVolVars());
236 auto elemFluxVarsCache =
localView(gridVariables_.gridFluxVarsCache());
238 for (
const auto& element : elements(problem_().gridGeometry().gridView()))
240 fvGeometry.bindElement(element);
241 elemVolVars.bindElement(element, fvGeometry, sol_);
242 elemFluxVarsCache.bindElement(element, fvGeometry, elemVolVars);
244 for (
const auto& scvf : scvfs(fvGeometry))
248 for (
auto& [name, surfaceData] : surfaces_)
250 if (considerScvf_(scvf, surfaceData))
252 const auto result = fluxType(element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
253 surfaceData.flux += result;
256 std::cout <<
"At element " << problem_().gridGeometry().elementMapper().index(element)
257 <<
": Flux at face " << scvf.ipGlobal() <<
": " << result <<
" (" << name <<
")" << std::endl;
265 bool considerScvf_(
const SubControlVolumeFace& scvf,
const SurfaceData& SurfaceData)
const
270 if (scvf.boundary() || !std::signbit(scvf.unitOuterNormal()[SurfaceData.normalDirectionIndex]))
276 void snapSurfaceToClosestFace_()
279 const auto gridView = problem_().gridGeometry().gridView();
281 for (
auto& [name, surfaceData] : surfaces_)
283 GeometriesEntitySet entitySet({surfaceData.surface});
286 problem_().gridGeometry().boundingBoxTree(), geometriesTree
289 if (intersectingElements.empty())
291 if (!nonIntersectingSurfaceIsError_)
294 std::cout <<
"surface boundaries: " << std::endl;
295 printSurfaceBoundaries_(surfaceData.surface);
297 DUNE_THROW(Dune::InvalidStateException,
"surface " << name <<
" does not intersect with any element");
300 std::vector<std::size_t> sortedResults;
301 sortedResults.reserve(gridView.size(0));
303 for (
const auto& i : intersectingElements)
304 sortedResults.push_back(i.first());
306 std::sort(sortedResults.begin(), sortedResults.end());
307 sortedResults.erase(std::unique(
308 sortedResults.begin(), sortedResults.end()
309 ), sortedResults.end());
312 GlobalPosition normalVector(0.0);
313 normalVector[surfaceData.normalDirectionIndex] = 1.0;
315 const auto& firstIntersectingElement = problem_().gridGeometry().element(sortedResults[0]);
316 Scalar
distance = std::numeric_limits<Scalar>::max();
317 bool snappingOcurred =
false;
319 GlobalPosition surfaceLowerLeft = surfaceData.surface.corner(0);
320 GlobalPosition surfaceUpperRight = surfaceData.surface.corner(3);
322 bool surfaceAlreadyOnFaces =
false;
323 for (
const auto& intersection : intersections(gridView, firstIntersectingElement))
325 if (surfaceAlreadyOnFaces)
329 if (abs(1.0 - abs(normalVector * intersection.centerUnitOuterNormal())) < 1e-8)
332 const auto getDistance = [](
const auto& p,
const auto& geo)
334 if constexpr (dim == 2)
340 const auto& geo = intersection.geometry();
341 if (
const Scalar d = getDistance(geo.center(), surfaceData.surface); d < 1e-8 *
diameter(geo))
344 surfaceAlreadyOnFaces =
true;
345 snappingOcurred =
false;
350 snappingOcurred =
true;
353 for (
int i = 0; i < surfaceData.surface.corners(); ++i)
355 const auto& faceCenter = geo.center();
356 surfaceLowerLeft[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
357 surfaceUpperRight[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
365 std::cout <<
"\n\nSurface '" << name <<
"' was automatically snapped to the closest faces" << std::endl;
366 std::cout <<
"Old surface boundaries: " << std::endl;
367 printSurfaceBoundaries_(surfaceData.surface);
370 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
371 inSurfaceAxes.set(surfaceData.normalDirectionIndex,
false);
372 surfaceData.surface =
Surface{surfaceLowerLeft, surfaceUpperRight, inSurfaceAxes};
374 std::cout <<
"New surface boundaries: " << std::endl;
375 printSurfaceBoundaries_(surfaceData.surface);
376 std::cout << std::endl;
383 for (
int i = 0; i <
surface.corners(); ++i)
384 std::cout <<
surface.corner(i) << std::endl;
387 const auto& problem_()
const {
return gridVariables_.curGridVolVars().problem(); }
389 std::map<std::string, SurfaceData> surfaces_;
390 const GridVariables& gridVariables_;
391 const SolutionVector& sol_;
392 const LocalResidual localResidual_;
394 bool nonIntersectingSurfaceIsError_;
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:56
Class used to calculate fluxes over axis-aligned surfaces.
Definition: fluxoveraxisalignedsurface.hh:37
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:135
void setNonIntersectingSurfaceIsError(bool isError=true)
Set if non-intersecting surfaces are treated as error.
Definition: fluxoveraxisalignedsurface.hh:226
void printAllFluxes() const
Prints all fluxes.
Definition: fluxoveraxisalignedsurface.hh:217
void addAxisAlignedSurface(const std::string &name, T &&surface)
Add an axis-aligned surface with a given name.
Definition: fluxoveraxisalignedsurface.hh:90
const std::map< std::string, SurfaceData > & surfaces() const
Provides access to all surfaces.
Definition: fluxoveraxisalignedsurface.hh:211
const auto & flux(const std::string &name) const
Return the flux over given surface.
Definition: fluxoveraxisalignedsurface.hh:203
void calculateAllFluxes()
Calculate the fluxes over all surfaces.
Definition: fluxoveraxisalignedsurface.hh:160
SurfaceT Surface
Definition: fluxoveraxisalignedsurface.hh:66
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:105
void calculateFluxes(const FluxType &fluxType)
Calculate the fluxes over all surfaces for a given flux type.
Definition: fluxoveraxisalignedsurface.hh:188
FluxOverAxisAlignedSurface(const GridVariables &gridVariables, const SolutionVector &sol, const LocalResidual &localResidual, bool nonIntersectingSurfaceIsError=false)
The constructor.
Definition: fluxoveraxisalignedsurface.hh:71
An interface for a set of geometric entities.
Definition: geometricentityset.hh:160
A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
Helper functions for distance queries.
An interface for a set of geometric entities.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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:28
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:282
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:256
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:135
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:102
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:26
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Algorithms that finds which geometric entities intersect.
Detect if a point intersects a geometry.
constexpr Surface surface
Definition: couplingmanager1d3d_surface.hh:37
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.