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

pointDataI.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/polyMesh.H>
00027 #include <OpenFOAM/transform.H>
00028 #include <meshTools/wallPoint.H>
00029 
00030 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00031 
00032 // Update this with w2 if w2 nearer to pt.
00033 inline bool Foam::pointData::update
00034 (
00035     const point& pt,
00036     const pointData& w2,
00037     const scalar tol
00038 )
00039 {
00040     scalar dist2 = magSqr(pt - w2.origin());
00041 
00042     if (!valid())
00043     {
00044         distSqr_ = dist2;
00045         origin_ = w2.origin();
00046         s_ = w2.s();
00047         v_ = w2.v();
00048 
00049         return true;
00050     }
00051 
00052 
00053 //    if (v_ != w2.v())
00054 //    {
00055 //        return false;
00056 //    }
00057 
00058 
00059     scalar diff = distSqr_ - dist2;
00060 
00061     if (diff < 0)
00062     {
00063         // already nearer to pt
00064         return false;
00065     }
00066 
00067     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
00068     {
00069         // don't propagate small changes
00070         return false;
00071     }
00072     else
00073     {
00074         // update with new values
00075         distSqr_ = dist2;
00076         origin_ = w2.origin();
00077         s_ = w2.s();
00078         v_ = w2.v();
00079 
00080         return true;
00081     }
00082 }
00083 
00084 
00085 // Update this with w2 (information on same point)
00086 inline bool Foam::pointData::update
00087 (
00088     const pointData& w2,
00089     const scalar tol
00090 )
00091 {
00092     if (!valid())
00093     {
00094         // current not yet set so use any value
00095         distSqr_ = w2.distSqr();
00096         origin_ = w2.origin();
00097         s_ = w2.s();
00098         v_ = w2.v();
00099         return true;
00100     }
00101 
00102 
00103 //    if (v_ != w2.v())
00104 //    {
00105 //        return false;
00106 //    }
00107 
00108 
00109     scalar diff = distSqr_ - w2.distSqr();
00110 
00111     if (diff < 0)
00112     {
00113         // already nearer to pt
00114         return false;
00115     }
00116 
00117     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
00118     {
00119         // don't propagate small changes
00120         return false;
00121     }
00122     else
00123     {
00124         // update with new values
00125         distSqr_ =  w2.distSqr();
00126         origin_ = w2.origin();
00127         s_ = w2.s();
00128         v_ = w2.v();
00129 
00130         return true;
00131     }
00132 }
00133 
00134 
00135 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00136 
00137 // Null constructor
00138 inline Foam::pointData::pointData()
00139 :
00140     origin_(wallPoint::greatPoint),
00141     distSqr_(GREAT),
00142     s_(GREAT),
00143     v_(wallPoint::greatPoint)
00144 {}
00145 
00146 
00147 // Construct from origin, distance
00148 inline Foam::pointData::pointData
00149 (
00150     const point& origin,
00151     const scalar distSqr,
00152     const scalar s,
00153     const vector& v
00154 )
00155 :
00156     origin_(origin),
00157     distSqr_(distSqr),
00158     s_(s),
00159     v_(v)
00160 {}
00161 
00162 
00163 // Construct as copy
00164 inline Foam::pointData::pointData(const pointData& wpt)
00165 :
00166     origin_(wpt.origin()),
00167     distSqr_(wpt.distSqr()),
00168     s_(wpt.s()),
00169     v_(wpt.v())
00170 {}
00171 
00172 
00173 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00174 
00175 inline const Foam::point& Foam::pointData::origin() const
00176 {
00177     return origin_;
00178 }
00179 
00180 
00181 inline Foam::scalar Foam::pointData::distSqr() const
00182 {
00183     return distSqr_;
00184 }
00185 
00186 
00187 inline Foam::scalar Foam::pointData::s() const
00188 {
00189     return s_;
00190 }
00191 
00192 
00193 inline const Foam::vector& Foam::pointData::v() const
00194 {
00195     return v_;
00196 }
00197 
00198 
00199 inline bool Foam::pointData::valid() const
00200 {
00201     return origin_ != wallPoint::greatPoint;
00202 }
00203 
00204 
00205 // Checks for cyclic points
00206 inline bool Foam::pointData::sameGeometry
00207 (
00208     const pointData& w2,
00209     const scalar tol
00210 ) const
00211 {
00212     scalar diff = Foam::mag(distSqr() - w2.distSqr());
00213 
00214     if (diff < SMALL)
00215     {
00216         return true;
00217     }
00218     else
00219     {
00220         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
00221         {
00222             return true;
00223         }
00224         else
00225         {
00226             return false;
00227         }
00228     }
00229 }
00230 
00231 
00232 inline void Foam::pointData::leaveDomain
00233 (
00234     const polyPatch& patch,
00235     const label patchPointI,
00236     const point& coord
00237 )
00238 {
00239     origin_ -= coord;
00240 }
00241 
00242 
00243 inline void Foam::pointData::transform(const tensor& rotTensor)
00244 {
00245     origin_ = Foam::transform(rotTensor, origin_);
00246 }
00247 
00248 
00249 // Update absolute geometric quantities. Note that distance (distSqr_)
00250 // is not affected by leaving/entering domain.
00251 inline void Foam::pointData::enterDomain
00252 (
00253     const polyPatch& patch,
00254     const label patchPointI,
00255     const point& coord
00256 )
00257 {
00258     // back to absolute form
00259     origin_ += coord;
00260 }
00261 
00262 
00263 // Update this with information from connected edge
00264 inline bool Foam::pointData::updatePoint
00265 (
00266     const polyMesh& mesh,
00267     const label pointI,
00268     const label edgeI,
00269     const pointData& edgeInfo,
00270     const scalar tol
00271 )
00272 {
00273     return
00274         update
00275         (
00276             mesh.points()[pointI],
00277             edgeInfo,
00278             tol
00279         );
00280     }
00281 
00282 
00283 // Update this with new information on same point
00284 inline bool Foam::pointData::updatePoint
00285 (
00286     const polyMesh& mesh,
00287     const label pointI,
00288     const pointData& newPointInfo,
00289     const scalar tol
00290 )
00291 {
00292     return
00293         update
00294         (
00295             mesh.points()[pointI],
00296             newPointInfo,
00297             tol
00298         );
00299 }
00300 
00301 
00302 // Update this with new information on same point. No extra information.
00303 inline bool Foam::pointData::updatePoint
00304 (
00305     const pointData& newPointInfo,
00306     const scalar tol
00307 )
00308 {
00309     return update(newPointInfo, tol);
00310 }
00311 
00312 
00313 // Update this with information from connected point
00314 inline bool Foam::pointData::updateEdge
00315 (
00316     const polyMesh& mesh,
00317     const label edgeI,
00318     const label pointI,
00319     const pointData& pointInfo,
00320     const scalar tol
00321 )
00322 {
00323     const pointField& points = mesh.points();
00324 
00325     const edge& e = mesh.edges()[edgeI];
00326 
00327     const point edgeMid(0.5*(points[e[0]] + points[e[1]]));
00328 
00329     return
00330         update
00331         (
00332             edgeMid,
00333             pointInfo,
00334             tol
00335         );
00336 }
00337 
00338 
00339 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00340 
00341 inline bool Foam::pointData::operator==(const pointData& rhs) const
00342 {
00343     return origin() == rhs.origin();
00344 }
00345 
00346 
00347 inline bool Foam::pointData::operator!=(const pointData& rhs) const
00348 {
00349     return !(*this == rhs);
00350 }
00351 
00352 
00353 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines