pointEdgePointI.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 | 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 "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::pointEdgePoint::update
34 (
35  const point& pt,
36  const pointEdgePoint& w2,
37  const scalar tol,
38  TrackingData& td
39 )
40 {
41  scalar dist2 = magSqr(pt - w2.origin());
42 
43  if (!valid(td))
44  {
45  // current not yet set so use any value
46  distSqr_ = dist2;
47  origin_ = w2.origin();
48 
49  return true;
50  }
51 
52  scalar diff = distSqr_ - dist2;
53 
54  if (diff < 0)
55  {
56  // already nearer to pt
57  return false;
58  }
59 
60  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
61  {
62  // don't propagate small changes
63  return false;
64  }
65  else
66  {
67  // update with new values
68  distSqr_ = dist2;
69  origin_ = w2.origin();
70 
71  return true;
72  }
73 }
74 
75 
76 // Update this with w2 (information on same point)
77 template<class TrackingData>
78 inline bool Foam::pointEdgePoint::update
79 (
80  const pointEdgePoint& w2,
81  const scalar tol,
82  TrackingData& td
83 )
84 {
85  if (!valid(td))
86  {
87  // current not yet set so use any value
88  distSqr_ = w2.distSqr();
89  origin_ = w2.origin();
90 
91  return true;
92  }
93 
94  scalar diff = distSqr_ - w2.distSqr();
95 
96  if (diff < 0)
97  {
98  // already nearer to pt
99  return false;
100  }
101 
102  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
103  {
104  // don't propagate small changes
105  return false;
106  }
107  else
108  {
109  // update with new values
110  distSqr_ = w2.distSqr();
111  origin_ = w2.origin();
112 
113  return true;
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
120 // Null constructor
122 :
123  origin_(point::max),
124  distSqr_(GREAT)
125 {}
126 
127 
128 // Construct from origin, distance
130 (
131  const point& origin,
132  const scalar distSqr
133 )
134 :
135  origin_(origin),
136  distSqr_(distSqr)
137 {}
138 
139 
140 // Construct as copy
142 :
143  origin_(wpt.origin()),
144  distSqr_(wpt.distSqr())
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  return origin_;
153 }
154 
155 
156 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
157 {
158  return distSqr_;
159 }
160 
161 
162 template<class TrackingData>
163 inline bool Foam::pointEdgePoint::valid(TrackingData& td) const
164 {
165  return origin_ != point::max;
166 }
167 
168 
169 // Checks for cyclic points
170 template<class TrackingData>
172 (
173  const pointEdgePoint& w2,
174  const scalar tol,
175  TrackingData& td
176 ) const
177 {
178  scalar diff = Foam::mag(distSqr() - w2.distSqr());
179 
180  if (diff < SMALL)
181  {
182  return true;
183  }
184  else
185  {
186  if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
187  {
188  return true;
189  }
190  else
191  {
192  return false;
193  }
194  }
195 }
196 
197 
198 template<class TrackingData>
200 (
201  const polyPatch& patch,
202  const label patchPointi,
203  const point& coord,
204  TrackingData& td
205 )
206 {
207  origin_ -= coord;
208 }
209 
210 
211 template<class TrackingData>
213 (
214  const tensor& rotTensor,
215  TrackingData& td
216 )
217 {
218  origin_ = Foam::transform(rotTensor, origin_);
219 }
220 
221 
222 // Update absolute geometric quantities. Note that distance (distSqr_)
223 // is not affected by leaving/entering domain.
224 template<class TrackingData>
226 (
227  const polyPatch& patch,
228  const label patchPointi,
229  const point& coord,
230  TrackingData& td
231 )
232 {
233  // back to absolute form
234  origin_ += coord;
235 }
236 
237 
238 // Update this with information from connected edge
239 template<class TrackingData>
241 (
242  const polyMesh& mesh,
243  const label pointi,
244  const label edgeI,
245  const pointEdgePoint& edgeInfo,
246  const scalar tol,
247  TrackingData& td
248 )
249 {
250  return update(mesh.points()[pointi], edgeInfo, tol, td);
251 }
252 
253 
254 // Update this with new information on same point
255 template<class TrackingData>
257 (
258  const polyMesh& mesh,
259  const label pointi,
260  const pointEdgePoint& newPointInfo,
261  const scalar tol,
262  TrackingData& td
263 )
264 {
265  return update(mesh.points()[pointi], newPointInfo, tol, td);
266 }
267 
268 
269 // Update this with new information on same point. No extra information.
270 template<class TrackingData>
272 (
273  const pointEdgePoint& newPointInfo,
274  const scalar tol,
275  TrackingData& td
276 )
277 {
278  return update(newPointInfo, tol, td);
279 }
280 
281 
282 // Update this with information from connected point
283 template<class TrackingData>
285 (
286  const polyMesh& mesh,
287  const label edgeI,
288  const label pointi,
289  const pointEdgePoint& pointInfo,
290  const scalar tol,
291  TrackingData& td
292 )
293 {
294  const edge& e = mesh.edges()[edgeI];
295  return update(e.centre(mesh.points()), pointInfo, tol, td);
296 }
297 
298 
299 template<class TrackingData>
300 inline bool Foam::pointEdgePoint::equal
301 (
302  const pointEdgePoint& rhs,
303  TrackingData& td
304 ) const
305 {
306  return operator==(rhs);
307 }
308 
309 
310 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
311 
313 const
314 {
315  return (origin() == rhs.origin()) && (distSqr() == rhs.distSqr());
316 }
317 
318 
320 const
321 {
322  return !(*this == rhs);
323 }
324 
325 
326 // ************************************************************************* //
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
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
point centre(const pointField &) const
Return centre (centroid)
Definition: edgeI.H:151
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
pointEdgePoint()
Construct null.
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
const point & origin() const
scalar distSqr() const
bool operator==(const pointEdgePoint &) const
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
3D tensor transformation operations.
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
bool equal(const pointEdgePoint &, TrackingData &td) const
Same (like operator==)
vector point
Point is a vector.
Definition: point.H:41
bool operator!=(const pointEdgePoint &) const
bool sameGeometry(const pointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
dimensioned< scalar > mag(const dimensioned< Type > &)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
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:465