fvMeshSubset.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  Map<label>& facesToSubset,
119  labelList& nCellsUsingFace
120  ) const;
121 
122  //- Subset of subset
123  static labelList subset
124  (
125  const label nElems,
126  const labelList& selectedElements,
127  const labelList& subsetMap
128  );
129 
130  //- Create zones for submesh
131  void subsetZones();
132 
133  //- Helper: extract cells-to-remove from cells-to-keep
134  labelList getCellsToRemove
135  (
136  const labelList& region,
137  const label currentRegion
138  ) const;
139 
140  //- Disallow default bitwise copy construct
141  fvMeshSubset(const fvMeshSubset&);
142 
143  //- Disallow default bitwise assignment
144  void operator=(const fvMeshSubset&);
145 
146 public:
147 
148  // Constructors
149 
150  //- Construct given a mesh to subset
151  explicit fvMeshSubset(const fvMesh&);
152 
153 
154  // Member Functions
155 
156  // Edit
157 
158  //- Set the subset. Create "oldInternalFaces" patch for exposed
159  // internal faces (patchID==-1) or use supplied patch.
160  // Does not handle coupled patches correctly if only one side
161  // gets deleted.
162  void setCellSubset
163  (
164  const labelHashSet& globalCellMap,
165  const label patchID = -1,
166  const bool syncPar = true
167  );
168 
169  //- Set the subset from all cells with region == currentRegion.
170  // Create "oldInternalFaces" patch for exposed
171  // internal faces (patchID==-1) or use supplied patch.
172  // Handles coupled patches by if necessary making coupled patch
173  // face part of patchID (so uncoupled)
174  void setLargeCellSubset
175  (
176  const labelList& region,
177  const label currentRegion,
178  const label patchID = -1,
179  const bool syncCouples = true
180  );
181 
182  //- setLargeCellSubset but with labelHashSet.
183  void setLargeCellSubset
184  (
185  const labelHashSet& globalCellMap,
186  const label patchID = -1,
187  const bool syncPar = true
188  );
189 
190 
191  //- Two step subsetting
192 
193  //- Get labels of exposed faces.
194  // These are
195  // - internal faces that become boundary faces
196  // - coupled faces that become uncoupled (since one of the
197  // sides gets deleted)
199  (
200  const labelList& region,
201  const label currentRegion,
202  const bool syncCouples = true
203  ) const;
204 
205  //- For every exposed face (from above getExposedFaces)
206  // used supplied (existing!) patch
207  void setLargeCellSubset
208  (
209  const labelList& region,
210  const label currentRegion,
211  const labelList& exposedFaces,
212  const labelList& patchIDs,
213  const bool syncCouples = true
214  );
215 
216 
217  // Access
218 
219  //- Original mesh
220  const fvMesh& baseMesh() const
221  {
222  return baseMesh_;
223  }
224 
225  //- Have subMesh?
226  bool hasSubMesh() const;
227 
228  //- Return reference to subset mesh
229  const fvMesh& subMesh() const;
230 
231  fvMesh& subMesh();
232 
233  //- Return point map
234  const labelList& pointMap() const;
235 
236  //- Return face map
237  const labelList& faceMap() const;
238 
239  //- Return face map with sign to encode flipped faces
240  const labelList& faceFlipMap() const;
241 
242  //- Return cell map
243  const labelList& cellMap() const;
244 
245  //- Return patch map
246  const labelList& patchMap() const;
247 
248 
249  // Field mapping
250 
251  //- Map volume field
252  template<class Type>
255  (
257  const fvMesh& sMesh,
258  const labelList& patchMap,
259  const labelList& cellMap,
260  const labelList& faceMap
261  );
262 
263  template<class Type>
266  (
268  ) const;
269 
270  //- Map surface field. Optionally negates value if flipping
271  // a face (from exposing an internal face)
272  template<class Type>
275  (
277  const fvMesh& sMesh,
278  const labelList& patchMap,
279  const labelList& cellMap,
280  const labelList& faceMap,
281  const bool negateIfFlipped = true
282  );
283 
284  template<class Type>
287  (
289  const bool negateIfFlipped = true
290  ) const;
291 
292  //- Map point field
293  template<class Type>
296  (
298  const pointMesh& sMesh,
299  const labelList& patchMap,
300  const labelList& pointMap
301  );
302 
303  template<class Type>
306  (
308  ) const;
309 
310  //- Map dimensioned field
311  template<class Type>
314  (
316  const fvMesh& sMesh,
317  const labelList& cellMap
318  );
319 
320  template<class Type>
323 };
324 
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace Foam
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
332 #ifdef NoRepository
333  #include "fvMeshSubsetInterpolate.C"
334 #endif
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
338 #endif
339 
340 // ************************************************************************* //
const labelList & patchMap() const
Return patch map.
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 & faceFlipMap() const
Return face map with sign to encode flipped faces.
Generic GeometricField class.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Two step subsetting.
const labelList & faceMap() const
Return face map.
void setCellSubset(const labelHashSet &globalCellMap, const label patchID=-1, const bool syncPar=true)
Set the subset. Create "oldInternalFaces" patch for exposed.
Definition: fvMeshSubset.C:483
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
const fvMesh & subMesh() const
Return reference to subset mesh.
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...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
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:885
const labelList & cellMap() const
Return cell map.
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:219
const labelList & pointMap() const
Return point map.
Namespace for OpenFOAM.
bool hasSubMesh() const
Have subMesh?