pointTopoDistanceData.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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::pointTopoDistanceData
29 
30 Description
31  For use with PointEdgeWave. Determines topological distance to
32  starting points. Templated on passive transported data.
33 
34 SourceFiles
35  pointTopoDistanceDataI.H
36  pointTopoDistanceData.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef pointTopoDistanceData_H
41 #define pointTopoDistanceData_H
42 
43 #include "point.H"
44 #include "tensor.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward Declarations
52 class polyPatch;
53 class polyMesh;
54 template<class Type> class pointTopoDistanceData;
55 template<class Type>
57 template<class Type>
58 Ostream& operator<<(Ostream&, const pointTopoDistanceData<Type>&);
59 
60 /*---------------------------------------------------------------------------*\
61  Class pointTopoDistanceData Declaration
62 \*---------------------------------------------------------------------------*/
63 
64 template<class Type>
66 {
67 
68  // Protected data
69 
70  //- Distance
71  label distance_;
72 
73  //- Starting data
74  Type data_;
75 
76 
77 public:
78 
79  typedef Type dataType;
80 
81 
82  // Constructors
83 
84  //- Construct null with invalid (-1) for distance, null constructor
85  // for data
86  inline pointTopoDistanceData();
87 
88  //- Construct from components
90  (
91  const label distance,
92  const Type& data
93  );
94 
95 
96  // Member Functions
97 
98  // Access
99 
100  inline label distance() const
101  {
102  return distance_;
103  }
105  inline const Type& data() const
106  {
107  return data_;
108  }
110 
111  // Needed by PointEdgeWave
112 
113  //- Changed or contains original (invalid) value
114  template<class TrackingData>
115  inline bool valid(TrackingData& td) const;
116 
117  //- Check for identical geometrical data (eg, cyclics checking)
118  template<class TrackingData>
119  inline bool sameGeometry
120  (
122  const scalar tol,
123  TrackingData& td
124  ) const;
125 
126  //- Convert origin to relative vector to leaving point
127  // (= point coordinate)
128  template<class TrackingData>
129  inline void leaveDomain
130  (
131  const polyPatch& patch,
132  const label patchPointi,
133  const point& pos,
134  TrackingData& td
135  );
136 
137  //- Convert relative origin to absolute by adding entering point
138  template<class TrackingData>
139  inline void enterDomain
140  (
141  const polyPatch& patch,
142  const label patchPointi,
143  const point& pos,
144  TrackingData& td
145  );
146 
147  //- Apply rotation matrix to origin
148  template<class TrackingData>
149  inline void transform
150  (
151  const tensor& rotTensor,
152  TrackingData& td
153  );
154 
155  //- Influence of edge on point
156  template<class TrackingData>
157  inline bool updatePoint
158  (
159  const polyMesh& mesh,
160  const label pointi,
161  const label edgeI,
162  const pointTopoDistanceData<Type>& edgeInfo,
163  const scalar tol,
164  TrackingData& td
165  );
166 
167  //- Influence of different value on same point.
168  // Merge new and old info.
169  template<class TrackingData>
170  inline bool updatePoint
171  (
172  const polyMesh& mesh,
173  const label pointi,
174  const pointTopoDistanceData<Type>& newPointInfo,
175  const scalar tol,
176  TrackingData& td
177  );
178 
179  //- Influence of different value on same point.
180  // No information about current position whatsoever.
181  template<class TrackingData>
182  inline bool updatePoint
183  (
184  const pointTopoDistanceData<Type>& newPointInfo,
185  const scalar tol,
186  TrackingData& td
187  );
188 
189  //- Influence of point on edge.
190  template<class TrackingData>
191  inline bool updateEdge
192  (
193  const polyMesh& mesh,
194  const label edgeI,
195  const label pointi,
196  const pointTopoDistanceData<Type>& pointInfo,
197  const scalar tol,
198  TrackingData& td
199  );
200 
201  //- Test for equality, with TrackingData
202  template<class TrackingData>
203  inline bool equal(const pointTopoDistanceData<Type>&, TrackingData&)
204  const;
205 
206 
207  // Member Operators
208 
209  // Needed for List IO
210  inline bool operator==(const pointTopoDistanceData<Type>&) const;
211  inline bool operator!=(const pointTopoDistanceData<Type>&) const;
212 
213 
214  // IOstream Operators
215 
216  friend Ostream& operator<< <Type>
217  (
218  Ostream&,
220  );
221  friend Istream& operator>> <Type>
222  (
223  Istream&,
225  );
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * Traits * * * * * * * * * * * * * * * * //
230 
231 //- Data are contiguous if data type is contiguous
232 template<class Type>
233 struct is_contiguous<pointTopoDistanceData<Type>> : is_contiguous<Type> {};
234 
235 //- Data are contiguous label if data type is label
236 template<class Type>
237 struct is_contiguous_label<pointTopoDistanceData<Type>> :
238  is_contiguous_label<Type> {};
239 
240 //- Data are contiguous scalar if data type is scalar
241 template<class Type>
242 struct is_contiguous_scalar<pointTopoDistanceData<Type>> :
243  is_contiguous_scalar<Type>{};
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 } // End namespace Foam
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #ifdef NoRepository
252  #include "pointTopoDistanceData.C"
253 #endif
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 #include "pointTopoDistanceDataI.H"
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 #endif
263 // ************************************************************************* //
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointTopoDistanceData< Type > &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool operator==(const pointTopoDistanceData< Type > &) const
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointTopoDistanceData< Type > &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
For use with PointEdgeWave. Determines topological distance to starting points. Templated on passive ...
bool operator!=(const pointTopoDistanceData< Type > &) const
dimensionedScalar pos(const dimensionedScalar &ds)
dynamicFvMesh & mesh
Istream & operator>>(Istream &, directionInfo &)
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
bool equal(const pointTopoDistanceData< Type > &, TrackingData &) const
Test for equality, with TrackingData.
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 sameGeometry(const pointTopoDistanceData< Type > &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
pointTopoDistanceData()
Construct null with invalid (-1) for distance, null constructor.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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>.
Namespace for OpenFOAM.