checkBoundaryFacesSharingTwoEdges.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "demandDrivenData.H"
30 #include "helperFunctionsPar.H"
31 #include "meshSurfaceEngine.H"
32 #include "decomposeFaces.H"
33 #include "decomposeCells.H"
34 
35 # ifdef USE_OMP
36 #include <omp.h>
37 # endif
38 
39 //#define DEBUGCheck
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * Private member functions * * * * * * * * * * * * * * * //
47 
49 {
51 }
52 
54 {
55  const meshSurfaceEngine& mse = meshSurface();
56 
57  const labelList& bp = mse.bp();
58  const edgeList& edges = mse.edges();
59  const VRWGraph& pointEdges = mse.boundaryPointEdges();
60 
61  const label nIntFaces = mesh_.nInternalFaces();
62  const faceListPMG& faces = mesh_.faces();
63 
64  //- find the internal faces attached to the boundary points
65  removeBndPoint_.setSize(pointEdges.size());
66  removeBndPoint_ = true;
67 
68  # ifdef USE_OMP
69  # pragma omp parallel for if( nIntFaces > 100 ) schedule(dynamic, 20)
70  # endif
71  for(label fI=0;fI<nIntFaces;++fI)
72  {
73  const face& f = faces[fI];
74 
75  forAll(f, pI)
76  {
77  const label bpI = bp[f[pI]];
78 
79  if( bpI < 0 )
80  continue;
81 
82  if( nBndFacesAtBndPoint_[bpI] == 2 )
83  {
84  const edge ePrev = f.faceEdge(f.rcIndex(pI));
85  const edge eNext = f.faceEdge(pI);
86 
87  bool foundNext(false), foundPrev(false);
88  forAllRow(pointEdges, bpI, peI)
89  {
90  const label beI = pointEdges(bpI, peI);
91 
92  if( edges[beI] == ePrev )
93  {
94  foundPrev = true;
95  }
96  else if( edges[beI] == eNext )
97  {
98  foundNext = true;
99  }
100  }
101 
102  if( !(foundPrev && foundNext) )
103  {
104  removeBndPoint_[bpI] = false;
105  }
106  }
107  else
108  {
109  removeBndPoint_[bpI] = false;
110  }
111  }
112  }
113 
114  if( Pstream::parRun() )
115  {
116  //- check processor faces
117  const PtrList<processorBoundaryPatch>& procBoundaries =
119  forAll(procBoundaries, patchI)
120  {
121  const label start = procBoundaries[patchI].patchStart();
122  const label end = start + procBoundaries[patchI].patchSize();
123 
124  # ifdef USE_OMP
125  # pragma omp parallel for schedule(dynamic, 10)
126  # endif
127  for(label faceI=start;faceI<end;++faceI)
128  {
129  const face& f = faces[faceI];
130 
131  forAll(f, pI)
132  {
133  const label bpI = bp[f[pI]];
134 
135  if( bpI < 0 )
136  continue;
137 
138  if( nBndFacesAtBndPoint_[bpI] == 2 )
139  {
140  const edge ePrev = f.faceEdge(f.rcIndex(pI));
141  const edge eNext = f.faceEdge(pI);
142 
143  bool foundNext(false), foundPrev(false);
144  forAllRow(pointEdges, bpI, peI)
145  {
146  const label beI = pointEdges(bpI, peI);
147 
148  if( edges[beI] == ePrev )
149  {
150  foundPrev = true;
151  }
152  else if( edges[beI] == eNext )
153  {
154  foundNext = true;
155  }
156  }
157 
158  if( !(foundPrev && foundNext) )
159  removeBndPoint_[bpI] = false;
160  }
161  else
162  {
163  removeBndPoint_[bpI] = false;
164  }
165  }
166  }
167  }
168 
169  //- make sure that all processors have the same information
170  const DynList<label>& bpNei = mse.bpNeiProcs();
171  const VRWGraph& bpAtProcs = mse.bpAtProcs();
172  const Map<label>& globalToLocal = mse.globalToLocalBndPointAddressing();
173 
174  std::map<label, labelLongList> exchangeData;
175  forAll(bpNei, i)
176  exchangeData.insert(std::make_pair(bpNei[i], labelLongList()));
177 
178  forAllConstIter(Map<label>, globalToLocal, it)
179  {
180  const label bpI = it();
181 
182  if( removeBndPoint_[bpI] )
183  continue;
184 
185  //- the point shall not be removed
186  forAllRow(bpAtProcs, bpI, i)
187  {
188  const label neiProc = bpAtProcs(bpI, i);
189  if( neiProc == Pstream::myProcNo() )
190  continue;
191 
192  exchangeData[neiProc].append(it.key());
193  }
194  }
195 
196  //- exchange data
197  labelLongList receivedData;
198  help::exchangeMap(exchangeData, receivedData);
199 
200  //- set remove flag to false
201  forAll(receivedData, i)
202  removeBndPoint_[globalToLocal[receivedData[i]]] = false;
203  }
204 }
205 
207 {
208  const meshSurfaceEngine& mse = meshSurface();
209  const VRWGraph& pointFaces = mse.pointFaces();
210 
211  nBndFacesAtBndPoint_.setSize(pointFaces.size());
213 
215  nBndFacesAtBndPoint_[bpI] = pointFaces.sizeOfRow(bpI);
216 
217  if( Pstream::parRun() )
218  {
219  const VRWGraph& bpAtProcs = mse.bpAtProcs();
220  const Map<label>& globalToLocal = mse.globalToLocalBndPointAddressing();
221  const DynList<label>& neiProcs = mse.bpNeiProcs();
222 
223  //- create data that shall be exhcnaged
224  std::map<label, labelLongList> exchangeData;
225  forAll(neiProcs, i)
226  exchangeData.insert(std::make_pair(neiProcs[i], labelLongList()));
227 
228  forAllConstIter(Map<label>, globalToLocal, it)
229  {
230  const label bpI = it();
231 
232  forAllRow(bpAtProcs, bpI, i)
233  {
234  const label neiProc = bpAtProcs(bpI, i);
235  if( neiProc == Pstream::myProcNo() )
236  continue;
237 
238  labelLongList& data = exchangeData[neiProc];
239  data.append(it.key());
240  data.append(nBndFacesAtBndPoint_[bpI]);
241  }
242  }
243 
244  //- exchange data with other processors
245  labelLongList receivedData;
246  help::exchangeMap(exchangeData, receivedData);
247 
248  label counter(0);
249  while( counter < receivedData.size() )
250  {
251  const label bpI = globalToLocal[receivedData[counter++]];
252  nBndFacesAtBndPoint_[bpI] += receivedData[counter++];
253  }
254  }
255 }
256 
258 {
259  const labelList& bp = meshSurface().bp();
260  const faceListPMG& faces = mesh_.faces();
261 
262  //- remove points which can be safely be removed
263  //- internal faces
264  const label nIntFaces = mesh_.nInternalFaces();
265 
266  # ifdef USE_OMP
267  # pragma omp parallel for if( nIntFaces > 100 ) schedule(dynamic, 10)
268  # endif
269  for(label faceI=0;faceI<nIntFaces;++faceI)
270  {
271  const face& f = faces[faceI];
272 
273  DynList<label> newF;
274  forAll(f, pI)
275  {
276  const label bpI = bp[f[pI]];
277 
278  if(
279  (bpI >= 0) && removeBndPoint_[bpI] &&
280  (nBndFacesAtBndPoint_[bpI] == 2)
281  )
282  continue;
283 
284  newF.append(f[pI]);
285  }
286 
287  if( newF.size() < f.size() )
288  {
289  face& mf = const_cast<face&>(f);
290  mf.setSize(newF.size());
291  forAll(mf, i)
292  mf[i] = newF[i];
293  }
294  }
295 
296  //- boundary faces
297  forAll(mesh_.boundaries(), patchI)
298  {
299  const label start = mesh_.boundaries()[patchI].patchStart();
300  const label end = start + mesh_.boundaries()[patchI].patchSize();
301 
302  # ifdef USE_OMP
303  # pragma omp parallel for if( end - start > 100 ) \
304  schedule(dynamic, 10)
305  # endif
306  for(label faceI=start;faceI<end;++faceI)
307  {
308  const face& f = faces[faceI];
309 
310  DynList<label> newF;
311  forAll(f, pI)
312  {
313  const label bpI = bp[f[pI]];
314 
315  if( removeBndPoint_[bpI] && (nBndFacesAtBndPoint_[bpI] == 2) )
316  continue;
317 
318  newF.append(f[pI]);
319  }
320 
321  if( newF.size() < f.size() )
322  {
323  face& mf = const_cast<face&>(f);
324  mf.setSize(newF.size());
325  forAll(mf, i)
326  mf[i] = newF[i];
327  }
328  }
329  }
330 
331  //- processor boundaries
332  forAll(mesh_.procBoundaries(), patchI)
333  {
334  const processorBoundaryPatch& patch = mesh_.procBoundaries()[patchI];
335  const label start = patch.patchStart();
336  const label end = start + patch.patchSize();
337 
338  # ifdef USE_OMP
339  # pragma omp parallel for if( patch.patchSize() > 100 ) \
340  schedule(dynamic, 10)
341  # endif
342  for(label faceI=start;faceI<end;++faceI)
343  {
344  const face& f = faces[faceI];
345 
346  DynList<label> newF;
347  forAll(f, pI)
348  {
349  const label bpI = bp[f[pI]];
350 
351  if(
352  (bpI >= 0) && removeBndPoint_[bpI] &&
353  (nBndFacesAtBndPoint_[bpI] == 2)
354  )
355  continue;
356 
357  newF.append(f[pI]);
358  }
359 
360  if( newF.size() < f.size() )
361  {
362  face& mf = const_cast<face&>(f);
363  mf.setSize(newF.size());
364 
365  if( !patch.owner() && (newF[0] != f[0]) )
366  {
367  forAll(mf, i)
368  mf[i] = newF[mf.rcIndex(i)];
369  }
370  else
371  {
372  forAll(mf, i)
373  mf[i] = newF[i];
374  }
375  }
376  }
377  }
378 }
379 
381 (
382  boolList& decomposeFace
383 )
384 {
385  const meshSurfaceEngine& mse = meshSurface();
386  const labelList& bp = mse.bp();
387  const faceList::subList& bFaces = mse.boundaryFaces();
388 
389  label nDecomposed(0);
390  const label nIntFaces = mesh_.nInternalFaces();
391 
392  # ifdef USE_OMP
393  # pragma omp parallel for if( bFaces.size() > 100 ) \
394  schedule(dynamic, 10) reduction(+ : nDecomposed)
395  # endif
396  forAll(bFaces, bfI)
397  {
398  const face& bf = bFaces[bfI];
399 
400  forAll(bf, pI)
401  {
402  const label bpI = bp[bf[pI]];
403 
404  if( nBndFacesAtBndPoint_[bpI] == 2 )
405  {
406  ++nDecomposed;
407  decomposeFace[nIntFaces+bfI] = true;
408  }
409  }
410  }
411 
412  reduce(nDecomposed, sumOp<label>());
413 
414  return nDecomposed;
415 }
416 
417 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
418 // Constructors
419 
421 (
423 )
424 :
425  mesh_(mesh),
426  meshSurfacePtr_(NULL),
427  nBndFacesAtBndPoint_(),
428  removeBndPoint_()
429 {
430 }
431 
432 // * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * * * //
433 
435 {
437 }
438 
439 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
440 
442 {
443  badPoints.clear();
444 
446 
447  const labelList& bPoints = meshSurface().boundaryPoints();
449  {
450  if( nBndFacesAtBndPoint_[bpI] != 2 )
451  continue;
452 
453  badPoints.insert(bPoints[bpI]);
454  }
455 }
456 
458 {
459  bool changed(false);
460 
462 
464 
466 
467  boolList decomposeFace(mesh_.faces().size(), false);
468  const label nDecomposed = findBndFacesForDecomposition(decomposeFace);
469 
470  Info << "Marked " << nDecomposed << " faces for decomposition" << endl;
471 
472  if( nDecomposed != 0 )
473  {
474  //- delete the mesh surface engine
476 
477  //- find cells which will be decomposed
478  boolList decomposeCell(mesh_.cells().size(), false);
479  const labelList& owner = mesh_.owner();
480  forAll(decomposeFace, faceI)
481  {
482  if( decomposeFace[faceI] )
483  decomposeCell[owner[faceI]];
484  }
485 
486  //- decompose marked faces
487  decomposeFaces(mesh_).decomposeMeshFaces(decomposeFace);
488 
489  //- decompose cells
490  VRWGraph pRegions(mesh_.points().size());
491  decomposeCells dc(mesh_);
492  dc.decomposeMesh(decomposeCell);
493 
494  changed = true;
495  }
496 
498 
499  return changed;
500 }
501 
502 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503 
504 } // End namespace Foam
505 
506 // ************************************************************************* //
Foam::meshSurfaceEngine::boundaryPointEdges
const VRWGraph & boundaryPointEdges() const
Definition: meshSurfaceEngineI.H:315
Foam::meshSurfaceEngine::bpAtProcs
const VRWGraph & bpAtProcs() const
processors which contain the vertex
Definition: meshSurfaceEngineI.H:451
Foam::polyMeshGenFaces::owner
const labelList & owner() const
owner and neighbour cells for faces
Definition: polyMeshGenFacesI.H:67
Foam::meshSurfaceEngine::bp
const labelList & bp() const
Definition: meshSurfaceEngineI.H:64
Foam::decomposeFaces
Definition: decomposeFaces.H:47
Foam::boundaryPatchBase::patchSize
label patchSize() const
Definition: boundaryPatchBase.H:161
helperFunctionsPar.H
Foam::decomposeCells::decomposeMesh
void decomposeMesh(const boolList &)
perform decomposition of selected cell into pyramids
Definition: decomposeCellsDecomposition.C:42
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
decomposeCells.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::meshSurfaceEngine::bpNeiProcs
const DynList< label > & bpNeiProcs() const
communication matrix for sending point data
Definition: meshSurfaceEngineI.H:470
Foam::meshSurfaceEngine::pointFaces
const VRWGraph & pointFaces() const
Definition: meshSurfaceEngineI.H:162
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
decomposeFaces.H
Foam::Map< label >
Foam::checkBoundaryFacesSharingTwoEdges::findBndFacesAtBndVertex
void findBndFacesAtBndVertex()
find the number of faces connected to the boundary vertex
Definition: checkBoundaryFacesSharingTwoEdges.C:206
Foam::polyMeshGen
Definition: polyMeshGen.H:46
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMeshGenPoints::points
const pointFieldPMG & points() const
access to points
Definition: polyMeshGenPointsI.H:44
Foam::checkBoundaryFacesSharingTwoEdges::improveTopology
bool improveTopology()
Definition: checkBoundaryFacesSharingTwoEdges.C:457
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::checkBoundaryFacesSharingTwoEdges::nBndFacesAtBndPoint_
labelList nBndFacesAtBndPoint_
number of boundary faces attached to a boundary vertex
Definition: checkBoundaryFacesSharingTwoEdges.H:66
Foam::HashSet< label, Hash< label > >
Foam::checkBoundaryFacesSharingTwoEdges::~checkBoundaryFacesSharingTwoEdges
~checkBoundaryFacesSharingTwoEdges()
Definition: checkBoundaryFacesSharingTwoEdges.C:434
Foam::polyMeshGenFaces::faces
const faceListPMG & faces() const
access to faces
Definition: polyMeshGenFacesI.H:43
Foam::polyMeshGenFaces::nInternalFaces
label nInternalFaces() const
return number of internal faces
Definition: polyMeshGenFacesI.H:48
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::LongList< label >
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::meshSurfaceEngine::boundaryFaces
const faceList::subList & boundaryFaces() const
Definition: meshSurfaceEngineI.H:103
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:59
Foam::Info
messageStream Info
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::VRWGraph::size
label size() const
Returns the number of rows.
Definition: VRWGraphI.H:122
Foam::checkBoundaryFacesSharingTwoEdges::findPoints
void findPoints(labelHashSet &badPoints)
find boundary points connected to two boundary faces, only
Definition: checkBoundaryFacesSharingTwoEdges.C:441
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
Foam::processorBoundaryPatch::owner
bool owner() const
check if the processor is the owner of the interface
Definition: processorBoundaryPatch.H:124
Foam::pointFieldPMG::size
label size() const
return the number of used elements
Definition: pointFieldPMGI.H:71
Foam::polyMeshGenModifier::removeUnusedVertices
void removeUnusedVertices()
remove unused vertices
Definition: polyMeshGenModifierRemoveUnusedVertices.C:37
Foam::checkBoundaryFacesSharingTwoEdges::meshSurface
const meshSurfaceEngine & meshSurface() const
returns mesh surface
Definition: checkBoundaryFacesSharingTwoEdges.H:76
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::checkBoundaryFacesSharingTwoEdges::mesh_
polyMeshGen & mesh_
Reference to polyMeshGen.
Definition: checkBoundaryFacesSharingTwoEdges.H:60
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
checkBoundaryFacesSharingTwoEdges.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
meshSurfaceEngine.H
Foam::processorBoundaryPatch
Definition: processorBoundaryPatch.H:47
Foam::checkBoundaryFacesSharingTwoEdges::findBndFacesForDecomposition
label findBndFacesForDecomposition(boolList &decomposeFace)
find boundary faces which shall be decomposed into triangles
Definition: checkBoundaryFacesSharingTwoEdges.C:381
Foam::DynList< label >
Foam::checkBoundaryFacesSharingTwoEdges::meshSurfacePtr_
meshSurfaceEngine * meshSurfacePtr_
pointer to meshSurfaceEngine
Definition: checkBoundaryFacesSharingTwoEdges.H:63
Foam::decomposeFaces::decomposeMeshFaces
void decomposeMeshFaces(const boolList &decomposeFace)
decompose selected faces into triangles using midnode subdivision
Definition: decomposeFaces.C:62
Foam::meshSurfaceEngine::edges
const edgeList & edges() const
Definition: meshSurfaceEngineI.H:296
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::checkBoundaryFacesSharingTwoEdges::createMeshSurface
void createMeshSurface() const
creates meshSurfaceEngine
Definition: checkBoundaryFacesSharingTwoEdges.C:48
Foam::checkBoundaryFacesSharingTwoEdges::findFacesAtBndEdge
void findFacesAtBndEdge()
Definition: checkBoundaryFacesSharingTwoEdges.C:53
Foam::decomposeCells
Definition: decomposeCells.H:48
Foam::sumOp
Definition: ops.H:162
Foam::checkBoundaryFacesSharingTwoEdges::removeExcessiveVertices
void removeExcessiveVertices()
remove vertices from the faces
Definition: checkBoundaryFacesSharingTwoEdges.C:257
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:473
f
labelList f(nPoints)
Foam::faceListPMG::size
label size() const
return the number of used elements
Definition: faceListPMGI.H:73
Foam::boundaryPatchBase::patchStart
label patchStart() const
Definition: boundaryPatchBase.H:151
Foam::polyMeshGenFaces::boundaries
const PtrList< boundaryPatch > & boundaries() const
ordinary boundaries
Definition: polyMeshGenFacesI.H:111
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::help::exchangeMap
void exchangeMap(const std::map< label, ListType > &m, LongList< T > &data, const Pstream::commsTypes commsType)
Definition: helperFunctionsPar.C:129
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::polyMeshGenFaces::procBoundaries
const PtrList< processorBoundaryPatch > & procBoundaries() const
inter-processor boundaries
Definition: polyMeshGenFacesI.H:106
Foam::meshSurfaceEngine::boundaryPoints
const labelList & boundaryPoints() const
Definition: meshSurfaceEngineI.H:84
Foam::checkBoundaryFacesSharingTwoEdges::removeBndPoint_
boolList removeBndPoint_
a list of boundary points which can be removed from the mesh
Definition: checkBoundaryFacesSharingTwoEdges.H:69
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::checkBoundaryFacesSharingTwoEdges::checkBoundaryFacesSharingTwoEdges
checkBoundaryFacesSharingTwoEdges(const checkBoundaryFacesSharingTwoEdges &)
Disallow default bitwise copy construct.
Foam::DynList::size
label size() const
Definition: DynListI.H:235
Foam::meshSurfaceEngine::globalToLocalBndPointAddressing
const Map< label > & globalToLocalBndPointAddressing() const
global point label to local label. Only for processors points
Definition: meshSurfaceEngineI.H:431
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::labelLongList
LongList< label > labelLongList
Definition: labelLongList.H:46
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
Foam::meshSurfaceEngine
Definition: meshSurfaceEngine.H:54
Foam::cellListPMG::size
label size() const
return the number of used elements
Definition: cellListPMGI.H:56
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304