wallDist.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallDist.H"
27 #include "wallPolyPatch.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(wallDist, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::wallDist::constructn() const
40 {
41  n_ = new volVectorField
42  (
43  IOobject
44  (
45  "n" & patchTypeName_,
46  mesh().time().timeName(),
47  mesh()
48  ),
49  mesh(),
51  patchDistMethod::patchTypes<vector>(mesh(), patchIDs_)
52  );
53 
54  const fvPatchList& patches = mesh().boundary();
55 
56  volVectorField::Boundary& nbf = n_->boundaryFieldRef();
57 
58  forAllConstIter(labelHashSet, patchIDs_, iter)
59  {
60  label patchi = iter.key();
61  nbf[patchi] == patches[patchi].nf();
62  }
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
68 Foam::wallDist::wallDist(const fvMesh& mesh, const word& patchTypeName)
69 :
71  patchIDs_(mesh.boundaryMesh().findPatchIDs<wallPolyPatch>()),
72  patchTypeName_(patchTypeName),
73  pdm_
74  (
76  (
77  static_cast<const fvSchemes&>(mesh)
78  .subDict(patchTypeName_ & "Dist"),
79  mesh,
80  patchIDs_
81  )
82  ),
83  y_
84  (
85  IOobject
86  (
87  "y" & patchTypeName_,
88  mesh.time().timeName(),
89  mesh
90  ),
91  mesh,
92  dimensionedScalar("y" & patchTypeName_, dimLength, small),
93  patchDistMethod::patchTypes<scalar>(mesh, patchIDs_)
94  ),
95  nRequired_
96  (
97  static_cast<const fvSchemes&>(mesh).subDict(patchTypeName_ & "Dist")
98  .lookupOrDefault<Switch>("nRequired", false)
99  )
100 {
101  if (nRequired_)
102  {
103  constructn();
104  }
105 
106  movePoints();
107 }
108 
109 
111 (
112  const fvMesh& mesh,
113  const labelHashSet& patchIDs,
114  const word& patchTypeName
115 )
116 :
118  patchIDs_(patchIDs),
119  patchTypeName_(patchTypeName),
120  pdm_
121  (
123  (
124  static_cast<const fvSchemes&>(mesh)
125  .subDict(patchTypeName_ & "Dist"),
126  mesh,
127  patchIDs_
128  )
129  ),
130  y_
131  (
132  IOobject
133  (
134  "y" & patchTypeName_,
135  mesh.time().timeName(),
136  mesh
137  ),
138  mesh,
139  dimensionedScalar("y" & patchTypeName_, dimLength, small),
140  patchDistMethod::patchTypes<scalar>(mesh, patchIDs_)
141  ),
142  nRequired_
143  (
144  static_cast<const fvSchemes&>(mesh).subDict(patchTypeName_ & "Dist")
145  .lookupOrDefault<Switch>("nRequired", false)
146  )
147 {
148  if (nRequired_)
149  {
150  constructn();
151  }
152 
153  movePoints();
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158 
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 {
167  if (!n_.valid())
168  {
170  << "n requested but 'nRequired' not specified in the "
171  << (patchTypeName_ & "Dist") << " dictionary" << nl
172  << " Recalculating y and n fields." << endl;
173 
174  nRequired_ = true;
175  constructn();
176  pdm_->correct(y_, n_());
177  }
178 
179  return n_();
180 }
181 
182 
184 {
185  if (pdm_->movePoints())
186  {
187  if (nRequired_)
188  {
189  return pdm_->correct(y_, n_());
190  }
191  else
192  {
193  return pdm_->correct(y_);
194  }
195  }
196  else
197  {
198  return false;
199  }
200 }
201 
202 
204 {
205  pdm_->topoChange(map);
206  movePoints();
207 }
208 
209 
211 {
212  pdm_->mapMesh(map);
213  movePoints();
214 }
215 
216 
218 {
219  // The y and n fields are registered and distributed automatically
220 }
221 
222 
223 // ************************************************************************* //
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static autoPtr< patchDistMethod > New(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:165
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
wallDist(const fvMesh &mesh, const word &patchTypeName="wall")
Construct from mesh and optional patch type name.
Definition: wallDist.C:68
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: wallDist.C:210
wordList patchTypes(nPatches)
virtual ~wallDist()
Destructor.
Definition: wallDist.C:159
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
Type of the boundary field.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
const dimensionSet dimLength
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: wallDist.C:217
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:87
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
virtual bool movePoints()
Update the y-field when the mesh moves.
Definition: wallDist.C:183
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
virtual void topoChange(const polyTopoChangeMap &)
Update the y-field when the mesh changes.
Definition: wallDist.C:203
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields...
Definition: wallDist.H:68
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Specialisation of patchDist for wall distance calculation.
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800