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 | 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 "Istream.H"
29 #include "cell.H"
30 #include "cellModel.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 :
36  m(nullptr)
37 {}
38 
39 
41 (
42  const cellModel& model,
43  const labelUList& labels,
44  const bool doCollapse
45 )
46 :
47  labelList(labels),
48  m(&model)
49 {
50  if (doCollapse)
51  {
52  collapse();
53  }
54 }
55 
56 
58 (
59  const cellModel& model,
60  labelUList&& labels,
61  const bool doCollapse
62 )
63 :
64  labelList(std::move(labels)),
65  m(&model)
66 {
67  if (doCollapse)
68  {
69  collapse();
70  }
71 }
72 
73 
75 (
76  const word& modelName,
77  const labelUList& labels,
78  const bool doCollapse
79 )
80 :
81  labelList(labels),
82  m(cellModel::ptr(modelName))
83 {
84  if (doCollapse)
85  {
86  collapse();
87  }
88 }
89 
90 
92 {
93  is >> *this;
94 }
95 
96 
98 {
99  return autoPtr<cellShape>::New(*this);
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 (
107  const UList<point>& meshPoints
108 ) const
109 {
110  // There are as many points as there labels for them
111  pointField p(size());
112 
113  // For each point in list, set it to the point in 'pnts' addressed
114  // by 'labs'
115  forAll(p, i)
116  {
117  p[i] = meshPoints[operator[](i)];
118  }
119 
120  // Return list
121  return p;
122 }
123 
124 
126 {
127  return *m;
128 }
129 
130 
132 (
133  const faceList& allFaces,
134  const cell& cFaces
135 ) const
136 {
137  // Faces in model order
138  faceList localFaces(faces());
139 
140  // Do linear match (usually cell shape is low complexity)
141 
142  labelList modelToMesh(localFaces.size(), -1);
143 
144  forAll(localFaces, i)
145  {
146  const face& localF = localFaces[i];
147 
148  forAll(cFaces, j)
149  {
150  label meshFacei = cFaces[j];
151 
152  if (allFaces[meshFacei] == localF)
153  {
154  modelToMesh[i] = meshFacei;
155 
156  break;
157  }
158  }
159  }
160 
161  return modelToMesh;
162 }
163 
164 
166 (
167  const edgeList& allEdges,
168  const labelList& cEdges
169 ) const
170 {
171  // Edges in model order
172  edgeList localEdges(edges());
173 
174  // Do linear match (usually cell shape is low complexity)
175 
176  labelList modelToMesh(localEdges.size(), -1);
177 
178  forAll(localEdges, i)
179  {
180  const edge& e = localEdges[i];
181 
182  forAll(cEdges, j)
183  {
184  label edgeI = cEdges[j];
185 
186  if (allEdges[edgeI] == e)
187  {
188  modelToMesh[i] = edgeI;
189 
190  break;
191  }
192  }
193  }
194 
195  return modelToMesh;
196 }
197 
198 
200 {
201  return m->faces(*this);
202 }
203 
204 
206 {
207  faceList oldFaces(faces());
208 
209  faceList newFaces(oldFaces.size());
210  label newFacei = 0;
211 
212  forAll(oldFaces, oldFacei)
213  {
214  const face& f = oldFaces[oldFacei];
215 
216  face& newF = newFaces[newFacei];
217 
218  newF.setSize(f.size());
219 
220  label newFp = 0;
221  label prevVertI = -1;
222 
223  forAll(f, fp)
224  {
225  label vertI = f[fp];
226 
227  if (vertI != prevVertI)
228  {
229  newF[newFp++] = vertI;
230 
231  prevVertI = vertI;
232  }
233  }
234 
235  if ((newFp > 1) && (newF[newFp-1] == newF[0]))
236  {
237  --newFp;
238  }
239 
240  if (newFp > 2)
241  {
242  // Size face and go to next one
243  newF.setSize(newFp);
244 
245  newFacei++;
246  }
247  }
248  newFaces.setSize(newFacei);
249 
250  return newFaces;
251 }
252 
253 
255 {
256  return m->nFaces();
257 }
258 
259 
261 {
262  return m->edges(*this);
263 }
264 
265 
267 {
268  return m->nEdges();
269 }
270 
271 
273 {
274  return size();
275 }
276 
277 
279 {
280  return m->centre(*this, points);
281 }
282 
283 
284 inline Foam::scalar Foam::cellShape::mag(const UList<point>& points) const
285 {
286  return m->mag(*this, points);
287 }
288 
289 
290 // ************************************************************************* //
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
cell.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::cellShape::centre
point centre(const UList< point > &points) const
Centroid of the cell.
Definition: cellShapeI.H:278
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::cellShape::meshFaces
labelList meshFaces(const faceList &allFaces, const cell &cFaces) const
Mesh face labels of this cell (in order of model)
Definition: cellShapeI.H:132
Foam::cellShape::faces
faceList faces() const
Faces of this cell.
Definition: cellShapeI.H:199
Foam::cellModel::ptr
static const cellModel * ptr(const modelType model)
Look up pointer to cellModel by enumeration, or nullptr on failure.
Definition: cellModels.C:120
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::cellShape::meshEdges
labelList meshEdges(const edgeList &allEdges, const labelList &cEdges) const
Mesh edge labels of this cell (in order of model)
Definition: cellShapeI.H:166
Foam::cellShape::points
pointField points(const UList< point > &meshPoints) const
Return the points corresponding to this cellShape.
Definition: cellShapeI.H:106
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::Field< vector >
Foam::cellShape::edges
edgeList edges() const
Edges of this cellShape.
Definition: cellShapeI.H:260
cellModel.H
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::cellShape::clone
autoPtr< cellShape > clone() const
Clone.
Definition: cellShapeI.H:97
Foam::Vector::centre
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:119
Istream.H
Foam::cellShape::cellShape
cellShape()
Construct null.
Definition: cellShapeI.H:34
Foam::cellShape::mag
scalar mag(const UList< point > &points) const
Scalar magnitude.
Definition: cellShapeI.H:284
Foam::cellShape::nFaces
label nFaces() const
Number of faces.
Definition: cellShapeI.H:254
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cellShape::nEdges
label nEdges() const
Number of edges.
Definition: cellShapeI.H:266
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::UList< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::cellShape::nPoints
label nPoints() const
Number of points.
Definition: cellShapeI.H:272
Foam::cellShape::collapsedFaces
faceList collapsedFaces() const
Collapsed faces of this cell.
Definition: cellShapeI.H:205
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::cellShape::model
const cellModel & model() const
Model reference.
Definition: cellShapeI.H:125