externalPointEdgePointI.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) 2013-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 template<class TrackingData>
32 inline bool Foam::externalPointEdgePoint::update
33 (
34  const point& pt,
35  const externalPointEdgePoint& w2,
36  const scalar tol,
37  TrackingData& td
38 )
39 {
40  scalar dist2 = magSqr(pt - w2.origin());
41 
42  if (!valid(td))
43  {
44  // current not yet set so use any value
45  distSqr_ = dist2;
46  origin_ = w2.origin();
47 
48  return true;
49  }
50 
51  scalar diff = distSqr_ - dist2;
52 
53  if (diff < 0)
54  {
55  // already nearer to pt
56  return false;
57  }
58 
59  if ((diff < small) || ((distSqr_ > small) && (diff/distSqr_ < tol)))
60  {
61  // don't propagate small changes
62  return false;
63  }
64  else
65  {
66  // update with new values
67  distSqr_ = dist2;
68  origin_ = w2.origin();
69 
70  return true;
71  }
72 }
73 
74 
75 template<class TrackingData>
76 inline bool Foam::externalPointEdgePoint::update
77 (
78  const externalPointEdgePoint& w2,
79  const scalar tol,
80  TrackingData& td
81 )
82 {
83  if (!valid(td))
84  {
85  // current not yet set so use any value
86  distSqr_ = w2.distSqr();
87  origin_ = w2.origin();
88 
89  return true;
90  }
91 
92  scalar diff = distSqr_ - w2.distSqr();
93 
94  if (diff < 0)
95  {
96  // already nearer to pt
97  return false;
98  }
99 
100  if ((diff < small) || ((distSqr_ > small) && (diff/distSqr_ < tol)))
101  {
102  // don't propagate small changes
103  return false;
104  }
105  else
106  {
107  // update with new values
108  distSqr_ = w2.distSqr();
109  origin_ = w2.origin();
110 
111  return true;
112  }
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
119 :
120  origin_(point::max),
121  distSqr_(great)
122 {}
123 
124 
126 (
127  const point& origin,
128  const scalar distSqr
129 )
130 :
131  origin_(origin),
132  distSqr_(distSqr)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  return origin_;
141 }
142 
143 
144 inline Foam::scalar Foam::externalPointEdgePoint::distSqr() const
145 {
146  return distSqr_;
147 }
148 
149 
150 template<class TrackingData>
151 inline bool Foam::externalPointEdgePoint::valid(TrackingData& td) const
152 {
153  return origin_ != point::max;
154 }
155 
156 
157 // Checks for cyclic points
158 template<class TrackingData>
160 (
161  const externalPointEdgePoint& w2,
162  const scalar tol,
163  TrackingData& td
164 ) const
165 {
166  scalar diff = Foam::mag(distSqr() - w2.distSqr());
167 
168  if (diff < small)
169  {
170  return true;
171  }
172  else
173  {
174  if ((distSqr() > small) && ((diff/distSqr()) < tol))
175  {
176  return true;
177  }
178  else
179  {
180  return false;
181  }
182  }
183 }
184 
185 
186 template<class TrackingData>
188 (
189  const polyPatch& patch,
190  const label patchFacei,
191  const transformer& transform,
192  TrackingData& td
193 )
194 {
195  origin_ = transform.transformPosition(origin_);
196 }
197 
198 
199 template<class TrackingData>
201 (
202  const polyMesh& mesh,
203  const label pointi,
204  const label edgeI,
205  const externalPointEdgePoint& edgeInfo,
206  const scalar tol,
207  TrackingData& td
208 )
209 {
210  return update(td.points_[pointi], edgeInfo, tol, td);
211 }
212 
213 
214 template<class TrackingData>
216 (
217  const polyMesh& mesh,
218  const label pointi,
219  const externalPointEdgePoint& newPointInfo,
220  const scalar tol,
221  TrackingData& td
222 )
223 {
224  return update(td.points_[pointi], newPointInfo, tol, td);
225 }
226 
227 
228 template<class TrackingData>
230 (
231  const externalPointEdgePoint& newPointInfo,
232  const scalar tol,
233  TrackingData& td
234 )
235 {
236  return update(newPointInfo, tol, td);
237 }
238 
239 
240 template<class TrackingData>
242 (
243  const polyMesh& mesh,
244  const label edgeI,
245  const label pointi,
246  const externalPointEdgePoint& pointInfo,
247  const scalar tol,
248  TrackingData& td
249 )
250 {
251  const edge& e = mesh.edges()[edgeI];
252  return update(e.centre(td.points_), pointInfo, tol, td);
253 }
254 
255 
256 template<class TrackingData>
258 (
259  const externalPointEdgePoint& rhs,
260  TrackingData& td
261 ) const
262 {
263  return operator==(rhs);
264 }
265 
266 
267 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
268 
269 inline bool Foam::externalPointEdgePoint::operator==
270 (
272 )
273 const
274 {
275  return (origin() == rhs.origin()) && (distSqr() == rhs.distSqr());
276 }
277 
278 
279 inline bool Foam::externalPointEdgePoint::operator!=
280 (
282 )
283 const
284 {
285  return !(*this == rhs);
286 }
287 
288 
289 // ************************************************************************* //
#define w2
Definition: blockCreate.C:32
static const Form max
Definition: VectorSpace.H:115
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool equal(const externalPointEdgePoint &, TrackingData &td) const
Equivalent to operator== with TrackingData.
bool sameGeometry(const externalPointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const externalPointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const externalPointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
bool valid(const PtrList< ModelType > &l)
const doubleScalar e
Definition: doubleScalar.H:105
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
vector point
Point is a vector.
Definition: point.H:41
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:403
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > magSqr(const dimensioned< Type > &)