pointEdgeStructuredWalkI.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-2016 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// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33
34template<class TrackingData>
36(
37 const pointEdgeStructuredWalk& w2,
38 const scalar tol,
39 TrackingData& td
40)
41{
42 if (!valid(td))
43 {
44 // current not yet set. Walked from w2 to here (=point0)
45 dist_ = w2.dist_ + mag(point0_-w2.previousPoint_);
46 previousPoint_ = point0_;
47 data_ = w2.data_;
48 index_ = w2.index_;
49
50 return true;
51 }
52
53 return false;
54}
55
56
57// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
60:
61 point0_(point::max),
62 previousPoint_(point::max),
63 dist_(0),
64 data_(Zero),
65 index_(-1)
66{}
67
68
70(
71 const point& point0,
72 const point& previousPoint,
73 const scalar dist,
74 const vector& data,
75 const label index
76)
77:
78 point0_(point0),
79 previousPoint_(previousPoint),
80 dist_(dist),
81 data_(data),
82 index_(index)
83{}
84
85
86// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87
89{
90 return point0_ != point::max;
91}
92
93
94template<class TrackingData>
95inline bool Foam::pointEdgeStructuredWalk::valid(TrackingData& td) const
96{
97 return previousPoint_ != point::max;
98}
99
100
101// Checks for cyclic points
102template<class TrackingData>
104(
106 const scalar tol,
107 TrackingData& td
108) const
109{
110 const scalar diff = Foam::mag(dist() - w2.dist());
111
112 if (diff < SMALL)
113 {
114 return true;
115 }
116 else
117 {
118 if ((dist() > SMALL) && ((diff/dist()) < tol))
119 {
120 return true;
121 }
122 else
123 {
124 return false;
125 }
126 }
127}
128
129
130template<class TrackingData>
132(
133 const polyPatch& patch,
134 const label patchPointi,
135 const point& coord,
136 TrackingData& td
137)
138{
139 previousPoint_ -= coord;
140}
141
142
143template<class TrackingData>
145(
146 const tensor& rotTensor,
147 TrackingData& td
148)
149{
150 previousPoint_ = Foam::transform(rotTensor, previousPoint_);
151}
152
153
154// Update absolute geometric quantities. Note that distance (dist_)
155// is not affected by leaving/entering domain.
156template<class TrackingData>
158(
159 const polyPatch& patch,
160 const label patchPointi,
161 const point& coord,
162 TrackingData& td
163)
164{
165 // back to absolute form
166 previousPoint_ += coord;
167}
168
169
170// Update this with information from connected edge
171template<class TrackingData>
173(
174 const polyMesh& mesh,
175 const label pointi,
176 const label edgeI,
177 const pointEdgeStructuredWalk& edgeInfo,
178 const scalar tol,
179 TrackingData& td
180)
181{
182 if (inZone())
183 {
184 return update(edgeInfo, tol, td);
185 }
186
187 return false;
188}
189
190
191// Update this with new information on same point
192template<class TrackingData>
194(
195 const polyMesh& mesh,
196 const label pointi,
197 const pointEdgeStructuredWalk& newPointInfo,
198 const scalar tol,
199 TrackingData& td
200)
201{
202 if (inZone())
203 {
204 return update(newPointInfo, tol, td);
205 }
206
207 return false;
208}
209
210
211// Update this with new information on same point. No extra information.
212template<class TrackingData>
214(
215 const pointEdgeStructuredWalk& newPointInfo,
216 const scalar tol,
217 TrackingData& td
218)
219{
220 return update(newPointInfo, tol, td);
221}
222
223
224// Update this with information from connected point
225template<class TrackingData>
227(
228 const polyMesh& mesh,
229 const label edgeI,
230 const label pointi,
231 const pointEdgeStructuredWalk& pointInfo,
232 const scalar tol,
233 TrackingData& td
234)
235{
236 if (inZone())
237 {
238 return update(pointInfo, tol, td);
239 }
240
241 return false;
242}
243
244
245template<class TrackingData>
247(
248 const pointEdgeStructuredWalk& rhs,
249 TrackingData& td
250) const
251{
252 return operator==(rhs);
253}
254
255
256// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
257
259(
260 const pointEdgeStructuredWalk& rhs
261) const
262{
263 return previousPoint_ == rhs.previousPoint_;
264}
265
266
268(
269 const pointEdgeStructuredWalk& rhs
270) const
271{
272 return !(*this == rhs);
273}
274
275
276// ************************************************************************* //
#define w2
Definition: blockCreate.C:35
bool valid() const
True if all internal ids are non-negative.
Database for solution data, solver performance and other reduced data.
Definition: data.H:58
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
virtual bool update()
Update the mesh for both mesh motion and topology change.
Default transformation behaviour.
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
Determines length of string of edges walked to point.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgeStructuredWalk &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool sameGeometry(const pointEdgeStructuredWalk &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
bool equal(const pointEdgeStructuredWalk &, TrackingData &) const
Test for equality, with TrackingData.
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
bool inZone() const
True if starting point is valid (ie, not point::max)
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgeStructuredWalk &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
mesh update()
dynamicFvMesh & mesh
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
3D tensor transformation operations.