patchEdgeFaceInfoI.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "polyMesh.H"
29 #include "transform.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 // Update this with w2 if w2 nearer to pt.
34 template<class TrackingData>
35 inline bool Foam::patchEdgeFaceInfo::update
36 (
37  const point& pt,
38  const patchEdgeFaceInfo& w2,
39  const scalar tol,
40  TrackingData& td
41 )
42 {
43  scalar dist2 = magSqr(pt - w2.origin());
44 
45  if (!valid(td))
46  {
47  // current not yet set so use any value
48  distSqr_ = dist2;
49  origin_ = w2.origin();
50 
51  return true;
52  }
53 
54  scalar diff = distSqr_ - dist2;
55 
56  if (diff < 0)
57  {
58  // already nearer to pt
59  return false;
60  }
61 
62  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
63  {
64  // don't propagate small changes
65  return false;
66  }
67  else
68  {
69  // update with new values
70  distSqr_ = dist2;
71  origin_ = w2.origin();
72 
73  return true;
74  }
75 }
76 
77 
78 // Update this with w2 (information on same edge)
79 template<class TrackingData>
80 inline bool Foam::patchEdgeFaceInfo::update
81 (
82  const patchEdgeFaceInfo& w2,
83  const scalar tol,
84  TrackingData& td
85 )
86 {
87  if (!valid(td))
88  {
89  // current not yet set so use any value
90  distSqr_ = w2.distSqr();
91  origin_ = w2.origin();
92 
93  return true;
94  }
95 
96  scalar diff = distSqr_ - w2.distSqr();
97 
98  if (diff < 0)
99  {
100  // already nearer to pt
101  return false;
102  }
103 
104  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
105  {
106  // don't propagate small changes
107  return false;
108  }
109  else
110  {
111  // update with new values
112  distSqr_ = w2.distSqr();
113  origin_ = w2.origin();
114 
115  return true;
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 // Null constructor
124 :
125  origin_(point::max),
126  distSqr_(sqr(GREAT))
127 {}
128 
129 
130 // Construct from origin, distance
132 (
133  const point& origin,
134  const scalar distSqr
135 )
136 :
137  origin_(origin),
138  distSqr_(distSqr)
139 {}
140 
141 
142 // Construct as copy
144 :
145  origin_(wpt.origin()),
146  distSqr_(wpt.distSqr())
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
154  return origin_;
155 }
156 
157 
158 inline Foam::scalar Foam::patchEdgeFaceInfo::distSqr() const
159 {
160  return distSqr_;
161 }
162 
163 
164 template<class TrackingData>
165 inline bool Foam::patchEdgeFaceInfo::valid(TrackingData& td) const
166 {
167  return origin_ != point::max;
168 }
169 
170 
171 template<class TrackingData>
173 (
174  const polyMesh& mesh,
175  const primitivePatch& patch,
176  const tensor& rotTensor,
177  const scalar tol,
178  TrackingData& td
179 )
180 {
181  origin_ = Foam::transform(rotTensor, origin_);
182 }
183 
184 
185 template<class TrackingData>
187 (
188  const polyMesh& mesh,
189  const primitivePatch& patch,
190  const label edgeI,
191  const label facei,
192  const patchEdgeFaceInfo& faceInfo,
193  const scalar tol,
194  TrackingData& td
195 )
196 {
197  const edge& e = patch.edges()[edgeI];
198  point eMid =
199  0.5
200  * (
201  patch.points()[patch.meshPoints()[e[0]]]
202  + patch.points()[patch.meshPoints()[e[1]]]
203  );
204  return update(eMid, faceInfo, tol, td);
205 }
206 
207 
208 template<class TrackingData>
210 (
211  const polyMesh& mesh,
212  const primitivePatch& patch,
213  const patchEdgeFaceInfo& edgeInfo,
214  const bool sameOrientation,
215  const scalar tol,
216  TrackingData& td
217 )
218 {
219  return update(edgeInfo, tol, td);
220 }
221 
222 
223 template<class TrackingData>
225 (
226  const polyMesh& mesh,
227  const primitivePatch& patch,
228  const label facei,
229  const label edgeI,
230  const patchEdgeFaceInfo& edgeInfo,
231  const scalar tol,
232  TrackingData& td
233 )
234 {
235  return update(patch.faceCentres()[facei], edgeInfo, tol, td);
236 }
237 
238 
239 template<class TrackingData>
241 (
242  const patchEdgeFaceInfo& rhs,
243  TrackingData& td
244 ) const
245 {
246  return operator==(rhs);
247 }
248 
249 
250 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
251 
252 inline bool Foam::patchEdgeFaceInfo::operator==
253 (
254  const Foam::patchEdgeFaceInfo& rhs
255 ) const
256 {
257  return origin() == rhs.origin();
258 }
259 
260 
261 inline bool Foam::patchEdgeFaceInfo::operator!=
262 (
263  const Foam::patchEdgeFaceInfo& rhs
264 ) const
265 {
266  return !(*this == rhs);
267 }
268 
269 
270 // ************************************************************************* //
Foam::Tensor< scalar >
update
mesh update()
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::patchEdgeFaceInfo::patchEdgeFaceInfo
patchEdgeFaceInfo()
Construct null.
Definition: patchEdgeFaceInfoI.H:123
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
Foam::patchEdgeFaceInfo::valid
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: patchEdgeFaceInfoI.H:165
Foam::label
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:62
Foam::patchEdgeFaceInfo::updateEdge
bool updateEdge(const polyMesh &mesh, const primitivePatch &patch, const label edgeI, const label facei, const patchEdgeFaceInfo &faceInfo, const scalar tol, TrackingData &td)
Influence of face on edge.
Definition: patchEdgeFaceInfoI.H:187
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::patchEdgeFaceInfo::transform
void transform(const polyMesh &mesh, const primitivePatch &patch, const tensor &rotTensor, const scalar tol, TrackingData &td)
Apply rotation matrix.
Definition: patchEdgeFaceInfoI.H:173
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
w2
#define w2
Definition: blockCreate.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::max
static const Vector< scalar > max
Definition: VectorSpace.H:117
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::patchEdgeFaceInfo::updateFace
bool updateFace(const polyMesh &mesh, const primitivePatch &patch, const label facei, const label edgeI, const patchEdgeFaceInfo &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on face.
Definition: patchEdgeFaceInfoI.H:225
Foam::Vector< scalar >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::patchEdgeFaceInfo::equal
bool equal(const patchEdgeFaceInfo &, TrackingData &td) const
Same (like operator==)
Definition: patchEdgeFaceInfoI.H:241
Foam::patchEdgeFaceInfo
Definition: patchEdgeFaceInfo.H:65
transform.H
3D tensor transformation operations.
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::patchEdgeFaceInfo::distSqr
scalar distSqr() const
Definition: patchEdgeFaceInfoI.H:158
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::patchEdgeFaceInfo::origin
const point & origin() const
Definition: patchEdgeFaceInfoI.H:152