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-------------------------------------------------------------------------------
11License
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
37namespace Foam
38{
41 (
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 // Compute an estimate of the centre as the average of the points
78 solveVector fCentre = p[f[0]];
79 for (label pi = 1; pi < nPoints; pi++)
80 {
81 fCentre += solveVector(p[f[pi]]);
82 }
83 fCentre /= nPoints;
84
85 // Compute the face area normal and unit normal by summing up the
86 // normals of the triangles formed by connecting each edge to the
87 // point average.
88 solveVector sumA = Zero;
89 for (label pi = 0; pi < nPoints; pi++)
90 {
91 const label nextPi(pi == nPoints-1 ? 0 : pi+1);
92 const solveVector nextPoint(p[f[nextPi]]);
93 const solveVector thisPoint(p[f[pi]]);
94
95 const solveVector a =
96 (nextPoint - thisPoint)^(fCentre - thisPoint);
97
98 sumA += a;
99 }
100 const solveVector sumAHat = normalised(sumA);
101
102 // Compute the area-weighted sum of the triangle centres. Note use
103 // the triangle area projected in the direction of the face normal
104 // as the weight, *not* the triangle area magnitude. Only the
105 // former makes the calculation independent of the initial estimate.
106 solveScalar sumAn = 0.0;
107 solveVector sumAnc = Zero;
108 for (label pi = 0; pi < nPoints; pi++)
109 {
110 const label nextPi(pi == nPoints-1 ? 0 : pi+1);
111 const solveVector nextPoint(p[f[nextPi]]);
112 const solveVector thisPoint(p[f[pi]]);
113
114 const solveVector c = thisPoint + nextPoint + fCentre;
115 const solveVector a =
116 (nextPoint - thisPoint)^(fCentre - thisPoint);
117
118 const scalar an = a & sumAHat;
119
120 sumAn += an;
121 sumAnc += an*c;
122 }
123
124 // Complete calculating centres and areas. If the face is too small
125 // for the sums to be reliably divided then just set the centre to
126 // the initial estimate.
127 if (sumAn > ROOTVSMALL)
128 {
129 fCtrs[facei] = (1.0/3.0)*sumAnc/sumAn;
130 }
131 else
132 {
133 fCtrs[facei] = fCentre;
134 }
135 fAreas[facei] = 0.5*sumA;
136 }
137 }
138}
139
140
141// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142
144(
145 const fvMesh& mesh,
146 const dictionary& dict
147)
148:
150{
151 // Force local calculation
152 movePoints();
153}
154
155
156// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157
159{
161
162 if (debug)
163 {
164 Pout<< "stabilisedFvGeometryScheme::movePoints() : "
165 << "recalculating primitiveMesh centres" << endl;
166 }
167
168 if
169 (
170 !mesh_.hasCellCentres()
171 && !mesh_.hasFaceCentres()
172 && !mesh_.hasCellVolumes()
173 && !mesh_.hasFaceAreas()
174 )
175 {
176 vectorField faceCentres(mesh_.nFaces());
177 vectorField faceAreas(mesh_.nFaces());
178
179 makeFaceCentresAndAreas
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
199 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
200 (
201 std::move(faceCentres),
202 std::move(faceAreas),
203 std::move(cellCentres),
204 std::move(cellVolumes)
205 );
206 }
207}
208
209
210// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
Default geometry calculation scheme. Slight stabilisation for bad meshes.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base class for geometry calculation schemes.
virtual void movePoints()
Update basic geometric properties from provided points.
const fvMesh & mesh() const
Return mesh reference.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
static void makeCellCentresAndVols(const primitiveMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, vectorField &cellCtrs, scalarField &cellVols)
Calculate cell centres and volumes from face properties.
void resetGeometry(pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes)
Reset the local geometry.
Geometry calculation scheme that implements face geometry calculation using normal-component-of-area ...
static void makeFaceCentresAndAreas(const polyMesh &mesh, const pointField &p, vectorField &fCtrs, vectorField &fAreas)
Calculate face area and centre weighted using pyramid height.
virtual void movePoints()
Do what is necessary if the mesh has moved.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
label nPoints
Namespace for OpenFOAM.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< solveScalar > solveVector
Definition: vector.H:64
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333