FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

lineI.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include <OpenFOAM/IOstreams.H>
00027 
00028 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00029 
00030 namespace Foam
00031 {
00032 
00033 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00034 
00035 template<class Point, class PointRef>
00036 inline line<Point, PointRef>::line(const Point& start, const Point& end)
00037 :
00038     a_(start),
00039     b_(end)
00040 {}
00041 
00042 
00043 template<class Point, class PointRef>
00044 inline line<Point, PointRef>::line(Istream& is)
00045 {
00046     // Read beginning of line point pair
00047     is.readBegin("line");
00048 
00049     is >> a_ >> b_;
00050 
00051     // Read end of line point pair
00052     is.readEnd("line");
00053 
00054     // Check state of Istream
00055     is.check("line::line(Istream& is)");
00056 }
00057 
00058 
00059 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00060 
00061 template<class Point, class PointRef>
00062 inline PointRef line<Point, PointRef>::start() const
00063 {
00064     return a_;
00065 }
00066 
00067 template<class Point, class PointRef>
00068 inline PointRef line<Point, PointRef>::end() const
00069 {
00070     return b_;
00071 }
00072 
00073 
00074 template<class Point, class PointRef>
00075 inline Point line<Point, PointRef>::centre() const
00076 {
00077     return 0.5*(a_ + b_);
00078 }
00079 
00080 
00081 template<class Point, class PointRef>
00082 inline scalar line<Point, PointRef>::mag() const
00083 {
00084     return ::Foam::mag(vec());
00085 }
00086 
00087 
00088 template<class Point, class PointRef>
00089 inline Point line<Point, PointRef>::vec() const
00090 {
00091     return b_ - a_;
00092 }
00093 
00094 
00095 template<class Point, class PointRef>
00096 PointHit<Point> line<Point, PointRef>::nearestDist(const Point& p) const
00097 {
00098     Point v = vec();
00099 
00100     Point w(p - a_);
00101 
00102     scalar c1 = v & w;
00103 
00104     if (c1 <= 0)
00105     {
00106         return PointHit<Point>(false, a_, Foam::mag(p - a_), true);
00107     }
00108 
00109     scalar c2 = v & v;
00110 
00111     if (c2 <= c1)
00112     {
00113         return PointHit<Point>(false, b_, Foam::mag(p - b_), true);
00114     }
00115 
00116     scalar b = c1/c2;
00117 
00118     Point pb(a_ + b*v);
00119 
00120     return PointHit<Point>(true, pb, Foam::mag(p - pb), false);
00121 }
00122 
00123 
00124 template<class Point, class PointRef>
00125 scalar line<Point, PointRef>::nearestDist
00126 (
00127     const line<Point, const Point&>& edge,
00128     Point& thisPt,
00129     Point& edgePt
00130 ) const
00131 {
00132     // From Mathworld Line-Line distance/(Gellert et al. 1989, p. 538).
00133     Point a(end() - start());
00134     Point b(edge.end() - edge.start());
00135     Point c(edge.start() - start());
00136 
00137     Point crossab = a ^ b;
00138     scalar magCrossSqr = magSqr(crossab);
00139 
00140     if (magCrossSqr > VSMALL)
00141     {
00142         scalar s = ((c ^ b) & crossab)/magCrossSqr;
00143         scalar t = ((c ^ a) & crossab)/magCrossSqr;
00144 
00145         // Check for end points outside of range 0..1
00146         if (s >= 0 && s <= 1 && t >= 0 && t <= 1)
00147         {
00148             // Both inside range 0..1
00149             thisPt = start() + a*s;
00150             edgePt = edge.start() + b*t;
00151         }
00152         else
00153         {
00154             // Do brute force. Distance of everything to everything.
00155             // Can quite possibly be improved!
00156 
00157             // From edge endpoints to *this
00158             PointHit<Point> this0(nearestDist(edge.start()));
00159             PointHit<Point> this1(nearestDist(edge.end()));
00160             scalar thisDist = min(this0.distance(), this1.distance());
00161 
00162             // From *this to edge
00163             PointHit<Point> edge0(edge.nearestDist(start()));
00164             PointHit<Point> edge1(edge.nearestDist(end()));
00165             scalar edgeDist = min(edge0.distance(), edge1.distance());
00166 
00167             if (thisDist < edgeDist)
00168             {
00169                 if (this0.distance() < this1.distance())
00170                 {
00171                     thisPt = this0.rawPoint();
00172                     edgePt = edge.start();
00173                 }
00174                 else
00175                 {
00176                     thisPt = this1.rawPoint();
00177                     edgePt = edge.end();
00178                 }
00179             }
00180             else
00181             {
00182                 if (edge0.distance() < edge1.distance())
00183                 {
00184                     thisPt = start();
00185                     edgePt = edge0.rawPoint();
00186                 }
00187                 else
00188                 {
00189                     thisPt = end();
00190                     edgePt = edge1.rawPoint();
00191                 }
00192             }
00193         }
00194     }
00195     else
00196     {
00197         // Parallel lines. Find overlap of both lines by projecting onto
00198         // direction vector (now equal for both lines).
00199 
00200         scalar edge0 = edge.start() & a;
00201         scalar edge1 = edge.end() & a;
00202         bool edgeOrder = edge0 < edge1;
00203 
00204         scalar minEdge = (edgeOrder ? edge0 : edge1);
00205         scalar maxEdge = (edgeOrder ? edge1 : edge0);
00206         const Point& minEdgePt = (edgeOrder ? edge.start() : edge.end());
00207         const Point& maxEdgePt = (edgeOrder ? edge.end() : edge.start());
00208 
00209         scalar this0 = start() & a;
00210         scalar this1 = end() & a;
00211         bool thisOrder = this0 < this1;
00212 
00213         scalar minThis = min(this0, this1);
00214         scalar maxThis = max(this1, this0);
00215         const Point& minThisPt = (thisOrder ? start() : end());
00216         const Point& maxThisPt = (thisOrder ? end() : start());
00217 
00218         if (maxEdge < minThis)
00219         {
00220             // edge completely below *this
00221             edgePt = maxEdgePt;
00222             thisPt = minThisPt;
00223         }
00224         else if (maxEdge < maxThis)
00225         {
00226             // maxEdge inside interval of *this
00227             edgePt = maxEdgePt;
00228             thisPt = nearestDist(edgePt).rawPoint();
00229         }
00230         else
00231         {
00232             // maxEdge outside. Check if minEdge inside.
00233             if (minEdge < minThis)
00234             {
00235                 // Edge completely envelops this. Take any this point and
00236                 // determine nearest on edge.
00237                 thisPt = minThisPt;
00238                 edgePt = edge.nearestDist(thisPt).rawPoint();
00239             }
00240             else if (minEdge < maxThis)
00241             {
00242                 // minEdge inside this interval.
00243                 edgePt = minEdgePt;
00244                 thisPt = nearestDist(edgePt).rawPoint();
00245             }
00246             else
00247             {
00248                 // minEdge outside this interval
00249                 edgePt = minEdgePt;
00250                 thisPt = maxThisPt;
00251             }
00252         }
00253     }
00254 
00255     return Foam::mag(thisPt - edgePt);
00256 }
00257 
00258 
00259 // * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * * //
00260 
00261 template<class Point, class PointRef>
00262 inline Istream& operator>>(Istream& is, line<Point, PointRef>& l)
00263 {
00264     // Read beginning of line point pair
00265     is.readBegin("line");
00266 
00267     is >> l.a_ >> l.b_;
00268 
00269     // Read end of line point pair
00270     is.readEnd("line");
00271 
00272     // Check state of Ostream
00273     is.check("Istream& operator>>(Istream&, line&)");
00274 
00275     return is;
00276 }
00277 
00278 
00279 template<class Point, class PointRef>
00280 inline Ostream& operator<<(Ostream& os, const line<Point, PointRef>& l)
00281 {
00282     os << token::BEGIN_LIST << l.a_ << token::SPACE << l.b_ << token::END_LIST;
00283     return os;
00284 }
00285 
00286 
00287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00288 
00289 } // End namespace Foam
00290 
00291 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines