Loading...
Searching...
No Matches
wallDistAddressing.H
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 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
26Class
27 Foam::wallDistAddressing
28
29Description
30 Variant of wallDist that uses meshWave and stores the addressing.
31
32 Can be used to transport other wall properties.
33 For example getting the nearest-wall normal into a volVectorField
34 \verbatim
35 // If not yet constructed construct from all wall patches
36 const auto& wDist = wallDistAddressing::New(mesh);
37
38 // Fill boundaryField with normal
39 for (const label patchi : wDist.patchIDs())
40 {
41 auto& pnf = n.boundaryFieldRef()[patchi];
42 pnf == pnf.patch().nf();
43 }
44
45 // Map data from nearest wall (using transformation if nearest is
46 // reached through a coupled patch with transformation)
47 wDist.map(n, mapDistribute::transform());
48 \endverbatim
49
50
51 Another use is e.g. to get the distance to cyclic patches. This enforces
52 registration under the supplied name.
53
54 \verbatim
55 // If not yet constructed construct from selected patches
56 const labelList& patchIDs = ..
57
58 const auto& wDist = wallDistAddressing::New
59 (
60 "myPatches", // registration name
61 mesh,
62 patchIDs
63 );
64
65 // Fill boundaryField with normal
66 for (const label patchi : wDist.patchIDs())
67 {
68 auto& pnf = n.boundaryFieldRef()[patchi];
69 pnf == pnf.patch().nf();
70 }
71
72 // Map data from nearest patch face (using transformation if nearest is
73 // reached through a coupled patch with transformation)
74 wDist.map(n, mapDistribute::transform());
75 \endverbatim
76
77Note
78 By default (correctWalls = true) all
79 cells point-connected to a wall explicitly search for the nearest
80 location on the point-surrouding wall faces (so override the wave
81 behaviour). This gets only done for processor-local wall faces.
82
83See also
84 Foam::wallDist
85
86SourceFiles
87 wallDistAddressing.C
88
89\*---------------------------------------------------------------------------*/
90
91#ifndef Foam_wallDistAddressing_H
92#define Foam_wallDistAddressing_H
93
94#include "MeshObject.H"
95#include "FaceCellWave.H"
96#include "cellDistFuncs.H"
97#include "volFields.H"
98
99// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100
101namespace Foam
102{
103
104// Forward Declarations
105class mapDistribute;
108/*---------------------------------------------------------------------------*\
109 Class wallDistAddressing Declaration
110\*---------------------------------------------------------------------------*/
111
113:
114 public MeshObject<fvMesh, UpdateableMeshObject, wallDistAddressing>,
115 public cellDistFuncs
116{
117 // Private Typedefs
118
119 typedef MeshObject
120 <
121 fvMesh,
124 > MeshObject_type;
125
126
127 //- Use fvMesh& instead of cellDistFuncs polyMesh&
129
130
131protected:
132
133 // Private Data
134
135 //- Set of patch IDs
136 const labelList patchIDs_;
137
138 //- Name for the patch set, e.g. "wall"
139 const word patchTypeName_;
141 //- Update wall distance every updateInterval_ steps
142 const label updateInterval_;
143
144 //- Do accurate distance calculation for near-wall cells.
145 const bool correctWalls_;
146
147 //- Flag to indicate whether the wall distance requires updating
148 bool requireUpdate_;
149
150 //- Distance-to-wall field
151 mutable volScalarField y_;
152
153 //- Number of wall faces
156 //- Map to pull wall face info to cell or boundary face
158
159 //- Indices of cells followed by boundary faces
161
162 //- Corresponding slot in mapPtr distribution result
164
165 //- Start of patches. Start of untransformedPatchStarts_[0] is end
166 // of cell information
168
169 // Corresponding data for transformed items
173
174
175 // Protected Member Functions
176
177 //- Extract FaceCellWave data
178 label getValues
179 (
184 ) const;
185
186 //- Store nearest-data to cell or boundary face
188 (
189 const label item,
190 const labelPair& data,
191 label& untransformi,
192 label& transformi,
193 labelPairList& transformedWallInfo
194 );
195
196 //- Extract nearest-patch distance data
197 void correct(volScalarField& y);
198
199
200 //- No copy construct
201 wallDistAddressing(const wallDistAddressing&) = delete;
202
203 //- No copy assignment
204 void operator=(const wallDistAddressing&) = delete;
205
206
207public:
208
209 // Declare name of the class and its debug switch
210 ClassName("wallDistAddressing");
211
212
213 // Constructors
214
215 //- Construct from mesh and near-wall behaviour
216 explicit wallDistAddressing
217 (
218 const fvMesh& mesh,
219 const bool correctWalls = true,
220 const label updateInterval = 1
221 );
222
223 //- Construct from patch type name (= registration name),
224 //- mesh and patch IDs.
226 (
227 const word& patchTypeName,
228 const fvMesh& mesh,
230 const bool correctWalls = true,
231 const label updateInterval = 1
232 );
233
235 //- Destructor
236 virtual ~wallDistAddressing();
237
238
239 // Member Functions
241 //- Return the patchIDs
242 const labelUList& patchIDs() const noexcept
243 {
244 return patchIDs_;
245 }
246
247 //- Return reference to cached distance-to-wall field. Unvisited
248 // elements set to GREAT
249 const volScalarField& y() const noexcept
250 {
251 return y_;
252 }
253
254 //- Update the y-field when the mesh moves
255 virtual bool movePoints();
256
257 //- Update the y-field when the mesh changes
258 virtual void updateMesh(const mapPolyMesh&);
259
260
261 // Get values from nearest patches
262
263 //- Collect patchFields from patchIDs into straight list
264 template<class Container, class Type>
265 tmp<Field<Type>> collectPatchFields(const Container& bfld) const;
266
267 //- Take collected/distributed patch field and fill volField
268 template<class VolField>
269 void extract
270 (
271 const Field<typename VolField::value_type>& patchFld,
272 VolField& fld
273 ) const;
274
275 //- Map nearest-patch information. Take wall patch values
276 // and transports these to the internal field and other patches.
277 // Use mapDistribute::transformPosition if transporting absolute
278 // coordinates.
279 template<class VolField, class TransformOp>
280 const VolField& map
281 (
282 VolField& fld,
283 const TransformOp& top = mapDistribute::transform()
284 ) const;
285};
286
287
288// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289
290} // End namespace Foam
291
292// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293
294#ifdef NoRepository
296#endif
297
298// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299
300#endif
301
302// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const polyMesh & mesh() const
Access mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Default transformation behaviour.
Class containing processor-to-processor mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition tmp.H:75
virtual bool movePoints()
Update the y-field when the mesh moves.
const word patchTypeName_
Name for the patch set, e.g. "wall".
ClassName("wallDistAddressing")
labelList untransformedPatchStarts_
Start of patches. Start of untransformedPatchStarts_[0] is end.
void extract(const Field< typename VolField::value_type > &patchFld, VolField &fld) const
Take collected/distributed patch field and fill volField.
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.
void operator=(const wallDistAddressing &)=delete
No copy assignment.
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.
tmp< Field< Type > > collectPatchFields(const Container &bfld) const
Collect patchFields from patchIDs into straight list.
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.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
thermo correct()
List< wallPoints > allFaceInfo(mesh_.nFaces())
List< wallPoints > allCellInfo(mesh_.nCells())
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< labelPair > labelPairList
List of labelPair.
Definition labelPair.H:33
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75