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-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 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
85  autoPtr<fvMesh> fvMeshSubsetPtr_;
86 
87  //- Point mapping array
88  labelList pointMap_;
89 
90  //- Face mapping array
91  labelList faceMap_;
92 
93  //- Cell mapping array
94  labelList cellMap_;
95 
96  //- Patch mapping array
97  labelList patchMap_;
98 
99  //- Optional face mapping array with flip encoded
100  mutable autoPtr<labelList> faceFlipMapPtr_;
101 
102 
103  // Private Member Functions
104 
105  //- Check if subset has been performed
106  bool checkCellSubset() const;
107 
108  //- Mark points in Map
109  static void markPoints(const labelList&, Map<label>&);
110 
111  //- Mark points (with 0) in labelList
112  static void markPoints(const labelList&, labelList&);
113 
114  //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
115  void doCoupledPatches
116  (
117  const bool syncPar,
118  labelList& nCellsUsingFace
119  ) const;
120 
121  //- Subset of subset
122  static labelList subset
123  (
124  const label nElems,
125  const labelList& selectedElements,
126  const labelList& subsetMap
127  );
128 
129  //- Create zones for submesh
130  void subsetZones();
131 
132  //- Helper: extract cells-to-remove from cells-to-keep
133  labelList getCellsToRemove
134  (
135  const labelList& region,
136  const label currentRegion
137  ) const;
138 
139  //- Disallow default bitwise copy construct
140  fvMeshSubset(const fvMeshSubset&);
141 
142  //- Disallow default bitwise assignment
143  void operator=(const fvMeshSubset&);
144 
145 public:
146 
147  // Constructors
148 
149  //- Construct given a mesh to subset
150  explicit fvMeshSubset(const fvMesh&);
151 
152 
153  // Member Functions
154 
155  // Edit
156 
157  //- Set the subset. Create "oldInternalFaces" patch for exposed
158  // internal faces (patchID==-1) or use supplied patch.
159  // Does not handle coupled patches correctly if only one side
160  // gets deleted.
161  void setCellSubset
162  (
163  const labelHashSet& globalCellMap,
164  const label patchID = -1
165  );
166 
167  //- Set the subset from all cells with region == currentRegion.
168  // Create "oldInternalFaces" patch for exposed
169  // internal faces (patchID==-1) or use supplied patch.
170  // Handles coupled patches by if necessary making coupled patch
171  // face part of patchID (so uncoupled)
172  void setLargeCellSubset
173  (
174  const labelList& region,
175  const label currentRegion,
176  const label patchID = -1,
177  const bool syncCouples = true
178  );
179 
180  //- setLargeCellSubset but with labelHashSet.
181  void setLargeCellSubset
182  (
183  const labelHashSet& globalCellMap,
184  const label patchID = -1,
185  const bool syncPar = true
186  );
187 
188 
189  //- Two step subsetting
190 
191  //- Get labels of exposed faces.
192  // These are
193  // - internal faces that become boundary faces
194  // - coupled faces that become uncoupled (since one of the
195  // sides gets deleted)
197  (
198  const labelList& region,
199  const label currentRegion,
200  const bool syncCouples = true
201  ) const;
202 
203  //- For every exposed face (from above getExposedFaces)
204  // used supplied (existing!) patch
205  void setLargeCellSubset
206  (
207  const labelList& region,
208  const label currentRegion,
209  const labelList& exposedFaces,
210  const labelList& patchIDs,
211  const bool syncCouples = true
212  );
213 
214 
215  // Access
216 
217  //- Original mesh
218  const fvMesh& baseMesh() const
219  {
220  return baseMesh_;
221  }
222 
223  //- Have subMesh?
224  bool hasSubMesh() const;
225 
226  //- Return reference to subset mesh
227  const fvMesh& subMesh() const;
228 
229  fvMesh& subMesh();
230 
231  //- Return point map
232  const labelList& pointMap() const;
233 
234  //- Return face map
235  const labelList& faceMap() const;
236 
237  //- Return face map with sign to encode flipped faces
238  const labelList& faceFlipMap() const;
239 
240  //- Return cell map
241  const labelList& cellMap() const;
242 
243  //- Return patch map
244  const labelList& patchMap() const;
245 
246 
247  // Field mapping
248 
249  //- Map volume field
250  template<class Type>
253  (
255  const fvMesh& sMesh,
256  const labelList& patchMap,
257  const labelList& cellMap,
258  const labelList& faceMap
259  );
260 
261  template<class Type>
264  (
266  ) const;
267 
268  //- Map surface field. Optionally negates value if flipping
269  // a face (from exposing an internal face)
270  template<class Type>
273  (
275  const fvMesh& sMesh,
276  const labelList& patchMap,
277  const labelList& cellMap,
278  const labelList& faceMap,
279  const bool negateIfFlipped = true
280  );
281 
282  template<class Type>
285  (
287  const bool negateIfFlipped = true
288  ) const;
289 
290  //- Map point field
291  template<class Type>
294  (
296  const pointMesh& sMesh,
297  const labelList& patchMap,
298  const labelList& pointMap
299  );
300 
301  template<class Type>
304  (
306  ) const;
307 
308  //- Map dimensioned field
309  template<class Type>
312  (
314  const fvMesh& sMesh,
315  const labelList& cellMap
316  );
317 
318  template<class Type>
321 };
322 
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 } // End namespace Foam
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 #ifdef NoRepository
331  #include "fvMeshSubsetInterpolate.C"
332 #endif
333 
334 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335 
336 #endif
337 
338 // ************************************************************************* //
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
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
const labelList & pointMap() const
Return point map.
Generic GeometricField class.
const fvMesh & subMesh() const
Return reference to subset mesh.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const labelList & faceMap() const
Return face map.
const labelList & patchMap() const
Return patch map.
bool hasSubMesh() const
Have subMesh?
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Two step subsetting.
const labelList & cellMap() const
Return cell map.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:217
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
A class for managing temporary objects.
Definition: PtrList.H:54
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.
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
void setCellSubset(const labelHashSet &globalCellMap, const label patchID=-1)
Set the subset. Create "oldInternalFaces" patch for exposed.
Definition: fvMeshSubset.C:410
Namespace for OpenFOAM.