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 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(stabilisedFvGeometryScheme, 0);
41  (
42  fvGeometryScheme,
43  stabilisedFvGeometryScheme,
44  dict
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 (
53  const polyMesh& mesh,
54  const pointField& p,
55  vectorField& fCtrs,
56  vectorField& fAreas
57 )
58 {
59  const faceList& fs = mesh.faces();
60 
61  forAll(fs, facei)
62  {
63  const labelList& f = fs[facei];
64  label nPoints = f.size();
65 
66  // If the face is a triangle, do a direct calculation for efficiency
67  // and to avoid round-off error-related problems
68  if (nPoints == 3)
69  {
70  fCtrs[facei] = (1.0/3.0)*(p[f[0]] + p[f[1]] + p[f[2]]);
71  fAreas[facei] = 0.5*((p[f[1]] - p[f[0]])^(p[f[2]] - p[f[0]]));
72  }
73 
74  // For more complex faces, decompose into triangles
75  else
76  {
77  typedef Vector<solveScalar> solveVector;
78 
79  // Compute an estimate of the centre as the average of the points
80  solveVector fCentre = p[f[0]];
81  for (label pi = 1; pi < nPoints; pi++)
82  {
83  fCentre += solveVector(p[f[pi]]);
84  }
85  fCentre /= nPoints;
86 
87  // Compute the face area normal and unit normal by summing up the
88  // normals of the triangles formed by connecting each edge to the
89  // point average.
90  solveVector sumA = Zero;
91  for (label pi = 0; pi < nPoints; pi++)
92  {
93  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
94  const solveVector nextPoint(p[f[nextPi]]);
95  const solveVector thisPoint(p[f[pi]]);
96 
97  const solveVector a =
98  (nextPoint - thisPoint)^(fCentre - thisPoint);
99 
100  sumA += a;
101  }
102  const solveVector sumAHat = normalised(sumA);
103 
104  // Compute the area-weighted sum of the triangle centres. Note use
105  // the triangle area projected in the direction of the face normal
106  // as the weight, *not* the triangle area magnitude. Only the
107  // former makes the calculation independent of the initial estimate.
108  solveScalar sumAn = 0.0;
109  solveVector sumAnc = Zero;
110  for (label pi = 0; pi < nPoints; pi++)
111  {
112  const label nextPi(pi == nPoints-1 ? 0 : pi+1);
113  const solveVector nextPoint(p[f[nextPi]]);
114  const solveVector thisPoint(p[f[pi]]);
115 
116  const solveVector c = thisPoint + nextPoint + fCentre;
117  const solveVector a =
118  (nextPoint - thisPoint)^(fCentre - thisPoint);
119 
120  const scalar an = a & sumAHat;
121 
122  sumAn += an;
123  sumAnc += an*c;
124  }
125 
126  // Complete calculating centres and areas. If the face is too small
127  // for the sums to be reliably divided then just set the centre to
128  // the initial estimate.
129  if (sumAn > ROOTVSMALL)
130  {
131  fCtrs[facei] = (1.0/3.0)*sumAnc/sumAn;
132  }
133  else
134  {
135  fCtrs[facei] = fCentre;
136  }
137  fAreas[facei] = 0.5*sumA;
138  }
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
144 
145 Foam::stabilisedFvGeometryScheme::stabilisedFvGeometryScheme
146 (
147  const fvMesh& mesh,
148  const dictionary& dict
149 )
150 :
152 {
153  // Force local calculation
154  movePoints();
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  if (debug)
163  {
164  Pout<< "stabilisedFvGeometryScheme::movePoints() : "
165  << "recalculating primitiveMesh centres" << endl;
166  }
167 
168  if
169  (
171  && !mesh_.hasFaceCentres()
172  && !mesh_.hasCellVolumes()
173  && !mesh_.hasFaceAreas()
174  )
175  {
176  vectorField faceCentres(mesh_.nFaces());
177  vectorField faceAreas(mesh_.nFaces());
178 
180  (
181  mesh_,
182  mesh_.points(),
183  faceCentres,
184  faceAreas
185  );
186 
187  vectorField cellCentres(mesh_.nCells());
188  scalarField cellVolumes(mesh_.nCells());
189 
191  (
192  mesh_,
193  faceCentres,
194  faceAreas,
195  cellCentres,
196  cellVolumes
197  );
198 
200  (
201  std::move(faceCentres),
202  std::move(faceAreas),
203  std::move(cellCentres),
204  std::move(cellVolumes)
205  );
206  }
207 }
208 
209 
210 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
p
volScalarField & p
Definition: createFieldRefs.H:8
primitiveMeshTools.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
stabilisedFvGeometryScheme.H
PrecisionAdaptor.H
Foam::basicFvGeometryScheme
Default geometry calculation scheme. Slight stabilisation for bad meshes.
Definition: basicFvGeometryScheme.H:50
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::primitiveMesh::hasCellVolumes
bool hasCellVolumes() const noexcept
Definition: primitiveMeshI.H:201
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::primitiveMeshTools::makeCellCentresAndVols
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
Definition: primitiveMeshTools.C:107
Foam::Field< vector >
Foam::primitiveMesh::hasCellCentres
bool hasCellCentres() const noexcept
Definition: primitiveMeshI.H:189
Foam::primitiveMesh::hasFaceAreas
bool hasFaceAreas() const noexcept
Definition: primitiveMeshI.H:207
Foam::stabilisedFvGeometryScheme::movePoints
virtual void movePoints()
Do what is necessary if the mesh has moved.
Definition: stabilisedFvGeometryScheme.C:160
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
Foam::stabilisedFvGeometryScheme::makeFaceCentresAndAreas
static void makeFaceCentresAndAreas(const polyMesh &mesh, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face area and centre weighted using pyramid height.
Definition: stabilisedFvGeometryScheme.C:52
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::primitiveMesh::resetGeometry
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
Definition: primitiveMesh.C:284
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::primitiveMesh::hasFaceCentres
bool hasFaceCentres() const noexcept
Definition: primitiveMeshI.H:195
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::fvGeometryScheme::mesh_
const fvMesh & mesh_
Hold reference to mesh.
Definition: fvGeometryScheme.H:72