wallDist.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) 2015-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "wallDist.H"
30 #include "wallPolyPatch.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(wallDist, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::wallDist::constructn() const
43 {
44  n_.reset
45  (
46  new volVectorField
47  (
48  IOobject
49  (
50  "n" & patchTypeName_,
51  mesh().time().timeName(),
52  mesh()
53  ),
54  mesh(),
56  patchDistMethod::patchTypes<vector>(mesh(), patchIDs_)
57  )
58  );
59 
60  const fvPatchList& patches = mesh().boundary();
61 
62  volVectorField::Boundary& nbf = n_.ref().boundaryFieldRef();
63 
64  for (const label patchi : patchIDs_)
65  {
66  nbf[patchi] == patches[patchi].nf();
67  }
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 Foam::wallDist::wallDist
74 (
75  const fvMesh& mesh,
76  const labelHashSet& patchIDs,
77  const word& patchTypeName
78 )
79 :
80  wallDist(mesh, word::null, patchIDs, patchTypeName)
81 {}
82 
83 
84 Foam::wallDist::wallDist
85 (
86  const fvMesh& mesh,
87  const word& defaultPatchDistMethod,
88  const labelHashSet& patchIDs,
89  const word& patchTypeName
90 )
91 :
93  patchIDs_(patchIDs),
94  patchTypeName_(patchTypeName),
95  dict_
96  (
97  static_cast<const fvSchemes&>(mesh).subOrEmptyDict
98  (
99  patchTypeName_ & "Dist"
100  )
101  ),
102  pdm_
103  (
105  (
106  dict_,
107  mesh,
108  patchIDs_,
109  defaultPatchDistMethod
110  )
111  ),
112  y_
113  (
114  IOobject
115  (
116  "y" & patchTypeName_,
117  mesh.time().timeName(),
118  mesh
119  ),
120  mesh,
121  dimensionedScalar("y" & patchTypeName_, dimLength, SMALL),
122  patchDistMethod::patchTypes<scalar>(mesh, patchIDs_)
123  ),
124  n_(volVectorField::null()),
125  updateInterval_(dict_.getOrDefault<label>("updateInterval", 1)),
126  nRequired_(dict_.getOrDefault("nRequired", false)),
127  requireUpdate_(true)
128 {
129  if (nRequired_)
130  {
131  constructn();
132  }
133 
134  movePoints();
135 }
136 
137 
138 Foam::wallDist::wallDist(const fvMesh& mesh, const word& patchTypeName)
139 :
140  wallDist
141  (
142  mesh,
143  mesh.boundaryMesh().findPatchIDs<wallPolyPatch>(),
144  patchTypeName
145  )
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152 {}
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
158 {
159  if (isNull(n_()))
160  {
162  << "n requested but 'nRequired' not specified in the "
163  << (patchTypeName_ & "Dist") << " dictionary" << nl
164  << " Recalculating y and n fields." << endl;
165 
166  nRequired_ = true;
167  constructn();
168  pdm_->correct(y_, n_.ref());
169  }
170 
171  return n_();
172 }
173 
174 
176 {
177  if
178  (
179  (updateInterval_ != 0)
180  && ((mesh_.time().timeIndex() % updateInterval_) == 0)
181  )
182  {
183  requireUpdate_ = true;
184  }
185 
186  if (requireUpdate_ && pdm_->movePoints())
187  {
188  DebugInfo<< "Updating wall distance" << endl;
189 
190  requireUpdate_ = false;
191 
192  if (nRequired_)
193  {
194  return pdm_->correct(y_, n_.ref());
195  }
196  else
197  {
198  return pdm_->correct(y_);
199  }
200  }
201 
202  return false;
203 }
204 
205 
206 void Foam::wallDist::updateMesh(const mapPolyMesh& mpm)
207 {
208  pdm_->updateMesh(mpm);
209 
210  // Force update if performing topology change
211  // Note: needed?
212  // - field would have been mapped, so if using updateInterval option (!= 1)
213  // live with error associated of not updating and use mapped values?
214  requireUpdate_ = true;
215  movePoints();
216 }
217 
218 
219 // ************************************************************************* //
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:150
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
wordList patchTypes(nPatches)
virtual ~wallDist()
Destructor.
Definition: wallDist.C:144
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
Type of boundary fields.
const dimensionSet dimless
Dimensionless.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:84
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
Definition: wallDist.C:199
virtual bool movePoints()
Update the y-field when the mesh moves.
Definition: wallDist.C:168
#define DebugInfo
Report an information message using Foam::Info.
bool isNull(const T *ptr)
True if ptr is a pointer (of type T) to the nullObject.
Definition: nullObject.H:227
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:447
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
IOobject(const IOobject &)=default
Copy construct.
const polyBoundaryMesh & patches
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields...
Definition: wallDist.H:71
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
Specialisation of patchDist for wall distance calculation.
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:57
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:705
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157