wallPointI.H
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 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 "polyMesh.H"
30 #include "transform.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 // Update this with w2 if w2 nearer to pt.
35 template<class TrackingData>
36 inline bool Foam::wallPoint::update
37 (
38  const point& pt,
39  const wallPoint& w2,
40  const scalar tol,
41  TrackingData& td
42 )
43 {
44  //Already done in calling algorithm
45  //if (w2.origin() == origin_)
46  //{
47  // // Shortcut. Same input so same distance.
48  // return false;
49  //}
50 
51  const scalar dist2 = magSqr(pt - w2.origin());
52 
53  if (!valid(td))
54  {
55  // current not yet set so use any value
56  distSqr_ = dist2;
57  origin_ = w2.origin();
58 
59  return true;
60  }
61 
62  const scalar diff = distSqr_ - dist2;
63 
64  if (diff < 0)
65  {
66  // already nearer to pt
67  return false;
68  }
69 
70  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
71  {
72  // don't propagate small changes
73  return false;
74  }
75  else
76  {
77  // update with new values
78  distSqr_ = dist2;
79  origin_ = w2.origin();
80 
81  return true;
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 :
90  origin_(point::max),
91  distSqr_(-GREAT)
92 {}
93 
94 
95 inline Foam::wallPoint::wallPoint(const point& origin, const scalar distSqr)
96 :
97  origin_(origin),
98  distSqr_(distSqr)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class TrackingData>
105 inline bool Foam::wallPoint::valid(TrackingData& td) const
106 {
107  return distSqr_ > -SMALL;
108 }
109 
110 
111 // Checks for cyclic faces
112 template<class TrackingData>
114 (
115  const polyMesh&,
116  const wallPoint& w2,
117  const scalar tol,
118  TrackingData& td
119 ) const
120 {
121  const scalar diff = mag(distSqr() - w2.distSqr());
122 
123  if (diff < SMALL)
124  {
125  return true;
126  }
127  else
128  {
129  if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
130  {
131  return true;
132  }
133  else
134  {
135  return false;
136  }
137  }
138 }
139 
140 
141 template<class TrackingData>
143 (
144  const polyMesh&,
145  const polyPatch&,
146  const label,
147  const point& faceCentre,
148  TrackingData& td
149 )
150 {
151  origin_ -= faceCentre;
152 }
153 
154 
155 template<class TrackingData>
156 inline void Foam::wallPoint::transform
157 (
158  const polyMesh&,
159  const tensor& rotTensor,
160  TrackingData& td
161 )
162 {
163  origin_ = Foam::transform(rotTensor, origin_);
164 }
165 
166 
167 // Update absolute geometric quantities. Note that distance (distSqr_)
168 // is not affected by leaving/entering domain.
169 template<class TrackingData>
171 (
172  const polyMesh&,
173  const polyPatch&,
174  const label,
175  const point& faceCentre,
176  TrackingData& td
177 )
178 {
179  // back to absolute form
180  origin_ += faceCentre;
181 }
182 
183 
184 // Update this with w2 if w2 nearer to pt.
185 template<class TrackingData>
186 inline bool Foam::wallPoint::updateCell
187 (
188  const polyMesh& mesh,
189  const label thisCelli,
190  const label neighbourFacei,
191  const wallPoint& neighbourWallInfo,
192  const scalar tol,
193  TrackingData& td
194 )
195 {
196  return
197  update
198  (
199  mesh.cellCentres()[thisCelli],
200  neighbourWallInfo,
201  tol,
202  td
203  );
204 }
205 
206 
207 // Update this with w2 if w2 nearer to pt.
208 template<class TrackingData>
209 inline bool Foam::wallPoint::updateFace
210 (
211  const polyMesh& mesh,
212  const label thisFacei,
213  const label neighbourCelli,
214  const wallPoint& neighbourWallInfo,
215  const scalar tol,
216  TrackingData& td
217 )
218 {
219  return
220  update
221  (
222  mesh.faceCentres()[thisFacei],
223  neighbourWallInfo,
224  tol,
225  td
226  );
227 }
228 
229 
230 // Update this with w2 if w2 nearer to pt.
231 template<class TrackingData>
232 inline bool Foam::wallPoint::updateFace
233 (
234  const polyMesh& mesh,
235  const label thisFacei,
236  const wallPoint& neighbourWallInfo,
237  const scalar tol,
238  TrackingData& td
239 )
240 {
241  return
242  update
243  (
244  mesh.faceCentres()[thisFacei],
245  neighbourWallInfo,
246  tol,
247  td
248  );
249 }
250 
251 
252 template<class TrackingData>
253 inline bool Foam::wallPoint::equal
254 (
255  const wallPoint& rhs,
256  TrackingData& td
257 ) const
258 {
259  return operator==(rhs);
260 }
261 
262 
263 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
264 
265 inline bool Foam::wallPoint::operator==
266 (
267  const wallPoint& rhs
268 ) const
269 {
270  return origin_ == rhs.origin_;
271 }
272 
273 
274 inline bool Foam::wallPoint::operator!=
275 (
276  const wallPoint& rhs
277 ) const
278 {
279  return !(*this == rhs);
280 }
281 
282 
283 // ************************************************************************* //
Foam::wallPoint::sameGeometry
bool sameGeometry(const polyMesh &, const wallPoint &, const scalar, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
Definition: wallPointI.H:114
Foam::wallPoint
Holds information regarding nearest wall point. Used in wall distance calculation.
Definition: wallPoint.H:65
Foam::Tensor< scalar >
Foam::wallPoint::valid
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Definition: wallPointI.H:105
update
mesh update()
Foam::wallPoint::updateFace
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const wallPoint &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: wallPointI.H:210
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
polyMesh.H
Foam::wallPoint::wallPoint
wallPoint()
Default construct.
Definition: wallPointI.H:88
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Foam::wallPoint::enterDomain
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
Definition: wallPointI.H:171
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::wallPoint::leaveDomain
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Definition: wallPointI.H:143
Foam::wallPoint::updateCell
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const wallPoint &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: wallPointI.H:187
w2
#define w2
Definition: blockCreate.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
Foam::wallPoint::equal
bool equal(const wallPoint &, TrackingData &td) const
Test for equality, with TrackingData.
Definition: wallPointI.H:254
transform.H
3D tensor transformation operations.
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::wallPoint::transform
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Definition: wallPointI.H:157