pointDataI.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) 2011-2019 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 #include "transform.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
34  s_(great),
35  v_(point::max)
36 {}
37 
38 
40 (
41  const point& origin,
42  const scalar distSqr,
43  const scalar s,
44  const vector& v
45 )
46 :
47  pointEdgePoint(origin, distSqr),
48  s_(s),
49  v_(v)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 inline Foam::scalar Foam::pointData::s() const
56 {
57  return s_;
58 }
59 
60 
61 inline const Foam::vector& Foam::pointData::v() const
62 {
63  return v_;
64 }
65 
66 
67 template<class TrackingData>
69 (
70  const tensor& rotTensor,
71  TrackingData& td
72 )
73 {
74  pointEdgePoint::transform(rotTensor, td);
75  v_ = Foam::transform(rotTensor, v_);
76 }
77 
78 
79 // Update this with information from connected edge
80 template<class TrackingData>
82 (
83  const polyMesh& mesh,
84  const label pointi,
85  const label edgeI,
86  const pointData& edgeInfo,
87  const scalar tol,
88  TrackingData& td
89 )
90 {
91  if
92  (
94  (
95  mesh,
96  pointi,
97  edgeI,
98  edgeInfo,
99  tol,
100  td
101  )
102  )
103  {
104  s_ = edgeInfo.s_;
105  v_ = edgeInfo.v_;
106  return true;
107  }
108  else
109  {
110  return false;
111  }
112 }
113 
114 // Update this with new information on same point
115 template<class TrackingData>
117 (
118  const polyMesh& mesh,
119  const label pointi,
120  const pointData& newPointInfo,
121  const scalar tol,
122  TrackingData& td
123 )
124 {
125  if
126  (
128  (
129  mesh,
130  pointi,
131  newPointInfo,
132  tol,
133  td
134  )
135  )
136  {
137  s_ = newPointInfo.s_;
138  v_ = newPointInfo.v_;
139  return true;
140  }
141  else
142  {
143  return false;
144  }
145 }
146 
147 
148 // Update this with new information on same point. No extra information.
149 template<class TrackingData>
151 (
152  const pointData& newPointInfo,
153  const scalar tol,
154  TrackingData& td
155 )
156 {
157  if (pointEdgePoint::updatePoint(newPointInfo, tol, td))
158  {
159  s_ = newPointInfo.s_;
160  v_ = newPointInfo.v_;
161  return true;
162  }
163  else
164  {
165  return false;
166  }
167 }
168 
169 
170 // Update this with information from connected point
171 template<class TrackingData>
172 inline bool Foam::pointData::updateEdge
173 (
174  const polyMesh& mesh,
175  const label edgeI,
176  const label pointi,
177  const pointData& pointInfo,
178  const scalar tol,
179  TrackingData& td
180 
181 )
182 {
183  if
184  (
186  (
187  mesh,
188  edgeI,
189  pointi,
190  pointInfo,
191  tol,
192  td
193  )
194  )
195  {
196  s_ = pointInfo.s_;
197  v_ = pointInfo.v_;
198  return true;
199  }
200  else
201  {
202  return false;
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
208 
210 const
211 {
212  return
214  && (s() == rhs.s())
215  && (v() == rhs.v());
216 }
217 
218 
220 const
221 {
222  return !(*this == rhs);
223 }
224 
225 
226 // ************************************************************************* //
bool operator!=(const pointData &) const
Definition: pointDataI.H:219
Variant of pointEdgePoint with some transported additional data. WIP - should be templated on data li...
Definition: pointData.H:61
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
Definition: pointDataI.H:69
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointData &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
Definition: pointDataI.H:82
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
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointData &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Definition: pointDataI.H:173
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
pointEdgePoint()
Construct null.
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
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.
pointData()
Construct null.
Definition: pointDataI.H:31
const vector & v() const
Definition: pointDataI.H:61
3D tensor transformation operations.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool operator==(const pointEdgePoint &) const
bool operator==(const pointData &) const
Definition: pointDataI.H:209
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
scalar s() const
Definition: pointDataI.H:55