ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / meshTools / PointEdgeWave / pointEdgePointI.H
blob370db9fb85402732dce2bda4693de7638fe49c6e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "polyMesh.H"
27 #include "transform.H"
29 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
31 // Update this with w2 if w2 nearer to pt.
32 template<class TrackingData>
33 inline bool Foam::pointEdgePoint::update
35     const point& pt,
36     const pointEdgePoint& w2,
37     const scalar tol,
38     TrackingData& td
41     scalar dist2 = magSqr(pt - w2.origin());
43     if (!valid(td))
44     {
45         // current not yet set so use any value
46         distSqr_ = dist2;
47         origin_ = w2.origin();
49         return true;
50     }
52     scalar diff = distSqr_ - dist2;
54     if (diff < 0)
55     {
56         // already nearer to pt
57         return false;
58     }
60     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
61     {
62         // don't propagate small changes
63         return false;
64     }
65     else
66     {
67         // update with new values
68         distSqr_ = dist2;
69         origin_ = w2.origin();
71         return true;
72     }
76 // Update this with w2 (information on same point)
77 template<class TrackingData>
78 inline bool Foam::pointEdgePoint::update
80     const pointEdgePoint& w2,
81     const scalar tol,
82     TrackingData& td
85     if (!valid(td))
86     {
87         // current not yet set so use any value
88         distSqr_ = w2.distSqr();
89         origin_ = w2.origin();
91         return true;
92     }
94     scalar diff = distSqr_ - w2.distSqr();
96     if (diff < 0)
97     {
98         // already nearer to pt
99         return false;
100     }
102     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
103     {
104         // don't propagate small changes
105         return false;
106     }
107     else
108     {
109         // update with new values
110         distSqr_ =  w2.distSqr();
111         origin_ = w2.origin();
113         return true;
114     }
118 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
120 // Null constructor
121 inline Foam::pointEdgePoint::pointEdgePoint()
123     origin_(point::max),
124     distSqr_(GREAT)
128 // Construct from origin, distance
129 inline Foam::pointEdgePoint::pointEdgePoint
131     const point& origin,
132     const scalar distSqr
135     origin_(origin),
136     distSqr_(distSqr)
140 // Construct as copy
141 inline Foam::pointEdgePoint::pointEdgePoint(const pointEdgePoint& wpt)
143     origin_(wpt.origin()),
144     distSqr_(wpt.distSqr())
148 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
150 inline const Foam::point& Foam::pointEdgePoint::origin() const
152     return origin_;
156 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
158     return distSqr_;
162 template<class TrackingData>
163 inline bool Foam::pointEdgePoint::valid(TrackingData& td) const
165     return origin_ != point::max;
169 // Checks for cyclic points
170 template<class TrackingData>
171 inline bool Foam::pointEdgePoint::sameGeometry
173     const pointEdgePoint& w2,
174     const scalar tol,
175     TrackingData& td
176 ) const
178     scalar diff = Foam::mag(distSqr() - w2.distSqr());
180     if (diff < SMALL)
181     {
182         return true;
183     }
184     else
185     {
186         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
187         {
188             return true;
189         }
190         else
191         {
192             return false;
193         }
194     }
198 template<class TrackingData>
199 inline void Foam::pointEdgePoint::leaveDomain
201     const polyPatch& patch,
202     const label patchPointI,
203     const point& coord,
204     TrackingData& td
207     origin_ -= coord;
211 template<class TrackingData>
212 inline void Foam::pointEdgePoint::transform
214     const tensor& rotTensor,
215     TrackingData& td
218     origin_ = Foam::transform(rotTensor, origin_);
222 // Update absolute geometric quantities. Note that distance (distSqr_)
223 // is not affected by leaving/entering domain.
224 template<class TrackingData>
225 inline void Foam::pointEdgePoint::enterDomain
227     const polyPatch& patch,
228     const label patchPointI,
229     const point& coord,
230     TrackingData& td
233     // back to absolute form
234     origin_ += coord;
238 // Update this with information from connected edge
239 template<class TrackingData>
240 inline bool Foam::pointEdgePoint::updatePoint
242     const polyMesh& mesh,
243     const label pointI,
244     const label edgeI,
245     const pointEdgePoint& edgeInfo,
246     const scalar tol,
247     TrackingData& td
250     return update(mesh.points()[pointI], edgeInfo, tol, td);
254 // Update this with new information on same point
255 template<class TrackingData>
256 inline bool Foam::pointEdgePoint::updatePoint
258     const polyMesh& mesh,
259     const label pointI,
260     const pointEdgePoint& newPointInfo,
261     const scalar tol,
262     TrackingData& td
265     return update(mesh.points()[pointI], newPointInfo, tol, td);
269 // Update this with new information on same point. No extra information.
270 template<class TrackingData>
271 inline bool Foam::pointEdgePoint::updatePoint
273     const pointEdgePoint& newPointInfo,
274     const scalar tol,
275     TrackingData& td
278     return update(newPointInfo, tol, td);
282 // Update this with information from connected point
283 template<class TrackingData>
284 inline bool Foam::pointEdgePoint::updateEdge
286     const polyMesh& mesh,
287     const label edgeI,
288     const label pointI,
289     const pointEdgePoint& pointInfo,
290     const scalar tol,
291     TrackingData& td
294     const edge& e = mesh.edges()[edgeI];
295     return update(e.centre(mesh.points()), pointInfo, tol, td);
299 template <class TrackingData>
300 inline bool Foam::pointEdgePoint::equal
302     const pointEdgePoint& rhs,
303     TrackingData& td
304 ) const
306     return operator==(rhs);
310 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
312 inline bool Foam::pointEdgePoint::operator==(const Foam::pointEdgePoint& rhs)
313  const
315     return origin() == rhs.origin();
319 inline bool Foam::pointEdgePoint::operator!=(const Foam::pointEdgePoint& rhs)
320  const
322     return !(*this == rhs);
326 // ************************************************************************* //