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

treeDataCell.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 "treeDataCell.H"
00027 #include "indexedOctree.H"
00028 #include <OpenFOAM/primitiveMesh.H>
00029 
00030 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00031 
00032 defineTypeNameAndDebug(Foam::treeDataCell, 0);
00033 
00034 
00035 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00036 
00037 Foam::treeBoundBox Foam::treeDataCell::calcCellBb(const label cellI) const
00038 {
00039     const cellList& cells = mesh_.cells();
00040     const faceList& faces = mesh_.faces();
00041     const pointField& points = mesh_.points();
00042 
00043     treeBoundBox cellBb
00044     (
00045         vector(GREAT, GREAT, GREAT),
00046         vector(-GREAT, -GREAT, -GREAT)
00047     );
00048 
00049     const cell& cFaces = cells[cellI];
00050 
00051     forAll(cFaces, cFaceI)
00052     {
00053         const face& f = faces[cFaces[cFaceI]];
00054 
00055         forAll(f, fp)
00056         {
00057             const point& p = points[f[fp]];
00058 
00059             cellBb.min() = min(cellBb.min(), p);
00060             cellBb.max() = max(cellBb.max(), p);
00061         }
00062     }
00063     return cellBb;
00064 }
00065 
00066 
00067 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00068 
00069 // Construct from components
00070 Foam::treeDataCell::treeDataCell
00071 (
00072     const bool cacheBb,
00073     const primitiveMesh& mesh,
00074     const labelList& cellLabels
00075 )
00076 :
00077     mesh_(mesh),
00078     cellLabels_(cellLabels),
00079     cacheBb_(cacheBb)
00080 {
00081     if (cacheBb_)
00082     {
00083         bbs_.setSize(cellLabels_.size());
00084 
00085         forAll(cellLabels_, i)
00086         {
00087             bbs_[i] = calcCellBb(cellLabels_[i]);
00088         }
00089     }
00090 }
00091 
00092 
00093 Foam::treeDataCell::treeDataCell
00094 (
00095     const bool cacheBb,
00096     const primitiveMesh& mesh
00097 )
00098 :
00099     mesh_(mesh),
00100     cellLabels_(identity(mesh_.nCells())),
00101     cacheBb_(cacheBb)
00102 {
00103     if (cacheBb_)
00104     {
00105         bbs_.setSize(cellLabels_.size());
00106 
00107         forAll(cellLabels_, i)
00108         {
00109             bbs_[i] = calcCellBb(cellLabels_[i]);
00110         }
00111     }
00112 }
00113 
00114 
00115 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00116 
00117 Foam::pointField Foam::treeDataCell::points() const
00118 {
00119     pointField cc(cellLabels_.size());
00120 
00121     forAll(cellLabels_, i)
00122     {
00123         cc[i] = mesh_.cellCentres()[cellLabels_[i]];
00124     }
00125 
00126     return cc;
00127 }
00128 
00129 
00130 // Check if any point on shape is inside cubeBb.
00131 bool Foam::treeDataCell::overlaps
00132 (
00133     const label index,
00134     const treeBoundBox& cubeBb
00135 ) const
00136 {
00137     if (cacheBb_)
00138     {
00139         return cubeBb.overlaps(bbs_[index]);
00140     }
00141     else
00142     {
00143         return cubeBb.overlaps(calcCellBb(cellLabels_[index]));
00144     }
00145 }
00146 
00147 
00148 // Calculate nearest point to sample. Updates (if any) nearestDistSqr, minIndex,
00149 // nearestPoint.
00150 void Foam::treeDataCell::findNearest
00151 (
00152     const labelList& indices,
00153     const point& sample,
00154 
00155     scalar& nearestDistSqr,
00156     label& minIndex,
00157     point& nearestPoint
00158 ) const
00159 {
00160     forAll(indices, i)
00161     {
00162         label index = indices[i];
00163         label cellI = cellLabels_[index];
00164         scalar distSqr = magSqr(sample - mesh_.cellCentres()[cellI]);
00165 
00166         if (distSqr < nearestDistSqr)
00167         {
00168             nearestDistSqr = distSqr;
00169             minIndex = index;
00170             nearestPoint = mesh_.cellCentres()[cellI];
00171         }
00172     }
00173 }
00174 
00175 
00176 bool Foam::treeDataCell::intersects
00177 (
00178     const label index,
00179     const point& start,
00180     const point& end,
00181     point& intersectionPoint
00182 ) const
00183 {
00184     // Do quick rejection test
00185     if (cacheBb_)
00186     {
00187         const treeBoundBox& cellBb = bbs_[index];
00188 
00189         if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
00190         {
00191             // start and end in same block outside of cellBb.
00192             return false;
00193         }
00194     }
00195     else
00196     {
00197         const treeBoundBox cellBb = calcCellBb(cellLabels_[index]);
00198 
00199         if ((cellBb.posBits(start) & cellBb.posBits(end)) != 0)
00200         {
00201             // start and end in same block outside of cellBb.
00202             return false;
00203         }
00204     }
00205 
00206 
00207     // Do intersection with all faces of cell
00208     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00209 
00210     // Disable picking up intersections behind us.
00211     scalar oldTol = intersection::setPlanarTol(0.0);
00212 
00213     const cell& cFaces = mesh_.cells()[cellLabels_[index]];
00214 
00215     const vector dir(end - start);
00216     scalar minDistSqr = magSqr(dir);
00217     bool hasMin = false;
00218 
00219     forAll(cFaces, i)
00220     {
00221         const face& f = mesh_.faces()[cFaces[i]];
00222 
00223         pointHit inter = f.ray
00224         (
00225             start,
00226             dir,
00227             mesh_.points(),
00228             intersection::HALF_RAY
00229         );
00230 
00231         if (inter.hit() && sqr(inter.distance()) <= minDistSqr)
00232         {
00233             // Note: no extra test on whether intersection is in front of us
00234             // since using half_ray AND zero tolerance. (note that tolerance
00235             // is used to look behind us)
00236             minDistSqr = sqr(inter.distance());
00237             intersectionPoint = inter.hitPoint();
00238             hasMin = true;
00239         }
00240     }
00241 
00242     // Restore picking tolerance
00243     intersection::setPlanarTol(oldTol);
00244 
00245     return hasMin;
00246 }
00247 
00248 
00249 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines