All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2019 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 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  return origin_;
104 }
105 
106 
108 {
109  return origin_;
110 }
111 
112 
113 inline Foam::scalar Foam::wallPoint::distSqr() const
114 {
115  return distSqr_;
116 }
117 
118 
119 inline Foam::scalar& Foam::wallPoint::distSqr()
120 {
121  return distSqr_;
122 }
123 
124 
125 template<class TrackingData>
126 inline bool Foam::wallPoint::valid(TrackingData& td) const
127 {
128  return distSqr_ > -small;
129 }
130 
131 
132 // Checks for cyclic faces
133 template<class TrackingData>
135 (
136  const polyMesh&,
137  const wallPoint& w2,
138  const scalar tol,
139  TrackingData& td
140 ) const
141 {
142  scalar diff = mag(distSqr() - w2.distSqr());
143 
144  if (diff < small)
145  {
146  return true;
147  }
148  else
149  {
150  if ((distSqr() > small) && ((diff/distSqr()) < tol))
151  {
152  return true;
153  }
154  else
155  {
156  return false;
157  }
158  }
159 }
160 
161 
162 template<class TrackingData>
164 (
165  const polyMesh&,
166  const polyPatch&,
167  const label,
168  const point& faceCentre,
169  TrackingData& td
170 )
171 {
172  origin_ -= faceCentre;
173 }
174 
175 
176 template<class TrackingData>
177 inline void Foam::wallPoint::transform
178 (
179  const polyMesh&,
180  const tensor& rotTensor,
181  TrackingData& td
182 )
183 {
184  origin_ = Foam::transform(rotTensor, origin_);
185 }
186 
187 
188 // Update absolute geometric quantities. Note that distance (distSqr_)
189 // is not affected by leaving/entering domain.
190 template<class TrackingData>
192 (
193  const polyMesh&,
194  const polyPatch&,
195  const label,
196  const point& faceCentre,
197  TrackingData& td
198 )
199 {
200  // back to absolute form
201  origin_ += faceCentre;
202 }
203 
204 
205 // Update this with w2 if w2 nearer to pt.
206 template<class TrackingData>
207 inline bool Foam::wallPoint::updateCell
208 (
209  const polyMesh& mesh,
210  const label thisCelli,
211  const label neighbourFacei,
212  const wallPoint& neighbourWallInfo,
213  const scalar tol,
214  TrackingData& td
215 )
216 {
217  return
218  update
219  (
220  mesh.cellCentres()[thisCelli],
221  neighbourWallInfo,
222  tol,
223  td
224  );
225  }
226 
227 
228 // Update this with w2 if w2 nearer to pt.
229 template<class TrackingData>
230 inline bool Foam::wallPoint::updateFace
231 (
232  const polyMesh& mesh,
233  const label thisFacei,
234  const label neighbourCelli,
235  const wallPoint& neighbourWallInfo,
236  const scalar tol,
237  TrackingData& td
238 )
239 {
240  return
241  update
242  (
243  mesh.faceCentres()[thisFacei],
244  neighbourWallInfo,
245  tol,
246  td
247  );
248 }
249 
250 // Update this with w2 if w2 nearer to pt.
251 template<class TrackingData>
252 inline bool Foam::wallPoint::updateFace
253 (
254  const polyMesh& mesh,
255  const label thisFacei,
256  const wallPoint& neighbourWallInfo,
257  const scalar tol,
258  TrackingData& td
259 )
260 {
261  return
262  update
263  (
264  mesh.faceCentres()[thisFacei],
265  neighbourWallInfo,
266  tol,
267  td
268  );
269 }
270 
271 
272 template<class TrackingData>
273 inline bool Foam::wallPoint::equal
274 (
275  const wallPoint& rhs,
276  TrackingData& td
277 ) const
278 {
279  return operator==(rhs);
280 }
281 
282 
283 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
284 
285 inline bool Foam::wallPoint::operator==(const Foam::wallPoint& rhs) const
286 {
287  return origin() == rhs.origin();
288 }
289 
290 
291 inline bool Foam::wallPoint::operator!=(const Foam::wallPoint& rhs) const
292 {
293  return !(*this == rhs);
294 }
295 
296 
297 // ************************************************************************* //
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
Definition: wallPointI.H:178
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:192
bool sameGeometry(const polyMesh &, const wallPoint &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Definition: wallPointI.H:135
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:164
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: wallPointI.H:126
bool operator!=(const wallPoint &) const
Definition: wallPointI.H:291
bool operator==(const wallPoint &) const
Definition: wallPointI.H:285
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:208
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:274
const point & origin() const
Definition: wallPointI.H:101
vector point
Point is a vector.
Definition: point.H:41
scalar distSqr() const
Definition: wallPointI.H:113
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:231
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