fvMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Class
25  Foam::fvMesh
26 
27 Description
28  Mesh data needed to do the Finite Volume discretisation.
29 
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.
40 
41 SourceFiles
42  fvMesh.C
43  fvMeshGeometry.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef fvMesh_H
48 #define fvMesh_H
49 
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"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 
71 class fvMeshLduAddressing;
72 class volMesh;
73 
74 
75 /*---------------------------------------------------------------------------*\
76  Class fvMesh Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 class fvMesh
80 :
81  public polyMesh,
82  public lduMesh,
83  public surfaceInterpolation,
84  public fvSchemes,
85  public fvSolution,
86  public data
87 {
88  // Private data
89 
90  //- Boundary mesh
91  fvBoundaryMesh boundary_;
92 
93 
94  // Demand-driven data
95 
96  mutable fvMeshLduAddressing* lduPtr_;
97 
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_;
103 
104  //- Cell volumes old time level
105  mutable void* VPtr_;
106 
107  //- Cell volumes old time level
108  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
109 
110  //- Cell volumes old-old time level
111  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
112 
113  //- Face area vectors
114  mutable slicedSurfaceVectorField* SfPtr_;
115 
116  //- Mag face area vectors
117  mutable surfaceScalarField* magSfPtr_;
118 
119  //- Cell centres
120  mutable slicedVolVectorField* CPtr_;
121 
122  //- Face centres
123  mutable slicedSurfaceVectorField* CfPtr_;
124 
125  //- Face motion fluxes
126  mutable surfaceScalarField* phiPtr_;
127 
128 
129  // Private Member Functions
130 
131  // Storage management
132 
133  //- Clear geometry but not the old-time cell volumes
134  void clearGeomNotOldVol();
135 
136  //- Clear geometry like clearGeomNotOldVol but recreate any
137  // geometric demand-driven data that was set
138  void updateGeomNotOldVol();
139 
140  //- Clear geometry
141  void clearGeom();
142 
143  //- Clear addressing
144  void clearAddressing(const bool isMeshUpdate = false);
145 
146  //- Preserve old volume(s)
147  void storeOldVol(const scalarField&);
148 
149 
150  // Make geometric data
151 
152  void makeSf() const;
153  void makeMagSf() const;
154 
155  void makeC() const;
156  void makeCf() const;
157 
158 
159  //- Disallow construct as copy
160  fvMesh(const fvMesh&);
161 
162  //- Disallow assignment
163  void operator=(const fvMesh&);
164 
165 
166 public:
167 
168  // Public typedefs
170  typedef fvMesh Mesh;
172 
173 
174  // Declare name of the class and its debug switch
175  ClassName("fvMesh");
176 
177 
178  // Constructors
179 
180  //- Construct from IOobject
181  explicit fvMesh(const IOobject& io);
182 
183  //- Construct from cellShapes with boundary.
184  fvMesh
185  (
186  const IOobject& io,
187  const Xfer<pointField>& points,
188  const cellShapeList& shapes,
189  const faceListList& boundaryFaces,
190  const wordList& boundaryPatchNames,
191  const PtrList<dictionary>& boundaryDicts,
192  const word& defaultBoundaryPatchName,
193  const word& defaultBoundaryPatchType,
194  const bool syncPar = true
195  );
196 
197  //- Construct from components without boundary.
198  // Boundary is added using addFvPatches() member function
199  fvMesh
200  (
201  const IOobject& io,
202  const Xfer<pointField>& points,
203  const Xfer<faceList>& faces,
204  const Xfer<labelList>& allOwner,
205  const Xfer<labelList>& allNeighbour,
206  const bool syncPar = true
207  );
208 
209  //- Construct without boundary from cells rather than owner/neighbour.
210  // Boundary is added using addPatches() member function
211  fvMesh
212  (
213  const IOobject& io,
214  const Xfer<pointField>& points,
215  const Xfer<faceList>& faces,
216  const Xfer<cellList>& cells,
217  const bool syncPar = true
218  );
219 
220 
221  //- Destructor
222  virtual ~fvMesh();
223 
224 
225  // Member Functions
226 
227  // Helpers
228 
229  //- Add boundary patches. Constructor helper
230  void addFvPatches
231  (
232  const List<polyPatch*>&,
233  const bool validBoundary = true
234  );
235 
236  //- Update the mesh based on the mesh files saved in time
237  // directories
238  virtual readUpdateState readUpdate();
239 
240 
241  // Access
242 
243  //- Return the top-level database
244  const Time& time() const
245  {
246  return polyMesh::time();
247  }
248 
249  //- Return the object registry - resolve conflict polyMesh/lduMesh
250  virtual const objectRegistry& thisDb() const
251  {
252  return polyMesh::thisDb();
253  }
254 
255  //- Return reference to name
256  // Note: name() is currently ambiguous due to derivation from
257  // surfaceInterpolation
258  const word& name() const
259  {
260  return polyMesh::name();
261  }
262 
263  //- Return reference to boundary mesh
264  const fvBoundaryMesh& boundary() const;
265 
266  //- Return ldu addressing
267  virtual const lduAddressing& lduAddr() const;
268 
269  //- Return a list of pointers for each patch
270  // with only those pointing to interfaces being set
271  virtual lduInterfacePtrsList interfaces() const
272  {
273  return boundary().interfaces();
274  }
275 
276  //- Return communicator used for parallel communication
277  virtual label comm() const
278  {
279  return polyMesh::comm();
280  }
281 
282  //- Internal face owner
283  const labelUList& owner() const
284  {
285  return lduAddr().lowerAddr();
286  }
287 
288  //- Internal face neighbour
289  const labelUList& neighbour() const
290  {
291  return lduAddr().upperAddr();
292  }
293 
294  //- Return cell volumes
295  const DimensionedField<scalar, volMesh>& V() const;
296 
297  //- Return old-time cell volumes
298  const DimensionedField<scalar, volMesh>& V0() const;
299 
300  //- Return old-old-time cell volumes
301  const DimensionedField<scalar, volMesh>& V00() const;
302 
303  //- Return sub-cycle cell volumes
305 
306  //- Return sub-cycl old-time cell volumes
308 
309  //- Return cell face area vectors
310  const surfaceVectorField& Sf() const;
311 
312  //- Return cell face area magnitudes
313  const surfaceScalarField& magSf() const;
314 
315  //- Return cell face motion fluxes
316  const surfaceScalarField& phi() const;
317 
318  //- Return cell centres as volVectorField
319  const volVectorField& C() const;
320 
321  //- Return face centres as surfaceVectorField
322  const surfaceVectorField& Cf() const;
323 
324  //- Return face deltas as surfaceVectorField
326 
327  //- Return a labelType of valid component indicators
328  // 1 : valid (solved)
329  // -1 : invalid (not solved)
330  template<class Type>
331  typename pTraits<Type>::labelType validComponents() const;
332 
333 
334  // Edit
335 
336  //- Clear all geometry and addressing
337  void clearOut();
338 
339  //- Update mesh corresponding to the given map
340  virtual void updateMesh(const mapPolyMesh& mpm);
341 
342  //- Move points, returns volumes swept by faces in motion
343  virtual tmp<scalarField> movePoints(const pointField&);
344 
345  //- Map all fields in time using given map.
346  virtual void mapFields(const mapPolyMesh& mpm);
347 
348  //- Remove boundary patches. Warning: fvPatchFields hold ref to
349  // these fvPatches.
350  void removeFvBoundary();
351 
352  //- Return cell face motion fluxes
354 
355  //- Return old-time cell volumes
357 
358 
359  // Write
360 
361  //- Write the underlying polyMesh and other data
362  virtual bool writeObject
363  (
367  const bool valid
368  ) const;
369 
370  //- Write mesh using IO settings from time
371  virtual bool write(const bool valid = true) const;
372 
373 
374  // Member Operators
375 
376  bool operator!=(const fvMesh&) const;
377  bool operator==(const fvMesh&) const;
378 };
379 
380 
381 template<>
383 fvMesh::validComponents<sphericalTensor>() const;
384 
385 
386 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
387 
388 } // End namespace Foam
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 #ifdef NoRepository
393  #include "fvMeshTemplates.C"
394  #include "fvPatchFvMeshTemplates.C"
395 #endif
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 #endif
400 
401 // ************************************************************************* //
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:230
const surfaceVectorField & Sf() const
Return cell face area vectors.
void clearAddressing()
Clear topological data.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:478
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:291
void addFvPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:460
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:877
const surfaceScalarField & phi() const
Return cell face motion fluxes.
Traits class for primitives.
Definition: pTraits.H:50
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:59
bool movePoints()
Do what is neccessary if the mesh has moved.
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:170
const cellList & cells() const
Cell to surface interpolation scheme. Included in fvMesh.
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:249
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:288
fvMesh Mesh
Definition: fvMesh.H:169
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
surfaceScalarField & setPhi()
Return cell face motion fluxes.
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:451
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:893
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:276
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:800
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:494
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:860
A class for handling words, derived from string.
Definition: word.H:59
virtual const labelUList & upperAddr() const =0
Return upper addressing.
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:551
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1036
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1200
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycl old-time cell volumes.
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Foam::fvMeshLduAddressing.
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
const Time & time() const
Return time.
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:282
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:47
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
bool operator==(const fvMesh &) const
Definition: fvMesh.C:899
Foam::fvBoundaryMesh.
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:50
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Version number type.
Definition: IOstream.H:96
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:487
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const volVectorField & C() const
Return cell centres as volVectorField.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:270
Namespace for OpenFOAM.
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:562
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
ClassName("fvMesh")