Loading...
Searching...
No Matches
edgeStats.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "edgeStats.H"
29#include "Time.H"
30#include "polyMesh.H"
31#include "Ostream.H"
32#include "twoDPointCorrector.H"
33#include "IOdictionary.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37const Foam::scalar Foam::edgeStats::edgeTol_ = 1e-3;
38
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42Foam::direction Foam::edgeStats::getNormalDir
43(
44 const twoDPointCorrector* correct2DPtr
45) const
46{
47 if (correct2DPtr)
48 {
49 const vector& normal = correct2DPtr->planeNormal();
50
51 if (mag(normal.x()) > 1-edgeTol_)
52 {
53 return vector::X;
54 }
55 else if (mag(normal.y()) > 1-edgeTol_)
56 {
57 return vector::Y;
58 }
59 else if (mag(normal.z()) > 1-edgeTol_)
60 {
61 return vector::Z;
62 }
63 }
64
65 return direction(3);
66}
67
68
69// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70
71Foam::edgeStats::edgeStats(const polyMesh& mesh)
72:
73 mesh_(mesh),
74 normalDir_(3)
75{
76 IOobject motionObj
77 (
78 "motionProperties",
79 mesh.time().constant(),
80 mesh,
81 IOobject::MUST_READ_IF_MODIFIED,
82 IOobject::NO_WRITE
83 );
84
85 if (motionObj.typeHeaderOk<IOdictionary>(true))
86 {
87 Info<< "Reading " << mesh.time().constant() / "motionProperties"
88 << endl << endl;
89
90 IOdictionary motionProperties(motionObj);
91
92 if (motionProperties.get<bool>("twoDMotion"))
93 {
94 Info<< "Correcting for 2D motion" << endl << endl;
95
96 twoDPointCorrector correct2D(mesh);
97
98 normalDir_ = getNormalDir(&correct2D);
99 }
100 }
101}
102
103
104// Construct from components
105Foam::edgeStats::edgeStats
106(
107 const polyMesh& mesh,
108 const twoDPointCorrector* correct2DPtr
109)
110:
111 mesh_(mesh),
112 normalDir_(getNormalDir(correct2DPtr))
113{}
114
115
116// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117
118Foam::scalar Foam::edgeStats::minLen(Ostream& os) const
119{
120 label nAny(0);
121 label nX(0);
122 label nY(0);
123 label nZ(0);
124
125 scalarMinMax limitsAny(GREAT, -GREAT);
126 scalarMinMax limitsX(limitsAny);
127 scalarMinMax limitsY(limitsAny);
128 scalarMinMax limitsZ(limitsAny);
129
130 const edgeList& edges = mesh_.edges();
131
132 for (const edge& e : edges)
133 {
134 vector eVec(e.vec(mesh_.points()));
135
136 scalar eMag = mag(eVec);
137
138 eVec /= eMag;
139
140 if (mag(eVec.x()) > 1-edgeTol_)
141 {
142 limitsX.add(eMag);
143 nX++;
144 }
145 else if (mag(eVec.y()) > 1-edgeTol_)
146 {
147 limitsY.add(eMag);
148 nY++;
149 }
150 else if (mag(eVec.z()) > 1-edgeTol_)
151 {
152 limitsZ.add(eMag);
153 nZ++;
154 }
155 else
156 {
157 limitsAny.add(eMag);
158 nAny++;
159 }
160 }
161
162 os << "Mesh bounding box:" << boundBox(mesh_.points()) << nl << nl
163 << "Mesh edge statistics:" << nl
164 << " x aligned : number:" << nX
165 << "\tminLen:" << limitsX.min() << "\tmaxLen:" << limitsX.max() << nl
166 << " y aligned : number:" << nY
167 << "\tminLen:" << limitsY.min() << "\tmaxLen:" << limitsY.max() << nl
168 << " z aligned : number:" << nZ
169 << "\tminLen:" << limitsZ.min() << "\tmaxLen:" << limitsZ.max() << nl
170 << " other : number:" << nAny
171 << "\tminLen:" << limitsAny.min()
172 << "\tmaxLen:" << limitsAny.max() << nl << endl;
173
174 if (normalDir_ == vector::X)
175 {
176 return Foam::min
177 (
178 limitsAny.min(),
179 Foam::min(limitsY.min(), limitsZ.min())
180 );
181 }
182 else if (normalDir_ == vector::Y)
183 {
184 return Foam::min
185 (
186 limitsAny.min(),
187 Foam::min(limitsX.min(), limitsZ.min())
188 );
189 }
190 else if (normalDir_ == vector::Z)
191 {
192 return Foam::min
193 (
194 limitsAny.min(),
195 Foam::min(limitsX.min(), limitsY.min())
196 );
197 }
198 else
199 {
200 return Foam::min
201 (
202 limitsAny.min(),
204 (
205 limitsX.min(),
206 Foam::min(limitsY.min(), limitsZ.min())
207 )
208 );
209 }
210}
211
212
213// ************************************************************************* //
static const scalar edgeTol_
Definition edgeStats.H:89
scalar minLen(Ostream &os) const
Calculate minimum edge length and print.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
List< edge > edgeList
List of edge.
Definition edgeList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & e