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 | Copyright (C) 2011-2016 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  //- Construct as copy
116  inline pointEdgePoint(const pointEdgePoint&);
117 
118 
119  // Member Functions
120 
121  // Access
122 
123  inline const point& origin() const;
124 
125  inline scalar distSqr() const;
126 
127 
128  // Needed by PointEdgeWave
129 
130  //- Check whether origin has been changed at all or
131  // still contains original (invalid) value.
132  template<class TrackingData>
133  inline bool valid(TrackingData& td) const;
134 
135  //- Check for identical geometrical data. Used for cyclics checking.
136  template<class TrackingData>
137  inline bool sameGeometry
138  (
139  const pointEdgePoint&,
140  const scalar tol,
141  TrackingData& td
142  ) const;
143 
144  //- Convert origin to relative vector to leaving point
145  // (= point coordinate)
146  template<class TrackingData>
147  inline void leaveDomain
148  (
149  const polyPatch& patch,
150  const label patchPointi,
151  const point& pos,
152  TrackingData& td
153  );
154 
155  //- Convert relative origin to absolute by adding entering point
156  template<class TrackingData>
157  inline void enterDomain
158  (
159  const polyPatch& patch,
160  const label patchPointi,
161  const point& pos,
162  TrackingData& td
163  );
164 
165  //- Apply rotation matrix to origin
166  template<class TrackingData>
167  inline void transform
168  (
169  const tensor& rotTensor,
170  TrackingData& td
171  );
172 
173  //- Influence of edge on point
174  template<class TrackingData>
175  inline bool updatePoint
176  (
177  const polyMesh& mesh,
178  const label pointi,
179  const label edgeI,
180  const pointEdgePoint& edgeInfo,
181  const scalar tol,
182  TrackingData& td
183  );
184 
185  //- Influence of different value on same point.
186  // Merge new and old info.
187  template<class TrackingData>
188  inline bool updatePoint
189  (
190  const polyMesh& mesh,
191  const label pointi,
192  const pointEdgePoint& newPointInfo,
193  const scalar tol,
194  TrackingData& td
195  );
196 
197  //- Influence of different value on same point.
198  // No information about current position whatsoever.
199  template<class TrackingData>
200  inline bool updatePoint
201  (
202  const pointEdgePoint& newPointInfo,
203  const scalar tol,
204  TrackingData& td
205  );
206 
207  //- Influence of point on edge.
208  template<class TrackingData>
209  inline bool updateEdge
210  (
211  const polyMesh& mesh,
212  const label edgeI,
213  const label pointi,
214  const pointEdgePoint& pointInfo,
215  const scalar tol,
216  TrackingData& td
217  );
218 
219  //- Same (like operator==)
220  template<class TrackingData>
221  inline bool equal(const pointEdgePoint&, TrackingData& td) const;
222 
223 
224  // Member Operators
225 
226  // Needed for List IO
227  inline bool operator==(const pointEdgePoint&) const;
228  inline bool operator!=(const pointEdgePoint&) const;
229 
230 
231  // IOstream Operators
232 
233  friend Ostream& operator<<(Ostream&, const pointEdgePoint&);
235 };
236 
237 
238 //- Data associated with pointEdgePoint type are contiguous
239 template<>
240 inline bool contiguous<pointEdgePoint>()
241 {
242  return true;
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace Foam
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #include "pointEdgePointI.H"
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #endif
257 
258 // ************************************************************************* //
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:53
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.