cellModel.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017 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 
29 #include "cellModel.H"
30 #include "pyramidPointFaceRef.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
35 (
36  const labelList& pointLabels,
37  const UList<point>& points
38 ) const
39 {
40  // Estimate cell centre by averaging the cell points
41  vector cEst = Zero;
42  for (const label pointi : pointLabels)
43  {
44  cEst += points[pointi];
45  }
46  cEst /= scalar(pointLabels.size());
47 
48 
49  // Calculate the centre by breaking the cell into pyramids and
50  // volume-weighted averaging their centres
51 
52  scalar sumV = 0;
53  vector sumVc = Zero;
54 
55  forAll(faces_, facei)
56  {
57  const Foam::face f(pointLabels, faces_[facei]);
58 
59  const scalar pyrVol = pyramidPointFaceRef(f, cEst).mag(points);
60 
61  if (pyrVol > SMALL)
62  {
64  << "zero or negative pyramid volume: " << -pyrVol
65  << " for face " << facei
66  << endl;
67  }
68 
69  sumV -= pyrVol;
70  sumVc -= pyrVol * pyramidPointFaceRef(f, cEst).centre(points);
71  }
72 
73  return sumVc/(sumV + VSMALL);
74 }
75 
76 
77 Foam::scalar Foam::cellModel::mag
78 (
79  const labelList& pointLabels,
80  const UList<point>& points
81 ) const
82 {
83  // Estimate cell centre by averaging the cell points
84  vector cEst = Zero;
85  for (const label pointi : pointLabels)
86  {
87  cEst += points[pointi];
88  }
89  cEst /= scalar(pointLabels.size());
90 
91 
92  // Calculate the magnitude by summing the -mags of the pyramids
93  // The sign change is because the faces point outwards
94  // and a pyramid is constructed from an inward pointing face
95  // and the base centre-apex vector
96 
97  scalar sumV = 0;
98 
99  forAll(faces_, facei)
100  {
101  const Foam::face f(pointLabels, faces_[facei]);
102 
103  const scalar pyrVol = pyramidPointFaceRef(f, cEst).mag(points);
104 
105  if (pyrVol > SMALL)
106  {
108  << "zero or negative pyramid volume: " << -pyrVol
109  << " for face " << facei
110  << endl;
111  }
112 
113  sumV -= pyrVol;
114  }
115 
116  return sumV;
117 }
118 
119 
120 // ************************************************************************* //
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::cellModel::centre
vector centre(const labelList &pointLabels, const UList< point > &points) const
Centroid of the cell.
Definition: cellModel.C:35
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
cellModel.H
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::cellModel::mag
scalar mag(const labelList &pointLabels, const UList< point > &points) const
Cell volume.
Definition: cellModel.C:78
Foam::pyramidPointFaceRef
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred points and faces.
Definition: pyramidPointFaceRef.H:48
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
pyramidPointFaceRef.H
pointLabels
labelList pointLabels(nPoints, -1)