wallPointI.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 "wallPoint.H"
27 #include "polyMesh.H"
28 #include "transformer.H"
29 #include "SubField.H"
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
33 template<class Derived>
34 template<class TrackingData>
36 (
37  const point& pt,
38  const WallPointBase<Derived>& w2,
39  const scalar tol,
40  TrackingData& td
41 )
42 {
43  const scalar dist2 = magSqr(pt - w2.origin());
44 
45  if (valid(td))
46  {
47  const scalar diff = distSqr() - dist2;
48 
49  if (diff < 0)
50  {
51  // already nearer to pt
52  return false;
53  }
54 
55  if ((diff < small) || ((distSqr() > small) && (diff/distSqr() < tol)))
56  {
57  // don't propagate small changes
58  return false;
59  }
60  }
61 
62  // Either *this is not yet valid or w2 is closer
63  distSqr() = dist2;
64  origin() = w2.origin();
65 
66  return true;
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 template<class Derived>
74 :
75  origin_(point::max),
76  distSqr_(-great)
77 {}
78 
79 
80 template<class Derived>
82 (
83  const point& origin,
84  const scalar distSqr
85 )
86 :
87  origin_(origin),
88  distSqr_(distSqr)
89 {}
90 
91 
92 template<class Derived>
94 (
95  const face& f,
96  const pointField& points,
97  const point& centre,
98  const scalar distSqr
99 )
100 :
101  origin_(centre),
102  distSqr_(distSqr)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
108 template<class Derived>
110 {
111  return origin_;
112 }
113 
114 
115 template<class Derived>
117 {
118  return origin_;
119 }
120 
121 
122 template<class Derived>
123 inline Foam::scalar Foam::WallPointBase<Derived>::distSqr() const
124 {
125  return distSqr_;
126 }
127 
128 
129 template<class Derived>
131 {
132  return distSqr_;
133 }
134 
135 
136 template<class Derived>
137 template<class TrackingData>
138 inline Foam::scalar Foam::WallPointBase<Derived>::dist(TrackingData& td) const
139 {
140  return valid(td) ? sqrt(distSqr_) : great;
141 }
142 
143 
144 template<class Derived>
145 template<class TrackingData>
146 inline bool Foam::WallPointBase<Derived>::valid(TrackingData& td) const
147 {
148  return distSqr_ > -small;
149 }
150 
151 
152 template<class Derived>
153 template<class TrackingData>
155 (
156  const polyMesh&,
157  const WallPointBase<Derived>& w2,
158  const scalar tol,
159  TrackingData& td
160 ) const
161 {
162  scalar diff = mag(distSqr() - w2.distSqr());
163 
164  if (diff < small)
165  {
166  return true;
167  }
168  else
169  {
170  if ((distSqr() > small) && ((diff/distSqr()) < tol))
171  {
172  return true;
173  }
174  else
175  {
176  return false;
177  }
178  }
179 }
180 
181 
182 template<class Derived>
183 template<class TrackingData>
185 (
186  const polyPatch& patch,
187  const label patchFacei,
188  const transformer& transform,
189  TrackingData& td
190 )
191 {
192  // Note that distSqr_ is not affected by crossing an interface
193  origin_ = transform.transformPosition(origin_);
194 }
195 
196 
197 template<class Derived>
198 template<class TrackingData>
200 (
201  const polyMesh& mesh,
202  const label thisCelli,
203  const label neighbourFacei,
204  const WallPointBase<Derived>& neighbourWallInfo,
205  const scalar tol,
206  TrackingData& td
207 )
208 {
209  return
210  static_cast<Derived&>(*this).update
211  (
212  mesh.cellCentres()[thisCelli],
213  neighbourWallInfo,
214  tol,
215  td
216  );
217 }
218 
219 
220 template<class Derived>
221 template<class TrackingData>
223 (
224  const polyMesh& mesh,
225  const label thisFacei,
226  const label neighbourCelli,
227  const WallPointBase<Derived>& neighbourWallInfo,
228  const scalar tol,
229  TrackingData& td
230 )
231 {
232  return
233  static_cast<Derived&>(*this).update
234  (
235  mesh.faceCentres()[thisFacei],
236  neighbourWallInfo,
237  tol,
238  td
239  );
240 }
241 
242 
243 template<class Derived>
244 template<class TrackingData>
246 (
247  const polyMesh& mesh,
248  const label thisFacei,
249  const WallPointBase<Derived>& neighbourWallInfo,
250  const scalar tol,
251  TrackingData& td
252 )
253 {
254  return
255  static_cast<Derived&>(*this).update
256  (
257  mesh.faceCentres()[thisFacei],
258  neighbourWallInfo,
259  tol,
260  td
261  );
262 }
263 
264 
265 template<class Derived>
266 template<class TrackingData>
268 (
269  const WallPointBase<Derived>& rhs,
270  TrackingData& td
271 ) const
272 {
273  return operator==(rhs);
274 }
275 
276 
277 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
278 
279 template<class Derived>
280 inline bool Foam::WallPointBase<Derived>::operator==
281 (
283 ) const
284 {
285  return origin() == rhs.origin();
286 }
287 
288 
289 template<class Derived>
290 inline bool Foam::WallPointBase<Derived>::operator!=
291 (
293 ) const
294 {
295  return !(*this == rhs);
296 }
297 
298 
299 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
300 
301 template<class Derived>
302 Foam::Ostream& Foam::operator<<(Ostream& os, const WallPointBase<Derived>& w)
303 {
304  return os << w.origin() << token::SPACE << w.distSqr();
305 }
306 
307 
308 template<class Derived>
310 {
311  return is >> w.origin() >> w.distSqr();
312 }
313 
314 
315 // ************************************************************************* //
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:403
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
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)
scalar distSqr() const
Definition: wallPointI.H:123
scalar dist(TrackingData &td) const
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
bool equal(const WallPointBase< Derived > &, TrackingData &td) const
Same (like operator==)
Definition: wallPointI.H:268
Istream & operator>>(Istream &, directionInfo &)
bool sameGeometry(const polyMesh &, const WallPointBase< Derived > &, const scalar, TrackingData &td) const
Check for identical geometrical data. Used for cyclics checking.
Definition: wallPointI.H:155
const point & origin() const
Definition: wallPointI.H:109
bool valid(TrackingData &td) const
Check whether origin has been changed at all or.
Definition: wallPointI.H:146
const vectorField & cellCentres() const
WallPointBase()
Construct null.
Definition: wallPointI.H:73
void transform(const polyPatch &patch, const label patchFacei, const transformer &transform, TrackingData &td)
Transform across an interface.
Definition: wallPointI.H:185
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
bool update(const point &, const WallPointBase< Derived > &w2, const scalar tol, TrackingData &td)
Evaluate distance to point. Update distSqr, origin from whomever.
Definition: wallPointI.H:36
const vectorField & faceCentres() const
bool update(const point &, const FvWallInfoDataBase< WallInfo, Type, Derived > &w2, const scalar tol, TrackingData &td)
Evaluate distance to point. Update distSqr, origin from whomever.
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const WallPointBase< Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: wallPointI.H:223
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const WallPointBase< Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: wallPointI.H:200