primitiveMeshFaceCentresAndAreas.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 Description
27  Calculate the face centres and areas.
28 
29  Calculate the centre by breaking the face into triangles using the face
30  centre and area-weighted averaging their centres. This method copes with
31  small face-concavity.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "primitiveMesh.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
40 {
41  if (debug)
42  {
43  Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
44  << "Calculating face centres and face areas"
45  << endl;
46  }
47 
48  // It is an error to attempt to recalculate faceCentres
49  // if the pointer is already set
50  if (faceCentresPtr_ || faceAreasPtr_)
51  {
53  << "Face centres or face areas already calculated"
54  << abort(FatalError);
55  }
56 
57  faceCentresPtr_ = new vectorField(nFaces());
58  vectorField& fCtrs = *faceCentresPtr_;
59 
60  faceAreasPtr_ = new vectorField(nFaces());
61  vectorField& fAreas = *faceAreasPtr_;
62 
63  makeFaceCentresAndAreas(points(), fCtrs, fAreas);
64 
65  if (debug)
66  {
67  Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
68  << "Finished calculating face centres and face areas"
69  << endl;
70  }
71 }
72 
73 
75 (
76  const pointField& p,
77  vectorField& fCtrs,
78  vectorField& fAreas
79 ) const
80 {
81  const faceList& fs = faces();
82 
83  forAll(fs, facei)
84  {
85  const labelList& f = fs[facei];
86  const label nPoints = f.size();
87 
88  // If the face is a triangle, do a direct calculation for efficiency
89  // and to avoid round-off error-related problems
90  if (nPoints == 3)
91  {
92  fCtrs[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
93  fAreas[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
94  }
95  else
96  {
97  typedef Vector<solveScalar> solveVector;
98 
99  solveVector sumN = Zero;
100  solveScalar sumA = 0.0;
101  solveVector sumAc = Zero;
102 
103  solveVector fCentre = p[f[0]];
104  for (label pi = 1; pi < nPoints; pi++)
105  {
106  fCentre += solveVector(p[f[pi]]);
107  }
108 
109  fCentre /= nPoints;
110 
111  for (label pi = 0; pi < nPoints; pi++)
112  {
113  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
114  const solveVector nextPoint(p[f[nextPi]]);
115  const solveVector thisPoint(p[f[pi]]);
116 
117  solveVector c = thisPoint + nextPoint + fCentre;
118  solveVector n = (nextPoint - thisPoint)^(fCentre - thisPoint);
119  solveScalar a = mag(n);
120  sumN += n;
121  sumA += a;
122  sumAc += a*c;
123  }
124 
125  // This is to deal with zero-area faces. Mark very small faces
126  // to be detected in e.g., processorPolyPatch.
127  if (sumA < ROOTVSMALL)
128  {
129  fCtrs[facei] = fCentre;
130  fAreas[facei] = Zero;
131  }
132  else
133  {
134  fCtrs[facei] = (1.0/3.0)*sumAc/sumA;
135  fAreas[facei] = 0.5*sumN;
136  }
137  }
138  }
139 }
140 
141 
142 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
143 
145 {
146  if (!faceCentresPtr_)
147  {
148  calcFaceCentresAndAreas();
149  }
150 
151  return *faceCentresPtr_;
152 }
153 
154 
156 {
157  if (!faceAreasPtr_)
158  {
159  calcFaceCentresAndAreas();
160  }
161 
162  return *faceAreasPtr_;
163 }
164 
165 
166 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::primitiveMesh::makeFaceCentresAndAreas
void makeFaceCentresAndAreas(const pointField &p, vectorField &fCtrs, vectorField &fAreas) const
Definition: primitiveMeshFaceCentresAndAreas.C:75
Foam::primitiveMesh::points
virtual const pointField & points() const =0
Return mesh points.
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::primitiveMesh::calcFaceCentresAndAreas
void calcFaceCentresAndAreas() const
Calculate face centres and areas.
Definition: primitiveMeshFaceCentresAndAreas.C:39
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
f
labelList f(nPoints)
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::List< face >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:144
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:155