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-2013 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 the mapping becomes a problem.
45  Do the new faces get the value of the internal face they came from?
46  What if e.g. the user supplied patch is a fixedValue 0? So for now
47  they get the face of existing patch face 0.
48 
49 SourceFiles
50  fvMeshSubset.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef fvMeshSubset_H
55 #define fvMeshSubset_H
56 
57 #include "fvMesh.H"
58 #include "pointMesh.H"
59 #include "GeometricField.H"
60 #include "HashSet.H"
61 #include "surfaceMesh.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class fvMeshSubset Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 class fvMeshSubset
73 {
74 
75 private:
76 
77  // Private data
78 
79  //- Mesh to subset from
80  const fvMesh& baseMesh_;
81 
82  //- Subset mesh pointer
83  autoPtr<fvMesh> fvMeshSubsetPtr_;
84 
85  //- Point mapping array
86  labelList pointMap_;
87 
88  //- Face mapping array
89  labelList faceMap_;
90 
91  //- Cell mapping array
92  labelList cellMap_;
93 
94  //- Patch mapping array
95  labelList patchMap_;
96 
97 
98  // Private Member Functions
99 
100  //- Check if subset has been performed
101  bool checkCellSubset() const;
102 
103  //- Mark points in Map
104  static void markPoints(const labelList&, Map<label>&);
105 
106  //- Mark points (with 0) in labelList
107  static void markPoints(const labelList&, labelList&);
108 
109  //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
110  void doCoupledPatches
111  (
112  const bool syncPar,
113  labelList& nCellsUsingFace
114  ) const;
115 
116  //- Subset of subset
117  static labelList subset
118  (
119  const label nElems,
120  const labelList& selectedElements,
121  const labelList& subsetMap
122  );
123 
124  //- Create zones for submesh
125  void subsetZones();
126 
127  //- Disallow default bitwise copy construct
128  fvMeshSubset(const fvMeshSubset&);
129 
130  //- Disallow default bitwise assignment
131  void operator=(const fvMeshSubset&);
132 
133 public:
134 
135  // Constructors
136 
137  //- Construct given a mesh to subset
138  explicit fvMeshSubset(const fvMesh&);
139 
140 
141  // Member Functions
142 
143  // Edit
144 
145  //- Set the subset. Create "oldInternalFaces" patch for exposed
146  // internal faces (patchID==-1) or use supplied patch.
147  // Does not handle coupled patches correctly if only one side
148  // gets deleted.
149  void setCellSubset
150  (
151  const labelHashSet& globalCellMap,
152  const label patchID = -1
153  );
154 
155  //- Set the subset from all cells with region == currentRegion.
156  // Create "oldInternalFaces" patch for exposed
157  // internal faces (patchID==-1) or use supplied patch.
158  // Handles coupled patches by if necessary making coupled patch
159  // face part of patchID (so uncoupled)
160  void setLargeCellSubset
161  (
162  const labelList& region,
163  const label currentRegion,
164  const label patchID = -1,
165  const bool syncCouples = true
166  );
167 
168  //- setLargeCellSubset but with labelHashSet.
169  void setLargeCellSubset
170  (
171  const labelHashSet& globalCellMap,
172  const label patchID = -1,
173  const bool syncPar = true
174  );
175 
176 
177  // Access
178 
179  //- Original mesh
180  const fvMesh& baseMesh() const
181  {
182  return baseMesh_;
183  }
184 
185  //- Have subMesh?
186  bool hasSubMesh() const;
187 
188  //- Return reference to subset mesh
189  const fvMesh& subMesh() const;
190 
191  fvMesh& subMesh();
192 
193  //- Return point map
194  const labelList& pointMap() const;
195 
196  //- Return face map
197  const labelList& faceMap() const;
198 
199  //- Return cell map
200  const labelList& cellMap() const;
201 
202  //- Return patch map
203  const labelList& patchMap() const;
204 
205 
206  // Field mapping
207 
208  //- Map volume field
209  template<class Type>
212  (
214  const fvMesh& sMesh,
215  const labelList& patchMap,
216  const labelList& cellMap,
217  const labelList& faceMap
218  );
219 
220  template<class Type>
223  (
225  ) const;
226 
227  //- Map surface field
228  template<class Type>
231  (
233  const fvMesh& sMesh,
234  const labelList& patchMap,
235  const labelList& faceMap
236  );
237 
238  template<class Type>
241  (
243  ) const;
244 
245  //- Map point field
246  template<class Type>
249  (
251  const pointMesh& sMesh,
252  const labelList& patchMap,
253  const labelList& pointMap
254  );
255 
256  template<class Type>
259  (
261  ) const;
262 };
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 } // End namespace Foam
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 #ifdef NoRepository
272 # include "fvMeshSubsetInterpolate.C"
273 #endif
274 
275 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276 
277 #endif
278 
279 // ************************************************************************* //
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
void setCellSubset(const labelHashSet &globalCellMap, const label patchID=-1)
Set the subset. Create "oldInternalFaces" patch for exposed.
Definition: fvMeshSubset.C:374
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:758
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 fvMesh & subMesh() const
Return reference to subset mesh.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const labelList & cellMap() const
Return cell map.
Namespace for OpenFOAM.
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:179
const labelList & pointMap() const
Return point map.
const labelList & patchMap() const
Return patch map.
Generic GeometricField class.
const labelList & faceMap() const
Return face map.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
A class for managing temporary objects.
Definition: PtrList.H:118
bool hasSubMesh() const
Have subMesh?
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:71