Fix tutorials: typo in tutorials/viscoelastic/viscoelasticFluidFoam/S-MDCPP/constant...
[OpenFOAM-1.6-ext.git] / src / finiteVolume / fvMesh / fvMeshSubset / fvMeshSubset.H
bloba7ba95c8742a31aab98a74eb0e2ad683e5963b89
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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
77     public regIOobject
80 public:
82     //- Patch-field subset interpolation class
83     class patchFieldSubset
84     :
85         public fvPatchFieldMapper
86     {
87         label sizeBeforeMapping_;
89         labelField directAddressing_;
91     public:
93         // Constructors
95             //- Construct from components
96             patchFieldSubset(const label sbm, const labelList& da)
97             :
98                 sizeBeforeMapping_(sbm),
99                 directAddressing_(da)
100             {}
102             //- Construct given addressing
103             patchFieldSubset(const fvMeshSubset& ms, const label patchNo)
104             :
105                 sizeBeforeMapping_
106                 (
107                     ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size()
108                 ),
109                 directAddressing_
110                 (
111                     static_cast<const labelField&>
112                     (
113                         labelField::subField
114                         (
115                             ms.faceMap(),
116                             ms.subMesh().boundary()[patchNo].size(),
117                             ms.subMesh().boundary()[patchNo].patch().start()
118                         )
119                     )
120                   - ms.baseMesh().boundary()
121                     [ms.patchMap()[patchNo]].patch().start()
122                 )
123             {
124                 // If patchNo supplied so exposed internal faces or uncoupled
125                 // patch faces get into existing patch what to do with
126                 // mapping? Truncate to 0 for now.
128                 const label sz =
129                     ms.baseMesh().boundary()[ms.patchMap()[patchNo]].size();
131                 forAll(directAddressing_, i)
132                 {
133                     label& address = directAddressing_[i];
135                     if (address < 0 || address >= sz)
136                     {
137                         address = 0;
138                     }
139                 }
140             }
143         // Destructor
145             virtual ~patchFieldSubset()
146             {}
149         // Member Functions
151             label size() const
152             {
153                 return directAddressing_.size();
154             }
156             label sizeBeforeMapping() const
157             {
158                 return directAddressing_.size();
159             }
161             bool direct() const
162             {
163                 return true;
164             }
166             const unallocLabelList& directAddressing() const
167             {
168                 return directAddressing_;
169             }
170     };
173     //- Patch-field subset interpolation class
174     class pointPatchFieldSubset
175     :
176         public PointPatchFieldMapper
177     {
178         //- Addressing
179         const labelList& directAddressing_;
182     public:
184         // Constructors
186             //- Construct given addressing
187             pointPatchFieldSubset(const labelList& directAddressing)
188             :
189                 directAddressing_(directAddressing)
190             {}
192         // Destructor
194             virtual ~pointPatchFieldSubset()
195             {}
198         // Member Functions
200             label size() const
201             {
202                 return directAddressing_.size();
203             }
205             label sizeBeforeMapping() const
206             {
207                 return directAddressing_.size();
208             }
210             bool direct() const
211             {
212                 return true;
213             }
215             const unallocLabelList& directAddressing() const
216             {
217                 return directAddressing_;
218             }
219     };
222 private:
224     // Private data
226         //- Reference to original mesh
227         const fvMesh& baseMesh_;
229         //- Subset fvMesh pointer
230         fvMesh* fvMeshSubsetPtr_;
232         //- Subset pointMesh pointer
233         mutable pointMesh* pointMeshSubsetPtr_;
235         //- Point mapping array
236         labelList pointMap_;
238         //- Face mapping array
239         labelList faceMap_;
241         //- Cell mapping array
242         labelList cellMap_;
244         //- Patch mapping array
245         labelList patchMap_;
248     // Private Member Functions
250         //- Disallow default bitwise copy construct
251         fvMeshSubset(const fvMeshSubset&);
253         //- Disallow default bitwise assignment
254         void operator=(const fvMeshSubset&);
257         //- Check if subset has been performed
258         bool checkCellSubset() const;
260         //- Mark points in Map
261         static void markPoints(const labelList&, Map<label>&);
263         //- Mark points (with 0) in labelList
264         static void markPoints(const labelList&, labelList&);
266         //- Adapt nCellsUsingFace for coupled faces becoming 'uncoupled'.
267         void doCoupledPatches
268         (
269             const bool syncPar,
270             labelList& nCellsUsingFace
271         ) const;
273         //- Subset of subset
274         static labelList subset
275         (
276             const label nElems,
277             const labelList& selectedElements,
278             const labelList& subsetMap
279         );
281         //- Create zones for submesh
282         void subsetZones();
284         //- Make dictionaries in case they are not present
285         //  Creates fvSchemes and fvSolution copying from base
286         void makeFvDictionaries();
289 public:
291     // Constructors
293         //- Construct given a mesh
294         explicit fvMeshSubset(const IOobject&, const fvMesh&);
297     // Destructor
299         ~fvMeshSubset();
302     // Member Functions
304         // Edit
306             //- Set the subset. Create "oldInternalFaces" patch for exposed
307             //  internal faces (patchID==-1) or use supplied patch.
308             //  Does not handle coupled patches correctly if only one side
309             //  gets deleted.
310             void setCellSubset
311             (
312                 const labelHashSet& globalCellMap,
313                 const label patchID = -1
314             );
316             //- Set the subset from all cells with region == currentRegion.
317             //  Create "oldInternalFaces" patch for exposed
318             //  internal faces (patchID==-1) or use supplied patch.
319             //  Handles coupled patches by if nessecary making coupled patch
320             //  face part of patchID (so uncoupled)
321             void setLargeCellSubset
322             (
323                 const labelList& region,
324                 const label currentRegion,
325                 const label patchID = -1,
326                 const bool syncCouples = true
327             );
329             //- setLargeCellSubset but with labelHashSet.
330             void setLargeCellSubset
331             (
332                 const labelHashSet& globalCellMap,
333                 const label patchID = -1,
334                 const bool syncPar = true
335             );
338         // Access
340             //- Return reference to base mesh
341             const fvMesh& baseMesh() const
342             {
343                 return baseMesh_;
344             }
346             //- Return reference to subset mesh
347             const fvMesh& subMesh() const;
349             fvMesh& subMesh();
351             //- Return reference to demand-driven subset pointMesh
352             const pointMesh& subPointMesh() const;
354             pointMesh& subPointMesh();
356             //- Return point map
357             const labelList& pointMap() const;
359             //- Return face map
360             const labelList& faceMap() const;
362             //- Return cell map
363             const labelList& cellMap() const;
365             //- Return patch map
366             const labelList& patchMap() const;
369         // Field mapping
371             //- Map volume field given map addressing
372             //  - it does not patch up the exposed patches with surface fields
373             //  - it does not handle patch type creation correctly
374             //  Please use only in dummy mapping operations
375             //  HJ, 31/Jul/2009
376             template<class Type>
377             static tmp<GeometricField<Type, fvPatchField, volMesh> >
378             meshToMesh
379             (
380                 const GeometricField<Type, fvPatchField, volMesh>&,
381                 const fvMesh& sMesh,
382                 const labelList& patchMap,
383                 const labelList& cellMap,
384                 const labelList& faceMap
385             );
387             //- Map volume field
388             template<class Type>
389             tmp<GeometricField<Type, fvPatchField, volMesh> >
390             interpolate
391             (
392                 const GeometricField<Type, fvPatchField, volMesh>&
393             ) const;
395             //- Map surface field
396             template<class Type>
397             static tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
398             meshToMesh
399             (
400                 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
401                 const fvMesh& sMesh,
402                 const labelList& patchMap,
403                 const labelList& faceMap
404             );
406             template<class Type>
407             tmp<GeometricField<Type, fvsPatchField, surfaceMesh> >
408             interpolate
409             (
410                 const GeometricField<Type, fvsPatchField, surfaceMesh>&
411             ) const;
413             //- Map point field
414             template<class Type>
415             static tmp<GeometricField<Type, pointPatchField, pointMesh> >
416             meshToMesh
417             (
418                 const GeometricField<Type, pointPatchField, pointMesh>&,
419                 const pointMesh& sMesh,
420                 const labelList& patchMap,
421                 const labelList& pointMap
422             );
424             template<class Type>
425             tmp<GeometricField<Type, pointPatchField, pointMesh> >
426             interpolate
427             (
428                 const GeometricField<Type, pointPatchField, pointMesh>&
429             ) const;
432         //- Write data
433         virtual bool writeData(Foam::Ostream&) const
434         {
435             return true;
436         }
441 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443 } // End namespace Foam
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447 #ifdef NoRepository
448 #   include "fvMeshSubsetInterpolate.C"
449 #endif
451 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453 #endif
455 // ************************************************************************* //