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