Loading...
Searching...
No Matches
wallDistAddressing.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) 2023-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "FaceCellWave.H"
29#include "wallDistAddressing.H"
30#include "wallPointAddressing.H"
32#include "surfaceFields.H"
33#include "wallPolyPatch.H"
35#include "OBJstream.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
48(
53) const
54{
55 // Extract volField. See patchWave::getValues
56
57 label nIllegal = 0;
58 forAll(allCellInfo, celli)
59 {
60 const scalar dist = allCellInfo[celli].distSqr();
61
62 if (allCellInfo[celli].valid(wave.data()))
63 {
64 y[celli] = Foam::sqrt(dist);
65 }
66 else
67 {
68 y[celli] = dist;
69 ++nIllegal;
70 }
71 }
72
73 // Copy boundary values
74 auto& bfld = y.boundaryFieldRef();
75
76 for (auto& pfld : bfld)
77 {
78 scalarField patchField(pfld.size(), Zero);
79
80 forAll(pfld, patchFacei)
81 {
82 const label meshFacei = pfld.patch().start() + patchFacei;
83 const scalar dist = allFaceInfo[meshFacei].distSqr();
84
85 if (allFaceInfo[meshFacei].valid(wave.data()))
86 {
87 // Adding SMALL to avoid problems with /0 in the turbulence
88 // models
89 patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
90 }
91 else
92 {
93 patchField[patchFacei] = dist;
94 ++nIllegal;
95 }
96 }
97 pfld = patchField;
98 }
99
100 return nIllegal;
101}
102
103
105(
106 const label item,
107 const labelPair& data,
108 label& untransformi,
109 label& transformi,
110 labelPairList& transformedWallInfo
111)
112{
113 const auto& gt = mesh_.globalData().globalTransforms();
114
115 if (gt.transformIndex(data) == gt.nullTransformIndex())
116 {
117 // Store item and global index of untransformed data
118 const label proci = gt.processor(data);
119 const label index = gt.index(data);
120 untransformedItems_[untransformi] = item;
121 untransformedSlots_[untransformi++] =
122 globalWallsPtr_().toGlobal(proci, index);
123 }
124 else
125 {
126 // Store item and transformed data
127 transformedItems_[transformi] = item;
128 transformedWallInfo[transformi++] = data;
129 }
130}
131
132
134{
135 y = dimensionedScalar("yWall", dimLength, GREAT);
136
137 const auto& C = mesh_.C();
138 const auto& gt = mesh_.globalData().globalTransforms();
139
140 label nWalls = 0;
141 for (const label patchi : patchIDs_)
142 {
143 nWalls += C.boundaryField()[patchi].size();
144 }
145
146 globalWallsPtr_.reset(new globalIndex(nWalls));
147 const globalIndex& globalWalls = globalWallsPtr_();
148
149 DynamicList<label> seedFaces(nWalls);
150 DynamicList<wallPointAddressing> seedInfo(nWalls);
151
152
153 nWalls = 0;
154 for (const label patchi : patchIDs_)
155 {
156 const auto& fc = C.boundaryField()[patchi];
157 const auto areaFraction(fc.patch().patch().areaFraction());
158
159 forAll(fc, patchFacei)
160 {
161 if
162 (
163 !areaFraction
164 || (areaFraction()[patchFacei] > 0.5) // mostly wall
165 )
166 {
167 seedFaces.append(fc.patch().start()+patchFacei);
168 seedInfo.append
169 (
171 (
172 fc[patchFacei],
173 gt.encode
174 (
176 nWalls,
177 gt.nullTransformIndex()
178 ),
179 scalar(0.0)
180 )
181 );
182 }
183 nWalls++;
184 }
185 }
186
188 // Initialise the passive data (since no mesh reference in the constructor
189 // of wallPointAddressing)
190 forAll(allCellInfo, celli)
191 {
192 allCellInfo[celli].data() = gt.encode
193 (
195 celli,
196 gt.nullTransformIndex()
197 );
198 }
200 forAll(allFaceInfo, facei)
201 {
202 allFaceInfo[facei].data() = gt.encode
203 (
205 facei,
206 gt.nullTransformIndex()
207 );
208 }
209
210 // Propagate information inwards
212 (
213 mesh_,
214 seedFaces,
215 seedInfo,
218 mesh_.globalData().nTotalCells()+1
219 );
220
221
222 // Extract wall distance
223 // ~~~~~~~~~~~~~~~~~~~~~
224
225 getValues(wave, allCellInfo, allFaceInfo, y);
226
227
228 const labelHashSet patchSet(patchIDs_);
229
230 // Correct wall cells for true distance
231 // - store for corrected cells the nearest wall face
232 // - update distance field (y) for these cells
233 Map<label> cellToWallFace;
234 if (correctWalls_)
235 {
236 cellToWallFace.reserve(nWalls);
237
239 {
240 // Correct across multiple patches
241 correctBoundaryCells
242 (
243 patchIDs_,
244 true, // do point-connected cells as well
245 y,
246 cellToWallFace
247 );
248 }
249 else
250 {
251 // Optional backwards compatibility
252 correctBoundaryFaceCells
253 (
254 patchSet,
255 y,
256 cellToWallFace
257 );
258 correctBoundaryPointCells
259 (
260 patchSet,
261 y,
262 cellToWallFace
263 );
264 }
265 }
266
267 // Make sure boundary values are up to date
268 y.correctBoundaryConditions();
269
270
271 // Extract all addressing
272 // ~~~~~~~~~~~~~~~~~~~~~~
273 // - untransformed : -local cell/boundary, -nearest (global) wall
274 // - transformed : -local cell/boundary, -nearest (global) wall
275
276 label untransformi = 0;
277 untransformedItems_.resize(mesh_.nCells()+mesh_.nBoundaryFaces());
278 untransformedSlots_.resize(untransformedItems_.size());
279
280 label transformi = 0;
281 transformedItems_.resize(mesh_.nCells()+mesh_.nBoundaryFaces());
282 labelPairList transformedWallInfo(transformedItems_.size());
283
284 for (label celli = 0; celli < mesh_.nCells(); celli++)
285 {
286 const auto cellFnd = cellToWallFace.find(celli);
287 if (cellFnd)
288 {
289 const label wallFacei = cellFnd();
290 const auto& data = allFaceInfo[wallFacei].data();
291 addItem(celli, data, untransformi, transformi, transformedWallInfo);
292 }
293 else if (allCellInfo[celli].valid(wave.data()))
294 {
295 const auto& data = allCellInfo[celli].data();
296 addItem(celli, data, untransformi, transformi, transformedWallInfo);
297 }
298 }
299 untransformedPatchStarts_.resize(mesh_.boundary().size()+1);
300 transformedPatchStarts_.resize(mesh_.boundary().size()+1);
301 for (const auto& pp : mesh_.boundary())
302 {
303 untransformedPatchStarts_[pp.index()] = untransformi;
304 transformedPatchStarts_[pp.index()] = transformi;
305
306 forAll(pp, patchFacei)
307 {
308 const label facei = pp.start()+patchFacei;
309 if (allFaceInfo[facei].valid(wave.data()))
310 {
311 const auto& data = allFaceInfo[facei].data();
312 addItem
313 (
314 facei,
315 data,
316 untransformi,
317 transformi,
318 transformedWallInfo
319 );
320 }
321 }
322 }
323
324 untransformedItems_.resize(untransformi);
325 untransformedSlots_.resize(untransformi);
326 untransformedPatchStarts_.back() = untransformi;
327 transformedItems_.resize(transformi);
328 transformedWallInfo.resize(transformi);
329 transformedPatchStarts_.back() = transformi;
330
331 if (debug)
332 {
333 Pout<< typeName
334 << " : untransformed:" << untransformi
335 << " transformed:" << transformi
336 << " out of nWalls:" << nWalls
337 << " untransformStart:" << flatOutput(untransformedPatchStarts_)
338 << " transformStart:" << flatOutput(transformedPatchStarts_)
339 << endl;
340 }
341
342 List<Map<label>> compactMap;
343 mapPtr_.reset
344 (
345 new mapDistribute
346 (
347 globalWalls, // globalIndex
348 untransformedSlots_, // untransformedElements
349
350 gt, // globalIndexAndTransform
351 transformedWallInfo, // transformedElements
352 transformedSlots_, // transformedIndices
353 compactMap
354 )
355 );
356
357 if (debug & 2)
358 {
359 // Use the obtained map to write nearest
360 {
361 OBJstream os(mesh_.polyMesh::path()/"nearest.obj");
362 Pout<< typeName
363 << " : writing line from cell/face to nearest wall face to "
364 << os.name() << endl;
365
366 // Seed all wall faces with centroid and override all other
367 // cells/faces with that of nearest wall
368 volVectorField wallCentre(mesh_.C());
369 this->map(wallCentre, mapDistribute::transformPosition());
370
371 // Seed all wall faces with normal and override all other
372 // cells/faces with that of nearest wall
374 (
376 (
377 "n" & patchTypeName_,
378 mesh_.time().timeName(),
379 mesh_,
381 ),
382 mesh_,
385 );
386 // Calculate normal
387 for (const label patchi : patchIDs_)
388 {
389 auto& pnf = wallNormal.boundaryFieldRef()[patchi];
390 pnf == pnf.patch().nf();
391 }
393
394 forAll(wallCentre, celli)
395 {
396 const point& cc = mesh_.C()[celli];
397 const point& wallC = wallCentre[celli];
398 const vector& wallN = wallNormal[celli];
399 os.write(linePointRef(cc, wallC), wallN, wallN);
400 }
402 }
403}
404
405
406// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
407
409(
410 const fvMesh& mesh,
411 const bool correctWalls,
412 const label updateInterval
413)
414:
415 // Register as "wallDistAddressing"
416 MeshObject_type(mesh),
418 patchIDs_(mesh.boundaryMesh().findPatchIDs<wallPolyPatch>().sortedToc()),
419 patchTypeName_("wall"),
420 updateInterval_(updateInterval),
421 correctWalls_(correctWalls),
422 requireUpdate_(true),
423 y_
424 (
426 (
427 "y" & patchTypeName_,
428 mesh.time().timeName(),
429 mesh.thisDb(),
430 IOobjectOption::NO_REGISTER
431 ),
432 mesh,
433 dimensionedScalar("y", dimLength, GREAT) // Same as wallDistData
434 )
435{
436 movePoints();
437}
438
439
441(
442 const word& patchTypeName,
443 const fvMesh& mesh,
444 const labelList& patchIDs,
445 const bool correctWalls,
446 const label updateInterval
447)
448:
449 MeshObject_type(patchTypeName, mesh),
451 patchIDs_(patchIDs),
452 patchTypeName_(patchTypeName),
453 updateInterval_(updateInterval),
454 correctWalls_(correctWalls),
455 requireUpdate_(true),
456 y_
457 (
459 (
460 "y" & patchTypeName_,
461 mesh.time().timeName(),
462 mesh.thisDb(),
463 IOobjectOption::NO_REGISTER
464 ),
465 mesh,
466 dimensionedScalar("y", dimLength, GREAT) // Same as wallDistData
467 )
469 movePoints();
470}
471
472
473// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
476{} // mapDistribute was forward declared
477
478
479// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
480
482{
483 if
484 (
485 (updateInterval_ > 0)
486 && ((mesh_.time().timeIndex() % updateInterval_) == 0)
487 )
488 {
489 requireUpdate_ = true;
490 }
491
492 if (requireUpdate_)
493 {
494 DebugInfo<< "Updating wall distance" << endl;
495
496 requireUpdate_ = false;
497
498 correct(y_);
499
500 return true;
501 }
502
503 return false;
504}
505
506
507void Foam::wallDistAddressing::updateMesh(const mapPolyMesh& mpm)
508{
509 // Force update if performing topology change
510 // Note: needed?
511 // - field would have been mapped, so if using updateInterval option (!= 1)
512 // live with error associated of not updating and use mapped values?
513 requireUpdate_ = true;
514 movePoints();
515}
516
517
518// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar y
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelList patchIDs
Graphite solid properties.
Definition C.H:49
C() noexcept
Default construct.
Definition C.C:36
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
const TrackingData & data() const noexcept
Additional data to be passed into container.
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
void resize(label newCapacity)
Rehash the hash table with new number of buckets. Currently identical to setCapacity().
Definition HashTable.C:722
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
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.
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.
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.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Default transformation behaviour for position.
Default transformation behaviour.
Class containing processor-to-processor mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static wordList patchTypes(const fvMesh &mesh, const labelHashSet &patchIDs)
Return the patch types for y and n.
Variant of wallDist that uses meshWave and stores the addressing.
virtual bool movePoints()
Update the y-field when the mesh moves.
const word patchTypeName_
Name for the patch set, e.g. "wall".
labelList untransformedPatchStarts_
Start of patches. Start of untransformedPatchStarts_[0] is end.
void correct(volScalarField &y)
Extract nearest-patch distance data.
volScalarField y_
Distance-to-wall field.
void addItem(const label item, const labelPair &data, label &untransformi, label &transformi, labelPairList &transformedWallInfo)
Store nearest-data to cell or boundary face.
const VolField & map(VolField &fld, const TransformOp &top=mapDistribute::transform()) const
Map nearest-patch information. Take wall patch values.
label getValues(const FaceCellWave< wallPointAddressing > &wave, const List< wallPointAddressing > &allCellInfo, const List< wallPointAddressing > &allFaceInfo, volScalarField &y) const
Extract FaceCellWave data.
const bool correctWalls_
Do accurate distance calculation for near-wall cells.
autoPtr< globalIndex > globalWallsPtr_
Number of wall faces.
labelList untransformedSlots_
Corresponding slot in mapPtr distribution result.
const labelUList & patchIDs() const noexcept
Return the patchIDs.
wallDistAddressing(const wallDistAddressing &)=delete
No copy construct.
virtual ~wallDistAddressing()
Destructor.
bool requireUpdate_
Flag to indicate whether the wall distance requires updating.
autoPtr< mapDistribute > mapPtr_
Map to pull wall face info to cell or boundary face.
labelList untransformedItems_
Indices of cells followed by boundary faces.
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
const label updateInterval_
Update wall distance every updateInterval_ steps.
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field. Unvisited.
const labelList patchIDs_
Set of patch IDs.
Holds information (coordinate and origin) regarding nearest wall point.
Foam::wallPolyPatch.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
thermo correct()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
vector wallNormal(Zero)
List< wallPoints > allFaceInfo(mesh_.nFaces())
List< wallPoints > allCellInfo(mesh_.nCells())
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to 'true' entries.
Definition BitOps.C:200
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.