cellShapeI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Istream.H"
27 #include "cell.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 :
33  m(NULL)
34 {}
35 
36 
38 (
39  const cellModel& M,
40  const labelList& l,
41  const bool doCollapse
42 )
43 :
44  labelList(l),
45  m(&M)
46 {
47  if (doCollapse)
48  {
49  collapse();
50  }
51 }
52 
53 
55 {
56  is >> *this;
57 }
58 
59 
61 {
62  return autoPtr<cellShape>(new cellShape(*this));
63 }
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 (
70  const pointField& meshPoints
71 ) const
72 {
73  // There are as many points as there labels for them
74  pointField p(size());
75 
76  // For each point in list, set it to the point in 'pnts' addressed
77  // by 'labs'
78  forAll(p, i)
79  {
80  p[i] = meshPoints[operator[](i)];
81  }
82 
83  // Return list
84  return p;
85 }
86 
87 
89 {
90  return *m;
91 }
92 
93 
95 (
96  const faceList& allFaces,
97  const cell& cFaces
98 ) const
99 {
100  // Faces in model order
101  faceList localFaces(faces());
102 
103  // Do linear match (usually cell shape is low complexity)
104 
105  labelList modelToMesh(localFaces.size(), -1);
106 
107  forAll(localFaces, i)
108  {
109  const face& localF = localFaces[i];
110 
111  forAll(cFaces, j)
112  {
113  label meshFacei = cFaces[j];
114 
115  if (allFaces[meshFacei] == localF)
116  {
117  modelToMesh[i] = meshFacei;
118 
119  break;
120  }
121  }
122  }
123 
124  return modelToMesh;
125 }
126 
127 
129 (
130  const edgeList& allEdges,
131  const labelList& cEdges
132 ) const
133 {
134  // Edges in model order
135  edgeList localEdges(edges());
136 
137  // Do linear match (usually cell shape is low complexity)
138 
139  labelList modelToMesh(localEdges.size(), -1);
140 
141  forAll(localEdges, i)
142  {
143  const edge& e = localEdges[i];
144 
145  forAll(cEdges, j)
146  {
147  label edgeI = cEdges[j];
148 
149  if (allEdges[edgeI] == e)
150  {
151  modelToMesh[i] = edgeI;
152 
153  break;
154  }
155  }
156  }
157 
158  return modelToMesh;
159 }
160 
161 
163 {
164  return m->faces(*this);
165 }
166 
167 
169 {
170  faceList oldFaces(faces());
171 
172  faceList newFaces(oldFaces.size());
173  label newFacei = 0;
174 
175  forAll(oldFaces, oldFacei)
176  {
177  const face& f = oldFaces[oldFacei];
178 
179  face& newF = newFaces[newFacei];
180 
181  newF.setSize(f.size());
182 
183  label newFp = 0;
184  label prevVertI = -1;
185 
186  forAll(f, fp)
187  {
188  label vertI = f[fp];
189 
190  if (vertI != prevVertI)
191  {
192  newF[newFp++] = vertI;
193 
194  prevVertI = vertI;
195  }
196  }
197 
198  if ((newFp > 1) && (newF[newFp-1] == newF[0]))
199  {
200  --newFp;
201  }
202 
203  if (newFp > 2)
204  {
205  // Size face and go to next one
206  newF.setSize(newFp);
207 
208  newFacei++;
209  }
210  }
211  newFaces.setSize(newFacei);
212 
213  return newFaces;
214 }
215 
216 
218 {
219  return m->nFaces();
220 }
221 
222 
224 {
225  return m->edges(*this);
226 }
227 
228 
230 {
231  return m->nEdges();
232 }
233 
234 
236 {
237  return size();
238 }
239 
240 
242 {
243  return m->centre(*this, points);
244 }
245 
246 
247 inline Foam::scalar Foam::cellShape::mag(const pointField& points) const
248 {
249  return m->mag(*this, points);
250 }
251 
252 
253 // ************************************************************************* //
label nFaces() const
Return number of faces.
Definition: cellModelI.H:62
const cellModel & model() const
Model reference.
Definition: cellShapeI.H:88
scalar mag(const pointField &) const
Scalar magnitude.
Definition: cellShapeI.H:247
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const double e
Elementary charge.
Definition: doubleFloat.H:78
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
void collapse()
Collapse shape to correct one after removing duplicate vertices.
Definition: cellShape.C:31
label nFaces() const
Number of faces.
Definition: cellShapeI.H:217
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
label nEdges() const
Number of edges.
Definition: cellShapeI.H:229
faceList faces(const labelList &pointLabels) const
Return list of faces.
Definition: cellModelI.H:97
scalar mag(const labelList &pointLabels, const pointField &points) const
Cell volume.
Definition: cellModel.C:89
vector centre(const labelList &pointLabels, const pointField &points) const
Vector centroid.
Definition: cellModel.C:32
faceList collapsedFaces() const
Collapsed faces of this cell.
Definition: cellShapeI.H:168
cellShape()
Construct null.
Definition: cellShapeI.H:31
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
autoPtr< cellShape > clone() const
Clone.
Definition: cellShapeI.H:60
edgeList edges(const labelList &pointLabels) const
Return list of edges.
Definition: cellModelI.H:70
pointField points(const pointField &meshPoints) const
Return the points corresponding to this cellShape.
Definition: cellShapeI.H:69
List< label > labelList
A List of labels.
Definition: labelList.H:56
edgeList edges() const
Edges of this cellShape.
Definition: cellShapeI.H:223
labelList meshFaces(const faceList &allFaces, const cell &) const
Mesh face labels of this cell (in order of model)
Definition: cellShapeI.H:95
labelList f(nPoints)
label size() const
Return the number of elements in the UList.
void setSize(const label)
Reset size of List.
Definition: List.C:295
faceList faces() const
Faces of this cell.
Definition: cellShapeI.H:162
point centre(const pointField &) const
Centroid of the cell.
Definition: cellShapeI.H:241
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Maps a geometry to a set of cell primitives, which enables geometric cell data to be calculated witho...
Definition: cellModel.H:64
label nEdges() const
Return number of edges.
Definition: cellModelI.H:56
label nPoints() const
Number of points.
Definition: cellShapeI.H:235
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
volScalarField & p
labelList meshEdges(const edgeList &allEdges, const labelList &) const
Mesh edge labels of this cell (in order of model)
Definition: cellShapeI.H:129