3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
25#define DUMUX_FREELOW_NAVIERSTOKES_FLUX_OVER_AXISALIGNED_SURFACE_HH
26
27#include <algorithm>
28#include <type_traits>
29#include <vector>
30
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/axisalignedcubegeometry.hh>
33
40
41namespace Dumux {
42
47template<class GridVariables, class SolutionVector, class LocalResidual>
49{
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;
58
59 static constexpr auto dim = GridView::dimension;
60 static constexpr auto dimWorld = GridView::dimensionworld;
61
62 static_assert(dim > 1, "Only implemented for dim > 1");
63
64 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
65
66 // in 2D, the surface is represented as a line
67 using SurfaceT = Dune::AxisAlignedCubeGeometry<Scalar, (dim == 2 ? 1 : 2), dimWorld>;
68
69 struct SurfaceData
70 {
71 SurfaceT surface;
72 std::size_t normalDirectionIndex;
73 NumEqVector flux;
74 };
75
76public:
77
78 using Surface = SurfaceT;
79
83 FluxOverAxisAlignedSurface(const GridVariables& gridVariables,
84 const SolutionVector& sol,
85 const LocalResidual& localResidual)
86 : gridVariables_(gridVariables)
87 , sol_(sol)
88 , localResidual_(localResidual)
89 {
90 verbose_ = getParamFromGroup<bool>(problem_().paramGroup(), "FluxOverAxisAlignedSurface.Verbose", false);
91 }
92
99 template<class T>
100 void addAxisAlignedSurface(const std::string& name, T&& surface)
101 {
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))
105 ));
106 }
107
115 void addAxisAlignedSurface(const std::string& name,
116 const GlobalPosition& lowerLeft,
117 const GlobalPosition& upperRight)
118 {
119 using std::abs;
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; });
122 if (it == v.end())
123 DUNE_THROW(Dune::InvalidStateException, "Surface is not axis-parallel!");
124
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);
128 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
129
130 surfaces_.emplace(std::make_pair(
131 name,
132 SurfaceData{
133 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
134 }
135 ));
136 }
137
145 void addAxisAlignedPlane(const std::string& name,
146 const GlobalPosition& center,
147 const std::size_t normalDirectionIndex)
148 {
149 GlobalPosition lowerLeft = gridVariables_.gridGeometry().bBoxMin();
150 GlobalPosition upperRight = gridVariables_.gridGeometry().bBoxMax();
151
152 lowerLeft[normalDirectionIndex] = center[normalDirectionIndex];
153 upperRight[normalDirectionIndex] = center[normalDirectionIndex];
154
155 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
156 inSurfaceAxes.set(normalDirectionIndex, false);
157 auto surface = Surface(lowerLeft, upperRight, inSurfaceAxes);
158
159 surfaces_.emplace(std::make_pair(
160 name,
161 SurfaceData{
162 std::move(surface), normalDirectionIndex, NumEqVector(0.0)
163 }
164 ));
165 }
166
171 {
172 auto fluxType = [this](const auto& element,
173 const auto& fvGeometry,
174 const auto& elemVolVars,
175 const auto& scvf,
176 const auto& elemFluxVarsCache)
177 {
178 return localResidual_.evalFlux(
179 problem_(), element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf
180 );
181 };
182
183 calculateFluxes(fluxType);
184 }
185
197 template<class FluxType>
198 void calculateFluxes(const FluxType& fluxType)
199 {
200 // make sure to reset all the values of the surfaces, in case this method has been called already before
201 for (auto& surface : surfaces_)
202 surface.second.flux = 0.0;
203
204 snapSurfaceToClosestFace_();
205 calculateFluxes_(fluxType);
206 }
207
213 const auto& flux(const std::string& name) const
214 {
215 return surfaces_.at(name).flux;
216 }
217
221 const std::map<std::string, SurfaceData>& surfaces() const
222 { return surfaces_; }
223
227 void printAllFluxes() const
228 {
229 for (const auto& [name, data] : surfaces_)
230 std::cout << "Flux over surface " << name << ": " << data.flux << std::endl;
231 }
232
233private:
234
235 template<class FluxType>
236 void calculateFluxes_(const FluxType& fluxType)
237 {
238 auto fvGeometry = localView(problem_().gridGeometry());
239 auto elemVolVars = localView(gridVariables_.curGridVolVars());
240 auto elemFluxVarsCache = localView(gridVariables_.gridFluxVarsCache());
241
242 for (const auto& element : elements(problem_().gridGeometry().gridView()))
243 {
244 fvGeometry.bindElement(element);
245 elemVolVars.bindElement(element, fvGeometry, sol_);
246 elemFluxVarsCache.bindElement(element, fvGeometry, elemVolVars);
247
248 for (const auto& scvf : scvfs(fvGeometry))
249 {
250 // iterate through all surfaces and check if the flux at the given position
251 // should be accounted for in the respective surface
252 for (auto& [name, surfaceData] : surfaces_)
253 {
254 if (considerScvf_(scvf, surfaceData))
255 {
256 const auto result = fluxType(element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
257 surfaceData.flux += result;
258
259 if (verbose_)
260 std::cout << "At element " << problem_().gridGeometry().elementMapper().index(element)
261 << ": Flux at face " << scvf.ipGlobal() << ": " << result << " (" << name << ")" << std::endl;
262 }
263 }
264 }
265 }
266 }
267
269 bool considerScvf_(const SubControlVolumeFace& scvf, const SurfaceData& SurfaceData) const
270 {
271 // In order to avoid considering scvfs at the same element intersection (and hence, the corresponding flux) twice,
272 // only use those with a unit outer normal pointing towards positive coordinate direction,
273 // unless the scvf lies on a boundary (then there is no second scvf).
274 if (scvf.boundary() || !std::signbit(scvf.unitOuterNormal()[SurfaceData.normalDirectionIndex]))
275 return intersectsPointGeometry(scvf.ipGlobal(), SurfaceData.surface);
276 else
277 return false;
278 }
279
280 void snapSurfaceToClosestFace_()
281 {
282 using GeometriesEntitySet = Dumux::GeometriesEntitySet<Surface>;
283 const auto gridView = problem_().gridGeometry().gridView();
284
285 for (auto& [name, surfaceData] : surfaces_)
286 {
287 GeometriesEntitySet entitySet({surfaceData.surface});
288 Dumux::BoundingBoxTree<GeometriesEntitySet> geometriesTree(std::make_shared<GeometriesEntitySet>(entitySet));
289 const auto intersectingElements = intersectingEntities(
290 problem_().gridGeometry().boundingBoxTree(), geometriesTree
291 );
292
293 if (intersectingElements.empty())
294 {
295 std::cout << "surface boundaries: " << std::endl;
296 printSurfaceBoundaries_(surfaceData.surface);
297
298 DUNE_THROW(Dune::InvalidStateException, "surface " << name << " does not intersect with any element");
299 }
300
301 std::vector<std::size_t> sortedResults;
302 sortedResults.reserve(gridView.size(0));
303
304 for (const auto& i : intersectingElements)
305 sortedResults.push_back(i.first());
306
307 std::sort(sortedResults.begin(), sortedResults.end());
308 sortedResults.erase(std::unique(
309 sortedResults.begin(), sortedResults.end()
310 ), sortedResults.end());
311
312 // pick the first intersecting element and make sure the surface snaps to the closest face with the same (or opposite facing) normal vector
313 GlobalPosition normalVector(0.0);
314 normalVector[surfaceData.normalDirectionIndex] = 1.0;
315
316 const auto& firstIntersectingElement = problem_().gridGeometry().element(sortedResults[0]);
317 Scalar distance = std::numeric_limits<Scalar>::max();
318 bool snappingOcurred = false;
319
320 GlobalPosition surfaceLowerLeft = surfaceData.surface.corner(0);
321 GlobalPosition surfaceUpperRight = surfaceData.surface.corner(3);
322
323 bool surfaceAlreadyOnFaces = false;
324 for (const auto& intersection : intersections(gridView, firstIntersectingElement))
325 {
326 if (surfaceAlreadyOnFaces)
327 continue;
328
329 using std::abs;
330 if (abs(1.0 - abs(normalVector * intersection.centerUnitOuterNormal())) < 1e-8)
331 {
332
333 const auto getDistance = [](const auto& p, const auto& geo)
334 {
335 if constexpr (dim == 2)
336 return distancePointSegment(p, geo);
337 else
338 return distancePointPolygon(p, geo);
339 };
340
341 const auto& geo = intersection.geometry();
342 if (const Scalar d = getDistance(geo.center(), surfaceData.surface); d < 1e-8 * diameter(geo))
343 {
344 // no snapping required, face already lies on surface
345 surfaceAlreadyOnFaces = true;
346 snappingOcurred = false;
347 }
348 else if (d < distance)
349 {
350 distance = d;
351 snappingOcurred = true;
352
353 // move the surface boundaries
354 for (int i = 0; i < surfaceData.surface.corners(); ++i)
355 {
356 const auto& faceCenter = geo.center();
357 surfaceLowerLeft[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
358 surfaceUpperRight[surfaceData.normalDirectionIndex] = faceCenter[surfaceData.normalDirectionIndex];
359 }
360 }
361 }
362 }
363
364 if (snappingOcurred)
365 {
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);
369
370 // overwrite the old surface with the new boundaries
371 auto inSurfaceAxes = std::move(std::bitset<dimWorld>{}.set());
372 inSurfaceAxes.set(surfaceData.normalDirectionIndex, false);
373 surfaceData.surface = Surface{surfaceLowerLeft, surfaceUpperRight, inSurfaceAxes};
374
375 std::cout << "New surface boundaries: " << std::endl;
376 printSurfaceBoundaries_(surfaceData.surface);
377 std::cout << std::endl;
378 }
379 }
380 }
381
382 void printSurfaceBoundaries_(const Surface& surface) const
383 {
384 for (int i = 0; i < surface.corners(); ++i)
385 std::cout << surface.corner(i) << std::endl;
386 }
387
388 const auto& problem_() const { return gridVariables_.curGridVolVars().problem(); }
389
390 std::map<std::string, SurfaceData> surfaces_;
391 const GridVariables& gridVariables_;
392 const SolutionVector& sol_;
393 const LocalResidual localResidual_; // store a copy of the local residual
394 bool verbose_;
395};
396
397} // end namespace Dumux
398
399#endif
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.
Helper functions for distance queries.
An interface for a set of geometric entities.
Algorithms that finds which geometric entities intersect.
Detect if a point intersects a geometry.
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:40
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:294
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:268
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:147
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:114
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:38
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:36
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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 &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: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:68
An interface for a set of geometric entities.
Definition: geometricentityset.hh:121