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 
71 
72 Foam::label Foam::patchWave::getValues(const MeshWave<wallPoint>& waveInfo)
73 {
74  const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
75  const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
76 
77  label nIllegal = 0;
78 
79  // Copy cell values
80  distance_.setSize(cellInfo.size());
81 
82  forAll(cellInfo, celli)
83  {
84  scalar dist = cellInfo[celli].distSqr();
85 
86  if (cellInfo[celli].valid(waveInfo.data()))
87  {
88  distance_[celli] = Foam::sqrt(dist);
89  }
90  else
91  {
92  distance_[celli] = dist;
93 
94  nIllegal++;
95  }
96  }
97 
98  // Copy boundary values
99  forAll(patchDistance_, patchi)
100  {
101  const polyPatch& patch = mesh().boundaryMesh()[patchi];
102 
103  // Allocate storage for patchDistance
104  scalarField* patchDistPtr = new scalarField(patch.size());
105 
106  patchDistance_.set(patchi, patchDistPtr);
107 
108  scalarField& patchField = *patchDistPtr;
109 
110  forAll(patchField, patchFacei)
111  {
112  label meshFacei = patch.start() + patchFacei;
113 
114  scalar dist = faceInfo[meshFacei].distSqr();
115 
116  if (faceInfo[meshFacei].valid(waveInfo.data()))
117  {
118  // Adding SMALL to avoid problems with /0 in the turbulence
119  // models
120  patchField[patchFacei] = Foam::sqrt(dist) + SMALL;
121  }
122  else
123  {
124  patchField[patchFacei] = dist;
125 
126  nIllegal++;
127  }
128  }
129  }
130  return nIllegal;
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135 
137 (
138  const polyMesh& mesh,
139  const labelHashSet& patchIDs,
140  const bool correctWalls
141 )
142 :
144  patchIDs_(patchIDs),
145  correctWalls_(correctWalls),
146  nUnset_(0),
147  distance_(mesh.nCells()),
148  patchDistance_(mesh.boundaryMesh().size())
149 {
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
155 
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
163 {
164  // Set initial changed faces: set wallPoint for wall faces to wall centre
165 
166  label nPatch = sumPatchSize(patchIDs_);
167 
168  List<wallPoint> faceDist(nPatch);
169  labelList changedFaces(nPatch);
170 
171  // Set to faceDist information to facecentre on walls.
172  setChangedFaces(patchIDs_, changedFaces, faceDist);
173 
174  // Do calculate wall distance by 'growing' from faces.
175  MeshWave<wallPoint> waveInfo
176  (
177  mesh(),
178  changedFaces,
179  faceDist,
180  mesh().globalData().nTotalCells()+1 // max iterations
181  );
182 
183  // Copy distance into return field
184  nUnset_ = getValues(waveInfo);
185 
186  // Correct wall cells for true distance
187  if (correctWalls_)
188  {
189  Map<label> nearestFace(2*nPatch);
190 
191  correctBoundaryFaceCells
192  (
193  patchIDs_,
194  distance_,
195  nearestFace
196  );
197 
198  correctBoundaryPointCells
199  (
200  patchIDs_,
201  distance_,
202  nearestFace
203  );
204  }
205 }
206 
207 
208 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
wallPoint.H
globalMeshData.H
Foam::Map< label >
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::MeshWave
FaceCellWave plus data.
Definition: MeshWave.H:59
Foam::patchWave::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: patchWave.C:162
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::cellDistFuncs::mesh
const polyMesh & mesh() const
Access mesh.
Definition: cellDistFuncs.H:92
Foam::cellDistFuncs
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:63
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
patchWave.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::patchWave::patchWave
patchWave(const polyMesh &mesh, const labelHashSet &patchIDs, bool correctWalls=true)
Construct from mesh and patches to initialize to 0 and flag.
Definition: patchWave.C:137
Foam::patchWave::~patchWave
virtual ~patchWave()
Destructor.
Definition: patchWave.C:156
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85