fvMeshSubset.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-2014 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 Class
25  Foam::fvMeshSubset
26 
27 Description
28  Post-processing mesh subset tool. Given the original mesh and the
29  list of selected cells, it creates the mesh consisting only of the
30  desired cells, with the mapping list for points, faces, and cells.
31 
32  Puts all exposed internal faces into either
33  - a user supplied patch
34  - a newly created patch "oldInternalFaces"
35 
36  - setCellSubset is for small subsets. Uses Maps to minimize memory.
37  - setLargeCellSubset is for largish subsets (>10% of mesh).
38  Uses labelLists instead.
39 
40  - setLargeCellSubset does coupled patch subsetting as well. If it detects
41  a face on a coupled patch 'losing' its neighbour it will move the
42  face into the oldInternalFaces patch.
43 
44  - if a user supplied patch is used it is up to the destination
45  patchField to handle exposed internal faces (mapping from face -1).
46  If not provided the default is to assign the internalField. All the
47  basic patch field types (e.g. fixedValue) will give a warning and
48  preferably derived patch field types should be used that know how to
49  handle exposed faces (e.g. use uniformFixedValue instead of fixedValue)
50 
51 SourceFiles
52  fvMeshSubset.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef fvMeshSubset_H
57 #define fvMeshSubset_H
58 
59 #include "fvMesh.H"
60 #include "pointMesh.H"
61 #include "GeometricField.H"
62 #include "HashSet.H"
63 #include "surfaceMesh.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 /*---------------------------------------------------------------------------*\
71  Class fvMeshSubset Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class fvMeshSubset
75 {
76 
77 private:
78 
79  // Private data
80 
81  //- Mesh to subset from
82  const fvMesh& baseMesh_;
83 
84  //- Subset mesh pointer
86 
87  //- Point mapping array
89 
90  //- Face mapping array
92 
93  //- Cell mapping array
95 
96  //- Patch mapping array
98 
100 
101 
102  // Private Member Functions
103 
104  //- Check if subset has been performed
105  bool checkCellSubset() const;
106 
107  //- Mark points in Map
108  static void markPoints(const labelList&, Map<label>&);
109 
110  //- Mark points (with 0) in labelList
111  static void markPoints(const labelList&, labelList&);
112 
113  //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
114  void doCoupledPatches
115  (
116  const bool syncPar,
117  labelList& nCellsUsingFace
118  ) const;
119 
120  //- Subset of subset
121  static labelList subset
122  (
123  const label nElems,
124  const labelList& selectedElements,
125  const labelList& subsetMap
126  );
127 
128  //- Create zones for submesh
129  void subsetZones();
130 
131  //- Helper: extract cells-to-remove from cells-to-keep
133  (
134  const labelList& region,
135  const label currentRegion
136  ) const;
137 
138  //- Disallow default bitwise copy construct
139  fvMeshSubset(const fvMeshSubset&);
140 
141  //- Disallow default bitwise assignment
142  void operator=(const fvMeshSubset&);
143 
144 public:
145 
146  // Constructors
147 
148  //- Construct given a mesh to subset
149  explicit fvMeshSubset(const fvMesh&);
150 
151 
152  // Member Functions
153 
154  // Edit
155 
156  //- Set the subset. Create "oldInternalFaces" patch for exposed
157  // internal faces (patchID==-1) or use supplied patch.
158  // Does not handle coupled patches correctly if only one side
159  // gets deleted.
160  void setCellSubset
161  (
162  const labelHashSet& globalCellMap,
163  const label patchID = -1
164  );
165 
166  //- Set the subset from all cells with region == currentRegion.
167  // Create "oldInternalFaces" patch for exposed
168  // internal faces (patchID==-1) or use supplied patch.
169  // Handles coupled patches by if necessary making coupled patch
170  // face part of patchID (so uncoupled)
171  void setLargeCellSubset
172  (
173  const labelList& region,
174  const label currentRegion,
175  const label patchID = -1,
176  const bool syncCouples = true
177  );
178 
179  //- setLargeCellSubset but with labelHashSet.
180  void setLargeCellSubset
181  (
182  const labelHashSet& globalCellMap,
183  const label patchID = -1,
184  const bool syncPar = true
185  );
186 
187 
188  //- Two step subsetting
189 
190  //- Get labels of exposed faces.
191  // These are
192  // - internal faces that become boundary faces
193  // - coupled faces that become uncoupled (since one of the
194  // sides gets deleted)
196  (
197  const labelList& region,
198  const label currentRegion,
199  const bool syncCouples = true
200  ) const;
201 
202  //- For every exposed face (from above getExposedFaces)
203  // used supplied (existing!) patch
204  void setLargeCellSubset
205  (
206  const labelList& region,
207  const label currentRegion,
208  const labelList& exposedFaces,
209  const labelList& patchIDs,
210  const bool syncCouples = true
211  );
212 
213 
214  // Access
215 
216  //- Original mesh
217  const fvMesh& baseMesh() const
218  {
219  return baseMesh_;
220  }
221 
222  //- Have subMesh?
223  bool hasSubMesh() const;
224 
225  //- Return reference to subset mesh
226  const fvMesh& subMesh() const;
227 
228  fvMesh& subMesh();
229 
230  //- Return point map
231  const labelList& pointMap() const;
232 
233  //- Return face map
234  const labelList& faceMap() const;
235 
236  //- Return face map with sign to encode flipped faces
237  const labelList& faceFlipMap() const;
238 
239  //- Return cell map
240  const labelList& cellMap() const;
241 
242  //- Return patch map
243  const labelList& patchMap() const;
244 
245 
246  // Field mapping
247 
248  //- Map volume field
249  template<class Type>
252  (
254  const fvMesh& sMesh,
255  const labelList& patchMap,
256  const labelList& cellMap,
257  const labelList& faceMap
258  );
259 
260  template<class Type>
263  (
265  ) const;
266 
267  //- Map surface field. Optionally negates value if flipping
268  // a face (from exposing an internal face)
269  template<class Type>
272  (
274  const fvMesh& sMesh,
275  const labelList& patchMap,
276  const labelList& cellMap,
277  const labelList& faceMap,
278  const bool negateIfFlipped = true
279  );
280 
281  template<class Type>
284  (
286  const bool negateIfFlipped = true
287  ) const;
288 
289  //- Map point field
290  template<class Type>
293  (
295  const pointMesh& sMesh,
296  const labelList& patchMap,
297  const labelList& pointMap
298  );
299 
300  template<class Type>
303  (
305  ) const;
306 
307  //- Map dimensioned fields
308  template<class Type>
311  (
313  const fvMesh& sMesh,
314  const labelList& cellMap
315  );
316 
317  template<class Type>
320  (
322  ) const;
323 
324 };
325 
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 } // End namespace Foam
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 #ifdef NoRepository
334 # include "fvMeshSubsetInterpolate.C"
335 #endif
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 #endif
340 
341 // ************************************************************************* //
Foam::fvMeshSubset::setLargeCellSubset
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:795
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Definition: fvMeshSubsetInterpolate.C:43
Foam::fvMeshSubset::patchMap_
labelList patchMap_
Patch mapping array.
Definition: fvMeshSubset.H:96
Foam::fvMeshSubset::fvMeshSubsetPtr_
autoPtr< fvMesh > fvMeshSubsetPtr_
Subset mesh pointer.
Definition: fvMeshSubset.H:84
Foam::fvMeshSubset::doCoupledPatches
void doCoupledPatches(const bool syncPar, labelList &nCellsUsingFace) const
Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
Definition: fvMeshSubset.C:96
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::fvMeshSubset::subsetZones
void subsetZones()
Create zones for submesh.
Definition: fvMeshSubset.C:243
Foam::fvMeshSubset::operator=
void operator=(const fvMeshSubset &)
Disallow default bitwise assignment.
Foam::Map< label >
Foam::fvMeshSubset::getCellsToRemove
labelList getCellsToRemove(const labelList &region, const label currentRegion) const
Helper: extract cells-to-remove from cells-to-keep.
Definition: fvMeshSubset.C:360
Foam::HashSet< label, Hash< label > >
Foam::fvMeshSubset::faceMap_
labelList faceMap_
Face mapping array.
Definition: fvMeshSubset.H:90
Foam::fvMeshSubset::cellMap
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubset.C:1542
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::fvMeshSubset::fvMeshSubset
fvMeshSubset(const fvMeshSubset &)
Disallow default bitwise copy construct.
Foam::fvMeshSubset::setCellSubset
void setCellSubset(const labelHashSet &globalCellMap, const label patchID=-1)
Set the subset. Create "oldInternalFaces" patch for exposed.
Definition: fvMeshSubset.C:410
fvMeshSubsetInterpolate.C
Foam::fvMeshSubset::baseMesh_
const fvMesh & baseMesh_
Mesh to subset from.
Definition: fvMeshSubset.H:81
Foam::fvMeshSubset::cellMap_
labelList cellMap_
Cell mapping array.
Definition: fvMeshSubset.H:93
Foam::fvMeshSubset::pointMap
const labelList & pointMap() const
Return point map.
Definition: fvMeshSubset.C:1488
Foam::fvMeshSubset::faceFlipMapPtr_
autoPtr< labelList > faceFlipMapPtr_
Definition: fvMeshSubset.H:98
Foam::fvMeshSubset::faceFlipMap
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
Definition: fvMeshSubset.C:1504
Foam::fvMeshSubset::getExposedFaces
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Two step subsetting.
Definition: fvMeshSubset.C:1405
HashSet.H
Foam::fvMeshSubset::subset
static labelList subset(const label nElems, const labelList &selectedElements, const labelList &subsetMap)
Subset of subset.
Definition: fvMeshSubset.C:204
surfaceMesh.H
Foam::fvMeshSubset::checkCellSubset
bool checkCellSubset() const
Check if subset has been performed.
Definition: fvMeshSubset.C:48
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::fvMeshSubset::patchMap
const labelList & patchMap() const
Return patch map.
Definition: fvMeshSubset.C:1550
GeometricField.H
Foam::fvMeshSubset::markPoints
static void markPoints(const labelList &, Map< label > &)
Mark points in Map.
Definition: fvMeshSubset.C:67
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
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::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubset.C:1472
Foam::fvMeshSubset::pointMap_
labelList pointMap_
Point mapping array.
Definition: fvMeshSubset.H:87
Foam::fvMeshSubset::baseMesh
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:216
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvMeshSubset::hasSubMesh
bool hasSubMesh() const
Have subMesh?
Definition: fvMeshSubset.C:1466
pointMesh.H
Foam::fvMeshSubset::faceMap
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubset.C:1496
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51