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-2024 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 {
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().name(),
47  mesh()
48  ),
49  mesh(),
51  patchDistMethod::patchTypes<vector>(mesh(), patchIndices_)
52  );
53 
54  const fvPatchList& patches = mesh().boundary();
55 
56  volVectorField::Boundary& nbf = n_->boundaryFieldRef();
57 
58  forAllConstIter(labelHashSet, patchIndices_, 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  patchIndices_(mesh.boundaryMesh().findIndices<wallPolyPatch>()),
72  patchTypeName_(patchTypeName),
73  pdm_
74  (
76  (
77  mesh.schemes().subDict(patchTypeName_ & "Dist"),
78  mesh,
79  patchIndices_
80  )
81  ),
82  y_
83  (
84  IOobject
85  (
86  "y" & patchTypeName_,
87  mesh.time().name(),
88  mesh
89  ),
90  mesh,
91  dimensionedScalar("y" & patchTypeName_, dimLength, small),
92  patchDistMethod::patchTypes<scalar>(mesh, patchIndices_)
93  ),
94  nRequired_
95  (
96  mesh
97  .schemes()
98  .subDict(patchTypeName_ & "Dist")
99  .lookupOrDefault<Switch>("nRequired", false)
100  )
101 {
102  if (nRequired_)
103  {
104  constructn();
105  pdm_->correct(y_, n_());
106  }
107  else
108  {
109  pdm_->correct(y_);
110  }
111 }
112 
113 
115 (
116  const fvMesh& mesh,
117  const labelHashSet& patchIDs,
118  const word& patchTypeName
119 )
120 :
122  patchIndices_(patchIDs),
123  patchTypeName_(patchTypeName),
124  pdm_
125  (
127  (
128  mesh.schemes().subDict(patchTypeName_ & "Dist"),
129  mesh,
130  patchIndices_
131  )
132  ),
133  y_
134  (
135  IOobject
136  (
137  "y" & patchTypeName_,
138  mesh.time().name(),
139  mesh
140  ),
141  mesh,
142  dimensionedScalar("y" & patchTypeName_, dimLength, small),
143  patchDistMethod::patchTypes<scalar>(mesh, patchIndices_)
144  ),
145  nRequired_
146  (
147  mesh
148  .schemes()
149  .subDict(patchTypeName_ & "Dist")
150  .lookupOrDefault<Switch>("nRequired", false)
151  )
152 {
153  if (nRequired_)
154  {
155  constructn();
156  pdm_->correct(y_, n_());
157  }
158  else
159  {
160  pdm_->correct(y_);
161  }
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
166 
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  if (!n_.valid())
176  {
178  << "n requested but 'nRequired' not specified in the "
179  << (patchTypeName_ & "Dist") << " dictionary" << nl
180  << " Recalculating y and n fields." << endl;
181 
182  nRequired_ = true;
183  constructn();
184  pdm_->correct(y_, n_());
185  }
186 
187  return n_();
188 }
189 
190 
191 // ************************************************************************* //
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
MeshObject types:
Definition: MeshObjects.H:67
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Generic GeometricField class.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:315
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:164
const word & name() const
Return name.
Definition: IOobject.H:307
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
Specialisation of patchDist for wall distance calculation.
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:71
wallDist(const fvMesh &mesh, const word &patchTypeName="wall")
Construct from mesh and optional patch type name.
Definition: wallDist.C:68
virtual ~wallDist()
Destructor.
Definition: wallDist.C:167
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition: wallDist.C:173
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:51
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
const dimensionSet dimless
const dimensionSet dimLength
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267
wordList patchTypes(nPatches)