BUGFIX: Seg-fault in multiphaseInterFoam. Author: Henrik Rusche. Merge: Hrvoje Jasak
[foam-extend-3.2.git] / src / finiteVolume / fvMesh / fvMeshSubset / fvMeshSubset.H
blobbf42952f6395510980d7d218c944e8f9833dde27
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::fvMeshSubset
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.
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.
49 SourceFiles
50     fvMeshSubset.C
52 \*---------------------------------------------------------------------------*/
54 #ifndef fvMeshSubset_H
55 #define fvMeshSubset_H
57 #include "fvMesh.H"
58 #include "pointMesh.H"
59 #include "fvPatchFieldMapper.H"
60 #include "PointPatchFieldMapper.H"
61 #include "GeometricField.H"
62 #include "HashSet.H"
63 #include "surfaceMesh.H"
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 namespace Foam
70 /*---------------------------------------------------------------------------*\
71                         Class fvMeshSubset Declaration
72 \*---------------------------------------------------------------------------*/
74 class fvMeshSubset
76     public regIOobject
79 public:
81     //- Patch-field subset interpolation class
82     class patchFieldSubset
83     :
84         public fvPatchFieldMapper
85     {
86         label sizeBeforeMapping_;
88         labelField directAddressing_;
90     public:
92         // Constructors
94             //- Construct from components
95             patchFieldSubset(const label sbm, const labelList& da)
96             :
97                 sizeBeforeMapping_(sbm),
98                 directAddressing_(da)
99             {}
101             //- Construct given addressing
102             patchFieldSubset(const fvMeshSubset& ms, const label patchNo)
103             :
104                 sizeBeforeMapping_
105                 (
106                     ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size()
107                 ),
108                 directAddressing_
109                 (
110                     static_cast<const labelField&>
111                     (
112                         labelField::subField
113                         (
114                             ms.faceMap(),
115                             ms.subMesh().boundary()[patchNo].size(),
116                             ms.subMesh().boundary()[patchNo].patch().start()
117                         )
118                     )
119                   - ms.baseMesh().boundary()
120                     [ms.patchMap()[patchNo]].patch().start()
121                 )
122             {
123                 // If patchNo supplied so exposed internal faces or uncoupled
124                 // patch faces get into existing patch what to do with
125                 // mapping? Truncate to 0 for now.
127                 const label sz =
128                     ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size();
130                 forAll(directAddressing_, i)
131                 {
132                     label& address = directAddressing_[i];
134                     if (address < 0 || address >= sz)
135                     {
136                         address = 0;
137                     }
138                 }
139             }
142         // Destructor
144             virtual ~patchFieldSubset()
145             {}
148         // Member Functions
150             label size() const
151             {
152                 return directAddressing_.size();
153             }
155             label sizeBeforeMapping() const
156             {
157                 return directAddressing_.size();
158             }
160             bool direct() const
161             {
162                 return true;
163             }
165             const unallocLabelList& directAddressing() const
166             {
167                 return directAddressing_;
168             }
169     };
172     //- Patch-field subset interpolation class
173     class pointPatchFieldSubset
174     :
175         public PointPatchFieldMapper
176     {
177         //- Addressing
178         const labelList& directAddressing_;
181     public:
183         // Constructors
185             //- Construct given addressing
186             pointPatchFieldSubset(const labelList& directAddressing)
187             :
188                 directAddressing_(directAddressing)
189             {}
191         // Destructor
193             virtual ~pointPatchFieldSubset()
194             {}
197         // Member Functions
199             label size() const
200             {
201                 return directAddressing_.size();
202             }
204             label sizeBeforeMapping() const
205             {
206                 return directAddressing_.size();
207             }
209             bool direct() const
210             {
211                 return true;
212             }
214             const unallocLabelList& directAddressing() const
215             {
216                 return directAddressing_;
217             }
218     };
221 private:
223     // Private data
225         //- Reference to original mesh
226         const fvMesh& baseMesh_;
228         //- Subset fvMesh pointer
229         fvMesh* fvMeshSubsetPtr_;
231         //- Subset pointMesh pointer
232         mutable pointMesh* pointMeshSubsetPtr_;
234         //- Point mapping array
235         labelList pointMap_;
237         //- Face mapping array
238         labelList faceMap_;
240         //- Cell mapping array
241         labelList cellMap_;
243         //- Patch mapping array
244         labelList patchMap_;
247     // Private Member Functions
249         //- Disallow default bitwise copy construct
250         fvMeshSubset(const fvMeshSubset&);
252         //- Disallow default bitwise assignment
253         void operator=(const fvMeshSubset&);
256         //- Check if subset has been performed
257         bool checkCellSubset() const;
259         //- Mark points in Map
260         static void markPoints(const labelList&, Map<label>&);
262         //- Mark points (with 0) in labelList
263         static void markPoints(const labelList&, labelList&);
265         //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
266         void doCoupledPatches
267         (
268             const bool syncPar,
269             labelList& nCellsUsingFace
270         ) const;
272         //- Subset of subset
273         static labelList subset
274         (
275             const label nElems,
276             const labelList& selectedElements,
277             const labelList& subsetMap
278         );
280         //- Create zones for submesh
281         void subsetZones();
283         //- Make dictionaries in case they are not present
284         //  Creates fvSchemes and fvSolution copying from base
285         void makeFvDictionaries();
288 public:
290     // Constructors
292         //- Construct given a mesh
293         explicit fvMeshSubset(const IOobject&, const fvMesh&);
296     // Destructor
298         ~fvMeshSubset();
301     // Member Functions
303         // Edit
305             //- Set the subset. Create "oldInternalFaces" patch for exposed
306             //  internal faces (patchID==-1) or use supplied patch.
307             //  Does not handle coupled patches correctly if only one side
308             //  gets deleted.
309             void setCellSubset
310             (
311                 const labelHashSet& globalCellMap,
312                 const label patchID = -1
313             );
315             //- Set the subset from all cells with region == currentRegion.
316             //  Create "oldInternalFaces" patch for exposed
317             //  internal faces (patchID==-1) or use supplied patch.
318             //  Handles coupled patches by if nessecary making coupled patch
319             //  face part of patchID (so uncoupled)
320             void setLargeCellSubset
321             (
322                 const labelList& region,
323                 const label currentRegion,
324                 const label patchID = -1,
325                 const bool syncCouples = true
326             );
328             //- setLargeCellSubset but with labelHashSet.
329             void setLargeCellSubset
330             (
331                 const labelHashSet& globalCellMap,
332                 const label patchID = -1,
333                 const bool syncPar = true
334             );
337         // Access
339             //- Return reference to base mesh
340             const fvMesh& baseMesh() const
341             {
342                 return baseMesh_;
343             }
345             //- Return reference to subset mesh
346             const fvMesh& subMesh() const;
348             fvMesh& subMesh();
350             //- Return reference to demand-driven subset pointMesh
351             const pointMesh& subPointMesh() const;
353             pointMesh& subPointMesh();
355             //- Return point map
356             const labelList& pointMap() const;
358             //- Return face map
359             const labelList& faceMap() const;
361             //- Return cell map
362             const labelList& cellMap() const;
364             //- Return patch map
365             const labelList& patchMap() const;
368         // Field mapping
370             //- Map volume field given map addressing
371             //  - it does not patch up the exposed patches with surface fields
372             //  - it does not handle patch type creation correctly
373             //  Please use only in dummy mapping operations
374             //  HJ, 31/Jul/2009
375             template<class Type>
376             static tmp<GeometricField<Type, fvPatchField, volMesh> >
377             meshToMesh
378             (
379                 const GeometricField<Type, fvPatchField, volMesh>&,
380                 const fvMesh& sMesh,
381                 const labelList& patchMap,
382                 const labelList& cellMap,
383                 const labelList& faceMap
384             );
386             //- Map volume field
387             template<class Type>
388             tmp<GeometricField<Type, fvPatchField, volMesh> >
389             interpolate
390             (
391                 const GeometricField<Type, fvPatchField, volMesh>&
392             ) const;
394             //- Map surface field
395             template<class Type>
396             static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
397             meshToMesh
398             (
399                 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
400                 const fvMesh& sMesh,
401                 const labelList& patchMap,
402                 const labelList& faceMap
403             );
405             template<class Type>
406             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
407             interpolate
408             (
409                 const GeometricField<Type, fvsPatchField, surfaceMesh>&
410             ) const;
412             //- Map point field
413             template<class Type>
414             static tmp<GeometricField<Type, pointPatchField, pointMesh> >
415             meshToMesh
416             (
417                 const GeometricField<Type, pointPatchField, pointMesh>&,
418                 const pointMesh& sMesh,
419                 const labelList& patchMap,
420                 const labelList& pointMap
421             );
423             template<class Type>
424             tmp<GeometricField<Type, pointPatchField, pointMesh> >
425             interpolate
426             (
427                 const GeometricField<Type, pointPatchField, pointMesh>&
428             ) const;
431         //- Write data
432         virtual bool writeData(Foam::Ostream&) const
433         {
434             return true;
435         }
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
442 } // End namespace Foam
444 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 #ifdef NoRepository
447 #   include "fvMeshSubsetInterpolate.C"
448 #endif
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
452 #endif
454 // ************************************************************************* //