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