initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / meshTools / PointEdgeWave / pointEdgePointI.H
blob360840f9106269304b70c11631f311988eb3c591
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "polyMesh.H"
28 #include "transform.H"
30 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
32 // Update this with w2 if w2 nearer to pt.
33 inline bool Foam::pointEdgePoint::update
35     const point& pt,
36     const pointEdgePoint& w2,
37     const scalar tol
40     scalar dist2 = magSqr(pt - w2.origin());
42     if (!valid())
43     {
44         // current not yet set so use any value
45         distSqr_ = dist2;
46         origin_ = w2.origin();
48         return true;
49     }        
51     scalar diff = distSqr_ - dist2;
53     if (diff < 0)
54     {
55         // already nearer to pt
56         return false;
57     }
59     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
60     {
61         // don't propagate small changes
62         return false;
63     }
64     else
65     {
66         // update with new values
67         distSqr_ = dist2;
68         origin_ = w2.origin();
70         return true;
71     }
74     
75 // Update this with w2 (information on same point)
76 inline bool Foam::pointEdgePoint::update
78     const pointEdgePoint& w2,
79     const scalar tol
82     if (!valid())
83     {
84         // current not yet set so use any value
85         distSqr_ = w2.distSqr();
86         origin_ = w2.origin();
88         return true;
89     }        
91     scalar diff = distSqr_ - w2.distSqr();
93     if (diff < 0)
94     {
95         // already nearer to pt
96         return false;
97     }
99     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
100     {
101         // don't propagate small changes
102         return false;
103     }
104     else
105     {
106         // update with new values
107         distSqr_ =  w2.distSqr();
108         origin_ = w2.origin();
110         return true;
111     }
114 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
116 // Null constructor
117 inline Foam::pointEdgePoint::pointEdgePoint()
119     origin_(greatPoint),
120     distSqr_(GREAT)
124 // Construct from origin, distance
125 inline Foam::pointEdgePoint::pointEdgePoint
127     const point& origin, 
128     const scalar distSqr
131     origin_(origin),
132     distSqr_(distSqr)
136 // Construct as copy
137 inline Foam::pointEdgePoint::pointEdgePoint(const pointEdgePoint& wpt)
139     origin_(wpt.origin()),
140     distSqr_(wpt.distSqr())
144 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
146 inline const Foam::point& Foam::pointEdgePoint::origin() const
148     return origin_;
152 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
154     return distSqr_;
158 inline bool Foam::pointEdgePoint::valid() const
160     return origin_ != greatPoint;
164 // Checks for cyclic points
165 inline bool Foam::pointEdgePoint::sameGeometry
167     const pointEdgePoint& w2,
168     const scalar tol
169 ) const
171     scalar diff = Foam::mag(distSqr() - w2.distSqr());
173     if (diff < SMALL)
174     {
175         return true;
176     }
177     else
178     {
179         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
180         {
181             return true;
182         }
183         else
184         {
185             return false;
186         }
187     }
191 inline void Foam::pointEdgePoint::leaveDomain
193     const polyPatch& patch,
194     const label patchPointI,
195     const point& coord
198     origin_ -= coord;
202 inline void Foam::pointEdgePoint::transform(const tensor& rotTensor)
204     origin_ = Foam::transform(rotTensor, origin_);
208 // Update absolute geometric quantities. Note that distance (distSqr_)
209 // is not affected by leaving/entering domain.
210 inline void Foam::pointEdgePoint::enterDomain
212     const polyPatch& patch,
213     const label patchPointI,
214     const point& coord
217     // back to absolute form
218     origin_ += coord;
222 // Update this with information from connected edge
223 inline bool Foam::pointEdgePoint::updatePoint
225     const polyMesh& mesh,
226     const label pointI,
227     const label edgeI,
228     const pointEdgePoint& edgeInfo,
229     const scalar tol
232     return 
233         update
234         (
235             mesh.points()[pointI],
236             edgeInfo,
237             tol
238         );
239     }    
242 // Update this with new information on same point
243 inline bool Foam::pointEdgePoint::updatePoint
245     const polyMesh& mesh,
246     const label pointI,
247     const pointEdgePoint& newPointInfo,
248     const scalar tol
251     return
252         update
253         (
254             mesh.points()[pointI],
255             newPointInfo,
256             tol
257         );
258 }    
261 // Update this with new information on same point. No extra information.
262 inline bool Foam::pointEdgePoint::updatePoint
264     const pointEdgePoint& newPointInfo,
265     const scalar tol
268     return update(newPointInfo, tol);
269 }    
272 // Update this with information from connected point
273 inline bool Foam::pointEdgePoint::updateEdge
275     const polyMesh& mesh,
276     const label edgeI,
277     const label pointI,
278     const pointEdgePoint& pointInfo,
279     const scalar tol
282     const pointField& points = mesh.points();
284     const edge& e = mesh.edges()[edgeI];
286     const point edgeMid(0.5*(points[e[0]] + points[e[1]]));
288     return
289         update
290         (
291             edgeMid,
292             pointInfo,
293             tol
294         );
295 }    
298 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
300 inline bool Foam::pointEdgePoint::operator==(const Foam::pointEdgePoint& rhs)
301  const
303     return origin() == rhs.origin();
307 inline bool Foam::pointEdgePoint::operator!=(const Foam::pointEdgePoint& rhs)
308  const
310     return !(*this == rhs);
314 // ************************************************************************* //