Loading...
Searching...
No Matches
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"
34#include "triangle.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
42 (
45 dict
46 );
47}
48
49
50// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
51
52namespace Foam
53{
54
55// Average of points
56// Note: use double precision to avoid overflows when summing
58(
59 const UList<point>& points,
61)
62{
63 doubleVector avg(Zero);
64
65 if (const auto n = pointLabels.size(); n)
66 {
67 for (const auto pointi : pointLabels)
68 {
69 avg += points[pointi];
70 }
71 avg /= n;
72 }
73
74 return avg;
75}
77} // End namespace Foam
78
79
80
81// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
82
84(
85 const polyMesh& mesh,
86 const pointField& p,
87 vectorField& fCtrs,
88 vectorField& fAreas
89)
90{
91 const faceList& fs = mesh.faces();
92
93 forAll(fs, facei)
94 {
95 const auto& f = fs[facei];
96 const label nPoints = f.size();
97
98 // If the face is a triangle, do a direct calculation for efficiency
99 // and to avoid round-off error-related problems
100 if (nPoints == 3)
101 {
102 const triPointRef tri(p, f[0], f[1], f[2]);
103 fCtrs[facei] = tri.centre();
104 fAreas[facei] = tri.areaNormal();
105 }
106
107 // For more complex faces, decompose into triangles
108 else
109 {
110 // Estimated centre by averaging the face points
111 const solveVector fCentre(pointsAverage(p, f));
112
113 // Compute the face area normal and unit normal by summing up the
114 // normals of the triangles formed by connecting each edge to the
115 // point average.
116 solveVector sumA(Zero);
117 for (label pi = 0; pi < nPoints; pi++)
118 {
119 const label nextPi(pi == nPoints-1 ? 0 : pi+1);
120 const solveVector nextPoint(p[f[nextPi]]);
121 const solveVector thisPoint(p[f[pi]]);
122
123 const solveVector a =
124 (nextPoint - thisPoint)^(fCentre - thisPoint);
125
126 sumA += a;
127 }
128 const solveVector sumAHat = normalised(sumA);
129
130 // Compute the area-weighted sum of the triangle centres. Note use
131 // the triangle area projected in the direction of the face normal
132 // as the weight, *not* the triangle area magnitude. Only the
133 // former makes the calculation independent of the initial estimate.
134 solveScalar sumAn(0);
135 solveVector sumAnc(Zero);
136 for (label pi = 0; pi < nPoints; pi++)
137 {
138 const label nextPi(pi == nPoints-1 ? 0 : pi+1);
139 const solveVector nextPoint(p[f[nextPi]]);
140 const solveVector thisPoint(p[f[pi]]);
141
142 const solveVector c = thisPoint + nextPoint + fCentre;
143 const solveVector a =
144 (nextPoint - thisPoint)^(fCentre - thisPoint);
145
146 const scalar an = a & sumAHat;
147
148 sumAn += an;
149 sumAnc += an*c;
150 }
151
152 // Complete calculating centres and areas. If the face is too small
153 // for the sums to be reliably divided then just set the centre to
154 // the initial estimate.
155 if (sumAn > ROOTVSMALL)
156 {
157 fCtrs[facei] = (1.0/3.0)*sumAnc/sumAn;
158 }
159 else
160 {
161 fCtrs[facei] = fCentre;
162 }
163 fAreas[facei] = 0.5*sumA;
165 }
166}
167
168
169// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170
171Foam::stabilisedFvGeometryScheme::stabilisedFvGeometryScheme
172(
173 const fvMesh& mesh,
174 const dictionary& dict
175)
176:
178{
179 // Force local calculation
180 movePoints();
181}
182
183
184// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185
187{
189
190 if (debug)
191 {
192 Pout<< "stabilisedFvGeometryScheme::movePoints() : "
193 << "recalculating primitiveMesh centres" << endl;
194 }
195
196 if
197 (
198 !mesh_.hasCellCentres()
199 && !mesh_.hasFaceCentres()
200 && !mesh_.hasCellVolumes()
201 && !mesh_.hasFaceAreas()
202 )
203 {
204 vectorField faceCentres(mesh_.nFaces());
205 vectorField faceAreas(mesh_.nFaces());
206
207 makeFaceCentresAndAreas
208 (
209 mesh_,
210 mesh_.points(),
211 faceCentres,
212 faceAreas
213 );
214
215 vectorField cellCentres(mesh_.nCells());
216 scalarField cellVolumes(mesh_.nCells());
217
219 (
220 mesh_,
221 faceCentres,
222 faceAreas,
223 cellCentres,
224 cellVolumes
225 );
226
227 const_cast<fvMesh&>(mesh_).primitiveMesh::resetGeometry
228 (
229 std::move(faceCentres),
230 std::move(faceAreas),
231 std::move(cellCentres),
232 std::move(cellVolumes)
233 );
234 }
235}
236
237
239(
240 const pointField& points,
241 const refPtr<pointField>& oldPoints,
242 pointField& faceCentres,
243 vectorField& faceAreas,
244 pointField& cellCentres,
245 scalarField& cellVolumes
246) const
247{
248 makeFaceCentresAndAreas
249 (
250 mesh_,
251 points,
252 faceCentres,
253 faceAreas
254 );
255
257 (
258 mesh_,
259 faceCentres,
260 faceAreas,
261 cellCentres,
262 cellVolumes
263 );
264
265 return true;
266}
267
268
269// ************************************************************************* //
constexpr scalar pi(M_PI)
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
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:133
Abstract base class for geometry calculation schemes.
const fvMesh & mesh_
Hold reference to mesh.
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:85
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
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.
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
Geometry calculation scheme that implements face geometry calculation using normal-component-of-area ...
virtual bool updateGeom(const pointField &points, const refPtr< pointField > &oldPoints, pointField &faceCentres, vectorField &faceAreas, pointField &cellCentres, scalarField &cellVolumes) const
Calculate geometry quantities using mesh topology and provided points. If oldPoints provided only doe...
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.
static vector areaNormal(const Point &p0, const Point &p1, const Point &p2)
The area normal for a triangle defined by three points (right-hand rule). Magnitude equal to area of ...
Definition triangleI.H:169
static Point centre(const Point &p0, const Point &p1, const Point &p2)
The centre (centroid) of three points.
Definition triangleI.H:144
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
const pointField & points
label nPoints
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
static Foam::doubleVector pointsAverage(const UList< point > &points, const labelUList &pointLabels)
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Vector< double > doubleVector
Definition vector.H:54
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition triangleFwd.H:46
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
Vector< solveScalar > solveVector
Definition vector.H:60
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299