wallFaceI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallFace.H"
27 #include "polyMesh.H"
28 #include "transformer.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class Derived>
33 template<class TrackingData>
35 (
36  const point& pt,
37  const WallFaceBase<Derived>& w2,
38  const scalar tol,
39  TrackingData& td
40 )
41 {
42  const scalar dist2 =
43  sqr
44  (
45  face(identity(w2.points().size()))
46  .nearestPoint(pt, w2.points())
47  .distance()
48  );
49 
50  if (valid(td))
51  {
52  scalar diff = distSqr() - dist2;
53 
54  if (diff < 0)
55  {
56  // already nearer to pt
57  return false;
58  }
59 
60  if ((diff < small) || ((distSqr() > small) && (diff/distSqr() < tol)))
61  {
62  // don't propagate small changes
63  return false;
64  }
65  }
66 
67  // Either *this is not yet valid or w2 is closer
68  distSqr() = dist2;
69  points() = w2.points();
70 
71  return true;
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
77 template<class Derived>
79 :
80  points_(0),
81  distSqr_(-great)
82 {}
83 
84 
85 template<class Derived>
87 (
88  const face& f,
89  const pointField& points,
90  const scalar distSqr
91 )
92 :
93  points_(f.points(points)),
94  distSqr_(distSqr)
95 {}
96 
97 
98 template<class Derived>
100 (
101  const face& f,
102  const pointField& points,
103  const point& centre,
104  const scalar distSqr
105 )
106 :
107  points_(f.points(points)),
108  distSqr_(distSqr)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class Derived>
116 {
117  return points_;
118 }
119 
120 
121 template<class Derived>
123 {
124  return points_;
125 }
126 
127 
128 template<class Derived>
129 inline Foam::scalar Foam::WallFaceBase<Derived>::distSqr() const
130 {
131  return distSqr_;
132 }
133 
134 
135 template<class Derived>
137 {
138  return distSqr_;
139 }
140 
141 
142 template<class Derived>
143 template<class TrackingData>
144 inline Foam::scalar Foam::WallFaceBase<Derived>::dist(TrackingData& td) const
145 {
146  return valid(td) ? sqrt(distSqr_) : great;
147 }
148 
149 
150 template<class Derived>
151 template<class TrackingData>
152 inline bool Foam::WallFaceBase<Derived>::valid(TrackingData& td) const
153 {
154  return distSqr_ > -small;
155 }
156 
157 
158 template<class Derived>
159 template<class TrackingData>
161 (
162  const polyMesh&,
163  const WallFaceBase<Derived>& w2,
164  const scalar tol,
165  TrackingData& td
166 ) const
167 {
168  scalar diff = mag(distSqr() - w2.distSqr());
169 
170  if (diff < small)
171  {
172  return true;
173  }
174  else
175  {
176  if ((distSqr() > small) && ((diff/distSqr()) < tol))
177  {
178  return true;
179  }
180  else
181  {
182  return false;
183  }
184  }
185 }
186 
187 
188 template<class Derived>
189 template<class TrackingData>
191 (
192  const polyPatch& patch,
193  const label patchFacei,
194  const transformer& transform,
195  TrackingData& td
196 )
197 {
198  transform.transformPosition(points_, points_);
199 }
200 
201 
202 template<class Derived>
203 template<class TrackingData>
205 (
206  const polyMesh& mesh,
207  const label thisCelli,
208  const label neighbourFacei,
209  const WallFaceBase<Derived>& neighbourWallInfo,
210  const scalar tol,
211  TrackingData& td
212 )
213 {
214  return
215  static_cast<Derived&>(*this).update
216  (
217  mesh.cellCentres()[thisCelli],
218  neighbourWallInfo,
219  tol,
220  td
221  );
222 }
223 
224 
225 template<class Derived>
226 template<class TrackingData>
228 (
229  const polyMesh& mesh,
230  const label thisFacei,
231  const label neighbourCelli,
232  const WallFaceBase<Derived>& neighbourWallInfo,
233  const scalar tol,
234  TrackingData& td
235 )
236 {
237  return
238  static_cast<Derived&>(*this).update
239  (
240  mesh.faceCentres()[thisFacei],
241  neighbourWallInfo,
242  tol,
243  td
244  );
245 }
246 
247 
248 template<class Derived>
249 template<class TrackingData>
251 (
252  const polyMesh& mesh,
253  const label thisFacei,
254  const WallFaceBase<Derived>& neighbourWallInfo,
255  const scalar tol,
256  TrackingData& td
257 )
258 {
259  return
260  static_cast<Derived&>(*this).update
261  (
262  mesh.faceCentres()[thisFacei],
263  neighbourWallInfo,
264  tol,
265  td
266  );
267 }
268 
269 
270 template<class Derived>
271 template<class TrackingData>
273 (
274  const WallFaceBase<Derived>& rhs,
275  TrackingData& td
276 ) const
277 {
278  return operator==(rhs);
279 }
280 
281 
282 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
283 
284 template<class Derived>
285 inline bool Foam::WallFaceBase<Derived>::operator==
286 (
287  const Foam::WallFaceBase<Derived>& rhs
288 ) const
289 {
290  return points() == rhs.points();
291 }
292 
293 
294 template<class Derived>
295 inline bool Foam::WallFaceBase<Derived>::operator!=
296 (
297  const Foam::WallFaceBase<Derived>& rhs
298 ) const
299 {
300  return !(*this == rhs);
301 }
302 
303 
304 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
305 
306 template<class Derived>
307 Foam::Ostream& Foam::operator<<(Ostream& os, const WallFaceBase<Derived>& w)
308 {
309  return os << w.points() << token::SPACE << w.distSqr();
310 }
311 
312 
313 template<class Derived>
315 {
316  return is >> w.points() >> w.distSqr();
317 }
318 
319 
320 // ************************************************************************* //
const pointField & points() const
Definition: wallFaceI.H:115
WallFaceBase()
Construct null.
Definition: wallFaceI.H:78
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:403
bool equal(const WallFaceBase< Derived > &, TrackingData &td) const
Test equality.
Definition: wallFaceI.H:273
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
dimensionedSymmTensor sqr(const dimensionedVector &dv)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
dimensionedScalar sqrt(const dimensionedScalar &ds)
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
bool sameGeometry(const polyMesh &, const WallFaceBase< Derived > &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for checking.
Definition: wallFaceI.H:161
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Definition: wallFaceI.H:191
const pointField & points
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
scalar distSqr() const
Definition: wallFaceI.H:129
Istream & operator>>(Istream &, directionInfo &)
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:64
const vectorField & cellCentres() const
scalar dist(TrackingData &td) const
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool valid(TrackingData &td) const
Check whether the WallFaceBase has been changed at all or still.
Definition: wallFaceI.H:152
const vectorField & faceCentres() const
pointHit nearestPoint(const point &p, const pointField &) const
Return nearest point to face.
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const WallFaceBase< Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: wallFaceI.H:228
bool update(const point &, const FvWallInfoDataBase< WallInfo, Type, Derived > &w2, const scalar tol, TrackingData &td)
Evaluate distance to point. Update distSqr, origin from whomever.
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const WallFaceBase< Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: wallFaceI.H:205
dimensioned< scalar > mag(const dimensioned< Type > &)
vector transformPosition(const vector &v) const
Transform the given position.
Definition: transformerI.H:153
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
bool update(const point &pt, const WallFaceBase< Derived > &w2, const scalar tol, TrackingData &td)
...
Definition: wallFaceI.H:35
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66