pointTopoDistanceDataI.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) 2013-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMesh.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 // Null constructor
35 :
36  data_(-1),
37  distance_(-1)
38 {}
39 
40 
41 // Construct from components
43 (
44  const label data,
45  const label distance
46 )
47 :
48  data_(data),
49  distance_(distance)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 template<class TrackingData>
56 inline bool Foam::pointTopoDistanceData::valid(TrackingData& td) const
57 {
58  return distance_ != -1;
59 }
60 
61 
62 // No geometric data so never any problem on cyclics
63 template<class TrackingData>
65 (
66  const pointTopoDistanceData&,
67  const scalar tol,
68  TrackingData& td
69 ) const
70 {
71  return true;
72 }
73 
74 
75 // No geometric data.
76 template<class TrackingData>
78 (
79  const polyPatch& patch,
80  const label patchFacei,
81  const transformer& transform,
82  TrackingData& td
83 )
84 {}
85 
86 
87 // Update this with information from connected edge
88 template<class TrackingData>
90 (
91  const polyMesh& mesh,
92  const label pointi,
93  const label edgeI,
94  const pointTopoDistanceData& edgeInfo,
95  const scalar tol,
96  TrackingData& td
97 )
98 {
99  if (distance_ == -1)
100  {
101  data_ = edgeInfo.data_;
102  distance_ = edgeInfo.distance_ + 1;
103  return true;
104  }
105  else
106  {
107  return false;
108  }
109 }
110 
111 
112 // Update this with new information on same point
113 template<class TrackingData>
115 (
116  const polyMesh& mesh,
117  const label pointi,
118  const pointTopoDistanceData& newPointInfo,
119  const scalar tol,
120  TrackingData& td
121 )
122 {
123  if (distance_ == -1)
124  {
125  operator=(newPointInfo);
126  return true;
127  }
128  else
129  {
130  return false;
131  }
132 }
133 
134 
135 // Update this with new information on same point. No extra information.
136 template<class TrackingData>
138 (
139  const pointTopoDistanceData& newPointInfo,
140  const scalar tol,
141  TrackingData& td
142 )
143 {
144  if (distance_ == -1)
145  {
146  operator=(newPointInfo);
147  return true;
148  }
149  else
150  {
151  return false;
152  }
153 }
154 
155 
156 // Update this with information from connected point
157 template<class TrackingData>
159 (
160  const polyMesh& mesh,
161  const label edgeI,
162  const label pointi,
163  const pointTopoDistanceData& pointInfo,
164  const scalar tol,
165  TrackingData& td
166 )
167 {
168  if (distance_ == -1)
169  {
170  operator=(pointInfo);
171  return true;
172  }
173  else
174  {
175  return false;
176  }
177 }
178 
179 
180 template<class TrackingData>
182 (
183  const pointTopoDistanceData& rhs,
184  TrackingData& td
185 ) const
186 {
187  return operator==(rhs);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
192 
193 inline bool Foam::pointTopoDistanceData::operator==
194 (
195  const Foam::pointTopoDistanceData& rhs
196 ) const
197 {
198  return data() == rhs.data() && distance() == rhs.distance();
199 }
200 
201 
202 inline bool Foam::pointTopoDistanceData::operator!=
203 (
204  const Foam::pointTopoDistanceData& rhs
205 ) const
206 {
207  return !(*this == rhs);
208 }
209 
210 
211 // ************************************************************************* //
Database for solution and other reduced data.
Definition: data.H:54
data(const objectRegistry &obr)
Construct for objectRegistry.
Definition: data.C:38
For use with PointEdgeWave. Determines topological distance to starting points.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointTopoDistanceData &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointTopoDistanceData &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool equal(const pointTopoDistanceData &, TrackingData &) const
Same (like operator==)
bool sameGeometry(const pointTopoDistanceData &, const scalar tol, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483