pointEdgePoint.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  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::pointEdgePoint
29 
30 Description
31  Holds information regarding nearest wall point. Used in PointEdgeWave.
32  (so not standard FaceCellWave)
33  To be used in wall distance calculation.
34 
35 SourceFiles
36  pointEdgePointI.H
37  pointEdgePoint.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef pointEdgePoint_H
42 #define pointEdgePoint_H
43 
44 #include "point.H"
45 #include "label.H"
46 #include "scalar.H"
47 #include "tensor.H"
48 #include "pTraits.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward Declarations
56 class polyPatch;
57 class polyMesh;
58 class pointEdgePoint;
59 
60 Istream& operator>>(Istream&, pointEdgePoint&);
61 Ostream& operator<<(Ostream&, const pointEdgePoint&);
62 
63 /*---------------------------------------------------------------------------*\
64  Class pointEdgePoint Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class pointEdgePoint
68 {
69  // Private Data
70 
71  //- Position of nearest wall center
72  point origin_;
73 
74  //- Normal distance (squared) from point to origin
75  scalar distSqr_;
76 
77 
78  // Private Member Functions
79 
80  //- Evaluate distance to point.
81  // Update distSqr, origin from whomever is nearer pt.
82  // \return true if w2 is closer to point, false otherwise.
83  template<class TrackingData>
84  inline bool update
85  (
86  const point&,
87  const pointEdgePoint& w2,
88  const scalar tol,
89  TrackingData& td
90  );
91 
92  //- Combine current with w2. Update distSqr, origin if w2 has smaller
93  // quantities and returns true.
94  template<class TrackingData>
95  inline bool update
96  (
97  const pointEdgePoint& w2,
98  const scalar tol,
99  TrackingData& td
100  );
101 
102 
103 public:
104 
105  // Constructors
106 
107  //- Default construct
108  inline pointEdgePoint();
109 
110  //- Construct from origin, distance
111  inline pointEdgePoint(const point& origin, const scalar distSqr);
112 
113 
114  // Member Functions
115 
116  // Access
117 
118  const point& origin() const
119  {
120  return origin_;
121  }
122  point& origin()
123  {
124  return origin_;
125  }
126 
127  scalar distSqr() const
128  {
129  return distSqr_;
130  }
131  scalar& distSqr()
132  {
133  return distSqr_;
134  }
135 
136 
137  // Needed by PointEdgeWave
138 
139  //- Changed or contains original (invalid) value
140  template<class TrackingData>
141  inline bool valid(TrackingData& td) const;
142 
143  //- Check for identical geometrical data (eg, cyclics checking)
144  template<class TrackingData>
145  inline bool sameGeometry
146  (
147  const pointEdgePoint&,
148  const scalar tol,
149  TrackingData& td
150  ) const;
151 
152  //- Convert origin to relative vector to leaving point
153  // (= point coordinate)
154  template<class TrackingData>
155  inline void leaveDomain
156  (
157  const polyPatch& patch,
158  const label patchPointi,
159  const point& pos,
160  TrackingData& td
161  );
162 
163  //- Convert relative origin to absolute by adding entering point
164  template<class TrackingData>
165  inline void enterDomain
166  (
167  const polyPatch& patch,
168  const label patchPointi,
169  const point& pos,
170  TrackingData& td
171  );
172 
173  //- Apply rotation matrix to origin
174  template<class TrackingData>
175  inline void transform
176  (
177  const tensor& rotTensor,
178  TrackingData& td
179  );
180 
181  //- Influence of edge on point
182  template<class TrackingData>
183  inline bool updatePoint
184  (
185  const polyMesh& mesh,
186  const label pointi,
187  const label edgeI,
188  const pointEdgePoint& edgeInfo,
189  const scalar tol,
190  TrackingData& td
191  );
192 
193  //- Influence of different value on same point.
194  // Merge new and old info.
195  template<class TrackingData>
196  inline bool updatePoint
197  (
198  const polyMesh& mesh,
199  const label pointi,
200  const pointEdgePoint& newPointInfo,
201  const scalar tol,
202  TrackingData& td
203  );
204 
205  //- Influence of different value on same point.
206  // No information about current position whatsoever.
207  template<class TrackingData>
208  inline bool updatePoint
209  (
210  const pointEdgePoint& newPointInfo,
211  const scalar tol,
212  TrackingData& td
213  );
214 
215  //- Influence of point on edge.
216  template<class TrackingData>
217  inline bool updateEdge
218  (
219  const polyMesh& mesh,
220  const label edgeI,
221  const label pointi,
222  const pointEdgePoint& pointInfo,
223  const scalar tol,
224  TrackingData& td
225  );
226 
227  //- Test for equality, with TrackingData
228  template<class TrackingData>
229  inline bool equal(const pointEdgePoint&, TrackingData& td) const;
230 
231 
232  // Member Operators
233 
234  //- Test for equality
235  inline bool operator==(const pointEdgePoint&) const;
236 
237  //- Test for inequality
238  inline bool operator!=(const pointEdgePoint&) const;
239 
240 
241  // IOstream Operators
242 
243  friend Ostream& operator<<(Ostream&, const pointEdgePoint&);
244  friend Istream& operator>>(Istream&, pointEdgePoint&);
245 };
246 
247 
248 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
249 
250 //- Contiguous data for pointEdgePoint
251 template<> struct is_contiguous<pointEdgePoint> : std::true_type {};
252 
253 //- Contiguous scalar data for pointEdgePoint
254 template<> struct is_contiguous_scalar<pointEdgePoint> : std::true_type {};
255 
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 } // End namespace Foam
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 #include "pointEdgePointI.H"
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 #endif
268 
269 // ************************************************************************* //
Foam::pointEdgePoint::pointEdgePoint
pointEdgePoint()
Default construct.
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:222
Foam::Tensor< scalar >
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:266
pointEdgePointI.H
Foam::pointEdgePoint::distSqr
scalar distSqr() const
Definition: pointEdgePoint.H:126
update
mesh update()
Foam::pointEdgePoint::valid
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Definition: pointEdgePointI.H:144
point.H
Foam::pointEdgePoint::equal
bool equal(const pointEdgePoint &, TrackingData &td) const
Test for equality, with TrackingData.
Definition: pointEdgePointI.H:282
Foam::pointEdgePoint::operator!=
bool operator!=(const pointEdgePoint &) const
Test for inequality.
Definition: pointEdgePointI.H:303
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:230
Foam::pointEdgePoint::operator>>
friend Istream & operator>>(Istream &, pointEdgePoint &)
Foam::pointEdgePoint::operator==
bool operator==(const pointEdgePoint &) const
Test for equality.
Definition: pointEdgePointI.H:294
tensor.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::pointEdgePoint::origin
const point & origin() const
Definition: pointEdgePoint.H:117
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::pointEdgePoint::operator<<
friend Ostream & operator<<(Ostream &, const pointEdgePoint &)
Foam::pointEdgePoint::distSqr
scalar & distSqr()
Definition: pointEdgePoint.H:130
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:207
Foam::pointEdgePoint::transform
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
Definition: pointEdgePointI.H:194
scalar.H
w2
#define w2
Definition: blockCreate.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
pTraits.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::is_contiguous_scalar
A template class to specify if a data type is composed solely of Foam::scalar elements.
Definition: contiguous.H:91
Foam::pointEdgePoint
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgePoint.H:66
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
label.H
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:181
Foam::pointEdgePoint::origin
point & origin()
Definition: pointEdgePoint.H:121
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
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 (eg, cyclics checking)
Definition: pointEdgePointI.H:153
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::is_contiguous
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:75
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177