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-------------------------------------------------------------------------------
11License
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
34template<class DataType>
36:
38{}
39
40
41template<class DataType>
43(
44 const point& origin,
45 const scalar distSqr,
46 const DataType& data
47)
48:
49 pointEdgePoint(origin, distSqr),
50 data_(data)
51{}
52
53
54// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55
56template<class DataType>
57template<class TrackingData>
59(
60 const tensor& rotTensor,
61 TrackingData& td
62)
63{
64 pointEdgePoint::transform(rotTensor, td);
65 data_ = Foam::transform(rotTensor, data_);
66}
67
68
69template<class DataType>
70template<class TrackingData>
72(
73 const polyMesh& mesh,
74 const label pointI,
75 const label edgeI,
76 const PointData<DataType>& edgeInfo,
77 const scalar tol,
78 TrackingData& td
79)
80{
81 if
82 (
84 (
85 mesh,
86 pointI,
87 edgeI,
88 edgeInfo,
89 tol,
90 td
91 )
92 )
93 {
94 data_ = edgeInfo.data_;
95
96 return true;
97 }
98
99 return false;
100}
101
102
103template<class DataType>
104template<class TrackingData>
106(
107 const polyMesh& mesh,
108 const label pointI,
109 const PointData<DataType>& newPointInfo,
110 const scalar tol,
111 TrackingData& td
112)
113{
114 if
115 (
117 (
118 mesh,
119 pointI,
120 newPointInfo,
121 tol,
122 td
123 )
124 )
125 {
126 data_ = newPointInfo.data_;
127
128 return true;
129 }
130
131 return false;
132}
133
134
135template<class DataType>
136template<class TrackingData>
138(
139 const PointData<DataType>& newPointInfo,
140 const scalar tol,
141 TrackingData& td
142)
143{
144 if (pointEdgePoint::updatePoint(newPointInfo, tol, td))
145 {
146 data_ = newPointInfo.data_;
147
148 return true;
149 }
150
151 return false;
152}
153
154
155template<class DataType>
156template<class TrackingData>
158(
159 const polyMesh& mesh,
160 const label edgeI,
161 const label pointI,
162 const PointData<DataType>& pointInfo,
163 const scalar tol,
164 TrackingData& td
165
166)
167{
168 if
169 (
171 (
172 mesh,
173 edgeI,
174 pointI,
175 pointInfo,
176 tol,
177 td
178 )
179 )
180 {
181 data_ = pointInfo.data_;
182
183 return true;
184 }
185
186 return false;
187}
188
189
190// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
191
192template<class DataType>
194(
195 const PointData<DataType>& rhs
196) const
197{
198 return pointEdgePoint::operator==(rhs) && (data() == rhs.data());
199}
200
201
202template<class DataType>
204(
205 const PointData<DataType>& rhs
206) const
207{
208 return !(*this == rhs);
209}
210
211
212// ************************************************************************* //
Variant of pointEdgePoint with some transported additional data. Templated on the transported data ty...
Definition: PointData.H:67
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:158
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:72
PointData()
Default construct.
Definition: PointDataI.H:35
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
Default transformation behaviour.
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
dynamicFvMesh & mesh
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
3D tensor transformation operations.