ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / meshTools / PatchEdgeFaceWave / patchEdgeFaceInfoI.H
blob881e586626a1f69760486c948a9a66ed0da341dd
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::patchEdgeFaceInfo::update
35     const point& pt,
36     const patchEdgeFaceInfo& 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 edge)
77 template<class TrackingData>
78 inline bool Foam::patchEdgeFaceInfo::update
80     const patchEdgeFaceInfo& 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::patchEdgeFaceInfo::patchEdgeFaceInfo()
123     origin_(point::max),
124     distSqr_(sqr(GREAT))
128 // Construct from origin, distance
129 inline Foam::patchEdgeFaceInfo::patchEdgeFaceInfo
131     const point& origin,
132     const scalar distSqr
135     origin_(origin),
136     distSqr_(distSqr)
140 // Construct as copy
141 inline Foam::patchEdgeFaceInfo::patchEdgeFaceInfo(const patchEdgeFaceInfo& wpt)
143     origin_(wpt.origin()),
144     distSqr_(wpt.distSqr())
148 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
150 inline const Foam::point& Foam::patchEdgeFaceInfo::origin() const
152     return origin_;
156 inline Foam::scalar Foam::patchEdgeFaceInfo::distSqr() const
158     return distSqr_;
162 template<class TrackingData>
163 inline bool Foam::patchEdgeFaceInfo::valid(TrackingData& td) const
165     return origin_ != point::max;
169 template<class TrackingData>
170 inline void Foam::patchEdgeFaceInfo::transform
172     const polyMesh& mesh,
173     const primitivePatch& patch,
174     const tensor& rotTensor,
175     const scalar tol,
176     TrackingData& td
179     origin_ = Foam::transform(rotTensor, origin_);
183 template<class TrackingData>
184 inline bool Foam::patchEdgeFaceInfo::updateEdge
186     const polyMesh& mesh,
187     const primitivePatch& patch,
188     const label edgeI,
189     const label faceI,
190     const patchEdgeFaceInfo& faceInfo,
191     const scalar tol,
192     TrackingData& td
195     const edge& e = patch.edges()[edgeI];
196     point eMid =
197         0.5
198       * (
199             patch.points()[patch.meshPoints()[e[0]]]
200           + patch.points()[patch.meshPoints()[e[1]]]
201         );
202     return update(eMid, faceInfo, tol, td);
206 template<class TrackingData>
207 inline bool Foam::patchEdgeFaceInfo::updateEdge
209     const polyMesh& mesh,
210     const primitivePatch& patch,
211     const patchEdgeFaceInfo& edgeInfo,
212     const bool sameOrientation,
213     const scalar tol,
214     TrackingData& td
217     return update(edgeInfo, tol, td);
221 template<class TrackingData>
222 inline bool Foam::patchEdgeFaceInfo::updateFace
224     const polyMesh& mesh,
225     const primitivePatch& patch,
226     const label faceI,
227     const label edgeI,
228     const patchEdgeFaceInfo& edgeInfo,
229     const scalar tol,
230     TrackingData& td
233     return update(patch.faceCentres()[faceI], edgeInfo, tol, td);
237 template <class TrackingData>
238 inline bool Foam::patchEdgeFaceInfo::equal
240     const patchEdgeFaceInfo& rhs,
241     TrackingData& td
242 ) const
244     return operator==(rhs);
248 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
250 inline bool Foam::patchEdgeFaceInfo::operator==
252     const Foam::patchEdgeFaceInfo& rhs
253 ) const
255     return origin() == rhs.origin();
259 inline bool Foam::patchEdgeFaceInfo::operator!=
261     const Foam::patchEdgeFaceInfo& rhs
262 ) const
264     return !(*this == rhs);
268 // ************************************************************************* //