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  Copyright (C) 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 // Update this with w2 if w2 nearer to pt.
35 template<class TrackingData>
36 inline bool Foam::patchEdgeFaceInfo::update
37 (
38  const point& pt,
39  const patchEdgeFaceInfo& w2,
40  const scalar tol,
41  TrackingData& td
42 )
43 {
44  const scalar dist2 = magSqr(pt - w2.origin());
45 
46  if (!valid(td))
47  {
48  // current not yet set so use any value
49  distSqr_ = dist2;
50  origin_ = w2.origin();
51 
52  return true;
53  }
54 
55  const scalar diff = distSqr_ - dist2;
56 
57  if (diff < 0)
58  {
59  // already nearer to pt
60  return false;
61  }
62 
63  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
64  {
65  // don't propagate small changes
66  return false;
67  }
68  else
69  {
70  // update with new values
71  distSqr_ = dist2;
72  origin_ = w2.origin();
73 
74  return true;
75  }
76 }
77 
78 
79 // Update this with w2 (information on same edge)
80 template<class TrackingData>
81 inline bool Foam::patchEdgeFaceInfo::update
82 (
83  const patchEdgeFaceInfo& w2,
84  const scalar tol,
85  TrackingData& td
86 )
87 {
88  if (!valid(td))
89  {
90  // current not yet set so use any value
91  distSqr_ = w2.distSqr();
92  origin_ = w2.origin();
93 
94  return true;
95  }
96 
97  const scalar diff = distSqr_ - w2.distSqr();
98 
99  if (diff < 0)
100  {
101  // already nearer to pt
102  return false;
103  }
104 
105  if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
106  {
107  // don't propagate small changes
108  return false;
109  }
110  else
111  {
112  // update with new values
113  distSqr_ = w2.distSqr();
114  origin_ = w2.origin();
115 
116  return true;
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
124 :
125  patchEdgeFaceInfo(point::max, sqr(GREAT))
126 {}
127 
128 
130 (
131  const point& origin,
132  const scalar distSqr
133 )
134 :
135  origin_(origin),
136  distSqr_(distSqr)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
142 template<class TrackingData>
143 inline bool Foam::patchEdgeFaceInfo::valid(TrackingData& td) const
144 {
145  return origin_ != point::max;
146 }
147 
148 
149 template<class TrackingData>
151 (
152  const polyMesh& mesh,
153  const primitivePatch& patch,
154  const tensor& rotTensor,
155  const scalar tol,
156  TrackingData& td
157 )
158 {
159  origin_ = Foam::transform(rotTensor, origin_);
160 }
161 
162 
163 template<class TrackingData>
165 (
166  const polyMesh& mesh,
167  const primitivePatch& patch,
168  const label edgeI,
169  const label facei,
170  const patchEdgeFaceInfo& faceInfo,
171  const scalar tol,
172  TrackingData& td
173 )
174 {
175  const edge& e = patch.edges()[edgeI];
176  const point eMid
177  (
178  0.5
179  * (
180  patch.points()[patch.meshPoints()[e[0]]]
181  + patch.points()[patch.meshPoints()[e[1]]]
182  )
183  );
184  return update(eMid, faceInfo, tol, td);
185 }
186 
187 
188 template<class TrackingData>
190 (
191  const polyMesh& mesh,
192  const primitivePatch& patch,
193  const patchEdgeFaceInfo& edgeInfo,
194  const bool sameOrientation,
195  const scalar tol,
196  TrackingData& td
197 )
198 {
199  return update(edgeInfo, tol, td);
200 }
201 
202 
203 template<class TrackingData>
205 (
206  const polyMesh& mesh,
207  const primitivePatch& patch,
208  const label facei,
209  const label edgeI,
210  const patchEdgeFaceInfo& edgeInfo,
211  const scalar tol,
212  TrackingData& td
213 )
214 {
215  return update(patch.faceCentres()[facei], edgeInfo, tol, td);
216 }
217 
218 
219 template<class TrackingData>
221 (
222  const patchEdgeFaceInfo& rhs,
223  TrackingData& td
224 ) const
225 {
226  return operator==(rhs);
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
231 
232 inline bool Foam::patchEdgeFaceInfo::operator==
233 (
234  const patchEdgeFaceInfo& rhs
235 ) const
236 {
237  return origin_ == rhs.origin_;
238 }
239 
240 
241 inline bool Foam::patchEdgeFaceInfo::operator!=
242 (
243  const Foam::patchEdgeFaceInfo& rhs
244 ) const
245 {
246  return !(*this == rhs);
247 }
248 
249 
250 // ************************************************************************* //
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()
Default construct.
Definition: patchEdgeFaceInfoI.H:123
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
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
Changed or contains original (invalid) value.
Definition: patchEdgeFaceInfoI.H:143
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:165
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:151
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< Cmpt >, Cmpt, 3 >::max
static const Vector< Cmpt > 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:205
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
Test for equality, with TrackingData.
Definition: patchEdgeFaceInfoI.H:221
Foam::patchEdgeFaceInfo
Definition: patchEdgeFaceInfo.H:64
transform.H
3D tensor transformation operations.
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79