version 3.8
fluxoveraxisalignedsurface.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
13#define DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
14
15#include <algorithm>
16#include <type_traits>
17#include <vector>
18
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/axisalignedcubegeometry.hh>
21
28
29namespace Dumux {
30
35template<class GridVariables, class SolutionVector, class LocalResidual>
37{
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;
46
47 static constexpr auto dim = GridView::dimension;
48 static constexpr auto dimWorld = GridView::dimensionworld;
49
50 static_assert(dim > 1, "Only implemented for dim > 1");
51
52 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
53
54 // in 2D, the surface is represented as a line
55 using SurfaceT = Dune::AxisAlignedCubeGeometry<Scalar, (dim == 2 ? 1 : 2), dimWorld>;
56
57 struct SurfaceData
58 {
59 SurfaceT surface;
60 std::size_t normalDirectionIndex;
61 NumEqVector flux;
62 };
63
64public:
65
66 using Surface = SurfaceT;
67
71 FluxOverAxisAlignedSurface(const GridVariables& gridVariables,
72 const SolutionVector& sol,
73 const LocalResidual& localResidual)
74 : gridVariables_(gridVariables)
75 , sol_(sol)
76 , localResidual_(localResidual)
77 {
78 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(), "FluxOverAxisAlignedSurface.Verbose", false);
79 }
80
87 template<class T>
88 void addAxisAlignedSurface(const std::string& name, T&& surface)
89 {
90 static_assert(std::is_same_v<std::decay_t<T>, Surface>);
91 surfaces_.emplace(std::make_pair(
92 name, std::make_pair(surface, NumEqVector(0.0))
93 ));
94 }
95
103 void addAxisAlignedSurface(const std::string& name,
104 const GlobalPosition& lowerLeft,
105 const GlobalPosition& upperRight)
106 {
107 using std::abs;
108 const GlobalPosition v = upperRight - lowerLeft;
109 const auto it = std::find_if(v.begin(), v.end(), [](const auto& x){ return abs(x) < 1e-20; });
110 if (it == v.end())
111 DUNE_THROW(Dune::InvalidStateException, "Surface is not axis-parallel!");
112
113 const std::size_t normalDirectionIndex = std::distance(v.begin(), it);
114 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
115 inSurfaceAxes.set(normalDirectionIndex, false);
116 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
117
118 surfaces_.emplace(std::make_pair(
119 name,
120 SurfaceData{
121 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
122 }
123 ));
124 }
125
133 void addAxisAlignedPlane(const std::string& name,
134 const GlobalPosition& center,
135 const std::size_t normalDirectionIndex)
136 {
137 GlobalPosition lowerLeft = gridVariables_.gridGeometry().bBoxMin();
138 GlobalPosition upperRight = gridVariables_.gridGeometry().bBoxMax();
139
140 lowerLeft[normalDirectionIndex] = center[normalDirectionIndex];
141 upperRight[normalDirectionIndex] = center[normalDirectionIndex];
142
143 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
144 inSurfaceAxes.set(normalDirectionIndex, false);
145 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
146
147 surfaces_.emplace(std::make_pair(
148 name,
149 SurfaceData{
150 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
151 }
152 ));
153 }
154
159 {
160 auto fluxType = [this](const auto& element,
161 const auto& fvGeometry,
162 const auto& elemVolVars,
163 const auto& scvf,
164 const auto& elemFluxVarsCache)
165 {
166 return localResidual_.evalFlux(
167 problem_(), element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf
168 );
169 };
170
171 calculateFluxes(fluxType);
172 }
173
185 template<class FluxType>
186 void calculateFluxes(const FluxType& fluxType)
187 {
188 // make sure to reset all the values of the surfaces, in case this method has been called already before
189 for (auto& surface : surfaces_)
190 surface.second.flux = 0.0;
191
192 snapSurfaceToClosestFace_();
193 calculateFluxes_(fluxType);
194 }
195
201 const auto& flux(const std::string& name) const
202 {
203 return surfaces_.at(name).flux;
204 }
205
209 const std::map<std::string, SurfaceData>& surfaces() const
210 { return surfaces_; }
211
215 void printAllFluxes() const
216 {
217 for (const auto& [name, data] : surfaces_)
218 std::cout << "Flux over surface " << name << ": " << data.flux << std::endl;
219 }
220
221private:
222
223 template<class FluxType>
224 void calculateFluxes_(const FluxType& fluxType)
225 {
226 auto fvGeometry = localView(problem_().gridGeometry());
227 auto elemVolVars = localView(gridVariables_.curGridVolVars());
228 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
229
230 for (const auto& element : elements(problem_().gridGeometry().gridView()))
231 {
232 fvGeometry.bindElement(element);
233 elemVolVars.bindElement(element, fvGeometry, sol_);
234 elemFluxVarsCache.bindElement(element, fvGeometry, elemVolVars);
235
236 for (const auto& scvf : scvfs(fvGeometry))
237 {
238 // iterate through all surfaces and check if the flux at the given position
239 // should be accounted for in the respective surface
240 for (auto& [name, surfaceData] : surfaces_)
241 {
242 if (considerScvf_(scvf, surfaceData))
243 {
244 const auto result = fluxType(element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
245 surfaceData.flux += result;
246
247 if (verbose_)
248 std::cout << "At element " << problem_().gridGeometry().elementMapper().index(element)
249 << ": Flux at face " << scvf.ipGlobal() << ": " << result << " (" << name << ")" << std::endl;
250 }
251 }
252 }
253 }
254 }
255
257 bool considerScvf_(const SubControlVolumeFace& scvf, const SurfaceData& SurfaceData) const
258 {
259 // In order to avoid considering scvfs at the same element intersection (and hence, the corresponding flux) twice,
260 // only use those with a unit outer normal pointing towards positive coordinate direction,
261 // unless the scvf lies on a boundary (then there is no second scvf).
262 if (scvf.boundary() || !std::signbit(scvf.unitOuterNormal()[SurfaceData.normalDirectionIndex]))
263 return intersectsPointGeometry(scvf.ipGlobal(), SurfaceData.surface);
264 else
265 return false;
266 }
267
268 void snapSurfaceToClosestFace_()
269 {
270 using GeometriesEntitySet = Dumux::GeometriesEntitySet<Surface>;
271 const auto gridView = problem_().gridGeometry().gridView();
272
273 for (auto& [name, surfaceData] : surfaces_)
274 {
275 GeometriesEntitySet entitySet({surfaceData.surface});
276 Dumux::BoundingBoxTree<GeometriesEntitySet> geometriesTree(std::make_shared<GeometriesEntitySet>(entitySet));
277 const auto intersectingElements = intersectingEntities(
278 problem_().gridGeometry().boundingBoxTree(), geometriesTree
279 );
280
281 if (intersectingElements.empty())
282 {
283 std::cout << "surface boundaries: " << std::endl;
284 printSurfaceBoundaries_(surfaceData.surface);
285
286 DUNE_THROW(Dune::InvalidStateException, "surface " << name << " does not intersect with any element");
287 }
288
289 std::vector<std::size_t> sortedResults;
290 sortedResults.reserve(gridView.size(0));
291
292 for (const auto& i : intersectingElements)
293 sortedResults.push_back(i.first());
294
295 std::sort(sortedResults.begin(), sortedResults.end());
296 sortedResults.erase(std::unique(
297 sortedResults.begin(), sortedResults.end()
298 ), sortedResults.end());
299
300 // pick the first intersecting element and make sure the surface snaps to the closest face with the same (or opposite facing) normal vector
301 GlobalPosition normalVector(0.0);
302 normalVector[surfaceData.normalDirectionIndex] = 1.0;
303
304 const auto& firstIntersectingElement = problem_().gridGeometry().element(sortedResults[0]);
305 Scalar distance = std::numeric_limits<Scalar>::max();
306 bool snappingOcurred = false;
307
308 GlobalPosition surfaceLowerLeft = surfaceData.surface.corner(0);
309 GlobalPosition surfaceUpperRight = surfaceData.surface.corner(3);
310
311 bool surfaceAlreadyOnFaces = false;
312 for (const auto& intersection : intersections(gridView, firstIntersectingElement))
313 {
314 if (surfaceAlreadyOnFaces)
315 continue;
316
317 using std::abs;
318 if (abs(1.0 - abs(normalVector * intersection.centerUnitOuterNormal())) < 1e-8)
319 {
320
321 const auto getDistance = [](const auto& p, const auto& geo)
322 {
323 if constexpr (dim == 2)
324 return distancePointSegment(p, geo);
325 else
326 return distancePointPolygon(p, geo);
327 };
328
329 const auto& geo = intersection.geometry();
330 if (const Scalar d = getDistance(geo.center(), surfaceData.surface); d < 1e-8 * diameter(geo))
331 {
332 // no snapping required, face already lies on surface
333 surfaceAlreadyOnFaces = true;
334 snappingOcurred = false;
335 }
336 else if (d < distance)
337 {
338 distance = d;
339 snappingOcurred = true;
340
341 // move the surface boundaries
342 for (int i = 0; i < surfaceData.surface.corners(); ++i)
343 {
344 const auto& faceCenter = geo.center();
345 surfaceLowerLeft[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
346 surfaceUpperRight[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
347 }
348 }
349 }
350 }
351
352 if (snappingOcurred)
353 {
354 std::cout << "\n\nSurface '" << name << "' was automatically snapped to the closest faces" << std::endl;
355 std::cout << "Old surface boundaries: " << std::endl;
356 printSurfaceBoundaries_(surfaceData.surface);
357
358 // overwrite the old surface with the new boundaries
359 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
360 inSurfaceAxes.set(surfaceData.normalDirectionIndex, false);
361 surfaceData.surface = Surface{surfaceLowerLeft, surfaceUpperRight, inSurfaceAxes};
362
363 std::cout << "New surface boundaries: " << std::endl;
364 printSurfaceBoundaries_(surfaceData.surface);
365 std::cout << std::endl;
366 }
367 }
368 }
369
370 void printSurfaceBoundaries_(const Surface& surface) const
371 {
372 for (int i = 0; i < surface.corners(); ++i)
373 std::cout << surface.corner(i) << std::endl;
374 }
375
376 const auto& problem_() const { return gridVariables_.curGridVolVars().problem(); }
377
378 std::map<std::string, SurfaceData> surfaces_;
379 const GridVariables& gridVariables_;
380 const SolutionVector& sol_;
381 const LocalResidual localResidual_; // store a copy of the local residual
382 bool verbose_;
383};
384
385} // end namespace Dumux
386
387#endif
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 &center, 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:133
FluxOverAxisAlignedSurface(const GridVariables &gridVariables, const SolutionVector &sol, const LocalResidual &localResidual)
The constructor.
Definition: fluxoveraxisalignedsurface.hh:71
void printAllFluxes() const
Prints all fluxes.
Definition: fluxoveraxisalignedsurface.hh:215
void addAxisAlignedSurface(const std::string &name, T &&surface)
Add an axis-aligned surface with a given name.
Definition: fluxoveraxisalignedsurface.hh:88
const std::map< std::string, SurfaceData > & surfaces() const
Provides access to all surfaces.
Definition: fluxoveraxisalignedsurface.hh:209
const auto & flux(const std::string &name) const
Return the flux over given surface.
Definition: fluxoveraxisalignedsurface.hh:201
void calculateAllFluxes()
Calculate the fluxes over all surfaces.
Definition: fluxoveraxisalignedsurface.hh:158
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:103
void calculateFluxes(const FluxType &fluxType)
Calculate the fluxes over all surfaces for a given flux type.
Definition: fluxoveraxisalignedsurface.hh:186
An interface for a set of geometric entities.
Definition: geometricentityset.hh:109
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
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.