wallPointI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "polyMesh.H"
27 #include "transform.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Update this with w2 if w2 nearer to pt.
32 template<class TrackingData>
33 inline bool Foam::wallPoint::update
34 (
35  const point& pt,
36  const wallPoint& w2,
37  const scalar tol,
38  TrackingData& td
39 )
40 {
41  // Already done in calling algorithm
42  // if (w2.origin() == origin_)
43  //{
44  // // Shortcut. Same input so same distance.
45  // return false;
46  //}
47 
48  scalar dist2 = magSqr(pt - w2.origin());
49 
50  if (!valid(td))
51  {
52  // current not yet set so use any value
53  distSqr_ = dist2;
54  origin_ = w2.origin();
55 
56  return true;
57  }
58 
59  scalar diff = distSqr_ - dist2;
60 
61  if (diff < 0)
62  {
63  // already nearer to pt
64  return false;
65  }
66 
67  if ((diff < small) || ((distSqr_ > small) && (diff/distSqr_ < tol)))
68  {
69  // don't propagate small changes
70  return false;
71  }
72  else
73  {
74  // update with new values
75  distSqr_ = dist2;
76  origin_ = w2.origin();
77 
78  return true;
79  }
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
86 :
87  origin_(point::max),
88  distSqr_(-great)
89 {}
90 
91 
92 inline Foam::wallPoint::wallPoint(const point& origin, const scalar distSqr)
93 :
94  origin_(origin),
95  distSqr_(distSqr)
96 {}
97 
98 
100 :
101  origin_(wpt.origin()),
102  distSqr_(wpt.distSqr())
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 {
110  return origin_;
111 }
112 
113 
115 {
116  return origin_;
117 }
118 
119 
120 inline Foam::scalar Foam::wallPoint::distSqr() const
121 {
122  return distSqr_;
123 }
124 
125 
126 inline Foam::scalar& Foam::wallPoint::distSqr()
127 {
128  return distSqr_;
129 }
130 
131 
132 template<class TrackingData>
133 inline bool Foam::wallPoint::valid(TrackingData& td) const
134 {
135  return distSqr_ > -small;
136 }
137 
138 
139 // Checks for cyclic faces
140 template<class TrackingData>
142 (
143  const polyMesh&,
144  const wallPoint& w2,
145  const scalar tol,
146  TrackingData& td
147 ) const
148 {
149  scalar diff = mag(distSqr() - w2.distSqr());
150 
151  if (diff < small)
152  {
153  return true;
154  }
155  else
156  {
157  if ((distSqr() > small) && ((diff/distSqr()) < tol))
158  {
159  return true;
160  }
161  else
162  {
163  return false;
164  }
165  }
166 }
167 
168 
169 template<class TrackingData>
171 (
172  const polyMesh&,
173  const polyPatch&,
174  const label,
175  const point& faceCentre,
176  TrackingData& td
177 )
178 {
179  origin_ -= faceCentre;
180 }
181 
182 
183 template<class TrackingData>
184 inline void Foam::wallPoint::transform
185 (
186  const polyMesh&,
187  const tensor& rotTensor,
188  TrackingData& td
189 )
190 {
191  origin_ = Foam::transform(rotTensor, origin_);
192 }
193 
194 
195 // Update absolute geometric quantities. Note that distance (distSqr_)
196 // is not affected by leaving/entering domain.
197 template<class TrackingData>
199 (
200  const polyMesh&,
201  const polyPatch&,
202  const label,
203  const point& faceCentre,
204  TrackingData& td
205 )
206 {
207  // back to absolute form
208  origin_ += faceCentre;
209 }
210 
211 
212 // Update this with w2 if w2 nearer to pt.
213 template<class TrackingData>
214 inline bool Foam::wallPoint::updateCell
215 (
216  const polyMesh& mesh,
217  const label thisCelli,
218  const label neighbourFacei,
219  const wallPoint& neighbourWallInfo,
220  const scalar tol,
221  TrackingData& td
222 )
223 {
224  return
225  update
226  (
227  mesh.cellCentres()[thisCelli],
228  neighbourWallInfo,
229  tol,
230  td
231  );
232  }
233 
234 
235 // Update this with w2 if w2 nearer to pt.
236 template<class TrackingData>
237 inline bool Foam::wallPoint::updateFace
238 (
239  const polyMesh& mesh,
240  const label thisFacei,
241  const label neighbourCelli,
242  const wallPoint& neighbourWallInfo,
243  const scalar tol,
244  TrackingData& td
245 )
246 {
247  return
248  update
249  (
250  mesh.faceCentres()[thisFacei],
251  neighbourWallInfo,
252  tol,
253  td
254  );
255 }
256 
257 // Update this with w2 if w2 nearer to pt.
258 template<class TrackingData>
259 inline bool Foam::wallPoint::updateFace
260 (
261  const polyMesh& mesh,
262  const label thisFacei,
263  const wallPoint& neighbourWallInfo,
264  const scalar tol,
265  TrackingData& td
266 )
267 {
268  return
269  update
270  (
271  mesh.faceCentres()[thisFacei],
272  neighbourWallInfo,
273  tol,
274  td
275  );
276 }
277 
278 
279 template<class TrackingData>
280 inline bool Foam::wallPoint::equal
281 (
282  const wallPoint& rhs,
283  TrackingData& td
284 ) const
285 {
286  return operator==(rhs);
287 }
288 
289 
290 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
291 
292 inline bool Foam::wallPoint::operator==(const Foam::wallPoint& rhs) const
293 {
294  return origin() == rhs.origin();
295 }
296 
297 
298 inline bool Foam::wallPoint::operator!=(const Foam::wallPoint& rhs) const
299 {
300  return !(*this == rhs);
301 }
302 
303 
304 // ************************************************************************* //
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Definition: wallPointI.H:185
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
Definition: wallPointI.H:199
bool sameGeometry(const polyMesh &, const wallPoint &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Definition: wallPointI.H:142
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face.
Definition: wallPointI.H:171
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: wallPointI.H:133
bool operator!=(const wallPoint &) const
Definition: wallPointI.H:298
bool operator==(const wallPoint &) const
Definition: wallPointI.H:292
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:215
3D tensor transformation operations.
const vectorField & cellCentres() const
wallPoint()
Construct null.
Definition: wallPointI.H:85
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Holds information regarding nearest wall point. Used in wall distance calculation.
Definition: wallPoint.H:63
const vectorField & faceCentres() const
bool equal(const wallPoint &, TrackingData &td) const
Same (like operator==)
Definition: wallPointI.H:281
const point & origin() const
Definition: wallPointI.H:108
vector point
Point is a vector.
Definition: point.H:41
scalar distSqr() const
Definition: wallPointI.H:120
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:238
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477