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 | Copyright (C) 2011-2016 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 "patchWave.H"
27 #include "polyMesh.H"
28 #include "wallPoint.H"
29 #include "globalMeshData.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::patchWave::setChangedFaces
34 (
35  const labelHashSet& patchIDs,
36  labelList& changedFaces,
37  List<wallPoint>& faceDist
38 ) const
39 {
40  const polyMesh& mesh = cellDistFuncs::mesh();
41 
42  label nChangedFaces = 0;
43 
44  forAll(mesh.boundaryMesh(), patchi)
45  {
46  if (patchIDs.found(patchi))
47  {
48  const polyPatch& patch = mesh.boundaryMesh()[patchi];
49 
50  forAll(patch.faceCentres(), patchFacei)
51  {
52  label meshFacei = patch.start() + patchFacei;
53 
54  changedFaces[nChangedFaces] = meshFacei;
55 
56  faceDist[nChangedFaces] =
57  wallPoint
58  (
59  patch.faceCentres()[patchFacei],
60  0.0
61  );
62 
63  nChangedFaces++;
64  }
65  }
66  }
67 }
68 
69 
70 Foam::label Foam::patchWave::getValues(const MeshWave<wallPoint>& waveInfo)
71 {
72  const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
73  const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
74 
75  label nIllegal = 0;
76 
77  // Copy cell values
78  distance_.setSize(cellInfo.size());
79 
80  forAll(cellInfo, celli)
81  {
82  scalar dist = cellInfo[celli].distSqr();
83 
84  if (cellInfo[celli].valid(waveInfo.data()))
85  {
86  distance_[celli] = Foam::sqrt(dist);
87  }
88  else
89  {
90  distance_[celli] = dist;
91 
92  nIllegal++;
93  }
94  }
95 
96  // Copy boundary values
97  forAll(patchDistance_, patchi)
98  {
99  const polyPatch& patch = mesh().boundaryMesh()[patchi];
100 
101  // Allocate storage for patchDistance
102  scalarField* patchDistPtr = new scalarField(patch.size());
103 
104  patchDistance_.set(patchi, patchDistPtr);
105 
106  scalarField& patchField = *patchDistPtr;
107 
108  forAll(patchField, patchFacei)
109  {
110  label meshFacei = patch.start() + patchFacei;
111 
112  scalar dist = faceInfo[meshFacei].distSqr();
113 
114  if (faceInfo[meshFacei].valid(waveInfo.data()))
115  {
116  // Adding SMALL to avoid problems with /0 in the turbulence
117  // models
118  patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
119  }
120  else
121  {
122  patchField[patchFacei] = dist;
123 
124  nIllegal++;
125  }
126  }
127  }
128  return nIllegal;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133 
135 (
136  const polyMesh& mesh,
137  const labelHashSet& patchIDs,
138  const bool correctWalls
139 )
140 :
141  cellDistFuncs(mesh),
142  patchIDs_(patchIDs),
143  correctWalls_(correctWalls),
144  nUnset_(0),
145  distance_(mesh.nCells()),
146  patchDistance_(mesh.boundaryMesh().size())
147 {
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  // Set initial changed faces: set wallPoint for wall faces to wall centre
163 
164  label nPatch = sumPatchSize(patchIDs_);
165 
166  List<wallPoint> faceDist(nPatch);
167  labelList changedFaces(nPatch);
168 
169  // Set to faceDist information to facecentre on walls.
170  setChangedFaces(patchIDs_, changedFaces, faceDist);
171 
172  // Do calculate wall distance by 'growing' from faces.
173  MeshWave<wallPoint> waveInfo
174  (
175  mesh(),
176  changedFaces,
177  faceDist,
178  mesh().globalData().nTotalCells()+1 // max iterations
179  );
180 
181  // Copy distance into return field
182  nUnset_ = getValues(waveInfo);
183 
184  // Correct wall cells for true distance
185  if (correctWalls_)
186  {
187  Map<label> nearestFace(2*nPatch);
188 
190  (
191  patchIDs_,
192  distance_,
193  nearestFace
194  );
195 
197  (
198  patchIDs_,
199  distance_,
200  nearestFace
201  );
202  }
203 }
204 
205 
206 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
label sumPatchSize(const labelHashSet &patchIDs) const
Sum of patch sizes (out of supplied subset of patches).
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
virtual ~patchWave()
Destructor.
Definition: patchWave.C:154
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
label nCells() const
dimensionedScalar sqrt(const dimensionedScalar &ds)
void correctBoundaryFaceCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to boundary (via face). Sets values in.
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:61
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
FaceCellWave plus data.
Definition: MeshWave.H:56
void correctBoundaryPointCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to wall (via point). Sets values in.
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual void correct()
Correct for mesh geom/topo changes.
Definition: patchWave.C:160
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const polyMesh & mesh() const
Access mesh.
Definition: cellDistFuncs.H:95
patchWave(const polyMesh &mesh, const labelHashSet &patchIDs, bool correctWalls=true)
Construct from mesh and patches to initialize to 0 and flag.
Definition: patchWave.C:135