cutCell.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) 2016-2017 OpenCFD Ltd.
9  Copyright (C) 2016-2017 DHI
10  Copyright (C) 2018-2019 Johan Roenby
11  Copyright (C) 2019-2020 DLR
12 -------------------------------------------------------------------------------
13 License
14  This file is part of OpenFOAM.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 \*---------------------------------------------------------------------------*/
30 
31 #include "cutCell.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 int Foam::cutCell::debug = 0;
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
47 (
48  const DynamicList<point>& cutFaceCentres,
49  const DynamicList<vector>& cutFaceAreas,
50  vector& subCellCentre, scalar& subCellVolume
51 )
52 {
53  // Clear the fields for accumulation
54  subCellCentre = Zero;
55  subCellVolume = Zero;
56 
57  // first estimate the approximate cell centre as the average of
58  // face centres
59 
60  vector cEst = average(cutFaceCentres);
61 
62  // Contribution to subcell centre and volume from cut faces
63  forAll(cutFaceCentres, facei)
64  {
65  // Calculate 3*face-pyramid volume
66  scalar pyr3Vol = max(
67  mag(cutFaceAreas[facei] & (cutFaceCentres[facei] - cEst)), VSMALL);
68 
69  // Calculate face-pyramid centre
70  vector pc = 0.75 * cutFaceCentres[facei] + 0.25 * cEst;
71 
72  // Accumulate volume-weighted face-pyramid centre
73  subCellCentre += pyr3Vol * pc;
74 
75  // Accumulate face-pyramid volume
76  subCellVolume += pyr3Vol;
77  }
78 
79  subCellCentre /= subCellVolume;
80  subCellVolume /= 3; // formula of pyramid
81 }
82 
83 
85 (
86  const DynamicList<DynamicList<point>>& faceEdges,
87  const vector& subCellCentre,
88  vector& faceArea,
89  vector& faceCentre
90 )
91 {
92  // Initial guess of face centre from edge points
93  point fCentre{Zero};
94  label nEdgePoints{0};
95  for (const DynamicList<point>& edgePoints : faceEdges)
96  {
97  for (const point& p : edgePoints)
98  {
99  fCentre += p;
100  nEdgePoints++;
101  }
102  }
103  if (nEdgePoints > 0)
104  {
105  fCentre /= nEdgePoints;
106  }
107 
108  vector sumN{Zero};
109  scalar sumA{0};
110  vector sumAc{Zero};
111 
112  // calculate area and centre
113  forAll(faceEdges, ei)
114  {
115  const DynamicList<point>& edgePoints = faceEdges[ei];
116  const label nPoints = edgePoints.size();
117  for (label pi = 0; pi < nPoints - 1; pi++)
118  {
119  const point& nextPoint = edgePoints[pi + 1];
120 
121  vector c = edgePoints[pi] + nextPoint + fCentre;
122  vector n =
123  (nextPoint - edgePoints[pi]) ^ (fCentre - edgePoints[pi]);
124  scalar a = mag(n);
125 
126  // Edges may have different orientation
127  sumN += Foam::sign(n & sumN) * n;
128  sumA += a;
129  sumAc += a * c;
130  }
131  }
132 
133  // This is to deal with zero-area faces. Mark very small faces
134  // to be detected in e.g., processorPolyPatch.
135  if (sumA < ROOTVSMALL)
136  {
137  faceCentre = fCentre;
138  faceArea = Zero;
139  }
140  else
141  {
142  faceCentre = (1.0/3.0)*sumAc/sumA;
143  faceArea = 0.5*sumN;
144  }
145 
146  // Check faceArea direction and change if not pointing in the subcell
147  if ((faceArea & (faceCentre - subCellCentre)) >= 0)
148  {
149  faceArea *= (-1);
150  }
151 }
152 
153 
155 (
156  const vector& faceArea,
157  const vector& faceCentre,
158  const DynamicList<DynamicList<point>>& faceEdges,
159  DynamicList<point>& facePoints
160 )
161 {
162  const vector zhat = normalised(faceArea);
163  vector xhat = faceEdges[0][0] - faceCentre;
164  xhat = (xhat - (xhat & zhat)*zhat);
165  xhat.normalise();
166  vector yhat = normalised(zhat ^ xhat);
167 
168  // Calculating all intersection points
169  DynamicList<point> unsortedFacePoints(3 * faceEdges.size());
170  DynamicList<scalar> unsortedFacePointAngles(3 * faceEdges.size());
171  for (const DynamicList<point>& edgePoints : faceEdges)
172  {
173  for (const point& p : edgePoints)
174  {
175  unsortedFacePoints.append(p);
176  unsortedFacePointAngles.append
177  (
179  (
180  ((p - faceCentre) & yhat),
181  ((p - faceCentre) & xhat)
182  )
183  );
184  }
185  }
186 
187  // Sorting face points by angle and inserting into facePoints
188  labelList order(Foam::sortedOrder(unsortedFacePointAngles));
189  facePoints.append(unsortedFacePoints[order[0]]);
190  for (label pi = 1; pi < order.size(); ++pi)
191  {
192  if
193  (
194  mag
195  (
196  unsortedFacePointAngles[order[pi]]
197  - unsortedFacePointAngles[order[pi - 1]]
198  ) > 1e-8)
199  {
200  facePoints.append(unsortedFacePoints[order[pi]]);
201  }
202  }
203 }
204 
205 
206 // ************************************************************************* //
cutCell.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< point >
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:312
Foam::cutCell::calcCellData
static void calcCellData(const DynamicList< point > &cutFaceCentres, const DynamicList< vector > &cutFaceAreas, vector &subCellCentre, scalar &subCellVolume)
Calculates volume and centre of the cutted cell.
Definition: cutCell.C:47
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cutCell::cutCell
cutCell(const fvMesh &unused)
Construct from fvMesh.
Definition: cutCell.C:40
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:474
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::cutCell::debug
static int debug
Definition: cutCell.H:95
Foam::cutCell::calcIsoFacePointsFromEdges
static void calcIsoFacePointsFromEdges(const vector &faceArea, const vector &faceCentre, const DynamicList< DynamicList< point >> &faceEdges, DynamicList< point > &facePoints)
Calculates the point of the cutting face.
Definition: cutCell.C:155
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:483
Foam::cutCell::calcGeomDataCutFace
static void calcGeomDataCutFace(const DynamicList< DynamicList< point >> &faceEdges, const vector &subCellCentre, vector &faceArea, vector &faceCentre)
Calculates area and centre of the cutting face.
Definition: cutCell.C:85
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328