Merge branch 'master' of ssh://git.code.sf.net/p/foam-extend/foam-extend-3.2
[foam-extend-3.2.git] / src / dbns / multigrid / mgMeshLevel / fineMgMeshLevel.H
blob1fb1186cf797f86d13f94dc9981f42cb9f2b73f7
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::fineMgMeshLevel
27 Description
28     Fine multigrid mesh level.  Pach fvMesh information into mgMeshLevel
29     interface.  Allow for access to fvMesh at the fine level.
31 SourceFiles
32     fineMgMeshLevel.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef fineMgMeshLevel_H
37 #define fineMgMeshLevel_H
39 #include "mgMeshLevel.H"
40 #include "fvMesh.H"
41 #include "volFields.H"
42 #include "surfaceFields.H"
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 namespace Foam
49 /*---------------------------------------------------------------------------*\
50                        Class fineMgMeshLevel Declaration
51 \*---------------------------------------------------------------------------*/
53 class fineMgMeshLevel
55     public mgMeshLevel
57     // Private data
59         //- Reference to the mesh
60         const fvMesh& mesh_;
63     // Private Member Functions
65         //- Disallow default bitwise copy construct
66         fineMgMeshLevel(const fineMgMeshLevel&);
68         //- Disallow default bitwise assignment
69         void operator=(const fineMgMeshLevel&);
72 public:
74     // Static data
76         //- Runtime type information
77         TypeName("fineMgMeshLevel");
80     // Constructors
82         //- Construct from mesh
83         fineMgMeshLevel(const fvMesh& mesh)
84         :
85             mesh_(mesh)
86         {}
89     //- Destructor
91         virtual ~fineMgMeshLevel()
92         {}
95     // Member Functions
97         //- Is this the finest level?
98         virtual bool finest() const
99         {
100             return true;
101         }
103         //- Return access to fine level
104         virtual const mgMeshLevel& fineLevel() const
105         {
106             return *this;
107         }
109         //- Return access to mesh
110         virtual const fvMesh& mesh() const
111         {
112             return mesh_;
113         }
116         // Sizes
118             //- Return number of cells
119             virtual label nCells() const
120             {
121                 return mesh_.nCells();
122             }
124             //- Return number of internal faces
125             virtual label nInternalFaces() const
126             {
127                 return mesh_.nInternalFaces();
128             }
130             //- Return number of patches
131             virtual label nPatches() const
132             {
133                 return mesh_.boundary().size();
134             }
137         // Addressing
139             //- Return owner addressing
140             virtual const unallocLabelList& owner() const
141             {
142                 return mesh_.owner();
143             }
145             //- Return neighbour addressing
146             const unallocLabelList& neighbour() const
147             {
148                 return mesh_.neighbour();
149             }
151             //- Return face-cell addressing for all patches
152             virtual const unallocLabelList& faceCells
153             (
154                 const label patchNo
155             ) const
156             {
157                 return mesh_.boundary()[patchNo].faceCells();
158             }
161         // Geometrical data
163             //- Return number of valid geometric dimensions in the mesh
164             virtual label nGeometricD() const
165             {
166                 return mesh_.nGeometricD();
167             }
169             //- Return cell centres
170             virtual const vectorField& cellCentres() const
171             {
172                 return mesh_.C().internalField();
173             }
175             //- Return face centres
176             virtual const vectorField& faceCentres() const
177             {
178                 return mesh_.Cf().internalField();
179             }
181             //- Return cell volumes
182             virtual const scalarField& cellVolumes() const
183             {
184                 return mesh_.V().field();
185             }
187             //- Return face area vectors
188             virtual const vectorField& faceAreas() const
189             {
190                 return mesh_.Sf().internalField();
191             }
193             //- Return face area magnitudes
194             virtual const scalarField& magFaceAreas() const
195             {
196                 return mesh_.magSf().internalField();
197             }
199             //- Return face area vectors for patches
200             virtual const vectorField& patchFaceAreas
201             (
202                 const label patchNo
203             ) const
204             {
205                 return mesh_.Sf().boundaryField()[patchNo];
206             }
208             //- Return face area magnitudes for patches
209             virtual const scalarField& magPatchFaceAreas
210             (
211                 const label patchNo
212             ) const
213             {
214                 return mesh_.magSf().boundaryField()[patchNo];
215             }
217             //- Return cell to face vectors for patches
218             //  Note: no clustering of boundary faces allowed:
219             //        using fine-level geometry
220             virtual tmp<vectorField> patchDeltaR
221             (
222                 const label patchNo
223             ) const
224             {
225                 return mesh_.Cf().boundaryField()[patchNo]
226                     - mesh_.C().boundaryField()[patchNo].patchInternalField();
227             }
229             //- Return neighbour cell to face vectors for coupled patches
230             //  Note: no clustering of boundary faces allowed:
231             //        using fine-level geometry
232             virtual tmp<vectorField> patchDeltaRNeighbour
233             (
234                 const label patchNo
235             ) const
236             {
237                 // Check.  HJ, 4/Jun/2013
238                 return patchDeltaR(patchNo)
239                     - mesh_.boundary()[patchNo].delta();
240             }
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 } // End namespace Foam
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 #endif
252 // ************************************************************************* //