Loading...
Searching...
No Matches
cellDistFuncs.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020,2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "cellDistFuncs.H"
30#include "polyMesh.H"
31#include "polyBoundaryMesh.H"
33#include "registerSwitch.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39defineTypeNameAndDebug(cellDistFuncs, 0);
40}
41
43
45(
46 "useCombinedWallPatch",
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
54Foam::cellDistFuncs::cellDistFuncs(const polyMesh& mesh)
56 mesh_(mesh)
57{}
58
59
60// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61
63(
65) const
66{
67 return mesh().boundaryMesh().patchSet(patchNames, false);
68}
69
70
71// size of largest patch (out of supplied subset of patches)
73(
75) const
76{
77 label maxSize = 0;
78
79 forAll(mesh().boundaryMesh(), patchi)
80 {
81 if (patchIDs.found(patchi))
82 {
83 const polyPatch& patch = mesh().boundaryMesh()[patchi];
84
85 maxSize = Foam::max(maxSize, patch.size());
86 }
87 }
88 return maxSize;
89}
90
91
92// sum of patch sizes (out of supplied subset of patches)
94(
96)
97const
98{
99 label sum = 0;
100
101 forAll(mesh().boundaryMesh(), patchi)
102 {
103 if (patchIDs.found(patchi))
104 {
105 const polyPatch& patch = mesh().boundaryMesh()[patchi];
106
107 sum += patch.size();
109 }
110 return sum;
111}
112
113
114// Gets nearest wall for cells next to wall
116(
117 const labelHashSet& patchIDs,
118 scalarField& wallDistCorrected,
119 Map<label>& nearestFace
120) const
121{
122 const auto& pbm = mesh().boundaryMesh();
123
124 // Size neighbours array for maximum possible (= size of largest patch)
125 DynamicList<label> neighbours(maxPatchSize(patchIDs));
126
127 // Correct all cells with face on wall
128 const vectorField& cellCentres = mesh().cellCentres();
129 const labelList& faceOwner = mesh().faceOwner();
130
131 forAll(pbm, patchi)
132 {
133 if (patchIDs.found(patchi))
134 {
135 const polyPatch& patch = pbm[patchi];
136 const auto areaFraction(patch.areaFraction());
137
138 // Check cells with face on wall
139 forAll(patch, patchFacei)
140 {
141 if (areaFraction && (areaFraction()[patchFacei] <= 0.5))
142 {
143 // For cyclicACMI: more cyclic than wall
144 continue;
145 }
146
147 getPointNeighbours(patch, patchFacei, neighbours);
148
149 label celli = faceOwner[patch.start() + patchFacei];
150
151 label minFacei = -1;
152
153 wallDistCorrected[celli] = smallestDist
154 (
155 cellCentres[celli],
156 patch,
157 neighbours,
158 minFacei
159 );
160
161 // Store wallCell and its nearest neighbour
162 nearestFace.insert(celli, patch.start()+minFacei);
164 }
165 }
166}
167
168
169// Correct all cells connected to wall (via point) and not in nearestFace
171(
172 const labelHashSet& patchIDs,
173 scalarField& wallDistCorrected,
174 Map<label>& nearestFace
175) const
176{
177 // Correct all (non-visited) cells with point on wall
178
179 const auto& pbm = mesh().boundaryMesh();
180 const vectorField& cellCentres = mesh().cellCentres();
181
182 forAll(pbm, patchi)
183 {
184 if (patchIDs.found(patchi))
185 {
186 const polyPatch& patch = pbm[patchi];
187 const auto& localFaces = patch.localFaces();
188 const labelList& meshPoints = patch.meshPoints();
189 const labelListList& pointFaces = patch.pointFaces();
190
191 bitSet isWallPoint(meshPoints.size(), true);
192 {
193 const auto areaFraction(patch.areaFraction());
194
195 // Check cells with face on wall
196 forAll(patch, patchFacei)
197 {
198 if (areaFraction && (areaFraction()[patchFacei] <= 0.5))
199 {
200 // For cyclicACMI: more cyclic than wall
201 isWallPoint.unset(localFaces[patchFacei]);
202 }
203 }
204 }
205
206
207 forAll(meshPoints, patchPointi)
208 {
209 const label vertI = meshPoints[patchPointi];
210
211 if (!isWallPoint[patchPointi])
212 {
213 continue;
214 }
215
216 const labelList& neighbours = mesh().pointCells(vertI);
217
218 for (const label celli : neighbours)
219 {
220 if (!nearestFace.found(celli))
221 {
222 const labelList& wallFaces = pointFaces[patchPointi];
223
224 label minFacei = -1;
225
226 wallDistCorrected[celli] = smallestDist
227 (
228 cellCentres[celli],
229 patch,
230 wallFaces,
231 minFacei
232 );
233
234 // Store wallCell and its nearest neighbour
235 nearestFace.insert(celli, patch.start()+minFacei);
236 }
238 }
239 }
240 }
241}
242
243
245(
246 const labelList& patchIDs,
247 const bool doPointCells,
248 scalarField& wallDistCorrected,
249 Map<label>& nearestFace
250) const
251{
252 label nWalls = 0;
253 {
254 for (const label patchi : patchIDs)
255 {
256 nWalls += mesh().boundaryMesh()[patchi].size();
257 }
258 }
259
260
262 {
263 for (const label patchi : patchIDs)
264 {
265 const auto& patch = mesh().boundaryMesh()[patchi];
266 forAll(patch, i)
267 {
268 faceLabels.append(patch.start()+i);
269 }
270 }
271 }
272
273 const uindirectPrimitivePatch wallPatch
274 (
276 mesh_.points()
277 );
278
279
280 // Correct all cells with face on wall
281 const vectorField& cellCentres = mesh().cellCentres();
282
283 DynamicList<label> neighbours;
284
285 nWalls = 0;
286 for (const label patchi : patchIDs)
287 {
288 const auto& patch = mesh().boundaryMesh()[patchi];
289 const auto areaFraction(patch.areaFraction());
290 const labelUList& faceCells = patch.faceCells();
291
292 // Check cells with face on wall
293 forAll(patch, patchFacei)
294 {
295 if (areaFraction && (areaFraction()[patchFacei] <= 0.5))
296 {
297 // For cyclicACMI: more cyclic than wall
298 }
299 else
300 {
301 getPointNeighbours(wallPatch, nWalls, neighbours);
302
303 const label celli = faceCells[patchFacei];
304
305 label minFacei = -1;
306 wallDistCorrected[celli] = smallestDist
307 (
308 cellCentres[celli],
309 wallPatch,
310 neighbours,
311 minFacei
312 );
313
314 // Store wallCell and its nearest neighbour
315 nearestFace.insert(celli, nWalls+minFacei);
316 }
317
318 nWalls++;
319 }
320 }
321
322 // Correct all cells with a point on the wall
323 if (doPointCells)
324 {
325 const auto& meshPoints = wallPatch.meshPoints();
326 const auto& localFaces = wallPatch.localFaces();
327
328 bitSet isWallPoint(meshPoints.size(), true);
329
330 nWalls = 0;
331 for (const label patchi : patchIDs)
332 {
333 const auto& patch = mesh().boundaryMesh()[patchi];
334 const auto areaFraction(patch.areaFraction());
335
336 // Check cells with face on wall
337 forAll(patch, patchFacei)
338 {
339 if (areaFraction && (areaFraction()[patchFacei] <= 0.5))
340 {
341 // For cyclicACMI: more cyclic than wall
342 isWallPoint.unset(localFaces[nWalls]);
343 }
344
345 nWalls++;
346 }
347 }
348
349 const auto& pointFaces = wallPatch.pointFaces();
350
351 for (const label patchPointi : isWallPoint)
352 {
353 const label verti = meshPoints[patchPointi];
354
355 const labelList& neighbours = mesh().pointCells(verti);
356
357 for (const label celli : neighbours)
358 {
359 if (!nearestFace.found(celli))
360 {
361 const labelList& wallFaces = pointFaces[patchPointi];
362
363 label minFacei = -1;
364
365 wallDistCorrected[celli] = smallestDist
366 (
367 cellCentres[celli],
368 wallPatch,
369 wallFaces,
370 minFacei
371 );
372
373 // Store wallCell and its nearest neighbour
374 nearestFace.insert(celli, wallPatch.addressing()[minFacei]);
375 }
376 }
377 }
378 }
379}
380
381
382// ************************************************************************* //
labelList faceLabels(nFaceLabels)
labelList patchIDs
const polyBoundaryMesh & pbm
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
const Addr & addressing() const noexcept
The addressing used for the list.
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const labelListList & pointFaces() const
Return point-face addressing.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A List with indirect addressing. Like IndirectList but does not store addressing.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
Definition bitSetI.H:540
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Collection of functions used in wall distance calculation.
label maxPatchSize(const labelHashSet &patchIDs) const
Size of largest patch (out of supplied subset of patches).
label sumPatchSize(const labelHashSet &patchIDs) const
Sum of patch sizes (out of supplied subset of patches).
labelHashSet getPatchIDs() const
Get patchIDs of/derived off certain type (e.g. 'processorPolyPatch').
const polyMesh & mesh() const
Access mesh.
void correctBoundaryFaceCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to boundary (via face). Sets values in.
scalar smallestDist(const point &p, const PatchType &patch, const labelUList &wallFaces, label &patchFacei) const
Calculate smallest true distance (and patch face index).
void getPointNeighbours(const PatchType &, const label patchFacei, DynamicList< label > &) const
Get faces sharing point with face on patch.
void correctBoundaryCells(const labelList &patchIDs, const bool doPointCells, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to any of the patches in patchIDs. Sets.
void correctBoundaryPointCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to wall (via point). Sets values in.
static bool useCombinedWallPatch
Use combined-wall-patches wall distance v.s. v2406 per-patch distance. Default is true.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const labelListList & pointCells() const
const vectorField & cellCentres() const
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
wordList patchNames(nPatches)
#define registerOptSwitch(Name, Type, SwitchVar)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299