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-2021 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 minimise 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 "HashSet.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 
67 class pointMesh;
68 class surfaceMesh;
69 
70 /*---------------------------------------------------------------------------*\
71  Class fvMeshSubset Declaration
72 \*---------------------------------------------------------------------------*/
73 
74 class fvMeshSubset
75 {
76  // Private Data
77 
78  //- Mesh to subset from
79  const fvMesh& baseMesh_;
80 
81  //- Subset mesh pointer
82  autoPtr<fvMesh> fvMeshSubsetPtr_;
83 
84  //- Point mapping array
85  labelList pointMap_;
86 
87  //- Face mapping array
88  labelList faceMap_;
89 
90  //- Cell mapping array
91  labelList cellMap_;
92 
93  //- Patch mapping array
94  labelList patchMap_;
95 
96  //- Optional face mapping array with flip encoded
97  mutable autoPtr<labelList> faceFlipMapPtr_;
98 
99 
100  // Private Member Functions
101 
102  //- Check if subset has been performed
103  bool checkCellSubset() const;
104 
105  //- Mark points in Map
106  static void markPoints(const labelList&, Map<label>&);
107 
108  //- Mark points (with 0) in labelList
109  static void markPoints(const labelList&, labelList&);
110 
111  //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
112  void doCoupledPatches
113  (
114  const bool syncPar,
115  Map<label>& facesToSubset,
116  labelList& nCellsUsingFace
117  ) const;
118 
119  //- Subset of subset
120  static labelList subset
121  (
122  const label nElems,
123  const labelList& selectedElements,
124  const labelList& subsetMap
125  );
126 
127  //- Create zones for submesh
128  void subsetZones();
129 
130  //- Helper: extract cells-to-remove from cells-to-keep
131  labelList getCellsToRemove
132  (
133  const labelList& region,
134  const label currentRegion
135  ) const;
136 
137 
138 public:
139 
140  // Constructors
141 
142  //- Construct given a mesh to subset
143  explicit fvMeshSubset(const fvMesh&);
144 
145  //- Disallow default bitwise copy construction
146  fvMeshSubset(const fvMeshSubset&) = delete;
147 
148 
149  // Member Functions
150 
151  // Edit
152 
153  //- Set the subset. Create "oldInternalFaces" patch for exposed
154  // internal faces (patchID==-1) or use supplied patch.
155  // Does not handle coupled patches correctly if only one side
156  // gets deleted.
157  void setCellSubset
158  (
159  const labelHashSet& globalCellMap,
160  const label patchID = -1,
161  const bool syncPar = true
162  );
163 
164  //- Set the subset from all cells with region == currentRegion.
165  // Create "oldInternalFaces" patch for exposed
166  // internal faces (patchID==-1) or use supplied patch.
167  // Handles coupled patches by if necessary making coupled patch
168  // face part of patchID (so uncoupled)
169  void setLargeCellSubset
170  (
171  const labelList& region,
172  const label currentRegion,
173  const label patchID = -1,
174  const bool syncCouples = true
175  );
176 
177  //- setLargeCellSubset but with labelHashSet.
178  void setLargeCellSubset
179  (
180  const labelHashSet& globalCellMap,
181  const label patchID = -1,
182  const bool syncPar = true
183  );
184 
185 
186  // Two step subsetting
187 
188  //- Get labels of exposed faces.
189  // These are
190  // - internal faces that become boundary faces
191  // - coupled faces that become uncoupled (since one of the
192  // sides gets deleted)
194  (
195  const labelList& region,
196  const label currentRegion,
197  const bool syncCouples = true
198  ) const;
199 
200  //- For every exposed face (from above getExposedFaces)
201  // used supplied (existing!) patch
202  void setLargeCellSubset
203  (
204  const labelList& region,
205  const label currentRegion,
206  const labelList& exposedFaces,
207  const labelList& patchIDs,
208  const bool syncCouples = true
209  );
210 
211 
212  // Access
213 
214  //- Original mesh
215  const fvMesh& baseMesh() const
216  {
217  return baseMesh_;
218  }
219 
220  //- Have subMesh?
221  bool hasSubMesh() const;
222 
223  //- Return reference to subset mesh
224  const fvMesh& subMesh() const;
225 
226  fvMesh& subMesh();
227 
228  //- Return point map
229  const labelList& pointMap() const;
230 
231  //- Return face map
232  const labelList& faceMap() const;
233 
234  //- Return face map with sign to encode flipped faces
235  const labelList& faceFlipMap() const;
236 
237  //- Return cell map
238  const labelList& cellMap() const;
239 
240  //- Return patch map
241  const labelList& patchMap() const;
242 
243 
244  // Field mapping
245 
246  //- Map volume field
247  template<class Type>
250  (
252  const fvMesh& sMesh,
253  const labelList& patchMap,
254  const labelList& cellMap,
255  const labelList& faceMap
256  );
257 
258  template<class Type>
261  (
263  ) const;
264 
265  //- Map surface field. Optionally negates value if flipping
266  // a face (from exposing an internal face)
267  template<class Type>
270  (
272  const fvMesh& sMesh,
273  const labelList& patchMap,
274  const labelList& cellMap,
275  const labelList& faceMap
276  );
277 
278  template<class Type>
281  (
283  ) const;
284 
285  //- Map point field
286  template<class Type>
289  (
291  const pointMesh& sMesh,
292  const labelList& patchMap,
293  const labelList& pointMap
294  );
295 
296  template<class Type>
299  (
301  ) const;
302 
303  //- Map dimensioned field
304  template<class Type>
307  (
309  const fvMesh& sMesh,
310  const labelList& cellMap
311  );
312 
313  template<class Type>
316 
317 
318  // Member Operators
319 
320  //- Disallow default bitwise assignment
321  void operator=(const fvMeshSubset&) = delete;
322 };
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 } // End namespace Foam
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
331 #ifdef NoRepository
332  #include "fvMeshSubsetInterpolate.C"
333 #endif
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 #endif
338 
339 // ************************************************************************* //
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.
fvMeshSubset(const fvMesh &)
Construct given a mesh to subset.
Definition: fvMeshSubset.C:464
const labelList & faceMap() const
Return face map.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const labelList & patchMap() const
Return patch map.
const fvMesh & subMesh() const
Return reference to subset mesh.
labelList getExposedFaces(const labelList &region, const label currentRegion, const bool syncCouples=true) const
Get labels of exposed faces.
const labelList & pointMap() const
Return point map.
const labelList & cellMap() const
Return cell 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:479
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
void operator=(const fvMeshSubset &)=delete
Disallow default bitwise assignment.
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.
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
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:881
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:214
Namespace for OpenFOAM.
bool hasSubMesh() const
Have subMesh?