version 3.10-dev
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 bool nonIntersectingSurfaceIsError = false)
75 : gridVariables_(gridVariables)
76 , sol_(sol)
77 , localResidual_(localResidual)
78 , nonIntersectingSurfaceIsError_(nonIntersectingSurfaceIsError)
79 {
80 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(), "FluxOverAxisAlignedSurface.Verbose", false);
81 }
82
89 template<class T>
90 void addAxisAlignedSurface(const std::string& name, T&& surface)
91 {
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))
95 ));
96 }
97
105 void addAxisAlignedSurface(const std::string& name,
106 const GlobalPosition& lowerLeft,
107 const GlobalPosition& upperRight)
108 {
109 using std::abs;
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; });
112 if (it == v.end())
113 DUNE_THROW(Dune::InvalidStateException, "Surface is not axis-parallel!");
114
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);
118 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
119
120 surfaces_.emplace(std::make_pair(
121 name,
122 SurfaceData{
123 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
124 }
125 ));
126 }
127
135 void addAxisAlignedPlane(const std::string& name,
136 const GlobalPosition& center,
137 const std::size_t normalDirectionIndex)
138 {
139 GlobalPosition lowerLeft = gridVariables_.gridGeometry().bBoxMin();
140 GlobalPosition upperRight = gridVariables_.gridGeometry().bBoxMax();
141
142 lowerLeft[normalDirectionIndex] = center[normalDirectionIndex];
143 upperRight[normalDirectionIndex] = center[normalDirectionIndex];
144
145 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
146 inSurfaceAxes.set(normalDirectionIndex, false);
147 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
148
149 surfaces_.emplace(std::make_pair(
150 name,
151 SurfaceData{
152 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
153 }
154 ));
155 }
156
161 {
162 auto fluxType = [this](const auto& element,
163 const auto& fvGeometry,
164 const auto& elemVolVars,
165 const auto& scvf,
166 const auto& elemFluxVarsCache)
167 {
168 return localResidual_.evalFlux(
169 problem_(), element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf
170 );
171 };
172
173 calculateFluxes(fluxType);
174 }
175
187 template<class FluxType>
188 void calculateFluxes(const FluxType& fluxType)
189 {
190 // make sure to reset all the values of the surfaces, in case this method has been called already before
191 for (auto& surface : surfaces_)
192 surface.second.flux = 0.0;
193
194 snapSurfaceToClosestFace_();
195 calculateFluxes_(fluxType);
196 }
197
203 const auto& flux(const std::string& name) const
204 {
205 return surfaces_.at(name).flux;
206 }
207
211 const std::map<std::string, SurfaceData>& surfaces() const
212 { return surfaces_; }
213
217 void printAllFluxes() const
218 {
219 for (const auto& [name, data] : surfaces_)
220 std::cout << "Flux over surface " << name << ": " << data.flux << std::endl;
221 }
222
226 void setNonIntersectingSurfaceIsError(bool isError = true)
227 { nonIntersectingSurfaceIsError_ = isError; }
228
229private:
230
231 template<class FluxType>
232 void calculateFluxes_(const FluxType& fluxType)
233 {
234 auto fvGeometry = localView(problem_().gridGeometry());
235 auto elemVolVars = localView(gridVariables_.curGridVolVars());
236 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
237
238 for (const auto& element : elements(problem_().gridGeometry().gridView()))
239 {
240 fvGeometry.bindElement(element);
241 elemVolVars.bindElement(element, fvGeometry, sol_);
242 elemFluxVarsCache.bindElement(element, fvGeometry, elemVolVars);
243
244 for (const auto& scvf : scvfs(fvGeometry))
245 {
246 // iterate through all surfaces and check if the flux at the given position
247 // should be accounted for in the respective surface
248 for (auto& [name, surfaceData] : surfaces_)
249 {
250 if (considerScvf_(scvf, surfaceData))
251 {
252 const auto result = fluxType(element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
253 surfaceData.flux += result;
254
255 if (verbose_)
256 std::cout << "At element " << problem_().gridGeometry().elementMapper().index(element)
257 << ": Flux at face " << scvf.ipGlobal() << ": " << result << " (" << name << ")" << std::endl;
258 }
259 }
260 }
261 }
262 }
263
265 bool considerScvf_(const SubControlVolumeFace& scvf, const SurfaceData& SurfaceData) const
266 {
267 // In order to avoid considering scvfs at the same element intersection (and hence, the corresponding flux) twice,
268 // only use those with a unit outer normal pointing towards positive coordinate direction,
269 // unless the scvf lies on a boundary (then there is no second scvf).
270 if (scvf.boundary() || !std::signbit(scvf.unitOuterNormal()[SurfaceData.normalDirectionIndex]))
271 return intersectsPointGeometry(scvf.ipGlobal(), SurfaceData.surface);
272 else
273 return false;
274 }
275
276 void snapSurfaceToClosestFace_()
277 {
278 using GeometriesEntitySet = Dumux::GeometriesEntitySet<Surface>;
279 const auto gridView = problem_().gridGeometry().gridView();
280
281 for (auto& [name, surfaceData] : surfaces_)
282 {
283 GeometriesEntitySet entitySet({surfaceData.surface});
284 Dumux::BoundingBoxTree<GeometriesEntitySet> geometriesTree(std::make_shared<GeometriesEntitySet>(entitySet));
285 const auto intersectingElements = intersectingEntities(
286 problem_().gridGeometry().boundingBoxTree(), geometriesTree
287 );
288
289 if (intersectingElements.empty())
290 {
291 if (!nonIntersectingSurfaceIsError_)
292 continue;
293
294 std::cout << "surface boundaries: " << std::endl;
295 printSurfaceBoundaries_(surfaceData.surface);
296
297 DUNE_THROW(Dune::InvalidStateException, "surface " << name << " does not intersect with any element");
298 }
299
300 std::vector<std::size_t> sortedResults;
301 sortedResults.reserve(gridView.size(0));
302
303 for (const auto& i : intersectingElements)
304 sortedResults.push_back(i.first());
305
306 std::sort(sortedResults.begin(), sortedResults.end());
307 sortedResults.erase(std::unique(
308 sortedResults.begin(), sortedResults.end()
309 ), sortedResults.end());
310
311 // pick the first intersecting element and make sure the surface snaps to the closest face with the same (or opposite facing) normal vector
312 GlobalPosition normalVector(0.0);
313 normalVector[surfaceData.normalDirectionIndex] = 1.0;
314
315 const auto& firstIntersectingElement = problem_().gridGeometry().element(sortedResults[0]);
316 Scalar distance = std::numeric_limits<Scalar>::max();
317 bool snappingOcurred = false;
318
319 GlobalPosition surfaceLowerLeft = surfaceData.surface.corner(0);
320 GlobalPosition surfaceUpperRight = surfaceData.surface.corner(3);
321
322 bool surfaceAlreadyOnFaces = false;
323 for (const auto& intersection : intersections(gridView, firstIntersectingElement))
324 {
325 if (surfaceAlreadyOnFaces)
326 continue;
327
328 using std::abs;
329 if (abs(1.0 - abs(normalVector * intersection.centerUnitOuterNormal())) < 1e-8)
330 {
331
332 const auto getDistance = [](const auto& p, const auto& geo)
333 {
334 if constexpr (dim == 2)
335 return distancePointSegment(p, geo);
336 else
337 return distancePointPolygon(p, geo);
338 };
339
340 const auto& geo = intersection.geometry();
341 if (const Scalar d = getDistance(geo.center(), surfaceData.surface); d < 1e-8 * diameter(geo))
342 {
343 // no snapping required, face already lies on surface
344 surfaceAlreadyOnFaces = true;
345 snappingOcurred = false;
346 }
347 else if (d < distance)
348 {
349 distance = d;
350 snappingOcurred = true;
351
352 // move the surface boundaries
353 for (int i = 0; i < surfaceData.surface.corners(); ++i)
354 {
355 const auto& faceCenter = geo.center();
356 surfaceLowerLeft[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
357 surfaceUpperRight[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
358 }
359 }
360 }
361 }
362
363 if (snappingOcurred)
364 {
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);
368
369 // overwrite the old surface with the new boundaries
370 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
371 inSurfaceAxes.set(surfaceData.normalDirectionIndex, false);
372 surfaceData.surface = Surface{surfaceLowerLeft, surfaceUpperRight, inSurfaceAxes};
373
374 std::cout << "New surface boundaries: " << std::endl;
375 printSurfaceBoundaries_(surfaceData.surface);
376 std::cout << std::endl;
377 }
378 }
379 }
380
381 void printSurfaceBoundaries_(const Surface& surface) const
382 {
383 for (int i = 0; i < surface.corners(); ++i)
384 std::cout << surface.corner(i) << std::endl;
385 }
386
387 const auto& problem_() const { return gridVariables_.curGridVolVars().problem(); }
388
389 std::map<std::string, SurfaceData> surfaces_;
390 const GridVariables& gridVariables_;
391 const SolutionVector& sol_;
392 const LocalResidual localResidual_; // store a copy of the local residual
393 bool verbose_;
394 bool nonIntersectingSurfaceIsError_;
395};
396
397} // end namespace Dumux
398
399#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: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
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.