initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / mesh / advanced / selectCells / edgeStats.C
blob2bc8c249169e64ba657a608b59a5df5041ca0f03
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "edgeStats.H"
28 #include "Time.H"
29 #include "polyMesh.H"
30 #include "Ostream.H"
31 #include "twoDPointCorrector.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 const Foam::scalar Foam::edgeStats::edgeTol_ = 1E-3;
39 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
41 Foam::direction Foam::edgeStats::getNormalDir
43     const twoDPointCorrector* correct2DPtr
44 ) const
46     direction dir = 3;
48     if (correct2DPtr)
49     {
50         const vector& normal = correct2DPtr->planeNormal();
52         if (mag(normal & vector(1, 0, 0)) > 1-edgeTol_)
53         {
54             dir = 0;
55         }
56         else if (mag(normal & vector(0, 1, 0)) > 1-edgeTol_)
57         {
58             dir = 1;
59         }
60         else if (mag(normal & vector(0, 0, 1)) > 1-edgeTol_)
61         {
62             dir = 2;
63         }
64     }
65     return dir;
69 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
71 // Construct from mesh
72 Foam::edgeStats::edgeStats(const polyMesh& mesh)
74     mesh_(mesh),
75     normalDir_(3)
77     IOobject motionObj
78     (
79         "motionProperties",
80         mesh.time().constant(),
81         mesh,
82         IOobject::MUST_READ,
83         IOobject::NO_WRITE
84     );
86     if (motionObj.headerOk())
87     {
88         Info<< "Reading " << mesh.time().constant() / "motionProperties"
89             << endl << endl;
91         IOdictionary motionProperties(motionObj);
93         Switch twoDMotion(motionProperties.lookup("twoDMotion"));
95         if (twoDMotion)
96         {
97             Info<< "Correcting for 2D motion" << endl << endl;
99             autoPtr<twoDPointCorrector> correct2DPtr
100             (
101                 new twoDPointCorrector(mesh)
102             );
104             normalDir_ = getNormalDir(&correct2DPtr());
105         }
106     }
110 // Construct from components
111 Foam::edgeStats::edgeStats
113     const polyMesh& mesh,
114     const twoDPointCorrector* correct2DPtr
117     mesh_(mesh),
118     normalDir_(getNormalDir(correct2DPtr))
122 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
124 Foam::scalar Foam::edgeStats::minLen(Ostream& os) const
126     label nX = 0;
127     label nY = 0;
128     label nZ = 0;
130     scalar minX = GREAT;
131     scalar maxX = -GREAT;
132     vector x(1, 0, 0);
134     scalar minY = GREAT;
135     scalar maxY = -GREAT;
136     vector y(0, 1, 0);
138     scalar minZ = GREAT;
139     scalar maxZ = -GREAT;
140     vector z(0, 0, 1);
142     scalar minOther = GREAT;
143     scalar maxOther = -GREAT;
145     const edgeList& edges = mesh_.edges();
147     forAll(edges, edgeI)
148     {
149         const edge& e = edges[edgeI];
151         vector eVec(e.vec(mesh_.points()));
153         scalar eMag = mag(eVec);
155         eVec /= eMag;
157         if (mag(eVec & x) > 1-edgeTol_)
158         {
159             minX = min(minX, eMag);
160             maxX = max(maxX, eMag);
161             nX++;
162         }
163         else if (mag(eVec & y) > 1-edgeTol_)
164         {
165             minY = min(minY, eMag);
166             maxY = max(maxY, eMag);
167             nY++;
168         }
169         else if (mag(eVec & z) > 1-edgeTol_)
170         {
171             minZ = min(minZ, eMag);
172             maxZ = max(maxZ, eMag);
173             nZ++;
174         }
175         else
176         {
177             minOther = min(minOther, eMag);
178             maxOther = max(maxOther, eMag);
179         }
180     }
182     os  << "Mesh bounding box:" << boundBox(mesh_.points()) << nl << nl
183         << "Mesh edge statistics:" << nl
184         << "    x aligned :  number:" << nX << "\tminLen:" << minX
185         << "\tmaxLen:" << maxX << nl
186         << "    y aligned :  number:" << nY << "\tminLen:" << minY
187         << "\tmaxLen:" << maxY << nl
188         << "    z aligned :  number:" << nZ << "\tminLen:" << minZ
189         << "\tmaxLen:" << maxZ << nl
190         << "    other     :  number:" << mesh_.nEdges() - nX - nY - nZ
191         << "\tminLen:" << minOther
192         << "\tmaxLen:" << maxOther << nl << endl;
194     if (normalDir_ == 0)
195     {
196         return min(minY, min(minZ, minOther));
197     }
198     else if (normalDir_ == 1)
199     {
200         return min(minX, min(minZ, minOther));
201     }
202     else if (normalDir_ == 2)
203     {
204         return min(minX, min(minY, minOther));
205     }
206     else
207     {
208         return min(minX, min(minY, min(minZ, minOther)));
209     }
213 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
216 // * * * * * * * * * * * * * * * Friend Functions  * * * * * * * * * * * * * //
219 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
222 // ************************************************************************* //