1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
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.
32 Puts all exposed internal faces into either
33 - a user supplied patch
34 - a newly created patch "oldInternalFaces"
36 - setCellSubset is for small subsets. Uses Maps to minimize memory.
37 - setLargeCellSubset is for largish subsets (>10% of mesh).
38 Uses labelLists instead.
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.
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.
52 \*---------------------------------------------------------------------------*/
54 #ifndef fvMeshSubset_H
55 #define fvMeshSubset_H
57 #include <finiteVolume/fvMesh.H>
58 #include <OpenFOAM/pointMesh.H>
59 #include <finiteVolume/fvPatchFieldMapper.H>
60 #include <OpenFOAM/pointPatchFieldMapper.H>
61 #include <OpenFOAM/GeometricField.H>
62 #include <OpenFOAM/HashSet.H>
63 #include <finiteVolume/surfaceMesh.H>
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 /*---------------------------------------------------------------------------*\
71 Class fvMeshSubset Declaration
72 \*---------------------------------------------------------------------------*/
79 //- Patch-field subset interpolation class
80 class patchFieldSubset
82 public fvPatchFieldMapper
84 const labelList& directAddressing_;
90 //- Construct given addressing
91 patchFieldSubset(const labelList& directAddressing)
93 directAddressing_(directAddressing)
98 virtual ~patchFieldSubset()
106 return directAddressing_.size();
114 const unallocLabelList& directAddressing() const
116 return directAddressing_;
121 //- Patch-field subset interpolation class
122 class pointPatchFieldSubset
124 public pointPatchFieldMapper
126 const labelList& directAddressing_;
132 //- Construct given addressing
133 pointPatchFieldSubset(const labelList& directAddressing)
135 directAddressing_(directAddressing)
140 virtual ~pointPatchFieldSubset()
148 return directAddressing_.size();
156 const unallocLabelList& directAddressing() const
158 return directAddressing_;
167 //- Mesh to subset from
168 const fvMesh& baseMesh_;
170 //- Subset mesh pointer
171 autoPtr<fvMesh> fvMeshSubsetPtr_;
173 //- Point mapping array
176 //- Face mapping array
179 //- Cell mapping array
182 //- Patch mapping array
186 // Private Member Functions
188 //- Check if subset has been performed
189 bool checkCellSubset() const;
191 //- Mark points in Map
192 static void markPoints(const labelList&, Map<label>&);
194 //- Mark points (with 0) in labelList
195 static void markPoints(const labelList&, labelList&);
197 //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
198 void doCoupledPatches
201 labelList& nCellsUsingFace
205 static labelList subset
208 const labelList& selectedElements,
209 const labelList& subsetMap
212 //- Create zones for submesh
215 //- Disallow default bitwise copy construct
216 fvMeshSubset(const fvMeshSubset&);
218 //- Disallow default bitwise assignment
219 void operator=(const fvMeshSubset&);
225 //- Construct given a mesh to subset
226 explicit fvMeshSubset(const fvMesh&);
233 //- Set the subset. Create "oldInternalFaces" patch for exposed
234 // internal faces (patchID==-1) or use supplied patch.
235 // Does not handle coupled patches correctly if only one side
239 const labelHashSet& globalCellMap,
240 const label patchID = -1
243 //- Set the subset from all cells with region == currentRegion.
244 // Create "oldInternalFaces" patch for exposed
245 // internal faces (patchID==-1) or use supplied patch.
246 // Handles coupled patches by if nessecary making coupled patch
247 // face part of patchID (so uncoupled)
248 void setLargeCellSubset
250 const labelList& region,
251 const label currentRegion,
252 const label patchID = -1,
253 const bool syncCouples = true
256 //- setLargeCellSubset but with labelHashSet.
257 void setLargeCellSubset
259 const labelHashSet& globalCellMap,
260 const label patchID = -1,
261 const bool syncPar = true
268 const fvMesh& baseMesh() const
273 //- Return reference to subset mesh
274 const fvMesh& subMesh() const;
279 const labelList& pointMap() const;
282 const labelList& faceMap() const;
285 const labelList& cellMap() const;
288 const labelList& patchMap() const;
295 static tmp<GeometricField<Type, fvPatchField, volMesh> >
298 const GeometricField<Type, fvPatchField, volMesh>&,
300 const labelList& patchMap,
301 const labelList& cellMap,
302 const labelList& faceMap
306 tmp<GeometricField<Type, fvPatchField, volMesh> >
309 const GeometricField<Type, fvPatchField, volMesh>&
312 //- Map surface field
314 static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
317 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
319 const labelList& patchMap,
320 const labelList& faceMap
324 tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
327 const GeometricField<Type, fvsPatchField, surfaceMesh>&
332 static tmp<GeometricField<Type, pointPatchField, pointMesh> >
335 const GeometricField<Type, pointPatchField, pointMesh>&,
336 const pointMesh& sMesh,
337 const labelList& patchMap,
338 const labelList& pointMap
342 tmp<GeometricField<Type, pointPatchField, pointMesh> >
345 const GeometricField<Type, pointPatchField, pointMesh>&
350 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
352 } // End namespace Foam
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 # include "fvMeshSubsetInterpolate.C"
360 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
364 // ************************ vim: set sw=4 sts=4 et: ************************ //