cellQuality.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-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "cellQuality.H"
29 #include "unitConversion.H"
30 #include "SubField.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 Foam::cellQuality::cellQuality(const polyMesh& mesh)
35 :
36  mesh_(mesh)
37 {}
38 
39 
40 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
41 
43 {
44  auto tresult = tmp<scalarField>::New(mesh_.nCells(), Zero);
45  auto& result = tresult.ref();
46 
47  scalarField sumArea(mesh_.nCells(), Zero);
48 
49  const vectorField& centres = mesh_.cellCentres();
50  const vectorField& areas = mesh_.faceAreas();
51 
52  const labelList& own = mesh_.faceOwner();
53  const labelList& nei = mesh_.faceNeighbour();
54 
55  forAll(nei, facei)
56  {
57  vector d = centres[nei[facei]] - centres[own[facei]];
58  vector s = areas[facei];
59  scalar magS = mag(s);
60 
61  scalar cosDDotS =
62  radToDeg(Foam::acos(min(1.0, (d & s)/(mag(d)*magS + VSMALL))));
63 
64  result[own[facei]] = max(cosDDotS, result[own[facei]]);
65 
66  result[nei[facei]] = max(cosDDotS, result[nei[facei]]);
67  }
68 
69  forAll(mesh_.boundaryMesh(), patchi)
70  {
71  const labelUList& faceCells =
72  mesh_.boundaryMesh()[patchi].faceCells();
73 
74  const vectorField::subField faceCentres =
75  mesh_.boundaryMesh()[patchi].faceCentres();
76 
77  const vectorField::subField faceAreas =
78  mesh_.boundaryMesh()[patchi].faceAreas();
79 
80  forAll(faceCentres, facei)
81  {
82  vector d = faceCentres[facei] - centres[faceCells[facei]];
83  vector s = faceAreas[facei];
84  scalar magS = mag(s);
85 
86  scalar cosDDotS =
87  radToDeg(Foam::acos(min(1.0, (d & s)/(mag(d)*magS + VSMALL))));
88 
89  result[faceCells[facei]] = max(cosDDotS, result[faceCells[facei]]);
90  }
91  }
92 
93  return tresult;
94 }
95 
96 
98 {
99  auto tresult = tmp<scalarField>::New(mesh_.nCells(), Zero);
100  auto& result = tresult.ref();
101 
102  scalarField sumArea(mesh_.nCells(), Zero);
103 
104  const vectorField& cellCtrs = mesh_.cellCentres();
105  const vectorField& faceCtrs = mesh_.faceCentres();
106  const vectorField& areas = mesh_.faceAreas();
107 
108  const labelList& own = mesh_.faceOwner();
109  const labelList& nei = mesh_.faceNeighbour();
110 
111  forAll(nei, facei)
112  {
113  scalar dOwn = mag
114  (
115  (faceCtrs[facei] - cellCtrs[own[facei]]) & areas[facei]
116  )/mag(areas[facei]);
117 
118  scalar dNei = mag
119  (
120  (cellCtrs[nei[facei]] - faceCtrs[facei]) & areas[facei]
121  )/mag(areas[facei]);
122 
123  point faceIntersection =
124  cellCtrs[own[facei]]
125  + (dOwn/(dOwn+dNei))*(cellCtrs[nei[facei]] - cellCtrs[own[facei]]);
126 
127  scalar skewness =
128  mag(faceCtrs[facei] - faceIntersection)
129  /(mag(cellCtrs[nei[facei]] - cellCtrs[own[facei]]) + VSMALL);
130 
131  result[own[facei]] = max(skewness, result[own[facei]]);
132 
133  result[nei[facei]] = max(skewness, result[nei[facei]]);
134  }
135 
136  forAll(mesh_.boundaryMesh(), patchi)
137  {
138  const labelUList& faceCells =
139  mesh_.boundaryMesh()[patchi].faceCells();
140 
141  const vectorField::subField faceCentres =
142  mesh_.boundaryMesh()[patchi].faceCentres();
143 
144  const vectorField::subField faceAreas =
145  mesh_.boundaryMesh()[patchi].faceAreas();
146 
147  forAll(faceCentres, facei)
148  {
149  vector n = faceAreas[facei]/mag(faceAreas[facei]);
150 
151  point faceIntersection =
152  cellCtrs[faceCells[facei]]
153  + ((faceCentres[facei] - cellCtrs[faceCells[facei]])&n)*n;
154 
155  scalar skewness =
156  mag(faceCentres[facei] - faceIntersection)
157  /(
158  mag(faceCentres[facei] - cellCtrs[faceCells[facei]])
159  + VSMALL
160  );
161 
162  result[faceCells[facei]] = max(skewness, result[faceCells[facei]]);
163  }
164  }
165 
166  return tresult;
167 }
168 
169 
171 {
172  auto tresult = tmp<scalarField>::New(mesh_.nFaces(), Zero);
173  auto& result = tresult.ref();
174 
175  const vectorField& centres = mesh_.cellCentres();
176  const vectorField& areas = mesh_.faceAreas();
177 
178  const labelList& own = mesh_.faceOwner();
179  const labelList& nei = mesh_.faceNeighbour();
180 
181  forAll(nei, facei)
182  {
183  vector d = centres[nei[facei]] - centres[own[facei]];
184  vector s = areas[facei];
185  scalar magS = mag(s);
186 
187  scalar cosDDotS =
188  radToDeg(Foam::acos(min(1.0, (d & s)/(mag(d)*magS + VSMALL))));
189 
190  result[facei] = cosDDotS;
191  }
192 
193  label globalFacei = mesh_.nInternalFaces();
194 
195  forAll(mesh_.boundaryMesh(), patchi)
196  {
197  const labelUList& faceCells =
198  mesh_.boundaryMesh()[patchi].faceCells();
199 
200  const vectorField::subField faceCentres =
201  mesh_.boundaryMesh()[patchi].faceCentres();
202 
203  const vectorField::subField faceAreas =
204  mesh_.boundaryMesh()[patchi].faceAreas();
205 
206  forAll(faceCentres, facei)
207  {
208  vector d = faceCentres[facei] - centres[faceCells[facei]];
209  vector s = faceAreas[facei];
210  scalar magS = mag(s);
211 
212  scalar cosDDotS =
213  radToDeg(Foam::acos(min(1.0, (d & s)/(mag(d)*magS + VSMALL))));
214 
215  result[globalFacei++] = cosDDotS;
216  }
217  }
218 
219  return tresult;
220 }
221 
222 
224 {
225  auto tresult = tmp<scalarField>::New(mesh_.nFaces(), Zero);
226  auto& result = tresult.ref();
227 
228  const vectorField& cellCtrs = mesh_.cellCentres();
229  const vectorField& faceCtrs = mesh_.faceCentres();
230  const vectorField& areas = mesh_.faceAreas();
231 
232  const labelList& own = mesh_.faceOwner();
233  const labelList& nei = mesh_.faceNeighbour();
234 
235  forAll(nei, facei)
236  {
237  scalar dOwn = mag
238  (
239  (faceCtrs[facei] - cellCtrs[own[facei]]) & areas[facei]
240  )/mag(areas[facei]);
241 
242  scalar dNei = mag
243  (
244  (cellCtrs[nei[facei]] - faceCtrs[facei]) & areas[facei]
245  )/mag(areas[facei]);
246 
247  point faceIntersection =
248  cellCtrs[own[facei]]
249  + (dOwn/(dOwn+dNei))*(cellCtrs[nei[facei]] - cellCtrs[own[facei]]);
250 
251  result[facei] =
252  mag(faceCtrs[facei] - faceIntersection)
253  /(mag(cellCtrs[nei[facei]] - cellCtrs[own[facei]]) + VSMALL);
254  }
255 
256 
257  label globalFacei = mesh_.nInternalFaces();
258 
259  forAll(mesh_.boundaryMesh(), patchi)
260  {
261  const labelUList& faceCells =
262  mesh_.boundaryMesh()[patchi].faceCells();
263 
264  const vectorField::subField faceCentres =
265  mesh_.boundaryMesh()[patchi].faceCentres();
266 
267  const vectorField::subField faceAreas =
268  mesh_.boundaryMesh()[patchi].faceAreas();
269 
270  forAll(faceCentres, facei)
271  {
272  vector n = faceAreas[facei]/mag(faceAreas[facei]);
273 
274  point faceIntersection =
275  cellCtrs[faceCells[facei]]
276  + ((faceCentres[facei] - cellCtrs[faceCells[facei]])&n)*n;
277 
278  result[globalFacei++] =
279  mag(faceCentres[facei] - faceIntersection)
280  /(
281  mag(faceCentres[facei] - cellCtrs[faceCells[facei]])
282  + VSMALL
283  );
284  }
285  }
286 
287  return tresult;
288 }
289 
290 
291 // ************************************************************************* //
dimensionedScalar acos(const dimensionedScalar &ds)
tmp< scalarField > nonOrthogonality() const
Return cell non-orthogonality.
Definition: cellQuality.C:35
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Unit conversion functions.
tmp< scalarField > faceSkewness() const
Return face skewness.
Definition: cellQuality.C:216
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
tmp< scalarField > faceNonOrthogonality() const
Return face non-orthogonality.
Definition: cellQuality.C:163
dynamicFvMesh & mesh
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vectorField & cellCentres() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
vector point
Point is a vector.
Definition: point.H:37
label nCells() const noexcept
Number of mesh cells.
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
const vectorField & faceAreas() const
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
tmp< scalarField > skewness() const
Return cell skewness.
Definition: cellQuality.C:90
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127