Loading...
Searching...
No Matches
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 Copyright (C) 2020-2022 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
29#include "cell.H"
30#include "pyramid.H"
31#include "primitiveMesh.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35const char* const Foam::cell::typeName = "cell";
36
37
38// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39
40Foam::labelList Foam::cell::labels(const faceUList& meshFaces) const
41{
42 const labelList& cFaces = *this;
43
44 label nVerts = 0;
45 for (const label facei : cFaces)
46 {
47 nVerts += meshFaces[facei].size();
48 }
49
50 labelList pointLabels(nVerts);
51
52 // The first face has no duplicates, can copy in values
53 const labelList& firstFace = meshFaces[cFaces[0]];
54
55 std::copy(firstFace.cbegin(), firstFace.cend(), pointLabels.begin());
56
57 // Now already contains some vertices
58 nVerts = firstFace.size();
59
60 // For the rest of the faces. For each vertex, check if the point is
61 // already inserted (up to nVerts, which now carries the number of real
62 // points. If not, add it at the end of the list.
63
64 for (label facei = 1; facei < cFaces.size(); ++facei)
65 {
66 for (const label curPoint : meshFaces[cFaces[facei]])
67 {
68 bool pointFound = false;
69
70 for (label checki = 0; checki < nVerts; ++checki)
71 {
72 if (curPoint == pointLabels[checki])
73 {
74 pointFound = true;
75 break;
76 }
77 }
78
79 if (!pointFound)
80 {
81 pointLabels[nVerts] = curPoint;
82 ++nVerts;
83 }
84 }
85 }
87 pointLabels.resize(nVerts);
88
89 return pointLabels;
90}
91
92
94(
95 const faceUList& meshFaces,
96 const UList<point>& meshPoints
97) const
98{
99 const labelList pointLabels = labels(meshFaces);
100
101 pointField allPoints(pointLabels.size());
102
103 forAll(allPoints, i)
104 {
105 allPoints[i] = meshPoints[pointLabels[i]];
106 }
107
108 return allPoints;
109}
110
111
112Foam::edgeList Foam::cell::edges(const faceUList& meshFaces) const
113{
114 const labelList& cFaces = *this;
115
116 label nEdges = 0;
117 for (const label facei : cFaces)
118 {
119 nEdges += meshFaces[facei].nEdges();
120 }
121
122 edgeList allEdges(nEdges);
123
124 nEdges = 0;
125
126 forAll(cFaces, facei)
127 {
128 for (const edge& curEdge : meshFaces[cFaces[facei]].edges())
129 {
130 bool edgeFound = false;
131
132 for (label checki = 0; checki < nEdges; ++checki)
133 {
134 if (curEdge == allEdges[checki])
135 {
136 edgeFound = true;
137 break;
138 }
139 }
140
141 if (!edgeFound)
142 {
143 allEdges[nEdges] = curEdge;
144 ++nEdges;
145 }
146 }
147 }
149 allEdges.resize(nEdges);
150
151 return allEdges;
152}
153
154
156(
157 const UList<point>& meshPoints,
158 const faceUList& meshFaces
159) const
160{
161 // When one wants to access the cell centre and magnitude, the
162 // functionality on the mesh level should be used in preference to the
163 // functions provided here. They do not rely to the functionality
164 // implemented here, provide additional checking and are more efficient.
165 // The cell::centre and cell::mag functions may be removed in the future.
166
167 // WARNING!
168 // In the old version of the code, it was possible to establish whether any
169 // of the pyramids had a negative volume, caused either by the poor
170 // estimate of the cell centre or by the fact that the cell was defined
171 // inside out. In the new description of the cell, this can only be
172 // established with the reference to the owner-neighbour face-cell
173 // relationship, which is not available on this level. Thus, all the
174 // pyramids are believed to be positive with no checking.
175
176 // Approximate cell centre as the area average of all face centres
177
178 vector ctr = Zero;
179 scalar sumArea = 0;
180
181 const labelList& cFaces = *this;
182
183 for (const label facei : cFaces)
184 {
185 const scalar magArea = meshFaces[facei].mag(meshPoints);
186 ctr += meshFaces[facei].centre(meshPoints)*magArea;
187 sumArea += magArea;
188 }
189
190 ctr /= sumArea + VSMALL;
191
192 // Calculate the centre by breaking the cell into pyramids and
193 // volume-weighted averaging their centres
194
195 scalar sumV = 0;
196 vector sumVc = Zero;
197
198 for (const label facei : cFaces)
199 {
200 const face& f = meshFaces[facei];
201
202 scalar pyrVol = pyramidPointFaceRef(f, ctr).mag(meshPoints);
203
204 // if pyramid inside-out because face points inwards invert
205 // N.B. pyramid remains unchanged
206 if (pyrVol < 0)
207 {
208 pyrVol = -pyrVol;
209 }
210
211 sumV += pyrVol;
212 sumVc += pyrVol * pyramidPointFaceRef(f, ctr).centre(meshPoints);
213 }
214
215 return sumVc/(sumV + VSMALL);
216}
217
218
219Foam::scalar Foam::cell::mag
220(
221 const UList<point>& meshPoints,
222 const faceUList& meshFaces
223) const
224{
225 // When one wants to access the cell centre and magnitude, the
226 // functionality on the mesh level should be used in preference to the
227 // functions provided here. They do not rely to the functionality
228 // implemented here, provide additional checking and are more efficient.
229 // The cell::centre and cell::mag functions may be removed in the future.
230
231 // WARNING! See cell::centre
232
233 const labelList& cFaces = *this;
234
235 // Approximate cell centre as the average of all face centres
236
237 vector ctr = Zero;
238 for (const label facei : cFaces)
239 {
240 ctr += meshFaces[facei].centre(meshPoints);
241 }
242 ctr /= cFaces.size();
243
244 // Calculate the magnitude by summing the mags of the pyramids
245 scalar sumV = 0;
246
247 for (const label facei : cFaces)
248 {
249 const face& f = meshFaces[facei];
250
251 sumV += ::Foam::mag(pyramidPointFaceRef(f, ctr).mag(meshPoints));
253
254 return sumV;
255}
256
257
260(
261 const UList<point>& meshPoints,
262 const faceUList& meshFaces
263) const
264{
265 Pair<point> bb(point::rootMax, point::rootMin);
266
267 for (const label facei : *this)
268 {
269 for (const label pointi : meshFaces[facei])
270 {
271 const point& p = meshPoints[pointi];
272
273 bb.first() = min(bb.first(), p);
274 bb.second() = max(bb.second(), p);
276 }
277
278 return bb;
279}
280
281
282Foam::Pair<Foam::point> Foam::cell::box(const primitiveMesh& mesh) const
284 return cell::box(mesh.points(), mesh.faces());
285}
286
287
288// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
289
290bool Foam::operator==(const cell& a, const cell& b)
291{
292 // Trivial reject: faces are different size
293 if (a.size() != b.size())
294 {
295 return false;
296 }
297
298 List<bool> fnd(a.size(), false);
299
300 for (const label curLabel : b)
301 {
302 bool found = false;
303
304 forAll(a, ai)
305 {
306 if (a[ai] == curLabel)
307 {
308 found = true;
309 fnd[ai] = true;
310 break;
311 }
312 }
313
314 if (!found)
315 {
316 return false;
317 }
318 }
319
320 // Any faces missed?
321 forAll(fnd, ai)
322 {
323 if (!fnd[ai])
324 {
325 return false;
326 }
327 }
328
329 return true;
330}
331
332
333// ************************************************************************* //
bool found
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
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
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
Definition UListI.H:468
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition UListI.H:424
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form rootMin
static const Form rootMax
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
labelList labels(const faceUList &meshFaces) const
Return unordered list of cell vertices given the list of faces.
Definition cell.C:33
edgeList edges(const faceUList &meshFaces) const
Return cell edges.
Definition cell.C:105
scalar mag(const UList< point > &meshPoints, const faceUList &meshFaces) const
Returns cell volume.
Definition cell.C:213
point centre(const UList< point > &meshPoints, const faceUList &meshFaces) const
Returns cell centre.
Definition cell.C:149
pointField points(const faceUList &meshFaces, const UList< point > &meshPoints) const
Return the cell vertices given the list of faces and mesh points.
Definition cell.C:87
Pair< point > box(const UList< point > &meshPoints, const faceUList &meshFaces) const
The bounding box for the cell.
Definition cell.C:253
static const char *const typeName
Definition cell.H:61
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Cell-face mesh analysis engine.
scalar mag(const UList< point > &points) const
Return scalar magnitude - returns volume of pyramid.
Definition pyramidI.H:69
Point centre(const UList< point > &points) const
Return centre (centroid).
Definition pyramidI.H:48
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
dynamicFvMesh & mesh
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred point and face.
Definition pyramid.H:70
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< face > faceUList
UList of faces.
Definition faceListFwd.H:43
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299