edgeStats.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "edgeStats.H"
27 #include "Time.H"
28 #include "polyMesh.H"
29 #include "Ostream.H"
30 #include "twoDPointCorrector.H"
31 #include "IOdictionary.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 const Foam::scalar Foam::edgeStats::edgeTol_ = 1e-3;
36 
37 
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 Foam::direction Foam::edgeStats::getNormalDir
42 (
43  const twoDPointCorrector* correct2DPtr
44 ) const
45 {
46  direction dir = 3;
47 
48  if (correct2DPtr)
49  {
50  const vector& normal = correct2DPtr->planeNormal();
51 
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;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 // Construct from mesh
72 Foam::edgeStats::edgeStats(const polyMesh& mesh)
73 :
74  mesh_(mesh),
75  normalDir_(3)
76 {
77  IOobject motionObj
78  (
79  "motionProperties",
80  mesh.time().constant(),
81  mesh,
84  );
85 
86  if (motionObj.typeHeaderOk<IOdictionary>(true))
87  {
88  Info<< "Reading " << mesh.time().constant() / "motionProperties"
89  << endl << endl;
90 
91  IOdictionary motionProperties(motionObj);
92 
93  Switch twoDMotion(motionProperties.lookup("twoDMotion"));
94 
95  if (twoDMotion)
96  {
97  Info<< "Correcting for 2D motion" << endl << endl;
98 
99  autoPtr<twoDPointCorrector> correct2DPtr
100  (
101  new twoDPointCorrector(mesh)
102  );
103 
104  normalDir_ = getNormalDir(&correct2DPtr());
105  }
106  }
107 }
108 
109 
110 // Construct from components
112 (
113  const polyMesh& mesh,
114  const twoDPointCorrector* correct2DPtr
115 )
116 :
117  mesh_(mesh),
118  normalDir_(getNormalDir(correct2DPtr))
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
124 Foam::scalar Foam::edgeStats::minLen(Ostream& os) const
125 {
126  label nX = 0;
127  label nY = 0;
128  label nZ = 0;
129 
130  scalar minX = great;
131  scalar maxX = -great;
132  vector x(1, 0, 0);
133 
134  scalar minY = great;
135  scalar maxY = -great;
136  vector y(0, 1, 0);
137 
138  scalar minZ = great;
139  scalar maxZ = -great;
140  vector z(0, 0, 1);
141 
142  scalar minOther = great;
143  scalar maxOther = -great;
144 
145  const edgeList& edges = mesh_.edges();
146 
147  forAll(edges, edgeI)
148  {
149  const edge& e = edges[edgeI];
150 
151  vector eVec(e.vec(mesh_.points()));
152 
153  scalar eMag = mag(eVec);
154 
155  eVec /= eMag;
156 
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  }
181 
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;
193 
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  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
214 
215 
216 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
217 
218 
219 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
220 
221 
222 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
edgeStats(const polyMesh &mesh)
Construct from mesh.
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static const scalar edgeTol_
Definition: edgeStats.H:75
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
scalar y
dynamicFvMesh & mesh
List< edge > edgeList
Definition: edgeList.H:38
static const char nl
Definition: Ostream.H:265
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
scalar minLen(Ostream &os) const
Calculate minimum edge length and print.