pointEdgePointI.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-2022 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 "transformer.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 
121 :
122  origin_(point::max),
123  distSqr_(great)
124 {}
125 
126 
128 (
129  const point& origin,
130  const scalar distSqr
131 )
132 :
133  origin_(origin),
134  distSqr_(distSqr)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
142  return origin_;
143 }
144 
145 
146 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
147 {
148  return distSqr_;
149 }
150 
151 
152 template<class TrackingData>
153 inline bool Foam::pointEdgePoint::valid(TrackingData& td) const
154 {
155  return origin_ != point::max;
156 }
157 
158 
159 // Checks for cyclic points
160 template<class TrackingData>
162 (
163  const pointEdgePoint& w2,
164  const scalar tol,
165  TrackingData& td
166 ) const
167 {
168  scalar diff = Foam::mag(distSqr() - w2.distSqr());
169 
170  if (diff < small)
171  {
172  return true;
173  }
174  else
175  {
176  if ((distSqr() > small) && ((diff/distSqr()) < tol))
177  {
178  return true;
179  }
180  else
181  {
182  return false;
183  }
184  }
185 }
186 
187 
188 template<class TrackingData>
190 (
191  const polyPatch& patch,
192  const label patchFacei,
193  const transformer& transform,
194  TrackingData& td
195 )
196 {
197  // Note that distSqr_ is not affected by crossing an interface
198  origin_ = transform.transformPosition(origin_);
199 }
200 
201 
202 // Update this with information from connected edge
203 template<class TrackingData>
205 (
206  const polyMesh& mesh,
207  const label pointi,
208  const label edgeI,
209  const pointEdgePoint& edgeInfo,
210  const scalar tol,
211  TrackingData& td
212 )
213 {
214  return update(mesh.points()[pointi], edgeInfo, tol, td);
215 }
216 
217 
218 // Update this with new information on same point
219 template<class TrackingData>
221 (
222  const polyMesh& mesh,
223  const label pointi,
224  const pointEdgePoint& newPointInfo,
225  const scalar tol,
226  TrackingData& td
227 )
228 {
229  return update(mesh.points()[pointi], newPointInfo, tol, td);
230 }
231 
232 
233 // Update this with new information on same point. No extra information.
234 template<class TrackingData>
236 (
237  const pointEdgePoint& newPointInfo,
238  const scalar tol,
239  TrackingData& td
240 )
241 {
242  return update(newPointInfo, tol, td);
243 }
244 
245 
246 // Update this with information from connected point
247 template<class TrackingData>
249 (
250  const polyMesh& mesh,
251  const label edgeI,
252  const label pointi,
253  const pointEdgePoint& pointInfo,
254  const scalar tol,
255  TrackingData& td
256 )
257 {
258  const edge& e = mesh.edges()[edgeI];
259  return update(e.centre(mesh.points()), pointInfo, tol, td);
260 }
261 
262 
263 template<class TrackingData>
264 inline bool Foam::pointEdgePoint::equal
265 (
266  const pointEdgePoint& rhs,
267  TrackingData& td
268 ) const
269 {
270  return operator==(rhs);
271 }
272 
273 
274 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
275 
277 const
278 {
279  return (origin() == rhs.origin()) && (distSqr() == rhs.distSqr());
280 }
281 
282 
284 const
285 {
286  return !(*this == rhs);
287 }
288 
289 
290 // ************************************************************************* //
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
pointEdgePoint()
Construct null.
const point & origin() const
point centre(const pointField &) const
Return centre (centroid)
Definition: edgeI.H:169
scalar distSqr() const
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
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.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
bool sameGeometry(const pointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
vector point
Point is a vector.
Definition: point.H:41
bool equal(const pointEdgePoint &, TrackingData &td) const
Same (like operator==)
dimensioned< scalar > mag(const dimensioned< Type > &)
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
bool operator!=(const pointEdgePoint &) const
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
bool operator==(const pointEdgePoint &) const