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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::pointEdgePoint
29
30Description
31 Holds information regarding nearest wall point. Used in PointEdgeWave.
32 (so not standard FaceCellWave)
33 To be used in wall distance calculation.
34
35SourceFiles
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
52namespace Foam
53{
54
55// Forward Declarations
56class polyPatch;
57class polyMesh;
58class pointEdgePoint;
59
60Istream& operator>>(Istream&, pointEdgePoint&);
61Ostream& operator<<(Ostream&, const pointEdgePoint&);
62
63/*---------------------------------------------------------------------------*\
64 Class pointEdgePoint Declaration
65\*---------------------------------------------------------------------------*/
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
103public:
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
118 const point& origin() const
119 {
120 return origin_;
122 point& origin()
123 {
124 return origin_;
125 }
127 scalar distSqr() const
128 {
129 return distSqr_;
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
245};
246
247
248// * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
249
250//- Contiguous data for pointEdgePoint
251template<> struct is_contiguous<pointEdgePoint> : std::true_type {};
252
253//- Contiguous scalar data for pointEdgePoint
254template<> 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// ************************************************************************* //
#define w2
Definition: blockCreate.C:35
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
const point & origin() const
bool sameGeometry(const pointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
bool equal(const pointEdgePoint &, TrackingData &td) const
Test for equality, with TrackingData.
pointEdgePoint()
Default construct.
bool operator==(const pointEdgePoint &) const
Test for equality.
bool operator!=(const pointEdgePoint &) const
Test for inequality.
friend Istream & operator>>(Istream &, pointEdgePoint &)
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
scalar distSqr() const
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
friend Ostream & operator<<(Ostream &, const pointEdgePoint &)
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Tensor of scalars, i.e. Tensor<scalar>.
mesh update()
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
vector point
Point is a vector.
Definition: point.H:43
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Istream & operator>>(Istream &, directionInfo &)
A template class to specify if a data type is composed solely of Foam::scalar elements.
Definition: contiguous.H:94
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:78