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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 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 \*---------------------------------------------------------------------------*/
28 
29 #include "polyMesh.H"
30 #include "transform.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class DataType>
36 (
37  const point& origin,
38  const scalar distSqr,
39  const DataType& data
40 )
41 :
42  pointEdgePoint(origin, distSqr),
43  data_(data)
44 {}
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
49 template<class DataType>
50 template<class TrackingData>
52 (
53  const tensor& rotTensor,
54  TrackingData& td
55 )
56 {
57  pointEdgePoint::transform(rotTensor, td);
58  data_ = Foam::transform(rotTensor, data_);
59 }
60 
61 
62 template<class DataType>
63 template<class TrackingData>
65 (
66  const polyMesh& mesh,
67  const label pointI,
68  const label edgeI,
69  const PointData<DataType>& edgeInfo,
70  const scalar tol,
71  TrackingData& td
72 )
73 {
74  if
75  (
77  (
78  mesh,
79  pointI,
80  edgeI,
81  edgeInfo,
82  tol,
83  td
84  )
85  )
86  {
87  data_ = edgeInfo.data_;
88 
89  return true;
90  }
91 
92  return false;
93 }
94 
95 
96 template<class DataType>
97 template<class TrackingData>
99 (
100  const polyMesh& mesh,
101  const label pointI,
102  const PointData<DataType>& newPointInfo,
103  const scalar tol,
104  TrackingData& td
105 )
106 {
107  if
108  (
110  (
111  mesh,
112  pointI,
113  newPointInfo,
114  tol,
115  td
116  )
117  )
118  {
119  data_ = newPointInfo.data_;
120 
121  return true;
122  }
123 
124  return false;
125 }
126 
127 
128 template<class DataType>
129 template<class TrackingData>
131 (
132  const PointData<DataType>& newPointInfo,
133  const scalar tol,
134  TrackingData& td
135 )
136 {
137  if (pointEdgePoint::updatePoint(newPointInfo, tol, td))
138  {
139  data_ = newPointInfo.data_;
140 
141  return true;
142  }
143 
144  return false;
145 }
146 
147 
148 template<class DataType>
149 template<class TrackingData>
151 (
152  const polyMesh& mesh,
153  const label edgeI,
154  const label pointI,
155  const PointData<DataType>& pointInfo,
156  const scalar tol,
157  TrackingData& td
158 
159 )
160 {
161  if
162  (
164  (
165  mesh,
166  edgeI,
167  pointI,
168  pointInfo,
169  tol,
170  td
171  )
172  )
173  {
174  data_ = pointInfo.data_;
175 
176  return true;
177  }
178 
179  return false;
180 }
181 
182 
183 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
184 
185 template<class DataType>
187 (
188  const PointData<DataType>& rhs
189 ) const
190 {
191  return pointEdgePoint::operator==(rhs) && (data() == rhs.data());
192 }
193 
194 
195 template<class DataType>
197 (
198  const PointData<DataType>& rhs
199 ) const
200 {
201  return !(*this == rhs);
202 }
203 
204 
205 // ************************************************************************* //
Variant of pointEdgePoint with some transported additional data. Templated on the transported data ty...
bool updatePoint(const polyMesh &mesh, const label pointI, const label edgeI, const PointData< DataType > &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
Definition: PointDataI.H:58
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointI, const PointData< DataType > &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Definition: PointDataI.H:144
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.
dynamicFvMesh & mesh
3D tensor transformation operations.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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 &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Tensor of scalars, i.e. Tensor<scalar>.
bool operator==(const pointEdgePoint &) const
Test for equality.
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to the data.
Definition: PointDataI.H:45
PointData()=default
Default construct.