initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / meshTools / cellDist / wallPoint / wallPointI.H
blobc9f76a35f2a1f519964ec450ed672b4df5ae34d0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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::wallPoint::update
35     const point& pt,
36     const wallPoint& w2,
37     const scalar tol
40     //Already done in calling algorithm
41     //if (w2.origin() == origin_)
42     //{
43     //    // Shortcut. Same input so same distance.
44     //    return false;
45     //}
47     scalar dist2 = magSqr(pt - w2.origin());
49     if (!valid())
50     {
51         // current not yet set so use any value
52         distSqr_ = dist2;
53         origin_ = w2.origin();
55         return true;
56     }        
58     scalar diff = distSqr_ - dist2;
60     if (diff < 0)
61     {
62         // already nearer to pt
63         return false;
64     }
66     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
67     {
68         // don't propagate small changes
69         return false;
70     }
71     else
72     {
73         // update with new values
74         distSqr_ = dist2;
75         origin_ = w2.origin();
77         return true;
78     }
80     
82 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
84 // Null constructor
85 inline Foam::wallPoint::wallPoint()
87     origin_(greatPoint),
88     distSqr_(GREAT)
92 // Construct from origin, distance
93 inline Foam::wallPoint::wallPoint(const point& origin, const scalar distSqr)
95     origin_(origin), distSqr_(distSqr)
99 // Construct as copy
100 inline Foam::wallPoint::wallPoint(const wallPoint& wpt)
102     origin_(wpt.origin()), distSqr_(wpt.distSqr())
106 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
108 inline const Foam::point& Foam::wallPoint::origin() const
110     return origin_;
114 inline Foam::point& Foam::wallPoint::origin()
116     return origin_;
120 inline Foam::scalar Foam::wallPoint::distSqr() const
122     return distSqr_;
126 inline Foam::scalar& Foam::wallPoint::distSqr()
128     return distSqr_;
132 inline bool Foam::wallPoint::valid() const
134     return origin_ != greatPoint;
138 // Checks for cyclic faces
139 inline bool Foam::wallPoint::sameGeometry
141     const polyMesh&,
142     const wallPoint& w2,
143     const scalar tol
145  const
147     scalar diff = mag(distSqr() - w2.distSqr());
149     if (diff < SMALL)
150     {
151         return true;
152     }
153     else
154     {
155         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
156         {
157             return true;
158         }
159         else
160         {
161             return false;
162         }
163     }
167 inline void Foam::wallPoint::leaveDomain
169     const polyMesh&,
170     const polyPatch&,
171     const label,
172     const point& faceCentre
175     origin_ -= faceCentre;
179 inline void Foam::wallPoint::transform
181     const polyMesh&,
182     const tensor& rotTensor
185     origin_ = Foam::transform(rotTensor, origin_);
189 // Update absolute geometric quantities. Note that distance (distSqr_)
190 // is not affected by leaving/entering domain.
191 inline void Foam::wallPoint::enterDomain
193     const polyMesh&,
194     const polyPatch&,
195     const label,
196     const point& faceCentre
199     // back to absolute form
200     origin_ += faceCentre;
204 // Update this with w2 if w2 nearer to pt.
205 inline bool Foam::wallPoint::updateCell
207     const polyMesh& mesh,
208     const label thisCellI,
209     const label neighbourFaceI,
210     const wallPoint& neighbourWallInfo,
211     const scalar tol
214     return 
215         update
216         (
217             mesh.cellCentres()[thisCellI],
218             neighbourWallInfo,
219             tol
220         );
221     }    
224 // Update this with w2 if w2 nearer to pt.
225 inline bool Foam::wallPoint::updateFace
227     const polyMesh& mesh,
228     const label thisFaceI,
229     const label neighbourCellI,
230     const wallPoint& neighbourWallInfo,
231     const scalar tol
234     return
235         update
236         (
237             mesh.faceCentres()[thisFaceI],
238             neighbourWallInfo,
239             tol
240         );
241 }    
243 // Update this with w2 if w2 nearer to pt.
244 inline bool Foam::wallPoint::updateFace
246     const polyMesh& mesh,
247     const label thisFaceI,
248     const wallPoint& neighbourWallInfo,
249     const scalar tol
252     return
253         update
254         (
255             mesh.faceCentres()[thisFaceI],
256             neighbourWallInfo,
257             tol
258         );
259 }    
262 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
264 inline bool Foam::wallPoint::operator==(const Foam::wallPoint& rhs) const
266     return origin() == rhs.origin();
270 inline bool Foam::wallPoint::operator!=(const Foam::wallPoint& rhs) const
272     return !(*this == rhs);
276 // ************************************************************************* //