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

edgeStats.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 "edgeStats.H"
00027 #include <OpenFOAM/Time.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <OpenFOAM/Ostream.H>
00030 #include <meshTools/twoDPointCorrector.H>
00031 
00032 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00033 
00034 const Foam::scalar Foam::edgeStats::edgeTol_ = 1E-3;
00035 
00036 
00037 
00038 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00039 
00040 Foam::direction Foam::edgeStats::getNormalDir
00041 (
00042     const twoDPointCorrector* correct2DPtr
00043 ) const
00044 {
00045     direction dir = 3;
00046 
00047     if (correct2DPtr)
00048     {
00049         const vector& normal = correct2DPtr->planeNormal();
00050 
00051         if (mag(normal & vector(1, 0, 0)) > 1-edgeTol_)
00052         {
00053             dir = 0;
00054         }
00055         else if (mag(normal & vector(0, 1, 0)) > 1-edgeTol_)
00056         {
00057             dir = 1;
00058         }
00059         else if (mag(normal & vector(0, 0, 1)) > 1-edgeTol_)
00060         {
00061             dir = 2;
00062         }
00063     }
00064     return dir;
00065 }
00066 
00067 
00068 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00069 
00070 // Construct from mesh
00071 Foam::edgeStats::edgeStats(const polyMesh& mesh)
00072 :
00073     mesh_(mesh),
00074     normalDir_(3)
00075 {
00076     IOobject motionObj
00077     (
00078         "motionProperties",
00079         mesh.time().constant(),
00080         mesh,
00081         IOobject::MUST_READ,
00082         IOobject::NO_WRITE
00083     );
00084 
00085     if (motionObj.headerOk())
00086     {
00087         Info<< "Reading " << mesh.time().constant() / "motionProperties"
00088             << endl << endl;
00089 
00090         IOdictionary motionProperties(motionObj);
00091 
00092         Switch twoDMotion(motionProperties.lookup("twoDMotion"));
00093 
00094         if (twoDMotion)
00095         {
00096             Info<< "Correcting for 2D motion" << endl << endl;
00097 
00098             autoPtr<twoDPointCorrector> correct2DPtr
00099             (
00100                 new twoDPointCorrector(mesh)
00101             );
00102 
00103             normalDir_ = getNormalDir(&correct2DPtr());
00104         }
00105     }
00106 }
00107 
00108 
00109 // Construct from components
00110 Foam::edgeStats::edgeStats
00111 (
00112     const polyMesh& mesh,
00113     const twoDPointCorrector* correct2DPtr
00114 )
00115 :
00116     mesh_(mesh),
00117     normalDir_(getNormalDir(correct2DPtr))
00118 {}
00119 
00120 
00121 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00122 
00123 Foam::scalar Foam::edgeStats::minLen(Ostream& os) const
00124 {
00125     label nX = 0;
00126     label nY = 0;
00127     label nZ = 0;
00128 
00129     scalar minX = GREAT;
00130     scalar maxX = -GREAT;
00131     vector x(1, 0, 0);
00132 
00133     scalar minY = GREAT;
00134     scalar maxY = -GREAT;
00135     vector y(0, 1, 0);
00136 
00137     scalar minZ = GREAT;
00138     scalar maxZ = -GREAT;
00139     vector z(0, 0, 1);
00140 
00141     scalar minOther = GREAT;
00142     scalar maxOther = -GREAT;
00143 
00144     const edgeList& edges = mesh_.edges();
00145 
00146     forAll(edges, edgeI)
00147     {
00148         const edge& e = edges[edgeI];
00149 
00150         vector eVec(e.vec(mesh_.points()));
00151 
00152         scalar eMag = mag(eVec);
00153 
00154         eVec /= eMag;
00155 
00156         if (mag(eVec & x) > 1-edgeTol_)
00157         {
00158             minX = min(minX, eMag);
00159             maxX = max(maxX, eMag);
00160             nX++;
00161         }
00162         else if (mag(eVec & y) > 1-edgeTol_)
00163         {
00164             minY = min(minY, eMag);
00165             maxY = max(maxY, eMag);
00166             nY++;
00167         }
00168         else if (mag(eVec & z) > 1-edgeTol_)
00169         {
00170             minZ = min(minZ, eMag);
00171             maxZ = max(maxZ, eMag);
00172             nZ++;
00173         }
00174         else
00175         {
00176             minOther = min(minOther, eMag);
00177             maxOther = max(maxOther, eMag);
00178         }
00179     }
00180 
00181     os  << "Mesh bounding box:" << boundBox(mesh_.points()) << nl << nl
00182         << "Mesh edge statistics:" << nl
00183         << "    x aligned :  number:" << nX << "\tminLen:" << minX
00184         << "\tmaxLen:" << maxX << nl
00185         << "    y aligned :  number:" << nY << "\tminLen:" << minY
00186         << "\tmaxLen:" << maxY << nl
00187         << "    z aligned :  number:" << nZ << "\tminLen:" << minZ
00188         << "\tmaxLen:" << maxZ << nl
00189         << "    other     :  number:" << mesh_.nEdges() - nX - nY - nZ
00190         << "\tminLen:" << minOther
00191         << "\tmaxLen:" << maxOther << nl << endl;
00192 
00193     if (normalDir_ == 0)
00194     {
00195         return min(minY, min(minZ, minOther));
00196     }
00197     else if (normalDir_ == 1)
00198     {
00199         return min(minX, min(minZ, minOther));
00200     }
00201     else if (normalDir_ == 2)
00202     {
00203         return min(minX, min(minY, minOther));
00204     }
00205     else
00206     {
00207         return min(minX, min(minY, min(minZ, minOther)));
00208     }
00209 }
00210 
00211 
00212 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00213 
00214 
00215 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
00216 
00217 
00218 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
00219 
00220 
00221 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines