cell.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "cell.H"
29 #include "pyramidPointFaceRef.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 const char* const Foam::cell::typeName = "cell";
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
38 {
39  // return the unordered list of vertex labels supporting the cell
40 
41  // count the maximum size of all vertices
42  label maxVert = 0;
43 
44  const labelList& faces = *this;
45 
46  forAll(faces, facei)
47  {
48  maxVert += f[faces[facei]].size();
49  }
50 
51  // set the fill-in list
52  labelList p(maxVert);
53 
54  // in the first face there is no duplicates
55  const labelList& first = f[faces[0]];
56 
57  forAll(first, pointi)
58  {
59  p[pointi] = first[pointi];
60  }
61 
62  // re-use maxVert to count the real vertices
63  maxVert = first.size();
64 
65  // go through the rest of the faces. For each vertex, check if the point is
66  // already inserted (up to maxVert, which now carries the number of real
67  // points. If not, add it at the end of the list.
68  for (label facei = 1; facei < faces.size(); facei++)
69  {
70  const labelList& curFace = f[faces[facei]];
71 
72  forAll(curFace, pointi)
73  {
74  const label curPoint = curFace[pointi];
75 
76  bool found = false;
77 
78  for (label checkI = 0; checkI < maxVert; checkI++)
79  {
80  if (curPoint == p[checkI])
81  {
82  found = true;
83  break;
84  }
85  }
86 
87  if (!found)
88  {
89  p[maxVert] = curPoint;
90 
91  maxVert++;
92  }
93  }
94  }
95 
96  // reset the size of the list
97  p.setSize(maxVert);
98 
99  return p;
100 }
101 
102 
104 (
105  const faceUList& f,
106  const UList<point>& meshPoints
107 ) const
108 {
109  const labelList pointLabels = labels(f);
110 
111  pointField p(pointLabels.size());
112 
113  forAll(p, i)
114  {
115  p[i] = meshPoints[pointLabels[i]];
116  }
117 
118  return p;
119 }
120 
121 
123 {
124  // return the lisf of cell edges
125 
126  const labelList& curFaces = *this;
127 
128  // create a list of edges
129  label maxNoEdges = 0;
130 
131  forAll(curFaces, facei)
132  {
133  maxNoEdges += f[curFaces[facei]].nEdges();
134  }
135 
136  edgeList allEdges(maxNoEdges);
137  label nEdges = 0;
138 
139  forAll(curFaces, facei)
140  {
141  const edgeList curFaceEdges = f[curFaces[facei]].edges();
142 
143  forAll(curFaceEdges, faceEdgeI)
144  {
145  const edge& curEdge = curFaceEdges[faceEdgeI];
146 
147  bool edgeFound = false;
148 
149  for (label addedEdgeI = 0; addedEdgeI < nEdges; addedEdgeI++)
150  {
151  if (allEdges[addedEdgeI] == curEdge)
152  {
153  edgeFound = true;
154 
155  break;
156  }
157  }
158 
159  if (!edgeFound)
160  {
161  // Add the new edge onto the list
162  allEdges[nEdges] = curEdge;
163  nEdges++;
164  }
165  }
166  }
167 
168  allEdges.setSize(nEdges);
169 
170  return allEdges;
171 }
172 
173 
175 (
176  const UList<point>& p,
177  const faceUList& f
178 ) const
179 {
180  // When one wants to access the cell centre and magnitude, the
181  // functionality on the mesh level should be used in preference to the
182  // functions provided here. They do not rely to the functionality
183  // implemented here, provide additional checking and are more efficient.
184  // The cell::centre and cell::mag functions may be removed in the future.
185 
186  // WARNING!
187  // In the old version of the code, it was possible to establish whether any
188  // of the pyramids had a negative volume, caused either by the poor
189  // estimate of the cell centre or by the fact that the cell was defined
190  // inside out. In the new description of the cell, this can only be
191  // established with the reference to the owner-neighbour face-cell
192  // relationship, which is not available on this level. Thus, all the
193  // pyramids are believed to be positive with no checking.
194 
195  // first calculate the approximate cell centre as the average of all
196  // face centres
197  vector cEst = Zero;
198  scalar sumArea = 0;
199 
200  const labelList& faces = *this;
201 
202  forAll(faces, facei)
203  {
204  scalar a = f[faces[facei]].mag(p);
205  cEst += f[faces[facei]].centre(p)*a;
206  sumArea += a;
207  }
208 
209  cEst /= sumArea + VSMALL;
210 
211  // Calculate the centre by breaking the cell into pyramids and
212  // volume-weighted averaging their centres
213  vector sumVc = Zero;
214 
215  scalar sumV = 0;
216 
217  forAll(faces, facei)
218  {
219  // calculate pyramid volume. If it is greater than zero, OK.
220  // If not, the pyramid is inside-out. Create a face with the opposite
221  // order and recalculate pyramid centre!
222  scalar pyrVol = pyramidPointFaceRef(f[faces[facei]], cEst).mag(p);
223  vector pyrCentre = pyramidPointFaceRef(f[faces[facei]], cEst).centre(p);
224 
225  // if pyramid inside-out because face points inwards invert
226  // N.B. pyramid remains unchanged
227  if (pyrVol < 0)
228  {
229  pyrVol = -pyrVol;
230  }
231 
232  sumVc += pyrVol*pyrCentre;
233  sumV += pyrVol;
234  }
235 
236  return sumVc/(sumV + VSMALL);
237 }
238 
239 
240 Foam::scalar Foam::cell::mag
241 (
242  const UList<point>& p,
243  const faceUList& f
244 ) const
245 {
246  // When one wants to access the cell centre and magnitude, the
247  // functionality on the mesh level should be used in preference to the
248  // functions provided here. They do not rely to the functionality
249  // implemented here, provide additional checking and are more efficient.
250  // The cell::centre and cell::mag functions may be removed in the future.
251 
252  // WARNING! See cell::centre
253 
254  // first calculate the approximate cell centre as the average of all
255  // face centres
256  vector cEst = Zero;
257  scalar nCellFaces = 0;
258 
259  const labelList& faces = *this;
260 
261  forAll(faces, facei)
262  {
263  cEst += f[faces[facei]].centre(p);
264  nCellFaces += 1;
265  }
266 
267  cEst /= nCellFaces;
268 
269  // Calculate the magnitude by summing the mags of the pyramids
270  scalar v = 0;
271 
272  forAll(faces, facei)
273  {
274  v += ::Foam::mag(pyramidPointFaceRef(f[faces[facei]], cEst).mag(p));
275  }
276 
277  return v;
278 }
279 
280 
281 // * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
282 
283 bool Foam::operator==(const cell& a, const cell& b)
284 {
285  // Trivial reject: faces are different size
286  if (a.size() != b.size())
287  {
288  return false;
289  }
290 
291  List<bool> fnd(a.size(), false);
292 
293  for (const label curLabel : b)
294  {
295  bool found = false;
296 
297  forAll(a, ai)
298  {
299  if (a[ai] == curLabel)
300  {
301  found = true;
302  fnd[ai] = true;
303  break;
304  }
305  }
306 
307  if (!found)
308  {
309  return false;
310  }
311  }
312 
313  // Any faces missed?
314  forAll(fnd, ai)
315  {
316  if (!fnd[ai])
317  {
318  return false;
319  }
320  }
321 
322  return true;
323 }
324 
325 
326 // ************************************************************************* //
cell.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::cell::typeName
static const char *const typeName
Definition: cell.H:63
Foam::cell::points
pointField points(const faceUList &f, const UList< point > &meshPoints) const
Return the cell vertices given the list of faces and mesh points.
Definition: cell.C:104
Foam::cell::centre
point centre(const UList< point > &p, const faceUList &f) const
Returns cell centre.
Definition: cell.C:175
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::cell::edges
edgeList edges(const faceUList &f) const
Return cell edges.
Definition: cell.C:122
Foam::cell::mag
scalar mag(const UList< point > &p, const faceUList &f) const
Returns cell volume.
Definition: cell.C:241
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::cell::labels
labelList labels(const faceUList &f) const
Return unordered list of cell vertices given the list of faces.
Definition: cell.C:37
found
bool found
Definition: TABSMDCalcMethod2.H:32
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::UList< face >
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::pyramidPointFaceRef
pyramid< point, const point &, const face & > pyramidPointFaceRef
Definition: pyramidPointFaceRef.H:47
pyramidPointFaceRef.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
pointLabels
labelList pointLabels(nPoints, -1)