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 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>
35
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39// Set initial set of changed faces (= all wall faces)
40template<class TransferType, class TrackingData>
42(
43 const labelHashSet& patchIDs,
44 labelList& changedFaces,
45 List<TransferType>& faceDist
46) const
47{
48 const polyMesh& mesh = cellDistFuncs::mesh();
49
50 label nChangedFaces = 0;
51
52 forAll(mesh.boundaryMesh(), patchi)
53 {
54 if (patchIDs.found(patchi))
55 {
56 const polyPatch& patch = mesh.boundaryMesh()[patchi];
57
58 const Field<Type>& patchField = initialPatchValuePtrs_[patchi];
59
60 forAll(patch.faceCentres(), patchFacei)
61 {
62 label meshFacei = patch.start() + patchFacei;
63
64 changedFaces[nChangedFaces] = meshFacei;
65
66 faceDist[nChangedFaces] =
67 TransferType
68 (
69 patch.faceCentres()[patchFacei],
70 patchField[patchFacei],
71 0.0
72 );
73
74 nChangedFaces++;
75 }
76 }
77 }
78}
79
80
81// Copy from MeshWave data into *this (distance) and field_ (transported data)
82template<class TransferType, class TrackingData>
84(
85 const MeshWave<TransferType, TrackingData>& waveInfo
86)
87{
88 const polyMesh& mesh = cellDistFuncs::mesh();
89
90 const List<TransferType>& cellInfo = waveInfo.allCellInfo();
91 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
92
93 label nIllegal = 0;
94
95 // Copy cell values
96 distance_.setSize(cellInfo.size());
97
98 forAll(cellInfo, celli)
99 {
100 const TransferType & wpn = cellInfo[celli];
101
102 scalar dist = wpn.distSqr();
103
104 if (cellInfo[celli].valid(waveInfo.data()))
105 {
106 distance_[celli] = Foam::sqrt(dist);
107
108 cellData_[celli] = cellInfo[celli].data();
109 }
110 else
111 {
112 // Illegal/unset value. What to do with data?
113 // Note: mag for now. Should maybe be member of TransferType?
114
115 distance_[celli] = mag(dist);
116
117 //cellData_[celli] = point::max;
118 cellData_[celli] = cellInfo[celli].data();
119
120 nIllegal++;
121 }
122 }
123
124 // Copy boundary values
125 forAll(patchDistance_, patchi)
126 {
127 const polyPatch& patch = mesh.boundaryMesh()[patchi];
128
129 // Allocate storage for patchDistance
130 scalarField* patchFieldPtr = new scalarField(patch.size());
131
132 patchDistance_.set(patchi, patchFieldPtr);
133
134 scalarField& patchField = *patchFieldPtr;
135
136 // Allocate storage for patchData
137 Field<Type>* patchDataFieldPtr = new Field<Type>(patch.size());
138
139 patchData_.set(patchi, patchDataFieldPtr);
140
141 Field<Type>& patchDataField = *patchDataFieldPtr;
142
143 // Copy distance and data
144 forAll(patchField, patchFacei)
145 {
146 label meshFacei = patch.start() + patchFacei;
147
148 scalar dist = faceInfo[meshFacei].distSqr();
149
150 if (faceInfo[meshFacei].valid(waveInfo.data()))
151 {
152 // Adding SMALL to avoid problems with /0 in the turbulence
153 // models
154 patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
155
156 patchDataField[patchFacei] = faceInfo[meshFacei].data();
157 }
158 else
159 {
160 // Illegal/unset value. What to do with data?
161
162 patchField[patchFacei] = mag(dist);
163
164 //patchDataField[patchFacei] = point::max;
165 patchDataField[patchFacei] = faceInfo[meshFacei].data();
166
167 nIllegal++;
168 }
169 }
170 }
171
172 return nIllegal;
173}
174
175
176// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
177
178// Construct from components
179template<class TransferType, class TrackingData>
181(
182 const polyMesh& mesh,
183 const labelHashSet& patchIDs,
184 const UPtrList<Field<Type>>& initialPatchValuePtrs,
185 const bool correctWalls,
186 TrackingData& td
187)
188:
190 patchIDs_(patchIDs),
191 initialPatchValuePtrs_(initialPatchValuePtrs),
192 correctWalls_(correctWalls),
193 td_(td),
194 nUnset_(0),
195 distance_(mesh.nCells()),
196 patchDistance_(mesh.boundaryMesh().size()),
197 cellData_(mesh.nCells()),
198 patchData_(mesh.boundaryMesh().size())
199{
201}
202
203
204// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
205
206template<class TransferType, class TrackingData>
208{}
209
210
211// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212
213// Correct for mesh geom/topo changes
214template<class TransferType, class TrackingData>
216{
217 //
218 // Set initial changed faces: set TransferType for wall faces
219 // to wall centre.
220 //
221
222 // Count walls
223 label nWalls = sumPatchSize(patchIDs_);
224
225 List<TransferType> faceDist(nWalls);
226 labelList changedFaces(nWalls);
227
228 setChangedFaces(patchIDs_, changedFaces, faceDist);
229
230 //
231 // Do calculate wall distance by 'growing' from faces.
232 //
233
235 (
236 mesh(),
237 changedFaces,
238 faceDist,
239 mesh().globalData().nTotalCells()+1, // max iterations
240 td_
241 );
242
243
244 //
245 // Copy distance into return field
246 //
247
248 nUnset_ = getValues(waveInfo);
249
250 //
251 // Correct wall cells for true distance
252 //
253
254 if (correctWalls_)
255 {
256 Map<label> nearestFace(2 * nWalls);
257
258 // Get distance and indices of nearest face
259 correctBoundaryFaceCells
260 (
261 patchIDs_,
262 distance_,
263 nearestFace
264 );
265
266 correctBoundaryPointCells
267 (
268 patchIDs_,
269 distance_,
270 nearestFace
271 );
272
273 // Transfer data from nearest face to cell
274 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
275
276 const labelList wallCells(nearestFace.toc());
277
278 forAll(wallCells, wallCelli)
279 {
280 label celli = wallCells[wallCelli];
281
282 label facei = nearestFace[celli];
283
284 cellData_[celli] = faceInfo[facei].data();
285 }
286 }
287}
288
289
290// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
FaceCellWave plus data.
Definition: MeshWave.H:62
const List< Type > & allFaceInfo() const
Get allFaceInfo.
Definition: MeshWave.H:125
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:237
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:64
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:69
virtual ~patchDataWave()
Destructor.
virtual void correct()
Correct for mesh geom/topo changes.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
dynamicFvMesh & mesh
const std::string patch
OpenFOAM patch number as a std::string.
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:333