primitiveMeshFaceCentresAndAreas.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-2023 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 Description
25  Calculate the face centres and areas.
26 
27  Calculate the centre by breaking the face into triangles using the face
28  centre and area-weighted averaging their centres. This method copes with
29  small face-concavity.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "primitiveMesh.H"
34 
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
39 {
40  if (debug)
41  {
42  Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
43  << "Calculating face centres and face areas"
44  << endl;
45  }
46 
47  // It is an error to attempt to recalculate faceCentres
48  // if the pointer is already set
49  if (faceCentresPtr_ || faceAreasPtr_ || magFaceAreasPtr_)
50  {
52  << "Face centres or face areas already calculated"
53  << abort(FatalError);
54  }
55 
56  faceCentresPtr_ = new vectorField(nFaces());
57  vectorField& fCtrs = *faceCentresPtr_;
58 
59  faceAreasPtr_ = new vectorField(nFaces());
60  vectorField& fAreas = *faceAreasPtr_;
61 
62  magFaceAreasPtr_ = new scalarField(nFaces());
63  scalarField& magfAreas = *magFaceAreasPtr_;
64 
65  makeFaceCentresAndAreas(points(), fCtrs, fAreas, magfAreas);
66 
67  if (debug)
68  {
69  Pout<< "primitiveMesh::calcFaceCentresAndAreas() : "
70  << "Finished calculating face centres and face areas"
71  << endl;
72  }
73 }
74 
75 
77 (
78  const pointField& p,
79  vectorField& fCtrs,
80  vectorField& fAreas,
81  scalarField& magfAreas
82 ) const
83 {
84  const faceList& fs = faces();
85 
86  forAll(fs, facei)
87  {
88  const Tuple2<vector, point> areaAndCentre =
90 
91  fCtrs[facei] = areaAndCentre.second();
92  fAreas[facei] = areaAndCentre.first();
93  magfAreas[facei] = max(mag(fAreas[facei]), rootVSmall);
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  if (!faceCentresPtr_)
103  {
104  calcFaceCentresAndAreas();
105  }
106 
107  return *faceCentresPtr_;
108 }
109 
110 
112 {
113  if (!faceAreasPtr_)
114  {
115  calcFaceCentresAndAreas();
116  }
117 
118  return *faceAreasPtr_;
119 }
120 
121 
123 {
124  if (!magFaceAreasPtr_)
125  {
126  calcFaceCentresAndAreas();
127  }
128 
129  return *magFaceAreasPtr_;
130 }
131 
132 
133 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
A List with indirect addressing.
Definition: UIndirectList.H:60
static Tuple2< vector, point > areaAndCentre(const PointField &)
Return vector area and centre point given face points.
const vectorField & faceCentres() const
label nFaces() const
const scalarField & magFaceAreas() const
virtual const pointField & points() const =0
Return mesh points.
const vectorField & faceAreas() const
void calcFaceCentresAndAreas() const
Calculate face centres and areas.
void makeFaceCentresAndAreas(const pointField &p, vectorField &fCtrs, vectorField &fAreas, scalarField &magfAreas) const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
volScalarField & p