Loading...
Searching...
No Matches
patchDataWave.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 "patchDataWave.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33template<class TransferType, class TrackingData>
34int Foam::patchDataWave<TransferType, TrackingData>::dummyTrackData_ = 12345;
35
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39// Set initial set of changed faces (= all wall faces)
40template<class TransferType, class TrackingData>
41void Foam::patchDataWave<TransferType, TrackingData>::setChangedFaces
42(
43 const labelHashSet& patchIDs,
44 DynamicList<label>& changedFaces,
45 DynamicList<TransferType>& faceDist
46) const
47{
48 const polyMesh& mesh = cellDistFuncs::mesh();
49
50 forAll(mesh.boundaryMesh(), patchi)
51 {
52 if (patchIDs.found(patchi))
53 {
54 const polyPatch& patch = mesh.boundaryMesh()[patchi];
55
56 const auto areaFraction(patch.areaFraction());
57
58 const auto faceCentres(patch.faceCentres());
59
60 const Field<Type>& patchField = initialPatchValuePtrs_[patchi];
61
62 forAll(faceCentres, patchFacei)
63 {
64 if
65 (
66 !areaFraction
67 || (areaFraction()[patchFacei] > 0.5) // mostly wall
68 )
69 {
70 label meshFacei = patch.start() + patchFacei;
71 changedFaces.append(meshFacei);
72
73 faceDist.append
74 (
75 TransferType
76 (
77 faceCentres[patchFacei],
78 patchField[patchFacei],
79 0.0
80 )
81 );
82 }
83 }
84 }
85 }
86}
87
88
89// Copy from MeshWave data into *this (distance) and field_ (transported data)
90template<class TransferType, class TrackingData>
91Foam::label Foam::patchDataWave<TransferType, TrackingData>::getValues
92(
93 const MeshWave<TransferType, TrackingData>& waveInfo
94)
95{
96 const polyMesh& mesh = cellDistFuncs::mesh();
97
98 const List<TransferType>& cellInfo = waveInfo.allCellInfo();
99 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
100
101 label nIllegal = 0;
102
103 // Copy cell values
104 distance_.setSize(cellInfo.size());
105
106 forAll(cellInfo, celli)
107 {
108 const TransferType & wpn = cellInfo[celli];
109
110 scalar dist = wpn.distSqr();
111
112 if (cellInfo[celli].valid(waveInfo.data()))
113 {
114 distance_[celli] = Foam::sqrt(dist);
115
116 cellData_[celli] = cellInfo[celli].data();
117 }
118 else
119 {
120 // Illegal/unset value. What to do with data?
121 // Note: mag for now. Should maybe be member of TransferType?
122
123 distance_[celli] = mag(dist);
124
125 //cellData_[celli] = point::max;
126 cellData_[celli] = cellInfo[celli].data();
127
128 nIllegal++;
129 }
130 }
131
132 // Copy boundary values
133 forAll(patchDistance_, patchi)
134 {
135 const polyPatch& patch = mesh.boundaryMesh()[patchi];
136
137 // Allocate storage for patchDistance
138 scalarField* patchFieldPtr = new scalarField(patch.size());
139
140 patchDistance_.set(patchi, patchFieldPtr);
141
142 scalarField& patchField = *patchFieldPtr;
143
144 // Allocate storage for patchData
145 Field<Type>* patchDataFieldPtr = new Field<Type>(patch.size());
146
147 patchData_.set(patchi, patchDataFieldPtr);
148
149 Field<Type>& patchDataField = *patchDataFieldPtr;
150
151 // Copy distance and data
152 forAll(patchField, patchFacei)
153 {
154 label meshFacei = patch.start() + patchFacei;
155
156 scalar dist = faceInfo[meshFacei].distSqr();
157
158 if (faceInfo[meshFacei].valid(waveInfo.data()))
159 {
160 // Adding SMALL to avoid problems with /0 in the turbulence
161 // models
162 patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
163
164 patchDataField[patchFacei] = faceInfo[meshFacei].data();
165 }
166 else
167 {
168 // Illegal/unset value. What to do with data?
169
170 patchField[patchFacei] = mag(dist);
171
172 //patchDataField[patchFacei] = point::max;
173 patchDataField[patchFacei] = faceInfo[meshFacei].data();
174
175 nIllegal++;
176 }
177 }
178 }
179
180 return nIllegal;
182
183
184// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
185
186// Construct from components
187template<class TransferType, class TrackingData>
189(
190 const polyMesh& mesh,
191 const labelHashSet& patchIDs,
192 const UPtrList<Field<Type>>& initialPatchValuePtrs,
193 const bool correctWalls,
194 TrackingData& td
195)
196:
198 patchIDs_(patchIDs),
199 initialPatchValuePtrs_(initialPatchValuePtrs),
200 correctWalls_(correctWalls),
201 td_(td),
202 nUnset_(0),
203 distance_(mesh.nCells()),
204 patchDistance_(mesh.boundaryMesh().size()),
205 cellData_(mesh.nCells()),
206 patchData_(mesh.boundaryMesh().size())
207{
210
211
212// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213
214template<class TransferType, class TrackingData>
217
218
219// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
220
221// Correct for mesh geom/topo changes
222template<class TransferType, class TrackingData>
224{
225 //
226 // Set initial changed faces: set TransferType for wall faces
227 // to wall centre.
228 //
229
230 // Count walls
231 label nWalls = sumPatchSize(patchIDs_);
232
233 DynamicList<TransferType> faceDist(nWalls);
234 DynamicList<label> changedFaces(nWalls);
235
236 setChangedFaces(patchIDs_, changedFaces, faceDist);
237
238 //
239 // Do calculate wall distance by 'growing' from faces.
240 //
241
243 (
244 mesh(),
245 changedFaces,
246 faceDist,
247 mesh().globalData().nTotalCells()+1, // max iterations
248 td_
249 );
250
251
252 //
253 // Copy distance into return field
254 //
255
256 nUnset_ = getValues(waveInfo);
257
258 //
259 // Correct wall cells for true distance
260 //
261
262 if (correctWalls_)
263 {
264 Map<label> nearestFace(2 * nWalls);
265
267 {
268 // Correct across multiple patches
269 correctBoundaryCells
270 (
271 patchIDs_.sortedToc(),
272 true, // do point-connected cells as well
273 distance_,
274 nearestFace
275 );
276 }
277 else
278 {
279 // Get distance and indices of nearest face
280 correctBoundaryFaceCells
281 (
282 patchIDs_,
283 distance_,
284 nearestFace
285 );
286
287 correctBoundaryPointCells
288 (
289 patchIDs_,
290 distance_,
291 nearestFace
292 );
293 }
294
295
296 // Transfer data from nearest face to cell
297 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
298
299 const labelList wallCells(nearestFace.toc());
300
301 forAll(wallCells, wallCelli)
302 {
303 label celli = wallCells[wallCelli];
304
305 label facei = nearestFace[celli];
306
307 cellData_[celli] = faceInfo[facei].data();
308 }
309 }
310}
311
312
313// ************************************************************************* //
labelList patchIDs
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
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
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
FaceCellWave plus data.
Definition MeshWave.H:58
const List< Type > & allFaceInfo() const noexcept
Get allFaceInfo.
Definition MeshWave.H:133
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition UListI.H:274
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
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 sumPatchSize(const labelHashSet &patchIDs) const
Sum of patch sizes (out of supplied subset of patches).
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.
virtual ~patchDataWave()
Destructor.
virtual void correct()
Correct for mesh geom/topo changes.
patchDataWave(const polyMesh &mesh, const labelHashSet &patchIDs, const UPtrList< Field< Type > > &initialPatchValuePtrs, const bool correctWalls=true, TrackingData &td=dummyTrackData_)
Construct from mesh, information on patches to initialize and flag.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const std::string patch
OpenFOAM patch number as a std::string.
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.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299