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 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 "pyramidPointFaceRef.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34const char* const Foam::cell::typeName = "cell";
35
36
37// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38
40{
41 const labelList& cFaces = *this;
42
43 label nVerts = 0;
44 for (const label facei : cFaces)
45 {
46 nVerts += meshFaces[facei].size();
47 }
48
49 labelList pointLabels(nVerts);
50
51 // The first face has no duplicates, can copy in values
52 const labelList& firstFace = meshFaces[cFaces[0]];
53
54 std::copy(firstFace.cbegin(), firstFace.cend(), pointLabels.begin());
55
56 // Now already contains some vertices
57 nVerts = firstFace.size();
58
59 // For the rest of the faces. For each vertex, check if the point is
60 // already inserted (up to nVerts, which now carries the number of real
61 // points. If not, add it at the end of the list.
62
63 for (label facei = 1; facei < cFaces.size(); ++facei)
64 {
65 for (const label curPoint : meshFaces[cFaces[facei]])
66 {
67 bool pointFound = false;
68
69 for (label checki = 0; checki < nVerts; ++checki)
70 {
71 if (curPoint == pointLabels[checki])
72 {
73 pointFound = true;
74 break;
75 }
76 }
77
78 if (!pointFound)
79 {
80 pointLabels[nVerts] = curPoint;
81 ++nVerts;
82 }
83 }
84 }
85
86 pointLabels.resize(nVerts);
87
88 return pointLabels;
89}
90
91
93(
94 const faceUList& meshFaces,
95 const UList<point>& meshPoints
96) const
97{
98 const labelList pointLabels = labels(meshFaces);
99
100 pointField allPoints(pointLabels.size());
101
102 forAll(allPoints, i)
103 {
104 allPoints[i] = meshPoints[pointLabels[i]];
105 }
106
107 return allPoints;
108}
109
110
112{
113 const labelList& cFaces = *this;
114
115 label nEdges = 0;
116 for (const label facei : cFaces)
117 {
118 nEdges += meshFaces[facei].nEdges();
119 }
120
121 edgeList allEdges(nEdges);
122
123 nEdges = 0;
124
125 forAll(cFaces, facei)
126 {
127 for (const edge& curEdge : meshFaces[cFaces[facei]].edges())
128 {
129 bool edgeFound = false;
130
131 for (label checki = 0; checki < nEdges; ++checki)
132 {
133 if (curEdge == allEdges[checki])
134 {
135 edgeFound = true;
136 break;
137 }
138 }
139
140 if (!edgeFound)
141 {
142 allEdges[nEdges] = curEdge;
143 ++nEdges;
144 }
145 }
146 }
147
148 allEdges.resize(nEdges);
149
150 return allEdges;
151}
152
153
155(
156 const UList<point>& meshPoints,
157 const faceUList& meshFaces
158) const
159{
160 // When one wants to access the cell centre and magnitude, the
161 // functionality on the mesh level should be used in preference to the
162 // functions provided here. They do not rely to the functionality
163 // implemented here, provide additional checking and are more efficient.
164 // The cell::centre and cell::mag functions may be removed in the future.
165
166 // WARNING!
167 // In the old version of the code, it was possible to establish whether any
168 // of the pyramids had a negative volume, caused either by the poor
169 // estimate of the cell centre or by the fact that the cell was defined
170 // inside out. In the new description of the cell, this can only be
171 // established with the reference to the owner-neighbour face-cell
172 // relationship, which is not available on this level. Thus, all the
173 // pyramids are believed to be positive with no checking.
174
175 // Approximate cell centre as the area average of all face centres
176
177 vector ctr = Zero;
178 scalar sumArea = 0;
179
180 const labelList& cFaces = *this;
181
182 for (const label facei : cFaces)
183 {
184 const scalar magArea = meshFaces[facei].mag(meshPoints);
185 ctr += meshFaces[facei].centre(meshPoints)*magArea;
186 sumArea += magArea;
187 }
188
189 ctr /= sumArea + VSMALL;
190
191 // Calculate the centre by breaking the cell into pyramids and
192 // volume-weighted averaging their centres
193
194 scalar sumV = 0;
195 vector sumVc = Zero;
196
197 for (const label facei : cFaces)
198 {
199 const face& f = meshFaces[facei];
200
201 scalar pyrVol = pyramidPointFaceRef(f, ctr).mag(meshPoints);
202
203 // if pyramid inside-out because face points inwards invert
204 // N.B. pyramid remains unchanged
205 if (pyrVol < 0)
206 {
207 pyrVol = -pyrVol;
208 }
209
210 sumV += pyrVol;
211 sumVc += pyrVol * pyramidPointFaceRef(f, ctr).centre(meshPoints);
212 }
213
214 return sumVc/(sumV + VSMALL);
215}
216
217
218Foam::scalar Foam::cell::mag
219(
220 const UList<point>& meshPoints,
221 const faceUList& meshFaces
222) const
223{
224 // When one wants to access the cell centre and magnitude, the
225 // functionality on the mesh level should be used in preference to the
226 // functions provided here. They do not rely to the functionality
227 // implemented here, provide additional checking and are more efficient.
228 // The cell::centre and cell::mag functions may be removed in the future.
229
230 // WARNING! See cell::centre
231
232 const labelList& cFaces = *this;
233
234 // Approximate cell centre as the average of all face centres
235
236 vector ctr = Zero;
237 for (const label facei : cFaces)
238 {
239 ctr += meshFaces[facei].centre(meshPoints);
240 }
241 ctr /= cFaces.size();
242
243 // Calculate the magnitude by summing the mags of the pyramids
244 scalar sumV = 0;
245
246 for (const label facei : cFaces)
247 {
248 const face& f = meshFaces[facei];
249
250 sumV += ::Foam::mag(pyramidPointFaceRef(f, ctr).mag(meshPoints));
251 }
252
253 return sumV;
254}
255
256
257// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
258
259bool Foam::operator==(const cell& a, const cell& b)
260{
261 // Trivial reject: faces are different size
262 if (a.size() != b.size())
263 {
264 return false;
265 }
266
267 List<bool> fnd(a.size(), false);
268
269 for (const label curLabel : b)
270 {
271 bool found = false;
272
273 forAll(a, ai)
274 {
275 if (a[ai] == curLabel)
276 {
277 found = true;
278 fnd[ai] = true;
279 break;
280 }
281 }
282
283 if (!found)
284 {
285 return false;
286 }
287 }
288
289 // Any faces missed?
290 forAll(fnd, ai)
291 {
292 if (!fnd[ai])
293 {
294 return false;
295 }
296 }
297
298 return true;
299}
300
301
302// ************************************************************************* //
bool found
DynamicList< label > & labels() const
Return temporary addressing.
Definition: Cloud.H:170
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:67
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
const_iterator cend() const noexcept
Return const_iterator to end traversing the constant UList.
Definition: UListI.H:364
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition: UListI.H:343
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
static const char *const typeName
Definition: cell.H:62
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Computes the magnitude of an input field.
Definition: mag.H:153
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:154
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
const pointField & points
pyramid< point, const point &, const face & > pyramidPointFaceRef
A pyramid using referred points and faces.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
labelList f(nPoints)
labelList pointLabels(nPoints, -1)
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333