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-2023 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::pointEdgePoint::update
33 (
34  const point& pt,
35  const pointEdgePoint& 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 :
79  origin_(point::max),
80  distSqr_(great)
81 {}
82 
83 
85 (
86  const point& origin,
87  const scalar distSqr
88 )
89 :
90  origin_(origin),
91  distSqr_(distSqr)
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
99  return origin_;
100 }
101 
102 
103 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
104 {
105  return distSqr_;
106 }
107 
108 
109 template<class TrackingData>
110 inline bool Foam::pointEdgePoint::valid(TrackingData& td) const
111 {
112  return origin_ != point::max;
113 }
114 
115 
116 template<class TrackingData>
118 (
119  const polyPatch& patch,
120  const label patchFacei,
121  const transformer& transform,
122  TrackingData& td
123 )
124 {
125  // Note that distSqr_ is not affected by crossing an interface
126  origin_ = transform.transformPosition(origin_);
127 }
128 
129 
130 // Update this with information from connected edge
131 template<class TrackingData>
133 (
134  const polyMesh& mesh,
135  const label pointi,
136  const label edgeI,
137  const pointEdgePoint& edgeInfo,
138  const scalar tol,
139  TrackingData& td
140 )
141 {
142  return update(mesh.points()[pointi], edgeInfo, tol, td);
143 }
144 
145 
146 // Update this with new information on same point
147 template<class TrackingData>
149 (
150  const polyMesh& mesh,
151  const label pointi,
152  const pointEdgePoint& newPointInfo,
153  const scalar tol,
154  TrackingData& td
155 )
156 {
157  return update(mesh.points()[pointi], newPointInfo, tol, td);
158 }
159 
160 
161 // Update this with information from connected point
162 template<class TrackingData>
164 (
165  const polyMesh& mesh,
166  const label edgeI,
167  const label pointi,
168  const pointEdgePoint& pointInfo,
169  const scalar tol,
170  TrackingData& td
171 )
172 {
173  const edge& e = mesh.edges()[edgeI];
174  return update(e.centre(mesh.points()), pointInfo, tol, td);
175 }
176 
177 
178 template<class TrackingData>
180 (
181  const pointEdgePoint& rhs,
182  TrackingData& td
183 ) const
184 {
185  return operator==(rhs);
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
190 
192 const
193 {
194  return (origin() == rhs.origin()) && (distSqr() == rhs.distSqr());
195 }
196 
197 
199 const
200 {
201  return !(*this == rhs);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
206 
207 inline Foam::Ostream& Foam::operator<<
208 (
209  Foam::Ostream& os,
210  const Foam::pointEdgePoint& wDist
211 )
212 {
213  return os << wDist.origin() << token::SPACE << wDist.distSqr();
214 }
215 
216 
217 inline Foam::Istream& Foam::operator>>
218 (
219  Foam::Istream& is,
220  Foam::pointEdgePoint& wDist
221 )
222 {
223  return is >> wDist.origin_ >> wDist.distSqr_;
224 }
225 
226 
227 // ************************************************************************* //
#define w2
Definition: blockCreate.C:32
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static const Form max
Definition: VectorSpace.H:120
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)...
scalar distSqr() const
bool equal(const pointEdgePoint &, TrackingData &td) const
Same (like operator==)
pointEdgePoint()
Construct null.
bool operator==(const pointEdgePoint &) const
bool operator!=(const pointEdgePoint &) const
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
const point & origin() const
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
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
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
const doubleScalar e
Definition: doubleScalar.H:106
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:409
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > magSqr(const dimensioned< Type > &)