patchWave.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "patchWave.H"
29 #include "polyMesh.H"
30 #include "wallPoint.H"
31 #include "globalMeshData.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::patchWave::setChangedFaces
36 (
37  const labelHashSet& patchIDs,
38  labelList& changedFaces,
39  List<wallPoint>& faceDist
40 ) const
41 {
42  const polyMesh& mesh = cellDistFuncs::mesh();
43 
44  label nChangedFaces = 0;
45 
46  forAll(mesh.boundaryMesh(), patchi)
47  {
48  if (patchIDs.found(patchi))
49  {
50  const polyPatch& patch = mesh.boundaryMesh()[patchi];
51 
52  forAll(patch.faceCentres(), patchFacei)
53  {
54  label meshFacei = patch.start() + patchFacei;
55 
56  changedFaces[nChangedFaces] = meshFacei;
57 
58  faceDist[nChangedFaces] =
59  wallPoint
60  (
61  patch.faceCentres()[patchFacei],
62  0.0
63  );
64 
65  nChangedFaces++;
66  }
67  }
68  }
69 
70  for (const label facei : sourceIDs_)
71  {
72  changedFaces[nChangedFaces] = facei;
73 
74  faceDist[nChangedFaces] =
75  wallPoint
76  (
77  mesh.faceCentres()[facei],
78  0.0
79  );
80 
81  nChangedFaces++;
82  }
83 }
84 
85 
86 Foam::label Foam::patchWave::getValues(const MeshWave<wallPoint>& waveInfo)
87 {
88  const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
89  const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
90 
91  label nIllegal = 0;
92 
93  // Copy cell values
94  distance_.setSize(cellInfo.size());
95 
96  forAll(cellInfo, celli)
97  {
98  scalar dist = cellInfo[celli].distSqr();
99 
100  if (cellInfo[celli].valid(waveInfo.data()))
101  {
102  distance_[celli] = Foam::sqrt(dist);
103  }
104  else
105  {
106  distance_[celli] = dist;
107 
108  nIllegal++;
109  }
110  }
111 
112  // Copy boundary values
113  forAll(patchDistance_, patchi)
114  {
115  const polyPatch& patch = mesh().boundaryMesh()[patchi];
116 
117  // Allocate storage for patchDistance
118  scalarField* patchDistPtr = new scalarField(patch.size());
119 
120  patchDistance_.set(patchi, patchDistPtr);
121 
122  scalarField& patchField = *patchDistPtr;
123 
124  forAll(patchField, patchFacei)
125  {
126  label meshFacei = patch.start() + patchFacei;
127 
128  scalar dist = faceInfo[meshFacei].distSqr();
129 
130  if (faceInfo[meshFacei].valid(waveInfo.data()))
131  {
132  // Adding SMALL to avoid problems with /0 in the turbulence
133  // models
134  patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
135  }
136  else
137  {
138  patchField[patchFacei] = dist;
139 
140  nIllegal++;
141  }
142  }
143  }
144  return nIllegal;
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
149 
151 (
152  const polyMesh& mesh,
153  const labelHashSet& patchIDs,
154  const bool correctWalls,
155  const labelList& sourceIDs
156 )
157 :
159  patchIDs_(patchIDs),
160  correctWalls_(correctWalls),
161  nUnset_(0),
162  distance_(mesh.nCells()),
163  patchDistance_(mesh.boundaryMesh().size()),
164  sourceIDs_(sourceIDs)
165 {
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
173 {}
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
179 {
180  // Set initial changed faces: set wallPoint for wall faces to wall centre
181 
182  label nPatch = sumPatchSize(patchIDs_) + sourceIDs_.size();
183 
184  List<wallPoint> faceDist(nPatch);
185  labelList changedFaces(nPatch);
186 
187  // Set to faceDist information to facecentre on walls.
188  setChangedFaces(patchIDs_, changedFaces, faceDist);
189 
190  // Do calculate wall distance by 'growing' from faces.
191  MeshWave<wallPoint> waveInfo
192  (
193  mesh(),
194  changedFaces,
195  faceDist,
196  mesh().globalData().nTotalCells()+1 // max iterations
197  );
198 
199  // Copy distance into return field
200  nUnset_ = getValues(waveInfo);
201 
202  // Correct wall cells for true distance
203  if (correctWalls_)
204  {
205  Map<label> nearestFace(2*nPatch);
206 
207  correctBoundaryFaceCells
208  (
209  patchIDs_,
210  distance_,
211  nearestFace
212  );
213 
214  correctBoundaryPointCells
215  (
216  patchIDs_,
217  distance_,
218  nearestFace
219  );
220  }
221 }
222 
223 
224 // ************************************************************************* //
const labelList patchIDs(pbm.indices(polyPatchNames, true))
virtual ~patchWave()
Destructor.
Definition: patchWave.C:165
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
dimensionedScalar sqrt(const dimensionedScalar &ds)
patchWave(const polyMesh &mesh, const labelHashSet &patchIDs, bool correctWalls=true, const labelList &sourceFaceIDs=labelList())
Construct from mesh and patches to initialize to 0 and flag.
Definition: patchWave.C:144
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:59
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
dynamicFvMesh & mesh
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
FaceCellWave plus data.
Definition: MeshWave.H:55
virtual void correct()
Correct for mesh geom/topo changes.
Definition: patchWave.C:171
const vectorField & faceCentres() const
const std::string patch
OpenFOAM patch number as a std::string.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
const polyMesh & mesh() const
Access mesh.
Definition: cellDistFuncs.H:98