polyMeshGenFacesI.H
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 \*---------------------------------------------------------------------------*/
28 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 #include "polyMeshGenFaces.H"
33 
34 # ifdef USE_OMP
35 #include <omp.h>
36 # endif
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 inline const faceListPMG& polyMeshGenFaces::faces() const
44 {
45  return faces_;
46 }
47 
49 {
50  if( !(ownerPtr_ && neighbourPtr_) )
51  {
52  # ifdef USE_OMP
53  if( omp_in_parallel() )
55  (
56  "inline label polyMeshGenFaces::nInternalFaces() const"
57  ) << "Calculating addressing inside a parallel region."
58  << " This is not thread safe" << exit(FatalError);
59  # endif
60 
62  }
63 
64  return nIntFaces_;
65 }
66 
67 inline const labelList& polyMeshGenFaces::owner() const
68 {
69  if( !ownerPtr_ )
70  {
71  # ifdef USE_OMP
72  if( omp_in_parallel() )
74  (
75  "inline label polyMeshGenFaces::owner() const"
76  ) << "Calculating addressing inside a parallel region."
77  << " This is not thread safe" << exit(FatalError);
78  # endif
79 
81  }
82 
83  return *ownerPtr_;
84 }
85 
87 {
88  if( !neighbourPtr_ )
89  {
90  # ifdef USE_OMP
91  if( omp_in_parallel() )
93  (
94  "inline label polyMeshGenFaces::neighbour() const"
95  ) << "Calculating addressing inside a parallel region."
96  << " This is not thread safe" << exit(FatalError);
97  # endif
98 
100  }
101 
102  return *neighbourPtr_;
103 }
104 
107 {
108  return procBoundaries_;
109 }
110 
112 {
113  return boundaries_;
114 }
115 
117 (
118  const label setID,
119  const label faceI
120 )
121 {
122  std::map<label, meshSubset>::iterator it = faceSubsets_.find(setID);
123  if( it == faceSubsets_.end() )
124  return;
125 
126  it->second.addElement(faceI);
127 }
128 
130 (
131  const label setI,
132  const label faceI)
133 {
134  std::map<label, meshSubset>::iterator it = faceSubsets_.find(setI);
135  if( it == faceSubsets_.end() )
136  return;
137 
138  it->second.removeElement(faceI);
139 }
140 
142 (
143  const label faceI,
144  DynList<label>& faceSubsets
145 ) const
146 {
147  faceSubsets.clear();
148 
149  std::map<label, meshSubset>::const_iterator it;
150  for
151  (
152  it=faceSubsets_.begin();
153  it!=faceSubsets_.end();
154  ++it
155  )
156  {
157  if( it->second.contains(faceI) )
158  faceSubsets.append(it->first);
159  }
160 }
161 
163 {
164  indices.clear();
165 
166  std::map<label, meshSubset>::const_iterator it;
167  for
168  (
169  it=faceSubsets_.begin();
170  it!=faceSubsets_.end();
171  ++it
172  )
173  indices.append(it->first);
174 }
175 
176 template<class ListType>
178 (
179  const label setI,
180  ListType& faceLabels
181 ) const
182 {
183  faceLabels.clear();
184 
185  std::map<label, meshSubset>::const_iterator it =
186  faceSubsets_.find(setI);
187  if( it == faceSubsets_.end() )
188  return;
189 
190  it->second.containedElements(faceLabels);
191 }
192 
193 template<class ListType>
194 inline void polyMeshGenFaces::updateFaceSubsets(const ListType& newFaceLabels)
195 {
196  for
197  (
198  std::map<label, meshSubset>::iterator it=faceSubsets_.begin();
199  it!=faceSubsets_.end();
200  ++it
201  )
202  it->second.updateSubset(newFaceLabels);
203 }
204 
205 inline void polyMeshGenFaces::updateFaceSubsets(const VRWGraph& newFacesForFace)
206 {
207  for
208  (
209  std::map<label, meshSubset>::iterator it=faceSubsets_.begin();
210  it!=faceSubsets_.end();
211  ++it
212  )
213  it->second.updateSubset(newFacesForFace);
214 }
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam::polyMeshGenFaces::neighbour
const labelList & neighbour() const
Definition: polyMeshGenFacesI.H:86
Foam::polyMeshGenFaces::owner
const labelList & owner() const
owner and neighbour cells for faces
Definition: polyMeshGenFacesI.H:67
Foam::polyMeshGenFaces::ownerPtr_
labelIOList * ownerPtr_
Definition: polyMeshGenFaces.H:73
Foam::polyMeshGenFaces::removeFaceFromSubset
void removeFaceFromSubset(const label, const label)
Definition: polyMeshGenFacesI.H:130
Foam::polyMeshGenFaces::faceSubsets_
std::map< label, meshSubset > faceSubsets_
face subsets
Definition: polyMeshGenFaces.H:68
Foam::polyMeshGenFaces::facesInSubset
void facesInSubset(const label, ListType &) const
Definition: polyMeshGenFacesI.H:178
Foam::polyMeshGenFaces::updateFaceSubsets
void updateFaceSubsets(const ListType &)
Definition: polyMeshGenFacesI.H:194
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
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::polyMeshGenFaces::procBoundaries_
PtrList< processorBoundaryPatch > procBoundaries_
Definition: polyMeshGenFaces.H:62
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
Foam::polyMeshGenFaces::neighbourPtr_
labelIOList * neighbourPtr_
Definition: polyMeshGenFaces.H:74
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::DynList< label >
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::polyMeshGenFaces::nIntFaces_
label nIntFaces_
number of internal faces, owner and neighbour
Definition: polyMeshGenFaces.H:72
Foam::polyMeshGenFaces::calculateOwnersAndNeighbours
virtual void calculateOwnersAndNeighbours() const =0
calculate owner and neighbour addressing
Foam::polyMeshGenFaces::faceInSubsets
void faceInSubsets(const label, DynList< label > &) const
Definition: polyMeshGenFacesI.H:142
Foam::polyMeshGenFaces::faces_
faceListPMG faces_
list of faces
Definition: polyMeshGenFaces.H:57
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::polyMeshGenFaces::procBoundaries
const PtrList< processorBoundaryPatch > & procBoundaries() const
inter-processor boundaries
Definition: polyMeshGenFacesI.H:106
Foam::polyMeshGenFaces::addFaceToSubset
void addFaceToSubset(const label, const label)
Definition: polyMeshGenFacesI.H:117
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::polyMeshGenFaces::faceSubsetIndices
void faceSubsetIndices(DynList< label > &) const
Definition: polyMeshGenFacesI.H:162
polyMeshGenFaces.H
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::polyMeshGenFaces::boundaries_
PtrList< boundaryPatch > boundaries_
boundary data
Definition: polyMeshGenFaces.H:65
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::DynList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: DynListI.H:279
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304