stabilisedFvGeometryScheme.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) 2020 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvMesh.H"
32 #include "PrecisionAdaptor.H"
33 #include "primitiveMeshTools.H"
34 #include "triangle.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(stabilisedFvGeometryScheme, 0);
42  (
43  fvGeometryScheme,
44  stabilisedFvGeometryScheme,
45  dict
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 (
54  const polyMesh& mesh,
55  const pointField& p,
56  vectorField& fCtrs,
57  vectorField& fAreas
58 )
59 {
60  const faceList& fs = mesh.faces();
61 
62  forAll(fs, facei)
63  {
64  const labelList& f = fs[facei];
65  label nPoints = f.size();
66 
67  // If the face is a triangle, do a direct calculation for efficiency
68  // and to avoid round-off error-related problems
69  if (nPoints == 3)
70  {
71  fCtrs[facei] = triPointRef::centre(p[f[0]], p[f[1]], p[f[2]]);
72  fAreas[facei] = triPointRef::areaNormal(p[f[0]], p[f[1]], p[f[2]]);
73  }
74 
75  // For more complex faces, decompose into triangles
76  else
77  {
78  // Compute an estimate of the centre as the average of the points
79  solveVector fCentre = p[f[0]];
80  for (label pi = 1; pi < nPoints; pi++)
81  {
82  fCentre += solveVector(p[f[pi]]);
83  }
84  fCentre /= nPoints;
85 
86  // Compute the face area normal and unit normal by summing up the
87  // normals of the triangles formed by connecting each edge to the
88  // point average.
89  solveVector sumA = Zero;
90  for (label pi = 0; pi < nPoints; pi++)
91  {
92  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
93  const solveVector nextPoint(p[f[nextPi]]);
94  const solveVector thisPoint(p[f[pi]]);
95 
96  const solveVector a =
97  (nextPoint - thisPoint)^(fCentre - thisPoint);
98 
99  sumA += a;
100  }
101  const solveVector sumAHat = normalised(sumA);
102 
103  // Compute the area-weighted sum of the triangle centres. Note use
104  // the triangle area projected in the direction of the face normal
105  // as the weight, *not* the triangle area magnitude. Only the
106  // former makes the calculation independent of the initial estimate.
107  solveScalar sumAn = 0.0;
108  solveVector sumAnc = Zero;
109  for (label pi = 0; pi < nPoints; pi++)
110  {
111  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
112  const solveVector nextPoint(p[f[nextPi]]);
113  const solveVector thisPoint(p[f[pi]]);
114 
115  const solveVector c = thisPoint + nextPoint + fCentre;
116  const solveVector a =
117  (nextPoint - thisPoint)^(fCentre - thisPoint);
118 
119  const scalar an = a & sumAHat;
120 
121  sumAn += an;
122  sumAnc += an*c;
123  }
124 
125  // Complete calculating centres and areas. If the face is too small
126  // for the sums to be reliably divided then just set the centre to
127  // the initial estimate.
128  if (sumAn > ROOTVSMALL)
129  {
130  fCtrs[facei] = (1.0/3.0)*sumAnc/sumAn;
131  }
132  else
133  {
134  fCtrs[facei] = fCentre;
135  }
136  fAreas[facei] = 0.5*sumA;
137  }
138  }
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143 
144 Foam::stabilisedFvGeometryScheme::stabilisedFvGeometryScheme
145 (
146  const fvMesh& mesh,
147  const dictionary& dict
148 )
149 :
150  basicFvGeometryScheme(mesh, dict)
151 {
152  // Force local calculation
153  movePoints();
154 }
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
162 
163  if (debug)
164  {
165  Pout<< "stabilisedFvGeometryScheme::movePoints() : "
166  << "recalculating primitiveMesh centres" << endl;
167  }
168 
169  if
170  (
171  !mesh_.hasCellCentres()
172  && !mesh_.hasFaceCentres()
173  && !mesh_.hasCellVolumes()
174  && !mesh_.hasFaceAreas()
175  )
176  {
177  vectorField faceCentres(mesh_.nFaces());
178  vectorField faceAreas(mesh_.nFaces());
179 
180  makeFaceCentresAndAreas
181  (
182  mesh_,
183  mesh_.points(),
184  faceCentres,
185  faceAreas
186  );
187 
188  vectorField cellCentres(mesh_.nCells());
189  scalarField cellVolumes(mesh_.nCells());
190 
192  (
193  mesh_,
194  faceCentres,
195  faceAreas,
196  cellCentres,
197  cellVolumes
198  );
199 
200  const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
201  (
202  std::move(faceCentres),
203  std::move(faceAreas),
204  std::move(cellCentres),
205  std::move(cellVolumes)
206  );
207  }
208 }
209 
210 
211 // ************************************************************************* //
dictionary dict
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Vector< solveScalar > solveVector
Definition: vector.H:60
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
dynamicFvMesh & mesh
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Point centre() const
Return centre (centroid)
Definition: triangleI.H:126
constexpr scalar pi(M_PI)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
virtual void movePoints()
Do what is necessary if the mesh has moved.
const fvMesh & mesh() const
Return mesh reference.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
const dimensionedScalar c
Speed of light in a vacuum.
virtual void movePoints()
Update basic geometric properties from provided points.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void makeFaceCentresAndAreas(const polyMesh &mesh, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face area and centre weighted using pyramid height.
volScalarField & p
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
vector areaNormal() const
The area normal - with magnitude equal to area of triangle.
Definition: triangleI.H:151
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127