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 | Copyright (C) 2015 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 
40 {
42  (
43  new volVectorField
44  (
45  IOobject
46  (
47  "n" & patchTypeName_,
48  mesh().time().timeName(),
49  mesh()
50  ),
51  mesh(),
53  patchDistMethod::patchTypes<vector>(mesh(), patchIDs_)
54  )
55  );
56 
57  const fvPatchList& patches = mesh().boundary();
58 
60  {
61  label patchi = iter.key();
62  n_().boundaryField()[patchi] == patches[patchi].nf();
63  }
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
69 Foam::wallDist::wallDist(const fvMesh& mesh, const word& patchTypeName)
70 :
72  patchIDs_(mesh.boundaryMesh().findPatchIDs<wallPolyPatch>()),
73  patchTypeName_(patchTypeName),
74  pdm_
75  (
77  (
78  static_cast<const fvSchemes&>(mesh)
79  .subDict(patchTypeName_ & "Dist"),
80  mesh,
81  patchIDs_
82  )
83  ),
84  y_
85  (
86  IOobject
87  (
88  "y" & patchTypeName_,
89  mesh.time().timeName(),
90  mesh
91  ),
92  mesh,
93  dimensionedScalar("y" & patchTypeName_, dimLength, SMALL),
94  patchDistMethod::patchTypes<scalar>(mesh, patchIDs_)
95  ),
96  nRequired_
97  (
98  static_cast<const fvSchemes&>(mesh).subDict(patchTypeName_ & "Dist")
99  .lookupOrDefault<Switch>("nRequired", false)
100  ),
101  n_(volVectorField::null())
102 {
103  if (nRequired_)
104  {
105  constructn();
106  }
107 
108  movePoints();
109 }
110 
111 
113 (
114  const fvMesh& mesh,
115  const labelHashSet& patchIDs,
116  const word& patchTypeName
117 )
118 :
120  patchIDs_(patchIDs),
121  patchTypeName_(patchTypeName),
122  pdm_
123  (
125  (
126  static_cast<const fvSchemes&>(mesh)
127  .subDict(patchTypeName_ & "Dist"),
128  mesh,
129  patchIDs_
130  )
131  ),
132  y_
133  (
134  IOobject
135  (
136  "y" & patchTypeName_,
137  mesh.time().timeName(),
138  mesh
139  ),
140  mesh,
141  dimensionedScalar("y" & patchTypeName_, dimLength, SMALL),
142  patchDistMethod::patchTypes<scalar>(mesh, patchIDs_)
143  ),
144  nRequired_
145  (
146  static_cast<const fvSchemes&>(mesh).subDict(patchTypeName_ & "Dist")
147  .lookupOrDefault<Switch>("nRequired", false)
148  ),
150 {
151  if (nRequired_)
152  {
153  constructn();
154  }
155 
156  movePoints();
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
161 
163 {}
164 
165 
166 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
167 
169 {
170  if (isNull(n_()))
171  {
173  << "n requested but 'nRequired' not specified in the "
174  << (patchTypeName_ & "Dist") << " dictionary" << nl
175  << " Recalculating y and n fields." << endl;
176 
177  nRequired_ = true;
178  constructn();
179  pdm_->correct(y_, n_());
180  }
181 
182  return n_();
183 }
184 
185 
187 {
188  if (pdm_->movePoints())
189  {
190  if (nRequired_)
191  {
192  return pdm_->correct(y_, n_());
193  }
194  else
195  {
196  return pdm_->correct(y_);
197  }
198  }
199  else
200  {
201  return false;
202  }
203 }
204 
205 
207 {
208  pdm_->updateMesh(mpm);
209  movePoints();
210 }
211 
212 
213 // ************************************************************************* //
wallDist.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::wallDist::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update the y-field when the mesh changes.
Definition: wallDist.C:206
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:50
wallPolyPatch.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::HashSet< label, Hash< label > >
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
Foam::UpdateableMeshObject
Definition: MeshObject.H:258
Foam::patchDistMethod::New
static autoPtr< patchDistMethod > New(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
Definition: patchDistMethod.C:53
patchTypes
wordList patchTypes(nPatches)
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::wallDist
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:68
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::mesh
const fvMesh & mesh() const
Definition: MeshObject.H:144
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::wallDist::movePoints
virtual bool movePoints()
Update the y-field when the mesh moves.
Definition: wallDist.C:186
Foam::wallDist::n
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:168
Foam::patchDistMethod
Specialisation of patchDist for wall distance calculation.
Definition: patchDistMethod.H:56
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::wallPolyPatch
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::wallDist::constructn
void constructn() const
Construct the normal-to-wall field as required.
Definition: wallDist.C:39
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::wallDist::n_
tmp< volVectorField > n_
Normal-to-wall field.
Definition: wallDist.H:90
patchi
label patchi
Definition: getPatchFieldScalar.H:1
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::isNull
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:40
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:59
Foam::wallDist::wallDist
wallDist(const wallDist &)
Disallow default bitwise copy construct.
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:81
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::wallDist::patchTypeName_
const word patchTypeName_
Name for the patch set, e.g. "wall".
Definition: wallDist.H:78
Foam::wallDist::patchIDs_
const labelHashSet patchIDs_
Set of patch IDs.
Definition: wallDist.H:75
Foam::GeometricField::null
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:30
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::wallDist::nRequired_
bool nRequired_
Flag to indicate if the distance-to-wall field is required.
Definition: wallDist.H:87
Foam::wallDist::~wallDist
virtual ~wallDist()
Destructor.
Definition: wallDist.C:162