Forward compatibility: flex
[foam-extend-3.2.git] / src / finiteVolume / fvMesh / singleCellFvMesh / singleCellFvMesh.H
blobeaf7230a531ed9f11796814d4e56dd8e935e23a9
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::singleCellFvMesh
27 Description
28     fvMesh as subset of other mesh. Consists of one cell and all original
29     bounday faces. Useful when manipulating boundary data. Single internal
30     cell only needed to be able to manipulate in a standard way.
32 SourceFiles
33     singleCellFvMesh.C
34     singleCellFvMeshInterpolate.C
36 \*---------------------------------------------------------------------------*/
38 #ifndef singleCellFvMesh_H
39 #define singleCellFvMesh_H
41 #include "fvPatchFieldMapper.H"
42 #include "fvMesh.H"
43 #include "labelIOList.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
50 /*---------------------------------------------------------------------------*\
51                       Class singleCellFvMesh Declaration
52 \*---------------------------------------------------------------------------*/
54 class singleCellFvMesh
56     public fvMesh
58     // Private data
60         const labelListIOList patchFaceAgglomeration_;
62         //- From patch faces back to agglomeration or fine mesh
63         labelListIOList patchFaceMap_;
65         //- From fine mesh faces to coarse mesh
66         labelIOList reverseFaceMap_;
68         //- From coarse points back to original mesh
69         labelIOList pointMap_;
71         //- From fine points to coarse mesh
72         labelIOList reversePointMap_;
75     // Private Member Functions
77         //- Calculate agglomerated mesh
78         void agglomerateMesh(const fvMesh&, const labelListList&);
81         //- Disallow default bitwise copy construct
82         singleCellFvMesh(const singleCellFvMesh&);
84         //- Disallow default bitwise assignment
85         void operator=(const singleCellFvMesh&);
88 public:
90         //- Patch field mapper class for non-agglomerated meshes
91         class directPatchFieldMapper
92         :
93             public fvPatchFieldMapper
94         {
95             // Private data
97                 //- Mapping
98                 const UList<label>& directAddressing_;
100                 //- Size before mapping
101                 const label sizeBeforeMapping_;
104         public:
106                 //- Construct given addressing
107                 directPatchFieldMapper
108                 (
109                     const UList<label>& directAddressing,
110                     const label sizeBeforeMapping
111                 )
112                 :
113                     directAddressing_(directAddressing),
114                     sizeBeforeMapping_(sizeBeforeMapping)
115                 {}
117                 virtual label size() const
118                 {
119                     return directAddressing_.size();
120                 }
122                 virtual label sizeBeforeMapping() const
123                 {
124                     return sizeBeforeMapping_;
125                 }
127                 virtual bool direct() const
128                 {
129                     return true;
130                 }
132                 virtual const UList<label>& directAddressing() const
133                 {
134                     return directAddressing_;
135                 }
136         };
138         //- Patch field mapper class for agglomerated meshes
139         class agglomPatchFieldMapper
140         :
141             public fvPatchFieldMapper
142         {
143             // Private data
145                 //- Addressing
146                 const labelListList& addressing_;
148                 //- Weights
149                 const scalarListList& weights_;
151                 //- Size before mapping
152                 label sizeBeforeMapping_;
155         public:
157                 //- Construct given addressing
158                 agglomPatchFieldMapper
159                 (
160                     const labelListList& addressing,
161                     const scalarListList& weights,
162                     const label sizeBeforeMapping
163                 )
164                 :
165                     addressing_(addressing),
166                     weights_(weights),
167                     sizeBeforeMapping_(sizeBeforeMapping)
168                 {}
170                 virtual label size() const
171                 {
172                     return addressing_.size();
173                 }
175                 virtual label sizeBeforeMapping() const
176                 {
177                     return sizeBeforeMapping_;
178                 }
180                 virtual bool direct() const
181                 {
182                     return false;
183                 }
185                 virtual const labelListList& addressing() const
186                 {
187                     return addressing_;
188                 }
190                 virtual const scalarListList& weights() const
191                 {
192                     return weights_;
193                 }
194         };
198     // Constructors
200         //- Construct from fvMesh and no agglomeration
201         singleCellFvMesh(const IOobject& io, const fvMesh&);
203         //- Construct from fvMesh and agglomeration of boundary faces.
204         //  agglomeration is per patch, per patch face index the agglomeration
205         //  the face goes into.
206         singleCellFvMesh
207         (
208             const IOobject& io,
209             const fvMesh&,
210             const labelListList& patchFaceAgglomeration
211         );
213         //- Read from IOobject
214         singleCellFvMesh(const IOobject& io);
216     // Member Functions
218         bool agglomerate() const
219         {
220             return patchFaceAgglomeration_.size() > 0;
221         }
223         //- From patchFace on this back to original mesh or agglomeration
224         const labelListList& patchFaceMap() const
225         {
226             return patchFaceMap_;
227         }
229         //- From point on this back to original mesh
230         const labelList& pointMap() const
231         {
232             return pointMap_;
233         }
235         //- From face on original mesh to face on this
236         const labelList& reverseFaceMap() const
237         {
238             return reverseFaceMap_;
239         }
241         //- From point on original mesh to point on this (or -1 for removed
242         //  points)
243         const labelList& reversePointMap() const
244         {
245             return reversePointMap_;
246         }
248         //- Map volField. Internal field set to average, patch fields straight
249         //  copies.
250         template<class Type>
251         tmp<GeometricField<Type, fvPatchField, volMesh> >
252         interpolate
253         (
254             const GeometricField<Type, fvPatchField, volMesh>&
255         ) const;
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 } // End namespace Foam
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 #ifdef NoRepository
267 #   include "singleCellFvMeshInterpolate.C"
268 #endif
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272 #endif
274 // ************************************************************************* //