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 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::wallDistAddressing
28 
29 Description
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 
77 Note
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 
83 See also
84  Foam::wallDist
85 
86 SourceFiles
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 
101 namespace Foam
102 {
103 
104 // Forward Declarations
105 class mapDistribute;
106 class wallPointAddressing;
108 /*---------------------------------------------------------------------------*\
109  Class wallDistAddressing Declaration
110 \*---------------------------------------------------------------------------*/
111 
112 class wallDistAddressing
113 :
114  public MeshObject<fvMesh, UpdateableMeshObject, wallDistAddressing>,
115  public cellDistFuncs
116 {
117  //- Use fvMesh& instead of cellDistFuncs polyMesh&
119 
120 protected:
121 
122  // Private Data
123 
124  //- Set of patch IDs
125  const labelList patchIDs_;
126 
127  //- Name for the patch set, e.g. "wall"
128  const word patchTypeName_;
130  //- Update wall distance every updateInterval_ steps
131  const label updateInterval_;
132 
133  //- Do accurate distance calculation for near-wall cells.
134  const bool correctWalls_;
135 
136  //- Flag to indicate whether the wall distance requires updating
137  bool requireUpdate_;
138 
139  //- Distance-to-wall field
140  mutable volScalarField y_;
141 
142  //- Number of wall faces
145  //- Map to pull wall face info to cell or boundary face
147 
148  //- Indices of cells followed by boundary faces
150 
151  //- Corresponding slot in mapPtr distribution result
153 
154  //- Start of patches. Start of untransformedPatchStarts_[0] is end
155  // of cell information
157 
158  // Corresponding data for transformed items
162 
163 
164  // Protected Member Functions
165 
166  //- Extract FaceCellWave data
167  label getValues
168  (
170  const List<wallPointAddressing>& allCellInfo,
171  const List<wallPointAddressing>& allFaceInfo,
173  ) const;
174 
175  //- Store nearest-data to cell or boundary face
176  void addItem
177  (
178  const label item,
179  const labelPair& data,
180  label& untransformi,
181  label& transformi,
182  labelPairList& transformedWallInfo
183  );
184 
185  //- Extract nearest-patch distance data
186  void correct(volScalarField& y);
187 
188 
189  //- No copy construct
190  wallDistAddressing(const wallDistAddressing&) = delete;
191 
192  //- No copy assignment
193  void operator=(const wallDistAddressing&) = delete;
194 
195 
196 public:
197 
198  // Declare name of the class and its debug switch
199  ClassName("wallDistAddressing");
200 
201 
202  // Constructors
203 
204  //- Construct from mesh and near-wall behaviour
205  explicit wallDistAddressing
206  (
207  const fvMesh& mesh,
208  const bool correctWalls = true,
209  const label updateInterval = 1
210  );
211 
212  //- Construct from patch type name (= registration name),
213  //- mesh and patch IDs.
215  (
216  const word& patchTypeName,
217  const fvMesh& mesh,
218  const labelList& patchIDs,
219  const bool correctWalls = true,
220  const label updateInterval = 1
221  );
222 
223 
224  //- Destructor
225  virtual ~wallDistAddressing();
226 
227 
228  // Member Functions
229 
230  //- Return the patchIDs
231  const labelUList& patchIDs() const noexcept
232  {
233  return patchIDs_;
234  }
235 
236  //- Return reference to cached distance-to-wall field. Unvisited
237  // elements set to GREAT
238  const volScalarField& y() const noexcept
239  {
240  return y_;
241  }
242 
243  //- Update the y-field when the mesh moves
244  virtual bool movePoints();
245 
246  //- Update the y-field when the mesh changes
247  virtual void updateMesh(const mapPolyMesh&);
248 
249 
250  // Get values from nearest patches
251 
252  //- Collect patchFields from patchIDs into straight list
253  template<class Container, class Type>
254  tmp<Field<Type>> collectPatchFields(const Container& bfld) const;
255 
256  //- Take collected/distributed patch field and fill volField
257  template<class VolField>
258  void extract
259  (
260  const Field<typename VolField::value_type>& patchFld,
261  VolField& fld
262  ) const;
263 
264  //- Map nearest-patch information. Take wall patch values
265  // and transports these to the internal field and other patches.
266  // Use mapDistribute::transformPosition if transporting absolute
267  // coordinates.
268  template<class VolField, class TransformOp>
269  const VolField& map
270  (
271  VolField& fld,
272  const TransformOp& top = mapDistribute::transform()
273  ) const;
274 };
275 
276 
277 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 
279 } // End namespace Foam
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 #ifdef NoRepository
285 #endif
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 #endif
290 
291 // ************************************************************************* //
const labelList patchIDs_
Set of patch IDs.
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field. Unvisited.
labelList untransformedItems_
Indices of cells followed by boundary faces.
const bool correctWalls_
Do accurate distance calculation for near-wall cells.
wallDistAddressing(const wallDistAddressing &)=delete
No copy construct.
const labelUList & patchIDs() const noexcept
Return the patchIDs.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:206
label getValues(const FaceCellWave< wallPointAddressing > &wave, const List< wallPointAddressing > &allCellInfo, const List< wallPointAddressing > &allFaceInfo, volScalarField &y) const
Extract FaceCellWave data.
void correct(volScalarField &y)
Extract nearest-patch distance data.
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:59
void extract(const Field< typename VolField::value_type > &patchFld, VolField &fld) const
Take collected/distributed patch field and fill volField.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
labelList untransformedSlots_
Corresponding slot in mapPtr distribution result.
tmp< Field< Type > > collectPatchFields(const Container &bfld) const
Collect patchFields from patchIDs into straight list.
Default transformation behaviour.
autoPtr< mapDistribute > mapPtr_
Map to pull wall face info to cell or boundary face.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList untransformedPatchStarts_
Start of patches. Start of untransformedPatchStarts_[0] is end.
virtual ~wallDistAddressing()
Destructor.
const label updateInterval_
Update wall distance every updateInterval_ steps.
const direction noexcept
Definition: Scalar.H:258
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
ClassName("wallDistAddressing")
autoPtr< globalIndex > globalWallsPtr_
Number of wall faces.
const VolField & map(VolField &fld, const TransformOp &top=mapDistribute::transform()) const
Map nearest-patch information. Take wall patch values.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const word patchTypeName_
Name for the patch set, e.g. "wall".
volScalarField y_
Distance-to-wall field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
void operator=(const wallDistAddressing &)=delete
No copy assignment.
virtual bool movePoints()
Update the y-field when the mesh moves.
bool requireUpdate_
Flag to indicate whether the wall distance requires updating.
const polyMesh & mesh() const
Access mesh.
Definition: cellDistFuncs.H:98
void addItem(const label item, const labelPair &data, label &untransformi, label &transformi, labelPairList &transformedWallInfo)
Store nearest-data to cell or boundary face.
Namespace for OpenFOAM.
Variant of wallDist that uses meshWave and stores the addressing.