ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / finiteVolume / fvMesh / fvMesh.H
blobc2d414dceaac2c0efaa276dfda595d51ce0da8f5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::fvMesh
27 Description
28     Mesh data needed to do the Finite Volume discretisation.
30     NOTE ON USAGE:
31     fvMesh contains all the topological and geometric information
32     related to the mesh.  It is also responsible for keeping the data
33     up-to-date.  This is done by deleting the cell volume, face area,
34     cell/face centre, addressing and other derived information as
35     required and recalculating it as necessary.  The fvMesh therefore
36     reserves the right to delete the derived information upon every
37     topological (mesh refinement/morphing) or geometric change (mesh
38     motion).  It is therefore unsafe to keep local references to the
39     derived data outside of the time loop.
41 SourceFiles
42     fvMesh.C
43     fvMeshGeometry.C
45 \*---------------------------------------------------------------------------*/
47 #ifndef fvMesh_H
48 #define fvMesh_H
50 #include "polyMesh.H"
51 #include "lduMesh.H"
52 #include "primitiveMesh.H"
53 #include "fvBoundaryMesh.H"
54 #include "surfaceInterpolation.H"
55 #include "fvSchemes.H"
56 #include "fvSolution.H"
57 #include "data.H"
58 #include "DimensionedField.H"
59 #include "volFieldsFwd.H"
60 #include "surfaceFieldsFwd.H"
61 #include "pointFieldsFwd.H"
62 #include "slicedVolFieldsFwd.H"
63 #include "slicedSurfaceFieldsFwd.H"
64 #include "className.H"
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 namespace Foam
71 class fvMeshLduAddressing;
72 class volMesh;
75 /*---------------------------------------------------------------------------*\
76                            Class fvMesh Declaration
77 \*---------------------------------------------------------------------------*/
79 class fvMesh
81     public polyMesh,
82     public lduMesh,
83     public surfaceInterpolation,
84     public fvSchemes,
85     public fvSolution,
86     public data
88     // Private data
90         //- Boundary mesh
91         fvBoundaryMesh boundary_;
94     // Demand-driven data
96         mutable fvMeshLduAddressing* lduPtr_;
98         //- Current time index for cell volumes
99         //  Note.  The whole mechanism will be replaced once the
100         //  dimensionedField is created and the dimensionedField
101         //  will take care of the old-time levels.
102         mutable label curTimeIndex_;
104         //- Cell volumes old time level
105         mutable void* VPtr_;
107         //- Cell volumes old time level
108         mutable DimensionedField<scalar, volMesh>* V0Ptr_;
110         //- Cell volumes old-old time level
111         mutable DimensionedField<scalar, volMesh>* V00Ptr_;
113         //- Face area vectors
114         mutable slicedSurfaceVectorField* SfPtr_;
116         //- Mag face area vectors
117         mutable surfaceScalarField* magSfPtr_;
119         //- Cell centres
120         mutable slicedVolVectorField* CPtr_;
122         //- Face centres
123         mutable slicedSurfaceVectorField* CfPtr_;
125         //- Face motion fluxes
126         mutable surfaceScalarField* phiPtr_;
129     // Private Member Functions
131         // Storage management
133             //- Clear geometry but not the old-time cell volumes
134             void clearGeomNotOldVol();
136             //- Clear geometry
137             void clearGeom();
139             //- Clear addressing
140             void clearAddressing();
143        // Make geometric data
145             void makeSf() const;
146             void makeMagSf() const;
148             void makeC() const;
149             void makeCf() const;
152         //- Disallow construct as copy
153         fvMesh(const fvMesh&);
155         //- Disallow assignment
156         void operator=(const fvMesh&);
159 public:
161     // Public typedefs
163         typedef fvMesh Mesh;
164         typedef fvBoundaryMesh BoundaryMesh;
167     // Declare name of the class and its debug switch
168     ClassName("fvMesh");
171     // Constructors
173         //- Construct from IOobject
174         explicit fvMesh(const IOobject& io);
176         //- Construct from components without boundary.
177         //  Boundary is added using addFvPatches() member function
178         fvMesh
179         (
180             const IOobject& io,
181             const Xfer<pointField>& points,
182             const Xfer<faceList>& faces,
183             const Xfer<labelList>& allOwner,
184             const Xfer<labelList>& allNeighbour,
185             const bool syncPar = true
186         );
188         //- Construct without boundary from cells rather than owner/neighbour.
189         //  Boundary is added using addPatches() member function
190         fvMesh
191         (
192             const IOobject& io,
193             const Xfer<pointField>& points,
194             const Xfer<faceList>& faces,
195             const Xfer<cellList>& cells,
196             const bool syncPar = true
197         );
200     //- Destructor
201     virtual ~fvMesh();
204     // Member Functions
206         // Helpers
208             //- Add boundary patches. Constructor helper
209             void addFvPatches
210             (
211                 const List<polyPatch*>&,
212                 const bool validBoundary = true
213             );
215             //- Update the mesh based on the mesh files saved in time
216             //  directories
217             virtual readUpdateState readUpdate();
220         // Access
222             //- Return the top-level database
223             const Time& time() const
224             {
225                 return polyMesh::time();
226             }
228             //- Return the object registry - resolve conflict polyMesh/lduMesh
229             virtual const objectRegistry& thisDb() const
230             {
231                 return polyMesh::thisDb();
232             }
234             //- Return reference to name
235             //  Note: name() is currently ambiguous due to derivation from
236             //  surfaceInterpolation
237             const word& name() const
238             {
239                 return polyMesh::name();
240             }
242             //- Return reference to boundary mesh
243             const fvBoundaryMesh& boundary() const;
245             //- Return ldu addressing
246             virtual const lduAddressing& lduAddr() const;
248             //- Return a list of pointers for each patch
249             //  with only those pointing to interfaces being set
250             virtual lduInterfacePtrsList interfaces() const
251             {
252                 return boundary().interfaces();
253             }
255             //- Internal face owner
256             const labelUList& owner() const
257             {
258                 return lduAddr().lowerAddr();
259             }
261             //- Internal face neighbour
262             const labelUList& neighbour() const
263             {
264                 return lduAddr().upperAddr();
265             }
267             //- Return cell volumes
268             const DimensionedField<scalar, volMesh>& V() const;
270             //- Return old-time cell volumes
271             const DimensionedField<scalar, volMesh>& V0() const;
273             //- Return old-old-time cell volumes
274             const DimensionedField<scalar, volMesh>& V00() const;
276             //- Return sub-cycle cell volumes
277             tmp<DimensionedField<scalar, volMesh> > Vsc() const;
279             //- Return sub-cycl old-time cell volumes
280             tmp<DimensionedField<scalar, volMesh> > Vsc0() const;
282             //- Return cell face area vectors
283             const surfaceVectorField& Sf() const;
285             //- Return cell face area magnitudes
286             const surfaceScalarField& magSf() const;
288             //- Return cell face motion fluxes
289             const surfaceScalarField& phi() const;
291             //- Return cell centres as volVectorField
292             const volVectorField& C() const;
294             //- Return face centres as surfaceVectorField
295             const surfaceVectorField& Cf() const;
298         // Edit
300             //- Clear all geometry and addressing
301             void clearOut();
303             //- Update mesh corresponding to the given map
304             virtual void updateMesh(const mapPolyMesh& mpm);
306             //- Move points, returns volumes swept by faces in motion
307             virtual tmp<scalarField> movePoints(const pointField&);
309             //- Map all fields in time using given map.
310             virtual void mapFields(const mapPolyMesh& mpm);
312             //- Remove boundary patches. Warning: fvPatchFields hold ref to
313             //  these fvPatches.
314             void removeFvBoundary();
316             //- Return cell face motion fluxes
317             surfaceScalarField& setPhi();
319             //- Return old-time cell volumes
320             DimensionedField<scalar, volMesh>& setV0();
323         // Write
325             //- Write the underlying polyMesh and other data
326             virtual bool writeObjects
327             (
328                 IOstream::streamFormat fmt,
329                 IOstream::versionNumber ver,
330                 IOstream::compressionType cmp
331             ) const;
333             //- Write mesh using IO settings from time
334             virtual bool write() const;
337     // Member Operators
339         bool operator!=(const fvMesh&) const;
340         bool operator==(const fvMesh&) const;
344 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
346 } // End namespace Foam
348 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350 #ifdef NoRepository
351 #   include "fvPatchFvMeshTemplates.C"
352 #endif
354 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356 #endif
358 // ************************************************************************* //