All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
pointEdgePoint.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-2019 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 Class
25  Foam::pointEdgePoint
26 
27 Description
28  Holds information regarding nearest wall point. Used in PointEdgeWave.
29  (so not standard FaceCellWave)
30  To be used in wall distance calculation.
31 
32 SourceFiles
33  pointEdgePointI.H
34  pointEdgePoint.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef pointEdgePoint_H
39 #define pointEdgePoint_H
40 
41 #include "point.H"
42 #include "label.H"
43 #include "scalar.H"
44 #include "tensor.H"
45 #include "pTraits.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of classes
53 class polyPatch;
54 class polyMesh;
55 
56 
57 // Forward declaration of friend functions and operators
58 
59 class pointEdgePoint;
60 
61 Istream& operator>>(Istream&, pointEdgePoint&);
62 Ostream& operator<<(Ostream&, const pointEdgePoint&);
63 
64 
65 /*---------------------------------------------------------------------------*\
66  Class pointEdgePoint Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 class pointEdgePoint
70 {
71  // Private Data
72 
73  //- Position of nearest wall center
74  point origin_;
75 
76  //- Normal distance (squared) from point to origin
77  scalar distSqr_;
78 
79 
80  // Private Member Functions
81 
82  //- Evaluate distance to point. Update distSqr, origin from whomever
83  // is nearer pt. Return true if w2 is closer to point,
84  // false otherwise.
85  template<class TrackingData>
86  inline bool update
87  (
88  const point&,
89  const pointEdgePoint& w2,
90  const scalar tol,
91  TrackingData& td
92  );
93 
94  //- Combine current with w2. Update distSqr, origin if w2 has smaller
95  // quantities and returns true.
96  template<class TrackingData>
97  inline bool update
98  (
99  const pointEdgePoint& w2,
100  const scalar tol,
101  TrackingData& td
102  );
103 
104 
105 public:
106 
107  // Constructors
108 
109  //- Construct null
110  inline pointEdgePoint();
111 
112  //- Construct from origin, distance
113  inline pointEdgePoint(const point&, const scalar);
114 
115 
116  // Member Functions
117 
118  // Access
119 
120  inline const point& origin() const;
121 
122  inline scalar distSqr() const;
123 
124 
125  // Needed by PointEdgeWave
126 
127  //- Check whether origin has been changed at all or
128  // still contains original (invalid) value.
129  template<class TrackingData>
130  inline bool valid(TrackingData& td) const;
131 
132  //- Check for identical geometrical data. Used for cyclics checking.
133  template<class TrackingData>
134  inline bool sameGeometry
135  (
136  const pointEdgePoint&,
137  const scalar tol,
138  TrackingData& td
139  ) const;
140 
141  //- Convert origin to relative vector to leaving point
142  // (= point coordinate)
143  template<class TrackingData>
144  inline void leaveDomain
145  (
146  const polyPatch& patch,
147  const label patchPointi,
148  const point& pos,
149  TrackingData& td
150  );
151 
152  //- Convert relative origin to absolute by adding entering point
153  template<class TrackingData>
154  inline void enterDomain
155  (
156  const polyPatch& patch,
157  const label patchPointi,
158  const point& pos,
159  TrackingData& td
160  );
161 
162  //- Apply rotation matrix to origin
163  template<class TrackingData>
164  inline void transform
165  (
166  const tensor& rotTensor,
167  TrackingData& td
168  );
169 
170  //- Influence of edge on point
171  template<class TrackingData>
172  inline bool updatePoint
173  (
174  const polyMesh& mesh,
175  const label pointi,
176  const label edgeI,
177  const pointEdgePoint& edgeInfo,
178  const scalar tol,
179  TrackingData& td
180  );
181 
182  //- Influence of different value on same point.
183  // Merge new and old info.
184  template<class TrackingData>
185  inline bool updatePoint
186  (
187  const polyMesh& mesh,
188  const label pointi,
189  const pointEdgePoint& newPointInfo,
190  const scalar tol,
191  TrackingData& td
192  );
193 
194  //- Influence of different value on same point.
195  // No information about current position whatsoever.
196  template<class TrackingData>
197  inline bool updatePoint
198  (
199  const pointEdgePoint& newPointInfo,
200  const scalar tol,
201  TrackingData& td
202  );
203 
204  //- Influence of point on edge.
205  template<class TrackingData>
206  inline bool updateEdge
207  (
208  const polyMesh& mesh,
209  const label edgeI,
210  const label pointi,
211  const pointEdgePoint& pointInfo,
212  const scalar tol,
213  TrackingData& td
214  );
215 
216  //- Same (like operator==)
217  template<class TrackingData>
218  inline bool equal(const pointEdgePoint&, TrackingData& td) const;
219 
220 
221  // Member Operators
222 
223  // Needed for List IO
224  inline bool operator==(const pointEdgePoint&) const;
225  inline bool operator!=(const pointEdgePoint&) const;
226 
227 
228  // IOstream Operators
229 
230  friend Ostream& operator<<(Ostream&, const pointEdgePoint&);
232 };
233 
234 
235 //- Data associated with pointEdgePoint type are contiguous
236 template<>
237 inline bool contiguous<pointEdgePoint>()
238 {
239  return true;
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace Foam
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #include "pointEdgePointI.H"
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 #endif
254 
255 // ************************************************************************* //
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
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()
Construct null.
#define w2
Definition: blockCreate.C:32
const point & origin() const
bool contiguous< pointEdgePoint >()
Data associated with pointEdgePoint type are contiguous.
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
scalar distSqr() const
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.
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. Used for cyclics checking.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
friend Istream & operator>>(Istream &, pointEdgePoint &)
Ostream & operator<<(Ostream &, const ensightPart &)
bool equal(const pointEdgePoint &, TrackingData &td) const
Same (like operator==)
friend Ostream & operator<<(Ostream &, const pointEdgePoint &)
bool operator!=(const pointEdgePoint &) const
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
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
Namespace for OpenFOAM.