1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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
29 Mesh data needed to do the Finite Volume discretisation.
32 fvMesh contains all the topological and geometric information
33 related to the mesh. It is also responsible for keeping the data
34 up-to-date. This is done by deleting the cell volume, face area,
35 cell/face centre, addressing and other derived information as
36 required and recalculating it as necessary. The fvMesh therefore
37 reserves the right to delete the derived information upon every
38 topological (mesh refinement/morphing) or geometric change (mesh
39 motion). It is therefore unsafe to keep local references to the
40 derived data outside of the time loop.
46 \*---------------------------------------------------------------------------*/
53 #include "primitiveMesh.H"
54 #include "fvBoundaryMesh.H"
55 #include "surfaceInterpolation.H"
56 #include "DimensionedField.H"
57 #include "volFieldsFwd.H"
58 #include "surfaceFieldsFwd.H"
59 #include "pointFieldsFwd.H"
60 #include "slicedVolFieldsFwd.H"
61 #include "slicedSurfaceFieldsFwd.H"
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 class fvMeshLduAddressing;
72 /*---------------------------------------------------------------------------*\
73 Class fvMesh Declaration
74 \*---------------------------------------------------------------------------*/
80 public surfaceInterpolation
85 fvBoundaryMesh boundary_;
90 mutable fvMeshLduAddressing* lduPtr_;
92 //- Current time index for cell volumes
93 // Note. The whole mechanism will be replaced once the
94 // dimensionedField is created and the dimensionedField
95 // will take care of the old-time levels.
96 mutable label curTimeIndex_;
99 mutable DimensionedField<scalar, volMesh>* VPtr_;
101 //- Cell volumes old time level
102 mutable DimensionedField<scalar, volMesh>* V0Ptr_;
104 //- Cell volumes old-old time level
105 mutable DimensionedField<scalar, volMesh>* V00Ptr_;
107 //- Face area vectors
108 mutable slicedSurfaceVectorField* SfPtr_;
110 //- Mag face area vectors
111 mutable surfaceScalarField* magSfPtr_;
114 mutable slicedVolVectorField* CPtr_;
117 mutable slicedSurfaceVectorField* CfPtr_;
119 //- Face motion fluxes
120 mutable surfaceScalarField* phiPtr_;
123 // Private Member Functions
125 // Make geometric data
129 //- Make face area magnitudes
130 void makeMagSf() const;
132 //- Make mesh motion fluxes
133 void makePhi() const;
135 //- Update mesh motion fluxes given swept volumes
136 void updatePhi(const scalarField& sweptVols) const;
138 //- Make cell centres
141 //- Make face centres
145 //- Disallow construct as copy
146 fvMesh(const fvMesh&);
148 //- Disallow assignment
149 void operator=(const fvMesh&);
152 // Storage management
154 //- Clear geometry but not the old-time cell volumes
155 void clearGeomNotOldVol();
161 void clearAddressing();
169 typedef fvBoundaryMesh BoundaryMesh;
172 // Declare name of the class and its debug switch
178 //- Construct from IOobject
179 explicit fvMesh(const IOobject& io);
181 //- Construct from components without boundary.
182 // Boundary is added using addFvPatches() member function
186 const Xfer<pointField>& points,
187 const Xfer<faceList>& faces,
188 const Xfer<labelList>& allOwner,
189 const Xfer<labelList>& allNeighbour,
190 const bool syncPar = true
193 //- Construct without boundary from cells rather than owner/neighbour.
194 // Boundary is added using addPatches() member function
198 const Xfer<pointField>& points,
199 const Xfer<faceList>& faces,
200 const Xfer<cellList>& cells,
201 const bool syncPar = true
204 //- Construct from cell shapes
208 const Xfer<pointField>& points,
209 const cellShapeList& shapes,
210 const faceListList& boundaryFaces,
211 const wordList& boundaryPatchNames,
212 const wordList& boundaryPatchTypes,
213 const word& defaultBoundaryPatchName,
214 const word& defaultBoundaryPatchType,
215 const wordList& boundaryPatchPhysicalTypes,
216 const bool syncPar = true
229 //- Add boundary patches. Constructor helper
232 const List<polyPatch*>&,
233 const bool validBoundary = true
236 //- Update the mesh based on the mesh files saved in time
238 virtual readUpdateState readUpdate();
243 //- Return the top-level database
244 const Time& time() const
246 return polyMesh::time();
249 //- Return the object registry - resolve conflict polyMesh/lduMesh
250 virtual const objectRegistry& thisDb() const
252 return polyMesh::thisDb();
255 //- Return reference to name
256 // Note: name() is currently ambiguous due to derivation from
257 // surfaceInterpolation
258 const word& name() const
260 return polyMesh::name();
263 //- Return reference to boundary mesh
264 const fvBoundaryMesh& boundary() const;
266 //- Return ldu addressing
267 virtual const lduAddressing& lduAddr() const;
269 //- Return a list of pointers for each patch
270 // with only those pointing to interfaces being set
271 virtual lduInterfacePtrsList interfaces() const
273 return boundary().interfaces();
276 //- Internal face owner
277 const unallocLabelList& owner() const
279 return lduAddr().lowerAddr();
282 //- Internal face neighbour
283 const unallocLabelList& neighbour() const
285 return lduAddr().upperAddr();
289 //- Return cell volumes
290 const DimensionedField<scalar, volMesh>& V() const;
292 //- Return old-time cell volumes
293 const DimensionedField<scalar, volMesh>& V0() const;
295 //- Return old-old-time cell volumes
296 const DimensionedField<scalar, volMesh>& V00() const;
298 //- Return sub-cycle cell volumes
299 tmp<DimensionedField<scalar, volMesh> > Vsc() const;
301 //- Return sub-cycl old-time cell volumes
302 tmp<DimensionedField<scalar, volMesh> > Vsc0() const;
304 //- Return cell face area vectors
305 const surfaceVectorField& Sf() const;
307 //- Return cell face area magnitudes
308 const surfaceScalarField& magSf() const;
310 //- Return cell face motion fluxes
311 const surfaceScalarField& phi() const;
313 //- Return cell centres as volVectorField
314 const volVectorField& C() const;
316 //- Return face centres as surfaceVectorField
317 const surfaceVectorField& Cf() const;
322 //- Clear all geometry and addressing
325 //- Update mesh corresponding to the given map
326 virtual void updateMesh(const mapPolyMesh& mpm);
328 //- Sync mesh update for topo change on other processors
329 // Locally, there is no topological change
330 virtual void syncUpdateMesh();
332 //- Move points, returns volumes swept by faces in motion
333 virtual tmp<scalarField> movePoints(const pointField&);
335 //- Map all fields in time using given map.
336 virtual void mapFields(const mapPolyMesh& mpm) const;
338 //- Map cell volumes in time using given map.
339 virtual void mapOldVolumes(const mapPolyMesh& mpm);
341 //- Remove boundary patches. Warning: fvPatchFields hold ref to
343 void removeFvBoundary();
345 //- Return cell face motion fluxes
346 surfaceScalarField& setPhi();
348 //- Return old-time cell volumes
349 DimensionedField<scalar, volMesh>& setV0();
354 //- Write the underlying polyMesh and other data
355 virtual bool writeObjects
357 IOstream::streamFormat fmt,
358 IOstream::versionNumber ver,
359 IOstream::compressionType cmp
362 //- Write mesh using IO settings from time
363 virtual bool write() const;
368 bool operator!=(const fvMesh&) const;
369 bool operator==(const fvMesh&) const;
373 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375 } // End namespace Foam
377 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
380 # include "fvPatchFvMeshTemplates.C"
383 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
387 // ************************************************************************* //