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

plane.C

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 "plane.H"
00027 #include <OpenFOAM/tensor.H>
00028 
00029 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00030 
00031 // Calculate base point and unit normal vector from plane equation
00032 void Foam::plane::calcPntAndVec(const scalarList& C)
00033 {
00034     if (mag(C[0]) > VSMALL)
00035     {
00036         basePoint_ = vector((-C[3]/C[0]), 0, 0);
00037     }
00038     else
00039     {
00040         if (mag(C[1]) > VSMALL)
00041         {
00042             basePoint_ = vector(0, (-C[3]/C[1]), 0);
00043         }
00044         else
00045         {
00046             if (mag(C[2]) > VSMALL)
00047             {
00048                 basePoint_ = vector(0, 0, (-C[3]/C[2]));
00049             }
00050             else
00051             {
00052                 FatalErrorIn("void plane::calcPntAndVec(const scalarList&)")
00053                     << "At least one plane coefficient must have a value"
00054                     << abort(FatalError);
00055             }
00056         }
00057     }
00058 
00059     unitVector_ = vector(C[0], C[1], C[2]);
00060     scalar magUnitVector(mag(unitVector_));
00061 
00062     if (magUnitVector < VSMALL)
00063     {
00064         FatalErrorIn("void plane::calcPntAndVec(const scalarList&)")
00065             << "Plane normal defined with zero length"
00066             << abort(FatalError);
00067     }
00068 
00069     unitVector_ /= magUnitVector;
00070 }
00071 
00072 
00073 void Foam::plane::calcPntAndVec
00074 (
00075     const point& point1,
00076     const point& point2,
00077     const point& point3
00078 )
00079 {
00080     basePoint_ = (point1 + point2 + point3)/3;
00081     vector line12 = point1 - point2;
00082     vector line23 = point2 - point3;
00083 
00084     if
00085     (
00086         mag(line12) < VSMALL
00087      || mag(line23) < VSMALL
00088      || mag(point3-point1) < VSMALL
00089     )
00090     {
00091         FatalErrorIn
00092         (
00093             "void plane::calcPntAndVec\n"
00094             "(\n"
00095             "    const point&,\n"
00096             "    const point&,\n"
00097             "    const point&\n"
00098             ")\n"
00099         ) << "Bad points." << abort(FatalError);
00100     }
00101 
00102     unitVector_ = line12 ^ line23;
00103     scalar magUnitVector(mag(unitVector_));
00104 
00105     if (magUnitVector < VSMALL)
00106     {
00107         FatalErrorIn
00108         (
00109             "void plane::calcPntAndVec\n"
00110             "(\n"
00111             "    const point&,\n"
00112             "    const point&,\n"
00113             "    const point&\n"
00114             ")\n"
00115         )   << "Plane normal defined with zero length"
00116             << abort(FatalError);
00117     }
00118 
00119     unitVector_ /= magUnitVector;
00120 }
00121 
00122 
00123 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00124 
00125 // Construct from normal vector through the origin
00126 Foam::plane::plane(const vector& normalVector)
00127 :
00128     unitVector_(normalVector),
00129     basePoint_(vector::zero)
00130 {}
00131 
00132 
00133 // Construct from point and normal vector
00134 Foam::plane::plane(const point& basePoint, const vector& normalVector)
00135 :
00136     unitVector_(normalVector),
00137     basePoint_(basePoint)
00138 {
00139     scalar magUnitVector(mag(unitVector_));
00140 
00141     if (magUnitVector > VSMALL)
00142     {
00143         unitVector_ /= magUnitVector;
00144     }
00145     else
00146     {
00147         FatalErrorIn("plane::plane(const point&, const vector&)")
00148             << "plane normal has zero length"
00149             << abort(FatalError);
00150     }
00151 }
00152 
00153 
00154 // Construct from plane equation
00155 Foam::plane::plane(const scalarList& C)
00156 {
00157     calcPntAndVec(C);
00158 }
00159 
00160 
00161 // Construct from three points
00162 Foam::plane::plane
00163 (
00164     const point& a,
00165     const point& b,
00166     const point& c
00167 )
00168 {
00169     calcPntAndVec(a, b, c);
00170 }
00171 
00172 
00173 // Construct from dictionary
00174 Foam::plane::plane(const dictionary& dict)
00175 :
00176     unitVector_(vector::zero),
00177     basePoint_(point::zero)
00178 {
00179     word planeType(dict.lookup("planeType"));
00180 
00181     if (planeType == "planeEquation")
00182     {
00183         const dictionary& subDict = dict.subDict("planeEquationDict");
00184         scalarList C(4);
00185 
00186         C[0] = readScalar(subDict.lookup("a"));
00187         C[1] = readScalar(subDict.lookup("b"));
00188         C[2] = readScalar(subDict.lookup("c"));
00189         C[3] = readScalar(subDict.lookup("d"));
00190 
00191         calcPntAndVec(C);
00192 
00193     }
00194     else if (planeType == "embeddedPoints")
00195     {
00196         const dictionary& subDict = dict.subDict("embeddedPoints");
00197 
00198         point point1(subDict.lookup("point1"));
00199         point point2(subDict.lookup("point2"));
00200         point point3(subDict.lookup("point3"));
00201 
00202         calcPntAndVec(point1, point2, point3);
00203     }
00204     else if (planeType == "pointAndNormal")
00205     {
00206         const dictionary& subDict = dict.subDict("pointAndNormalDict");
00207 
00208         basePoint_ = subDict.lookup("basePoint");
00209         unitVector_ = subDict.lookup("normalVector");
00210         unitVector_ /= mag(unitVector_);
00211     }
00212     else
00213     {
00214         FatalIOErrorIn
00215         (
00216             "plane::plane(const dictionary&)",
00217             dict
00218         )
00219             << "Invalid plane type: " << planeType
00220             << abort(FatalIOError);
00221     }
00222 }
00223 
00224 
00225 // Construct from Istream. Assumes point and normal vector.
00226 Foam::plane::plane(Istream& is)
00227 :
00228     unitVector_(is),
00229     basePoint_(is)
00230 {
00231     scalar magUnitVector(mag(unitVector_));
00232 
00233     if (magUnitVector > VSMALL)
00234     {
00235         unitVector_ /= magUnitVector;
00236     }
00237     else
00238     {
00239         FatalErrorIn("plane::plane(Istream& is)")
00240             << "plane normal has zero length"
00241             << abort(FatalError);
00242     }
00243 }
00244 
00245 
00246 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00247 
00248 // Return plane normal vector
00249 const Foam::vector& Foam::plane::normal() const
00250 {
00251     return unitVector_;
00252 }
00253 
00254 
00255 // Return plane base point
00256 const Foam::point& Foam::plane::refPoint() const
00257 {
00258     return basePoint_;
00259 }
00260 
00261 
00262 // Return coefficients for plane equation: ax + by + cz + d = 0
00263 Foam::FixedList<Foam::scalar, 4> Foam::plane::planeCoeffs() const
00264 {
00265     FixedList<scalar, 4> C(4);
00266 
00267     scalar magX = mag(unitVector_.x());
00268     scalar magY = mag(unitVector_.y());
00269     scalar magZ = mag(unitVector_.z());
00270 
00271     if (magX > magY)
00272     {
00273         if (magX > magZ)
00274         {
00275             C[0] = 1;
00276             C[1] = unitVector_.y()/unitVector_.x();
00277             C[2] = unitVector_.z()/unitVector_.x();
00278         }
00279         else
00280         {
00281             C[0] = unitVector_.x()/unitVector_.z();
00282             C[1] = unitVector_.y()/unitVector_.z();
00283             C[2] = 1;
00284         }
00285     }
00286     else
00287     {
00288         if (magY > magZ)
00289         {
00290             C[0] = unitVector_.x()/unitVector_.y();
00291             C[1] = 1;
00292             C[2] = unitVector_.z()/unitVector_.y();
00293         }
00294         else
00295         {
00296             C[0] = unitVector_.x()/unitVector_.z();
00297             C[1] = unitVector_.y()/unitVector_.z();
00298             C[2] = 1;
00299         }
00300     }
00301 
00302     C[3] = - C[0] * basePoint_.x()
00303            - C[1] * basePoint_.y()
00304            - C[2] * basePoint_.z();
00305 
00306     return C;
00307 }
00308 
00309 
00310 // Return nearest point in the plane for the given point
00311 Foam::point Foam::plane::nearestPoint(const point& p) const
00312 {
00313     return p - unitVector_*((p - basePoint_) & unitVector_);
00314 }
00315 
00316 
00317 // Return distance from the given point to the plane
00318 Foam::scalar Foam::plane::distance(const point& p) const
00319 {
00320     return mag((p - basePoint_) & unitVector_);
00321 }
00322 
00323 
00324 // Cutting point for plane and line defined by origin and direction
00325 Foam::scalar Foam::plane::normalIntersect
00326 (
00327     const point& pnt0,
00328     const vector& dir
00329 ) const
00330 {
00331     scalar denom = stabilise((dir & unitVector_), VSMALL);
00332 
00333     return ((basePoint_ - pnt0) & unitVector_)/denom;
00334 }
00335 
00336 
00337 // Cutting line of two planes
00338 Foam::plane::ray Foam::plane::planeIntersect(const plane& plane2) const
00339 {
00340     // Mathworld plane-plane intersection. Assume there is a point on the
00341     // intersection line with z=0 and solve the two plane equations
00342     // for that (now 2x2 equation in x and y)
00343     // Better: use either z=0 or x=0 or y=0.
00344 
00345     const vector& n1 = normal();
00346     const vector& n2 = plane2.normal();
00347 
00348     const point& p1 = refPoint();
00349     const point& p2 = plane2.refPoint();
00350 
00351     scalar n1p1 = n1&p1;
00352     scalar n2p2 = n2&p2;
00353 
00354     vector dir = n1 ^ n2;
00355 
00356     // Determine zeroed out direction (can be x,y or z) by looking at which
00357     // has the largest component in dir.
00358     scalar magX = mag(dir.x());
00359     scalar magY = mag(dir.y());
00360     scalar magZ = mag(dir.z());
00361 
00362     direction iZero, i1, i2;
00363 
00364     if (magX > magY)
00365     {
00366         if (magX > magZ)
00367         {
00368             iZero = 0;
00369             i1 = 1;
00370             i2 = 2;
00371         }
00372         else
00373         {
00374             iZero = 2;
00375             i1 = 0;
00376             i2 = 1;
00377         }
00378     }
00379     else
00380     {
00381         if (magY > magZ)
00382         {
00383             iZero = 1;
00384             i1 = 2;
00385             i2 = 0;
00386         }
00387         else
00388         {
00389             iZero = 2;
00390             i1 = 0;
00391             i2 = 1;
00392         }
00393     }
00394 
00395     vector pt;
00396 
00397     pt[iZero] = 0;
00398     pt[i1] = (n2[i2]*n1p1 - n1[i2]*n2p2) / (n1[i1]*n2[i2] - n2[i1]*n1[i2]);
00399     pt[i2] = (n2[i1]*n1p1 - n1[i1]*n2p2) / (n1[i2]*n2[i1] - n1[i1]*n2[i2]);
00400 
00401     return ray(pt, dir);
00402 }
00403 
00404 
00405 // Cutting point of three planes
00406 Foam::point Foam::plane::planePlaneIntersect
00407 (
00408     const plane& plane2,
00409     const plane& plane3
00410 ) const
00411 {
00412     FixedList<scalar, 4> coeffs1(planeCoeffs());
00413     FixedList<scalar, 4> coeffs2(plane2.planeCoeffs());
00414     FixedList<scalar, 4> coeffs3(plane3.planeCoeffs());
00415 
00416     tensor a
00417     (
00418         coeffs1[0],coeffs1[1],coeffs1[2],
00419         coeffs2[0],coeffs2[1],coeffs2[2],
00420         coeffs3[0],coeffs3[1],coeffs3[2]
00421     );
00422 
00423     vector b(coeffs1[3],coeffs2[3],coeffs3[3]);
00424 
00425     return (inv(a) & (-b));
00426 }
00427 
00428 
00429 void Foam::plane::writeDict(Ostream& os) const
00430 {
00431     os.writeKeyword("planeType") << "pointAndNormal"
00432         << token::END_STATEMENT << nl;
00433     os  << indent << "pointAndNormalDict" << nl
00434         << indent << token::BEGIN_BLOCK << incrIndent << nl;
00435     os.writeKeyword("basePoint") << basePoint_ << token::END_STATEMENT << nl;
00436     os.writeKeyword("normalVector") << unitVector_ << token::END_STATEMENT
00437         << nl;
00438     os << decrIndent << indent << token::END_BLOCK << endl;
00439 }
00440 
00441 
00442 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
00443 
00444 bool Foam::operator==(const plane& a, const plane& b)
00445 {
00446     if (a.basePoint_ == b.basePoint_ && a.unitVector_ == b.unitVector_)
00447     {
00448         return true;
00449     }
00450     else
00451     {
00452         return false;
00453     }
00454 }
00455 
00456 bool Foam::operator!=(const plane& a, const plane& b)
00457 {
00458     return !(a == b);
00459 }
00460 
00461 
00462 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
00463 
00464 Foam::Ostream& Foam::operator<<(Ostream& os, const plane& a)
00465 {
00466     os  << a.unitVector_ << token::SPACE << a.basePoint_;
00467 
00468     return os;
00469 }
00470 
00471 
00472 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines