initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fvMesh / fvMeshSubset / fvMeshSubset.H
blob11c1af1bc9102021f8ecf5d26e090afff3f57f2f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "labelHashSet.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         mutable autoPtr<pointMesh> pointMeshSubsetPtr_;
176         //- Point mapping array
177         labelList pointMap_;
179         //- Face mapping array
180         labelList faceMap_;
182         //- Cell mapping array
183         labelList cellMap_;
185         //- Patch mapping array
186         labelList patchMap_;
189     // Private Member Functions
191         //- Check if subset has been performed
192         bool checkCellSubset() const;
194         //- Mark points in Map
195         static void markPoints(const labelList&, Map<label>&); 
197         //- Mark points (with 0) in labelList
198         static void markPoints(const labelList&, labelList&); 
200         //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
201         void doCoupledPatches
202         (
203             const bool syncPar,
204             labelList& nCellsUsingFace
205         ) const;
207         //- Subset of subset
208         static labelList subset
209         (
210             const label nElems,
211             const labelList& selectedElements,
212             const labelList& subsetMap
213         );
215         //- Create zones for submesh
216         void subsetZones();
218         //- Disallow default bitwise copy construct
219         fvMeshSubset(const fvMeshSubset&);
221         //- Disallow default bitwise assignment
222         void operator=(const fvMeshSubset&);
224 public:
226     // Constructors
228         //- Construct given a mesh to subset
229         explicit fvMeshSubset(const fvMesh&);
232     // Member Functions
234         // Edit
236             //- Set the subset. Create "oldInternalFaces" patch for exposed
237             //  internal faces (patchID==-1) or use supplied patch.
238             //  Does not handle coupled patches correctly if only one side
239             //  gets deleted.
240             void setCellSubset
241             (
242                 const labelHashSet& globalCellMap,
243                 const label patchID = -1
244             );
246             //- Set the subset from all cells with region == currentRegion.
247             //  Create "oldInternalFaces" patch for exposed
248             //  internal faces (patchID==-1) or use supplied patch.
249             //  Handles coupled patches by if nessecary making coupled patch
250             //  face part of patchID (so uncoupled)
251             void setLargeCellSubset
252             (
253                 const labelList& region,
254                 const label currentRegion,
255                 const label patchID = -1,
256                 const bool syncCouples = true
257             );
259             //- setLargeCellSubset but with labelHashSet.
260             void setLargeCellSubset
261             (
262                 const labelHashSet& globalCellMap,
263                 const label patchID = -1,
264                 const bool syncPar = true
265             );
268         // Access
270             //- Original mesh
271             const fvMesh& baseMesh() const
272             {
273                 return baseMesh_;
274             }
276             //- Return reference to subset mesh
277             const fvMesh& subMesh() const;
279             fvMesh& subMesh();
281             //- Return reference to demand-driven subset pointMesh
282             const pointMesh& subPointMesh() const;
284             pointMesh& subPointMesh();
286             //- Return point map
287             const labelList& pointMap() const;
289             //- Return face map
290             const labelList& faceMap() const;
292             //- Return cell map
293             const labelList& cellMap() const;
295             //- Return patch map
296             const labelList& patchMap() const;
299         // Field mapping
301             //- Map volume field
302             template<class Type>
303             static tmp<GeometricField<Type, fvPatchField, volMesh> >
304             interpolate
305             (
306                 const GeometricField<Type, fvPatchField, volMesh>&,
307                 const fvMesh& sMesh,
308                 const labelList& patchMap,
309                 const labelList& cellMap,
310                 const labelList& faceMap
311             );
313             template<class Type>
314             tmp<GeometricField<Type, fvPatchField, volMesh> >
315             interpolate
316             (
317                 const GeometricField<Type, fvPatchField, volMesh>&
318             ) const;
320             //- Map surface field
321             template<class Type>
322             static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
323             interpolate
324             (
325                 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
326                 const fvMesh& sMesh,
327                 const labelList& patchMap,
328                 const labelList& faceMap
329             );
331             template<class Type>
332             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
333             interpolate
334             (
335                 const GeometricField<Type, fvsPatchField, surfaceMesh>&
336             ) const;
338             //- Map point field
339             template<class Type>
340             static tmp<GeometricField<Type, pointPatchField, pointMesh> >
341             interpolate
342             (
343                 const GeometricField<Type, pointPatchField, pointMesh>&,
344                 const pointMesh& sMesh,
345                 const objectRegistry& reg,
346                 const labelList& patchMap,
347                 const labelList& pointMap
348             );
350             template<class Type>
351             tmp<GeometricField<Type, pointPatchField, pointMesh> >
352             interpolate
353             (
354                 const GeometricField<Type, pointPatchField, pointMesh>&
355             ) const;
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361 } // End namespace Foam
363 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
365 #ifdef NoRepository
366 #   include "fvMeshSubsetInterpolate.C"
367 #endif
369 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 #endif
373 // ************************************************************************* //