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-2015 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  Class pointEdgePoint Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class pointEdgePoint
61 {
62  // Private data
63 
64  //- Position of nearest wall center
65  point origin_;
66 
67  //- Normal distance (squared) from point to origin
68  scalar distSqr_;
69 
70 
71  // Private Member Functions
72 
73  //- Evaluate distance to point. Update distSqr, origin from whomever
74  // is nearer pt. Return true if w2 is closer to point,
75  // false otherwise.
76  template<class TrackingData>
77  inline bool update
78  (
79  const point&,
80  const pointEdgePoint& w2,
81  const scalar tol,
82  TrackingData& td
83  );
84 
85  //- Combine current with w2. Update distSqr, origin if w2 has smaller
86  // quantities and returns true.
87  template<class TrackingData>
88  inline bool update
89  (
90  const pointEdgePoint& w2,
91  const scalar tol,
92  TrackingData& td
93  );
94 
95 
96 public:
97 
98  // Constructors
99 
100  //- Construct null
101  inline pointEdgePoint();
102 
103  //- Construct from origin, distance
104  inline pointEdgePoint(const point&, const scalar);
105 
106  //- Construct as copy
107  inline pointEdgePoint(const pointEdgePoint&);
108 
109 
110  // Member Functions
111 
112  // Access
113 
114  inline const point& origin() const;
115 
116  inline scalar distSqr() const;
117 
118 
119  // Needed by PointEdgeWave
120 
121  //- Check whether origin has been changed at all or
122  // still contains original (invalid) value.
123  template<class TrackingData>
124  inline bool valid(TrackingData& td) const;
125 
126  //- Check for identical geometrical data. Used for cyclics checking.
127  template<class TrackingData>
128  inline bool sameGeometry
129  (
130  const pointEdgePoint&,
131  const scalar tol,
132  TrackingData& td
133  ) const;
134 
135  //- Convert origin to relative vector to leaving point
136  // (= point coordinate)
137  template<class TrackingData>
138  inline void leaveDomain
139  (
140  const polyPatch& patch,
141  const label patchPointI,
142  const point& pos,
143  TrackingData& td
144  );
145 
146  //- Convert relative origin to absolute by adding entering point
147  template<class TrackingData>
148  inline void enterDomain
149  (
150  const polyPatch& patch,
151  const label patchPointI,
152  const point& pos,
153  TrackingData& td
154  );
155 
156  //- Apply rotation matrix to origin
157  template<class TrackingData>
158  inline void transform
159  (
160  const tensor& rotTensor,
161  TrackingData& td
162  );
163 
164  //- Influence of edge on point
165  template<class TrackingData>
166  inline bool updatePoint
167  (
168  const polyMesh& mesh,
169  const label pointI,
170  const label edgeI,
171  const pointEdgePoint& edgeInfo,
172  const scalar tol,
173  TrackingData& td
174  );
175 
176  //- Influence of different value on same point.
177  // Merge new and old info.
178  template<class TrackingData>
179  inline bool updatePoint
180  (
181  const polyMesh& mesh,
182  const label pointI,
183  const pointEdgePoint& newPointInfo,
184  const scalar tol,
185  TrackingData& td
186  );
187 
188  //- Influence of different value on same point.
189  // No information about current position whatsoever.
190  template<class TrackingData>
191  inline bool updatePoint
192  (
193  const pointEdgePoint& newPointInfo,
194  const scalar tol,
195  TrackingData& td
196  );
197 
198  //- Influence of point on edge.
199  template<class TrackingData>
200  inline bool updateEdge
201  (
202  const polyMesh& mesh,
203  const label edgeI,
204  const label pointI,
205  const pointEdgePoint& pointInfo,
206  const scalar tol,
207  TrackingData& td
208  );
209 
210  //- Same (like operator==)
211  template<class TrackingData>
212  inline bool equal(const pointEdgePoint&, TrackingData& td) const;
213 
214 
215  // Member Operators
216 
217  // Needed for List IO
218  inline bool operator==(const pointEdgePoint&) const;
219  inline bool operator!=(const pointEdgePoint&) const;
220 
221 
222  // IOstream Operators
223 
224  friend Ostream& operator<<(Ostream&, const pointEdgePoint&);
226 };
227 
228 
229 //- Data associated with pointEdgePoint type are contiguous
230 template<>
231 inline bool contiguous<pointEdgePoint>()
232 {
233  return true;
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 } // End namespace Foam
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 #include "pointEdgePointI.H"
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #endif
248 
249 // ************************************************************************* //
Foam::pointEdgePoint::pointEdgePoint
pointEdgePoint()
Construct null.
Definition: pointEdgePointI.H:121
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::contiguous< pointEdgePoint >
bool contiguous< pointEdgePoint >()
Data associated with pointEdgePoint type are contiguous.
Definition: pointEdgePoint.H:230
pointEdgePointI.H
Foam::pointEdgePoint::origin_
point origin_
Position of nearest wall center.
Definition: pointEdgePoint.H:64
Foam::pointEdgePoint::valid
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: pointEdgePointI.H:163
point.H
Foam::pointEdgePoint::equal
bool equal(const pointEdgePoint &, TrackingData &td) const
Same (like operator==)
Definition: pointEdgePointI.H:301
Foam::pointEdgePoint::operator!=
bool operator!=(const pointEdgePoint &) const
Definition: pointEdgePointI.H:319
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:285
Foam::pointEdgePoint::operator>>
friend Istream & operator>>(Istream &, pointEdgePoint &)
Foam::pointEdgePoint::operator==
bool operator==(const pointEdgePoint &) const
Definition: pointEdgePointI.H:312
tensor.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::label
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
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::pointEdgePoint::operator<<
friend Ostream & operator<<(Ostream &, const pointEdgePoint &)
Foam::pointEdgePoint::transform
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
Definition: pointEdgePointI.H:213
scalar.H
Foam::pointEdgePoint::update
bool update(const point &, const pointEdgePoint &w2, const scalar tol, TrackingData &td)
Evaluate distance to point. Update distSqr, origin from whomever.
Definition: pointEdgePointI.H:34
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
pTraits.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:200
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointEdgePoint::origin
const point & origin() const
Definition: pointEdgePointI.H:150
Foam::pointEdgePoint
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgePoint.H:59
Foam::pointEdgePoint::distSqr
scalar distSqr() const
Definition: pointEdgePointI.H:156
Foam::Vector< scalar >
label.H
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:226
Foam::pointEdgePoint::distSqr_
scalar distSqr_
Normal distance (squared) from point to origin.
Definition: pointEdgePoint.H:67
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::pointEdgePoint::sameGeometry
bool sameGeometry(const pointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Definition: pointEdgePointI.H:172
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:241
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190