initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / fvMesh / fvMeshSubset / fvMeshSubset.H
blob9ad70db6473e14811f34450c8904704009c00480
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::fvMeshSubset
28 Description
29     Post-processing mesh subset tool.  Given the original mesh and the
30     list of selected cells, it creates the mesh consisting only of the
31     desired cells, with the mapping list for points, faces, and cells.
33     Puts all exposed internal faces into either
34     - a user supplied patch
35     - a newly created patch "oldInternalFaces"
37     - setCellSubset is for small subsets. Uses Maps to minimize memory.
38     - setLargeCellSubset is for largish subsets (>10% of mesh).
39       Uses labelLists instead.
41     - setLargeCellSubset does coupled patch subsetting as well. If it detects
42       a face on a coupled patch 'losing' its neighbour it will move the
43       face into the oldInternalFaces patch.
45     - if a user supplied patch is used the mapping becomes a problem.
46     Do the new faces get the value of the internal face they came from?
47     What if e.g. the user supplied patch is a fixedValue 0? So for now
48     they get the face of existing patch face 0.
50 SourceFiles
51     fvMeshSubset.C
53 \*---------------------------------------------------------------------------*/
55 #ifndef fvMeshSubset_H
56 #define fvMeshSubset_H
58 #include "fvMesh.H"
59 #include "pointMesh.H"
60 #include "fvPatchFieldMapper.H"
61 #include "pointPatchFieldMapper.H"
62 #include "GeometricField.H"
63 #include "HashSet.H"
64 #include "surfaceMesh.H"
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 namespace Foam
71 /*---------------------------------------------------------------------------*\
72                         Class fvMeshSubset Declaration
73 \*---------------------------------------------------------------------------*/
75 class fvMeshSubset
78 public:
80     //- Patch-field subset interpolation class
81     class patchFieldSubset
82     :
83         public fvPatchFieldMapper
84     {
85         const labelList& directAddressing_;
87     public:
89         // Constructors
91             //- Construct given addressing
92             patchFieldSubset(const labelList& directAddressing)
93             :
94                 directAddressing_(directAddressing)
95             {}
97         // Destructor
99             virtual ~patchFieldSubset()
100             {}
103         // Member Functions
105             label size() const
106             {
107                 return directAddressing_.size();
108             }
110             bool direct() const
111             {
112                 return true;
113             }
115             const unallocLabelList& directAddressing() const
116             {
117                 return directAddressing_;
118             }
119     };
122     //- Patch-field subset interpolation class
123     class pointPatchFieldSubset
124     :
125         public pointPatchFieldMapper
126     {
127         const labelList& directAddressing_;
129     public:
131         // Constructors
133             //- Construct given addressing
134             pointPatchFieldSubset(const labelList& directAddressing)
135             :
136                 directAddressing_(directAddressing)
137             {}
139         // Destructor
141             virtual ~pointPatchFieldSubset()
142             {}
145         // Member Functions
147             label size() const
148             {
149                 return directAddressing_.size();
150             }
152             bool direct() const
153             {
154                 return true;
155             }
157             const unallocLabelList& directAddressing() const
158             {
159                 return directAddressing_;
160             }
161     };
164 private:
166     // Private data
168         //- Mesh to subset from
169         const fvMesh& baseMesh_;
171         //- Subset mesh pointer
172         autoPtr<fvMesh> fvMeshSubsetPtr_;
174         //- Point mapping array
175         labelList pointMap_;
177         //- Face mapping array
178         labelList faceMap_;
180         //- Cell mapping array
181         labelList cellMap_;
183         //- Patch mapping array
184         labelList patchMap_;
187     // Private Member Functions
189         //- Check if subset has been performed
190         bool checkCellSubset() const;
192         //- Mark points in Map
193         static void markPoints(const labelList&, Map<label>&); 
195         //- Mark points (with 0) in labelList
196         static void markPoints(const labelList&, labelList&); 
198         //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
199         void doCoupledPatches
200         (
201             const bool syncPar,
202             labelList& nCellsUsingFace
203         ) const;
205         //- Subset of subset
206         static labelList subset
207         (
208             const label nElems,
209             const labelList& selectedElements,
210             const labelList& subsetMap
211         );
213         //- Create zones for submesh
214         void subsetZones();
216         //- Disallow default bitwise copy construct
217         fvMeshSubset(const fvMeshSubset&);
219         //- Disallow default bitwise assignment
220         void operator=(const fvMeshSubset&);
222 public:
224     // Constructors
226         //- Construct given a mesh to subset
227         explicit fvMeshSubset(const fvMesh&);
230     // Member Functions
232         // Edit
234             //- Set the subset. Create "oldInternalFaces" patch for exposed
235             //  internal faces (patchID==-1) or use supplied patch.
236             //  Does not handle coupled patches correctly if only one side
237             //  gets deleted.
238             void setCellSubset
239             (
240                 const labelHashSet& globalCellMap,
241                 const label patchID = -1
242             );
244             //- Set the subset from all cells with region == currentRegion.
245             //  Create "oldInternalFaces" patch for exposed
246             //  internal faces (patchID==-1) or use supplied patch.
247             //  Handles coupled patches by if nessecary making coupled patch
248             //  face part of patchID (so uncoupled)
249             void setLargeCellSubset
250             (
251                 const labelList& region,
252                 const label currentRegion,
253                 const label patchID = -1,
254                 const bool syncCouples = true
255             );
257             //- setLargeCellSubset but with labelHashSet.
258             void setLargeCellSubset
259             (
260                 const labelHashSet& globalCellMap,
261                 const label patchID = -1,
262                 const bool syncPar = true
263             );
266         // Access
268             //- Original mesh
269             const fvMesh& baseMesh() const
270             {
271                 return baseMesh_;
272             }
274             //- Return reference to subset mesh
275             const fvMesh& subMesh() const;
277             fvMesh& subMesh();
279             //- Return point map
280             const labelList& pointMap() const;
282             //- Return face map
283             const labelList& faceMap() const;
285             //- Return cell map
286             const labelList& cellMap() const;
288             //- Return patch map
289             const labelList& patchMap() const;
292         // Field mapping
294             //- Map volume field
295             template<class Type>
296             static tmp<GeometricField<Type, fvPatchField, volMesh> >
297             interpolate
298             (
299                 const GeometricField<Type, fvPatchField, volMesh>&,
300                 const fvMesh& sMesh,
301                 const labelList& patchMap,
302                 const labelList& cellMap,
303                 const labelList& faceMap
304             );
306             template<class Type>
307             tmp<GeometricField<Type, fvPatchField, volMesh> >
308             interpolate
309             (
310                 const GeometricField<Type, fvPatchField, volMesh>&
311             ) const;
313             //- Map surface field
314             template<class Type>
315             static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
316             interpolate
317             (
318                 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
319                 const fvMesh& sMesh,
320                 const labelList& patchMap,
321                 const labelList& faceMap
322             );
324             template<class Type>
325             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
326             interpolate
327             (
328                 const GeometricField<Type, fvsPatchField, surfaceMesh>&
329             ) const;
331             //- Map point field
332             template<class Type>
333             static tmp<GeometricField<Type, pointPatchField, pointMesh> >
334             interpolate
335             (
336                 const GeometricField<Type, pointPatchField, pointMesh>&,
337                 const pointMesh& sMesh,
338                 const labelList& patchMap,
339                 const labelList& pointMap
340             );
342             template<class Type>
343             tmp<GeometricField<Type, pointPatchField, pointMesh> >
344             interpolate
345             (
346                 const GeometricField<Type, pointPatchField, pointMesh>&
347             ) const;
351 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353 } // End namespace Foam
355 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357 #ifdef NoRepository
358 #   include "fvMeshSubsetInterpolate.C"
359 #endif
361 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
363 #endif
365 // ************************************************************************* //