BUGFIX: Seg-fault in multiphaseInterFoam. Author: Henrik Rusche. Merge: Hrvoje Jasak
[foam-extend-3.2.git] / src / finiteVolume / fvMesh / fvMesh.H
blob3bb3a115f452398ecb766e17d62373d1a68a3f2f
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::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 "DimensionedField.H"
56 #include "volFieldsFwd.H"
57 #include "surfaceFieldsFwd.H"
58 #include "pointFieldsFwd.H"
59 #include "slicedVolFieldsFwd.H"
60 #include "slicedSurfaceFieldsFwd.H"
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 namespace Foam
67 class fvMeshLduAddressing;
68 class volMesh;
71 /*---------------------------------------------------------------------------*\
72                            Class fvMesh Declaration
73 \*---------------------------------------------------------------------------*/
75 class fvMesh
77     public polyMesh,
78     public lduMesh,
79     public surfaceInterpolation
81     // Private data
83         //- Boundary mesh
84         fvBoundaryMesh boundary_;
87     // Demand-driven data
89         mutable fvMeshLduAddressing* lduPtr_;
91         //- Current time index for cell volumes
92         //  Note.  The whole mechanism will be replaced once the
93         //  dimensionedField is created and the dimensionedField
94         //  will take care of the old-time levels.
95         mutable label curTimeIndex_;
97         //- Cell volumes
98         mutable DimensionedField<scalar, volMesh>* VPtr_;
100         //- Cell volumes old time level
101         mutable DimensionedField<scalar, volMesh>* V0Ptr_;
103         //- Cell volumes old-old time level
104         mutable DimensionedField<scalar, volMesh>* V00Ptr_;
106         //- Face area vectors
107         mutable slicedSurfaceVectorField* SfPtr_;
109         //- Mag face area vectors
110         mutable surfaceScalarField* magSfPtr_;
112         //- Cell centres
113         mutable slicedVolVectorField* CPtr_;
115         //- Face centres
116         mutable slicedSurfaceVectorField* CfPtr_;
118         //- Face motion fluxes
119         mutable surfaceScalarField* phiPtr_;
122     // Private Member Functions
124        // Make geometric data
126             void makeSf() const;
128             //- Make face area magnitudes
129             void makeMagSf() const;
131             //- Make mesh motion fluxes
132             void makePhi() const;
134             //- Update mesh motion fluxes given swept volumes
135             void updatePhi(const scalarField& sweptVols) const;
137             //- Make cell centres
138             void makeC() const;
140             //- Make face centres
141             void makeCf() const;
144         //- Disallow construct as copy
145         fvMesh(const fvMesh&);
147         //- Disallow assignment
148         void operator=(const fvMesh&);
151         // Storage management
153             //- Clear geometry but not the old-time cell volumes
154             void clearGeomNotOldVol();
156             //- Clear geometry
157             void clearGeom();
159             //- Clear addressing
160             void clearAddressing();
163 public:
165     // Public typedefs
167         typedef fvMesh Mesh;
168         typedef fvBoundaryMesh BoundaryMesh;
171     // Declare name of the class and its debug switch
172     TypeName("fvMesh");
175     // Constructors
177         //- Construct from IOobject
178         explicit fvMesh(const IOobject& io);
180         //- Construct from components without boundary.
181         //  Boundary is added using addFvPatches() member function
182         fvMesh
183         (
184             const IOobject& io,
185             const Xfer<pointField>& points,
186             const Xfer<faceList>& faces,
187             const Xfer<labelList>& allOwner,
188             const Xfer<labelList>& allNeighbour,
189             const bool syncPar = true
190         );
192         //- Construct without boundary from cells rather than owner/neighbour.
193         //  Boundary is added using addPatches() member function
194         fvMesh
195         (
196             const IOobject& io,
197             const Xfer<pointField>& points,
198             const Xfer<faceList>& faces,
199             const Xfer<cellList>& cells,
200             const bool syncPar = true
201         );
203         //- Construct from cell shapes
204         fvMesh
205         (
206             const IOobject& io,
207             const Xfer<pointField>& points,
208             const cellShapeList& shapes,
209             const faceListList& boundaryFaces,
210             const wordList& boundaryPatchNames,
211             const wordList& boundaryPatchTypes,
212             const word& defaultBoundaryPatchName,
213             const word& defaultBoundaryPatchType,
214             const wordList& boundaryPatchPhysicalTypes,
215             const bool syncPar = true
216         );
219     // Destructor
221         virtual ~fvMesh();
224     // Member Functions
226         // Helpers
228             //- Add boundary patches. Constructor helper
229             void addFvPatches
230             (
231                 const List<polyPatch*>&,
232                 const bool validBoundary = true
233             );
235             //- Update the mesh based on the mesh files saved in time
236             //  directories
237             virtual readUpdateState readUpdate();
240         // Access
242             //- Return the top-level database
243             const Time& time() const
244             {
245                 return polyMesh::time();
246             }
248             //- Return the object registry - resolve conflict polyMesh/lduMesh
249             virtual const objectRegistry& thisDb() const
250             {
251                 return polyMesh::thisDb();
252             }
254             //- Return reference to name
255             //  Note: name() is currently ambiguous due to derivation from
256             //  surfaceInterpolation
257             const word& name() const
258             {
259                 return polyMesh::name();
260             }
262             //- Return reference to boundary mesh
263             const fvBoundaryMesh& boundary() const;
265             //- Return ldu addressing
266             virtual const lduAddressing& lduAddr() const;
268             //- Return a list of pointers for each patch
269             //  with only those pointing to interfaces being set
270             virtual lduInterfacePtrsList interfaces() const
271             {
272                 return boundary().interfaces();
273             }
275             //- Internal face owner
276             const unallocLabelList& owner() const
277             {
278                 return lduAddr().lowerAddr();
279             }
281             //- Internal face neighbour
282             const unallocLabelList& neighbour() const
283             {
284                 return lduAddr().upperAddr();
285             }
288             //- Return cell volumes
289             const DimensionedField<scalar, volMesh>& V() const;
291             //- Return old-time cell volumes
292             const DimensionedField<scalar, volMesh>& V0() const;
294             //- Return old-old-time cell volumes
295             const DimensionedField<scalar, volMesh>& V00() const;
297             //- Return sub-cycle cell volumes
298             tmp<DimensionedField<scalar, volMesh> > Vsc() const;
300             //- Return sub-cycl old-time cell volumes
301             tmp<DimensionedField<scalar, volMesh> > Vsc0() const;
303             //- Return cell face area vectors
304             const surfaceVectorField& Sf() const;
306             //- Return cell face area magnitudes
307             const surfaceScalarField& magSf() const;
309             //- Return cell face motion fluxes
310             const surfaceScalarField& phi() const;
312             //- Return cell centres as volVectorField
313             const volVectorField& C() const;
315             //- Return face centres as surfaceVectorField
316             const surfaceVectorField& Cf() const;
319         // Edit
321             //- Clear all geometry and addressing
322             void clearOut();
324             //- Update mesh corresponding to the given map
325             virtual void updateMesh(const mapPolyMesh& mpm);
327             //- Sync mesh update for topo change on other processors
328             //  Locally, there is no topological change
329             virtual void syncUpdateMesh();
331             //- Move points, returns volumes swept by faces in motion
332             virtual tmp<scalarField> movePoints(const pointField&);
334             //- Map all fields in time using given map.
335             virtual void mapFields(const mapPolyMesh& mpm) const;
337             //- Map cell volumes in time using given map.
338             virtual void mapOldVolumes(const mapPolyMesh& mpm);
340             //- Remove boundary patches. Warning: fvPatchFields hold ref to
341             //  these fvPatches.
342             void removeFvBoundary();
344             //- Return cell face motion fluxes
345             surfaceScalarField& setPhi();
347             //- Return old-time cell volumes
348             DimensionedField<scalar, volMesh>& setV0();
351         // Write
353             //- Write the underlying polyMesh and other data
354             virtual bool writeObjects
355             (
356                 IOstream::streamFormat fmt,
357                 IOstream::versionNumber ver,
358                 IOstream::compressionType cmp
359             ) const;
361             //- Write mesh using IO settings from time
362             virtual bool write() const;
365     // Member Operators
367         bool operator!=(const fvMesh&) const;
368         bool operator==(const fvMesh&) const;
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
374 } // End namespace Foam
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
378 #ifdef NoRepository
379 #   include "fvPatchFvMeshTemplates.C"
380 #endif
382 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 #endif
386 // ************************************************************************* //