BUGFIX: Seg-fault in multiphaseInterFoam. Author: Henrik Rusche. Merge: Hrvoje Jasak
[foam-extend-3.2.git] / src / meshTools / PointEdgeWave / pointEdgePointI.H
bloba4058291bd3462a2bad590b355b7962288b63f6f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  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 inline bool Foam::pointEdgePoint::update
34     const point& pt,
35     const pointEdgePoint& w2,
36     const scalar tol
39     scalar dist2 = magSqr(pt - w2.origin());
41     if (!valid())
42     {
43         // current not yet set so use any value
44         distSqr_ = dist2;
45         origin_ = w2.origin();
47         return true;
48     }
50     scalar diff = distSqr_ - dist2;
52     if (diff < 0)
53     {
54         // already nearer to pt
55         return false;
56     }
58     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
59     {
60         // don't propagate small changes
61         return false;
62     }
63     else
64     {
65         // update with new values
66         distSqr_ = dist2;
67         origin_ = w2.origin();
69         return true;
70     }
74 // Update this with w2 (information on same point)
75 inline bool Foam::pointEdgePoint::update
77     const pointEdgePoint& w2,
78     const scalar tol
81     if (!valid())
82     {
83         // current not yet set so use any value
84         distSqr_ = w2.distSqr();
85         origin_ = w2.origin();
87         return true;
88     }
90     scalar diff = distSqr_ - w2.distSqr();
92     if (diff < 0)
93     {
94         // already nearer to pt
95         return false;
96     }
98     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
99     {
100         // don't propagate small changes
101         return false;
102     }
103     else
104     {
105         // update with new values
106         distSqr_ =  w2.distSqr();
107         origin_ = w2.origin();
109         return true;
110     }
113 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
115 // Null constructor
116 inline Foam::pointEdgePoint::pointEdgePoint()
118     origin_(greatPoint),
119     distSqr_(GREAT)
123 // Construct from origin, distance
124 inline Foam::pointEdgePoint::pointEdgePoint
126     const point& origin,
127     const scalar distSqr
130     origin_(origin),
131     distSqr_(distSqr)
135 // Construct as copy
136 inline Foam::pointEdgePoint::pointEdgePoint(const pointEdgePoint& wpt)
138     origin_(wpt.origin()),
139     distSqr_(wpt.distSqr())
143 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
145 inline const Foam::point& Foam::pointEdgePoint::origin() const
147     return origin_;
151 inline Foam::scalar Foam::pointEdgePoint::distSqr() const
153     return distSqr_;
157 inline bool Foam::pointEdgePoint::valid() const
159     return origin_ != greatPoint;
163 // Checks for cyclic points
164 inline bool Foam::pointEdgePoint::sameGeometry
166     const pointEdgePoint& w2,
167     const scalar tol
168 ) const
170     scalar diff = Foam::mag(distSqr() - w2.distSqr());
172     if (diff < SMALL)
173     {
174         return true;
175     }
176     else
177     {
178         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
179         {
180             return true;
181         }
182         else
183         {
184             return false;
185         }
186     }
190 inline void Foam::pointEdgePoint::leaveDomain
192     const polyPatch& patch,
193     const label patchPointI,
194     const point& coord
197     origin_ -= coord;
201 inline void Foam::pointEdgePoint::transform(const tensor& rotTensor)
203     origin_ = Foam::transform(rotTensor, origin_);
207 // Update absolute geometric quantities. Note that distance (distSqr_)
208 // is not affected by leaving/entering domain.
209 inline void Foam::pointEdgePoint::enterDomain
211     const polyPatch& patch,
212     const label patchPointI,
213     const point& coord
216     // back to absolute form
217     origin_ += coord;
221 // Update this with information from connected edge
222 inline bool Foam::pointEdgePoint::updatePoint
224     const polyMesh& mesh,
225     const label pointI,
226     const label edgeI,
227     const pointEdgePoint& edgeInfo,
228     const scalar tol
231     return
232         update
233         (
234             mesh.allPoints()[pointI],
235             edgeInfo,
236             tol
237         );
238     }
241 // Update this with new information on same point
242 inline bool Foam::pointEdgePoint::updatePoint
244     const polyMesh& mesh,
245     const label pointI,
246     const pointEdgePoint& newPointInfo,
247     const scalar tol
250     return
251         update
252         (
253             mesh.allPoints()[pointI],
254             newPointInfo,
255             tol
256         );
260 // Update this with new information on same point. No extra information.
261 inline bool Foam::pointEdgePoint::updatePoint
263     const pointEdgePoint& newPointInfo,
264     const scalar tol
267     return update(newPointInfo, tol);
271 // Update this with information from connected point
272 inline bool Foam::pointEdgePoint::updateEdge
274     const polyMesh& mesh,
275     const label edgeI,
276     const label pointI,
277     const pointEdgePoint& pointInfo,
278     const scalar tol
281     const pointField& points = mesh.points();
283     const edge& e = mesh.edges()[edgeI];
285     const point edgeMid(0.5*(points[e[0]] + points[e[1]]));
287     return
288         update
289         (
290             edgeMid,
291             pointInfo,
292             tol
293         );
297 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
299 inline bool Foam::pointEdgePoint::operator==(const Foam::pointEdgePoint& rhs)
300  const
302     return origin() == rhs.origin();
306 inline bool Foam::pointEdgePoint::operator!=(const Foam::pointEdgePoint& rhs)
307  const
309     return !(*this == rhs);
313 // ************************************************************************* //