nearWallDist.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2020,2024 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 "nearWallDist.H"
30 #include "cellDistFuncs.H"
31 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::nearWallDist::calculate()
36 {
37  const cellDistFuncs wallUtils(mesh_);
38 
39  const volVectorField& cellCentres = mesh_.C();
40 
41 
43  {
44  // Collect indices of wall patches
45 
46  DynamicList<label> wallPatchIDs(mesh_.boundary().size());
47  label nWalls = 0;
48 
49  forAll(mesh_.boundary(), patchi)
50  {
51  if (isA<wallFvPatch>(mesh_.boundary()[patchi]))
52  {
53  wallPatchIDs.append(patchi);
54  nWalls += mesh_.boundary()[patchi].size();
55  }
56  else
57  {
58  // Set distance to 0
59  operator[](patchi) = 0.0;
60  }
61  }
62 
63 
64  // Collect all mesh faces of wall patches
65 
66  DynamicList<label> faceLabels(nWalls);
67 
68  for (const label patchi : wallPatchIDs)
69  {
70  const fvPatch& patch = mesh_.boundary()[patchi];
71  const auto& pp = patch.patch();
72 
73  forAll(patch, i)
74  {
75  faceLabels.append(pp.start()+i);
76  }
77  }
78 
79  const uindirectPrimitivePatch wallPatch
80  (
81  UIndirectList<face>(mesh_.faces(), faceLabels),
82  mesh_.points()
83  );
84 
85  DynamicList<label> neighbours;
86 
87  nWalls = 0;
88  for (const label patchi : wallPatchIDs)
89  {
90  const fvPatch& patch = mesh_.boundary()[patchi];
91  const labelUList& faceCells = patch.patch().faceCells();
92 
93  fvPatchScalarField& ypatch = operator[](patchi);
94 
95  forAll(patch, patchFacei)
96  {
97  // Get point connected neighbours (in wallPatch indices!)
98  wallUtils.getPointNeighbours(wallPatch, nWalls, neighbours);
99 
100  label minFacei = -1;
101  ypatch[patchFacei] = wallUtils.smallestDist
102  (
103  cellCentres[faceCells[patchFacei]],
104  wallPatch,
105  neighbours,
106  minFacei
107  );
108 
109  nWalls++;
110  }
111  }
112  }
113  else
114  {
115  // Get patch ids of walls
116  const labelHashSet wallPatchIDs(wallUtils.getPatchIDs<wallPolyPatch>());
117 
118  // Size neighbours array for maximum possible
119  DynamicList<label> neighbours(wallUtils.maxPatchSize(wallPatchIDs));
120 
121 
122  // Correct all cells with face on wall
123 
124  forAll(mesh_.boundary(), patchi)
125  {
126  fvPatchScalarField& ypatch = operator[](patchi);
127 
128  const fvPatch& patch = mesh_.boundary()[patchi];
129 
130  if (isA<wallFvPatch>(patch))
131  {
132  const polyPatch& pPatch = patch.patch();
133 
134  const labelUList& faceCells = patch.faceCells();
135 
136  // Check cells with face on wall
137  forAll(patch, patchFacei)
138  {
139  wallUtils.getPointNeighbours
140  (
141  pPatch,
142  patchFacei,
143  neighbours
144  );
145 
146  label minFacei = -1;
147 
148  ypatch[patchFacei] = wallUtils.smallestDist
149  (
150  cellCentres[faceCells[patchFacei]],
151  pPatch,
152  neighbours,
153  minFacei
154  );
155  }
156  }
157  else
158  {
159  ypatch = 0.0;
160  }
161  }
162  }
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167 
168 Foam::nearWallDist::nearWallDist(const Foam::fvMesh& mesh)
169 :
170  volScalarField::Boundary
171  (
172  mesh.boundary(),
173  mesh.V(), // Dummy internal field,
174  fvPatchFieldBase::calculatedType()
175  ),
176  mesh_(mesh)
177 {
178  calculate();
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 {}
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
191 {
192  if (mesh_.topoChanging())
193  {
194  const DimensionedField<scalar, volMesh>& V = mesh_.V();
195  const fvBoundaryMesh& bnd = mesh_.boundary();
196 
197  this->setSize(bnd.size());
198  forAll(*this, patchi)
199  {
200  this->set
201  (
202  patchi,
204  (
206  bnd[patchi],
207  V
208  )
209  );
210  }
211  }
212 
213  calculate();
214 }
215 
216 
217 // ************************************************************************* //
faceListList boundary
const Type & operator[](const labelPair &index) const
Const access to a single field element via (fieldi, elemi)
Definition: FieldField.C:351
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
labelList faceLabels(nFaceLabels)
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static bool useCombinedWallPatch
Use combined-wall-patches wall distance v.s. v2406 per-patch distance. Default is true...
Definition: cellDistFuncs.H:92
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
points setSize(newPointi)
const word calculatedType
A calculated patch field type.
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
virtual ~nearWallDist()
Destructor.
Definition: nearWallDist.C:177
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
const std::string patch
OpenFOAM patch number as a std::string.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
const volVectorField & C() const
Return cell centres as volVectorField.
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:183
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.