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-2022 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 Foam_pointEdgePoint_H
42 #define Foam_pointEdgePoint_H
43 
44 #include "label.H"
45 #include "scalar.H"
46 #include "point.H"
47 #include "tensor.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // Forward Declarations
55 class polyPatch;
56 class polyMesh;
57 class pointEdgePoint;
58 
59 Istream& operator>>(Istream&, pointEdgePoint&);
60 Ostream& operator<<(Ostream&, const pointEdgePoint&);
61 
62 /*---------------------------------------------------------------------------*\
63  Class pointEdgePoint Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class pointEdgePoint
67 {
68  // Private Data
69 
70  //- Position of nearest wall center
71  point origin_;
72 
73  //- Normal distance (squared) from point to origin
74  scalar distSqr_;
75 
76 
77  // Private Member Functions
78 
79  //- Evaluate distance to point.
80  // Update distSqr, origin from whomever is nearer pt.
81  // \return true if w2 is closer to point, false otherwise.
82  template<class TrackingData>
83  inline bool update
84  (
85  const point&,
86  const pointEdgePoint& w2,
87  const scalar tol,
88  TrackingData& td
89  );
90 
91  //- Combine current with w2. Update distSqr, origin if w2 has smaller
92  // quantities and returns true.
93  template<class TrackingData>
94  inline bool update
95  (
96  const pointEdgePoint& w2,
97  const scalar tol,
98  TrackingData& td
99  );
100 
101 
102 public:
103 
104  // Constructors
105 
106  //- Default construct. Max point
107  inline pointEdgePoint();
108 
109  //- Construct from origin, distance
110  inline pointEdgePoint(const point& origin, const scalar distSqr);
111 
112 
113  // Member Functions
114 
115  // Access
116 
117  const point& origin() const noexcept { return origin_; }
118 
119  point& origin() noexcept { return origin_; }
120 
121  scalar distSqr() const noexcept { return distSqr_; }
122 
123  scalar& distSqr() noexcept { return distSqr_; }
124 
125 
126  // Member Operators
127 
128  //- Test for equality
129  inline bool operator==(const pointEdgePoint&) const;
131  //- Test for inequality
132  inline bool operator!=(const pointEdgePoint&) const;
133 
134 
135  // IOstream Operators
136 
137  friend Ostream& operator<<(Ostream&, const pointEdgePoint&);
139 
140 
141  // Wave Methods
142 
143  // Needed by PointEdgeWave
144 
145  //- Changed or contains original (invalid) value
146  template<class TrackingData>
147  inline bool valid(TrackingData& td) const;
148 
149  //- Check for identical geometrical data (eg, cyclics checking)
150  template<class TrackingData>
151  inline bool sameGeometry
152  (
153  const pointEdgePoint&,
154  const scalar tol,
155  TrackingData& td
156  ) const;
157 
158  //- Convert origin to relative vector to leaving point
159  // (= point coordinate)
160  template<class TrackingData>
161  inline void leaveDomain
162  (
163  const polyPatch& patch,
164  const label patchPointi,
165  const point& pos,
166  TrackingData& td
167  );
168 
169  //- Convert relative origin to absolute by adding entering point
170  template<class TrackingData>
171  inline void enterDomain
172  (
173  const polyPatch& patch,
174  const label patchPointi,
175  const point& pos,
176  TrackingData& td
177  );
178 
179  //- Apply rotation matrix to origin
180  template<class TrackingData>
181  inline void transform
182  (
183  const tensor& rotTensor,
184  TrackingData& td
185  );
186 
187  //- Influence of edge on point
188  template<class TrackingData>
189  inline bool updatePoint
190  (
191  const polyMesh& mesh,
192  const label pointi,
193  const label edgeI,
194  const pointEdgePoint& edgeInfo,
195  const scalar tol,
196  TrackingData& td
197  );
198 
199  //- Influence of different value on same point.
200  // Merge new and old info.
201  template<class TrackingData>
202  inline bool updatePoint
203  (
204  const polyMesh& mesh,
205  const label pointi,
206  const pointEdgePoint& newPointInfo,
207  const scalar tol,
208  TrackingData& td
209  );
210 
211  //- Influence of different value on same point.
212  // No information about current position whatsoever.
213  template<class TrackingData>
214  inline bool updatePoint
215  (
216  const pointEdgePoint& newPointInfo,
217  const scalar tol,
218  TrackingData& td
219  );
220 
221  //- Influence of point on edge.
222  template<class TrackingData>
223  inline bool updateEdge
224  (
225  const polyMesh& mesh,
226  const label edgeI,
227  const label pointi,
228  const pointEdgePoint& pointInfo,
229  const scalar tol,
230  TrackingData& td
231  );
232 
233  //- Test for equality, with TrackingData
234  template<class TrackingData>
235  inline bool equal(const pointEdgePoint&, TrackingData& td) const;
236 };
237 
238 
239 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
240 
241 //- Contiguous data for pointEdgePoint
242 template<> struct is_contiguous<pointEdgePoint> : std::true_type {};
243 
244 //- Contiguous scalar data for pointEdgePoint
245 template<> struct is_contiguous_scalar<pointEdgePoint> : std::true_type {};
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace Foam
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 #include "pointEdgePointI.H"
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
pointEdgePoint()
Default construct. Max point.
#define w2
Definition: blockCreate.C:28
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
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.
const point & origin() const noexcept
dimensionedScalar pos(const dimensionedScalar &ds)
dynamicFvMesh & mesh
Istream & operator>>(Istream &, directionInfo &)
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
bool sameGeometry(const pointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
friend Istream & operator>>(Istream &, pointEdgePoint &)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
mesh update()
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
A template class to specify that a data type can be considered as being contiguous in memory...
Definition: contiguous.H:70
const std::string patch
OpenFOAM patch number as a std::string.
bool equal(const pointEdgePoint &, TrackingData &td) const
Test for equality, with TrackingData.
friend Ostream & operator<<(Ostream &, const pointEdgePoint &)
bool operator!=(const pointEdgePoint &) const
Test for inequality.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Tensor of scalars, i.e. Tensor<scalar>.
scalar distSqr() const noexcept
bool operator==(const pointEdgePoint &) const
Test for equality.
Namespace for OpenFOAM.