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

wallPointI.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 
00029 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00030 
00031 // Update this with w2 if w2 nearer to pt.
00032 inline bool Foam::wallPoint::update
00033 (
00034     const point& pt,
00035     const wallPoint& w2,
00036     const scalar tol
00037 )
00038 {
00039     //Already done in calling algorithm
00040     //if (w2.origin() == origin_)
00041     //{
00042     //    // Shortcut. Same input so same distance.
00043     //    return false;
00044     //}
00045 
00046     scalar dist2 = magSqr(pt - w2.origin());
00047 
00048     if (!valid())
00049     {
00050         // current not yet set so use any value
00051         distSqr_ = dist2;
00052         origin_ = w2.origin();
00053 
00054         return true;
00055     }        
00056 
00057     scalar diff = distSqr_ - dist2;
00058 
00059     if (diff < 0)
00060     {
00061         // already nearer to pt
00062         return false;
00063     }
00064 
00065     if ((diff < SMALL) || ((distSqr_ > SMALL) && (diff/distSqr_ < tol)))
00066     {
00067         // don't propagate small changes
00068         return false;
00069     }
00070     else
00071     {
00072         // update with new values
00073         distSqr_ = dist2;
00074         origin_ = w2.origin();
00075 
00076         return true;
00077     }
00078 }
00079     
00080 
00081 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00082 
00083 // Null constructor
00084 inline Foam::wallPoint::wallPoint()
00085 :
00086     origin_(greatPoint),
00087     distSqr_(GREAT)
00088 {}
00089 
00090 
00091 // Construct from origin, distance
00092 inline Foam::wallPoint::wallPoint(const point& origin, const scalar distSqr)
00093 :
00094     origin_(origin), distSqr_(distSqr)
00095 {}
00096 
00097 
00098 // Construct as copy
00099 inline Foam::wallPoint::wallPoint(const wallPoint& wpt)
00100 :
00101     origin_(wpt.origin()), distSqr_(wpt.distSqr())
00102 {}
00103 
00104 
00105 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00106 
00107 inline const Foam::point& Foam::wallPoint::origin() const
00108 {
00109     return origin_;
00110 }
00111 
00112 
00113 inline Foam::point& Foam::wallPoint::origin()
00114 {
00115     return origin_;
00116 }
00117 
00118 
00119 inline Foam::scalar Foam::wallPoint::distSqr() const
00120 {
00121     return distSqr_;
00122 }
00123 
00124 
00125 inline Foam::scalar& Foam::wallPoint::distSqr()
00126 {
00127     return distSqr_;
00128 }
00129 
00130 
00131 inline bool Foam::wallPoint::valid() const
00132 {
00133     return origin_ != greatPoint;
00134 }
00135 
00136 
00137 // Checks for cyclic faces
00138 inline bool Foam::wallPoint::sameGeometry
00139 (
00140     const polyMesh&,
00141     const wallPoint& w2,
00142     const scalar tol
00143 )
00144  const
00145 {
00146     scalar diff = mag(distSqr() - w2.distSqr());
00147 
00148     if (diff < SMALL)
00149     {
00150         return true;
00151     }
00152     else
00153     {
00154         if ((distSqr() > SMALL) && ((diff/distSqr()) < tol))
00155         {
00156             return true;
00157         }
00158         else
00159         {
00160             return false;
00161         }
00162     }
00163 }
00164 
00165 
00166 inline void Foam::wallPoint::leaveDomain
00167 (
00168     const polyMesh&,
00169     const polyPatch&,
00170     const label,
00171     const point& faceCentre
00172 )
00173 {
00174     origin_ -= faceCentre;
00175 }
00176 
00177 
00178 inline void Foam::wallPoint::transform
00179 (
00180     const polyMesh&,
00181     const tensor& rotTensor
00182 )
00183 {
00184     origin_ = Foam::transform(rotTensor, origin_);
00185 }
00186 
00187 
00188 // Update absolute geometric quantities. Note that distance (distSqr_)
00189 // is not affected by leaving/entering domain.
00190 inline void Foam::wallPoint::enterDomain
00191 (
00192     const polyMesh&,
00193     const polyPatch&,
00194     const label,
00195     const point& faceCentre
00196 )
00197 {
00198     // back to absolute form
00199     origin_ += faceCentre;
00200 }
00201 
00202 
00203 // Update this with w2 if w2 nearer to pt.
00204 inline bool Foam::wallPoint::updateCell
00205 (
00206     const polyMesh& mesh,
00207     const label thisCellI,
00208     const label neighbourFaceI,
00209     const wallPoint& neighbourWallInfo,
00210     const scalar tol
00211 )
00212 {
00213     return 
00214         update
00215         (
00216             mesh.cellCentres()[thisCellI],
00217             neighbourWallInfo,
00218             tol
00219         );
00220     }    
00221 
00222 
00223 // Update this with w2 if w2 nearer to pt.
00224 inline bool Foam::wallPoint::updateFace
00225 (
00226     const polyMesh& mesh,
00227     const label thisFaceI,
00228     const label neighbourCellI,
00229     const wallPoint& neighbourWallInfo,
00230     const scalar tol
00231 )
00232 {
00233     return
00234         update
00235         (
00236             mesh.faceCentres()[thisFaceI],
00237             neighbourWallInfo,
00238             tol
00239         );
00240 }    
00241 
00242 // Update this with w2 if w2 nearer to pt.
00243 inline bool Foam::wallPoint::updateFace
00244 (
00245     const polyMesh& mesh,
00246     const label thisFaceI,
00247     const wallPoint& neighbourWallInfo,
00248     const scalar tol
00249 )
00250 {
00251     return
00252         update
00253         (
00254             mesh.faceCentres()[thisFaceI],
00255             neighbourWallInfo,
00256             tol
00257         );
00258 }    
00259 
00260 
00261 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00262 
00263 inline bool Foam::wallPoint::operator==(const Foam::wallPoint& rhs) const
00264 {
00265     return origin() == rhs.origin();
00266 }
00267 
00268 
00269 inline bool Foam::wallPoint::operator!=(const Foam::wallPoint& rhs) const
00270 {
00271     return !(*this == rhs);
00272 }
00273 
00274 
00275 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines