version 3.8
walldistance.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_DISCRETIZATION_WALL_DISTANCE_HH
13#define DUMUX_DISCRETIZATION_WALL_DISTANCE_HH
14
15#include <vector>
16
17#include <dune/common/parallel/mpihelper.hh>
18#include <dune/common/shared_ptr.hh>
19#include <dune/common/reservedvector.hh>
20#include <dune/grid/common/partitionset.hh>
21
22#include <dumux/common/tag.hh>
26
27namespace Dumux {
28
35template<class GridGeometry, template<class> class DistanceField = AABBDistanceField>
37{
38 using GridView = typename GridGeometry::GridView;
39 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
40 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
41 using FVElementGeometry = typename GridGeometry::LocalView;
42 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
43 using Scalar = typename GridView::Grid::ctype;
44 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
45
46 static constexpr int dim = GridView::dimension;
47 static constexpr int dimWorld = GridView::dimensionworld;
48 static_assert(dim == dimWorld, "Wall distances cannot be computed for embedded surface or network domains.");
49
50 using CornerStorage = Dune::ReservedVector<GlobalPosition, (1<<(GridView::dimension-1))>;
51
52 // We use a simplified geometry type here which allows much easier MPI communication
53 // for parallel runs compared to the Dune geometry types (due to fixed-size storage).
54 // This type only implements a minimal part of the Geometry interface.
55 struct SimpleGeometry
56 {
57 SimpleGeometry() = default;
58 SimpleGeometry(CornerStorage&& corners)
59 : corners_(std::move(corners))
60 , center_(0.0)
61 {
62 for (int i = 0; i < corners_.size(); ++i)
63 center_ += corners_[i];
64 center_ /= corners_.size();
65 }
66
67 using GlobalCoordinate = GlobalPosition;
68 using ctype = typename GlobalCoordinate::value_type;
69 static constexpr int coorddimension = GridView::dimensionworld;
70 static constexpr int mydimension = GridView::dimension-1;
71
72 std::size_t corners() const
73 { return corners_.size(); }
74
75 const auto& corner(int i) const
76 { return corners_[i]; }
77
78 const auto& center() const
79 { return center_; }
80
81 private:
82 CornerStorage corners_;
83 GlobalCoordinate center_;
84 };
85
86 struct WallDataImpl
87 {
88 GridIndexType eIdx;
89 GridIndexType scvfIdx;
90 GlobalPosition scvfOuterNormal;
91 int rank; // for parallel runs
92 };
93
94 struct ElementCenters {};
95 struct VertexCenters {};
96
97public:
100 using WallData = WallDataImpl;
101
108 template<class LocationTag, class ScvfSelectionFunctor>
109 WallDistance(std::shared_ptr<const GridGeometry> gridGeometry, LocationTag tag, const ScvfSelectionFunctor& select)
110 : gridGeometry_(gridGeometry)
111 {
112 initializeWallDistance_(select, tag);
113 }
114
121 template<class LocationTag>
122 WallDistance(std::shared_ptr<const GridGeometry> gridGeometry, LocationTag tag)
123 : WallDistance(gridGeometry, tag, [](const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) { return true; }) {}
124
126 template<class LocationTag, class ScvfSelectionFunctor>
127 WallDistance(const GridGeometry& gridGeometry, LocationTag tag, const ScvfSelectionFunctor& select)
128 : WallDistance(Dune::stackobject_to_shared_ptr(gridGeometry), tag, select) {}
129
131 template<class LocationTag>
132 WallDistance(const GridGeometry& gridGeometry, LocationTag tag)
133 : WallDistance(Dune::stackobject_to_shared_ptr(gridGeometry), tag) {}
134
140 const std::vector<Scalar>& wallDistance() const
141 { return distance_; }
142
148 const std::vector<WallData>& wallData() const
149 { return wallData_; }
150
151private:
158 template<class ConsiderFaceFunction, class LocationTag>
159 void initializeWallDistance_(const ConsiderFaceFunction& considerFace, LocationTag loc)
160 {
161 const std::size_t numSamplingPoints = (loc == atElementCenters)
162 ? gridGeometry_->gridView().size(0)
163 : gridGeometry_->gridView().size(dim);
164 // Reset the containers.
165 wallData_.resize(numSamplingPoints);
166 distance_.resize(numSamplingPoints, std::numeric_limits<Scalar>::max());
167
168 std::vector<SimpleGeometry> wallGeometries;
169 wallGeometries.reserve(gridGeometry_->numBoundaryScvf());
170
171 std::vector<WallData> tempWallData;
172 tempWallData.reserve(gridGeometry_->numBoundaryScvf());
173
174 // Loop over all elements: find all wall scvfs.
175 auto fvGeometry = localView(*gridGeometry_);
176 for (const auto& element : elements(gridGeometry_->gridView(), Dune::Partitions::interior))
177 {
178 fvGeometry.bindElement(element);
179 if (!fvGeometry.hasBoundaryScvf())
180 continue;
181
182 const auto eIdx = gridGeometry_->elementMapper().index(element);
183
184 for (const auto& scvf : scvfs(fvGeometry))
185 {
186 if (scvf.boundary() && considerFace(fvGeometry, scvf))
187 {
188 const auto& geo = fvGeometry.geometry(scvf);
189 CornerStorage corners;
190 for (int i = 0; i < geo.corners(); ++i)
191 corners.push_back(geo.corner(i));
192
193 wallGeometries.emplace_back(std::move(corners));
194 tempWallData.push_back(WallData{
195 eIdx, scvf.index(), scvf.unitOuterNormal(), gridGeometry_->gridView().comm().rank()
196 });
197 }
198 }
199 }
200
201#if HAVE_MPI
202 // Handle parallel runs. We need to prepare a global vector of wall geometries,
203 // containing the wall geometries of each process in order to get a correct distance field.
204 const bool isParallel = gridGeometry_->gridView().comm().size() > 1;
205 std::vector<SimpleGeometry> globalWallGeometries;
206 std::vector<WallData> globalTempWallData;
207 const auto distanceField = [&]
208 {
209 if (isParallel)
210 {
211 const auto& communication = gridGeometry_->gridView().comm();
212 const int totalNumberOfBoundaryGeometries = communication.sum(wallGeometries.size());
213 globalWallGeometries.resize(totalNumberOfBoundaryGeometries);
214 globalTempWallData.resize(totalNumberOfBoundaryGeometries);
215
216 // prepare a displacement vector
217 std::vector<int> numGeosPerProcLocal{static_cast<int>(wallGeometries.size())};
218 std::vector<int> numGeosPerProcGlobal(communication.size());
219 communication.allgather(numGeosPerProcLocal.data(), 1, numGeosPerProcGlobal.data());
220
221 std::vector<int> disp(communication.size(), 0);
222 disp[1] = numGeosPerProcGlobal[0];
223 for (int i = 2; i < numGeosPerProcGlobal.size(); ++i)
224 disp[i] = disp[i-1] + numGeosPerProcGlobal[i-1];
225
226 // concatenate the wall geometries and temp scvf data of each process into a global vector
227 communication.allgatherv(
228 wallGeometries.data(),
229 wallGeometries.size(),
230 globalWallGeometries.data(),
231 numGeosPerProcGlobal.data(),
232 disp.data()
233 );
234
235 communication.allgatherv(
236 tempWallData.data(),
237 tempWallData.size(),
238 globalTempWallData.data(),
239 numGeosPerProcGlobal.data(),
240 disp.data()
241 );
242
243 // pass the global vector of wall geometries to the distance field
244 return DistanceField<SimpleGeometry>(globalWallGeometries);
245 }
246 else
247 return DistanceField<SimpleGeometry>(wallGeometries);
248 }();
249#else
250 const DistanceField<SimpleGeometry> distanceField(wallGeometries);
251#endif
252
253 // compute sampling points
254 std::vector<GlobalPosition> points(numSamplingPoints);
255 if (loc == atElementCenters)
256 for (const auto& element : elements(gridGeometry_->gridView()))
257 points[gridGeometry_->elementMapper().index(element)] = element.geometry().center();
258 else
259 for (const auto& vertex : vertices(gridGeometry_->gridView()))
260 points[gridGeometry_->vertexMapper().index(vertex)] = vertex.geometry().corner(0);
261
262 // get the actual distances (this is the most expensive part)
263 if (loc == atElementCenters)
264 {
265 const auto kernel = [&](std::size_t eIdx){
266 const auto [d, idx] = distanceField.distanceAndIndex(points[eIdx]);
267 distance_[eIdx] = d;
268#if HAVE_MPI
269 wallData_[eIdx] = isParallel ? globalTempWallData[idx] : tempWallData[idx];
270#else
271 wallData_[eIdx] = tempWallData[idx];
272#endif
273 };
274
275 runKernel_(numSamplingPoints, kernel);
276 }
277 else
278 {
279 const auto kernel = [&](std::size_t vIdx){
280 const auto [d, idx] = distanceField.distanceAndIndex(points[vIdx]);
281 distance_[vIdx] = d;
282#if HAVE_MPI
283 wallData_[vIdx] = isParallel ? globalTempWallData[idx] : tempWallData[idx];
284#else
285 wallData_[vIdx] = tempWallData[idx];
286#endif
287 };
288
289 runKernel_(numSamplingPoints, kernel);
290 }
291 }
292
293 template<class Kernel>
294 void runKernel_(std::size_t size, const Kernel& kernel)
295 {
296 // parallelize, if we have enough work (enough evaluation points)
297 if (size > 10000)
298 Dumux::parallelFor(size, [&](const std::size_t i){ kernel(i); });
299 else
300 for (std::size_t i = 0; i < size; ++i) kernel(i);
301 }
302
303 std::vector<Scalar> distance_;
304 std::vector<WallData> wallData_;
305 std::shared_ptr<const GridGeometry> gridGeometry_;
306};
307
308} // end namespace Dumux
309
310#endif
Class to calculate the wall distance at every element or vertex of a grid.
Definition: walldistance.hh:37
static constexpr auto atVertexCenters
Definition: walldistance.hh:99
const std::vector< Scalar > & wallDistance() const
Returns a vector storing the distance from each DOF location to the nearest wall. For the atElementCe...
Definition: walldistance.hh:140
WallDistance(std::shared_ptr< const GridGeometry > gridGeometry, LocationTag tag)
Constructs a new wall distance object.
Definition: walldistance.hh:122
WallDistance(const GridGeometry &gridGeometry, LocationTag tag, const ScvfSelectionFunctor &select)
caller has to make sure the lifetime of grid geometry exceeds the lifetime of wall distance
Definition: walldistance.hh:127
WallDataImpl WallData
Definition: walldistance.hh:100
static constexpr auto atElementCenters
Definition: walldistance.hh:98
WallDistance(std::shared_ptr< const GridGeometry > gridGeometry, LocationTag tag, const ScvfSelectionFunctor &select)
Constructs a new wall distance object.
Definition: walldistance.hh:109
const std::vector< WallData > & wallData() const
Returns a vector storing additional information about the nearest scvf on the wall (element index and...
Definition: walldistance.hh:148
WallDistance(const GridGeometry &gridGeometry, LocationTag tag)
caller has to make sure the lifetime of grid geometry exceeds the lifetime of wall distance
Definition: walldistance.hh:132
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
AABBDistanceField< Geometry > DistanceField
Class to calculate the closest distance from a point to a given set of geometries describing the doma...
Definition: distancefield.hh:117
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
Defines the index types used for grid and local indices.
constexpr Kernel kernel
Definition: couplingmanager1d3d_kernel.hh:38
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Parallel for loop (multithreading)
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
Helper class to create (named and comparable) tagged types Tags any given type. The tagged type is eq...
Definition: tag.hh:30
Helper class to create (named and comparable) tagged types.