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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017,2022 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::fvMesh
29 
30 Description
31  Mesh data needed to do the Finite Volume discretisation.
32 
33  NOTE ON USAGE:
34  fvMesh contains all the topological and geometric information
35  related to the mesh. It is also responsible for keeping the data
36  up-to-date. This is done by deleting the cell volume, face area,
37  cell/face centre, addressing and other derived information as
38  required and recalculating it as necessary. The fvMesh therefore
39  reserves the right to delete the derived information upon every
40  topological (mesh refinement/morphing) or geometric change (mesh
41  motion). It is therefore unsafe to keep local references to the
42  derived data outside of the time loop.
43 
44 SourceFiles
45  fvMesh.C
46  fvMeshGeometry.C
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #ifndef Foam_fvMesh_H
51 #define Foam_fvMesh_H
52 
53 #include "polyMesh.H"
54 #include "lduMesh.H"
55 #include "primitiveMesh.H"
56 #include "fvBoundaryMesh.H"
57 #include "surfaceInterpolation.H"
58 #include "fvSchemes.H"
59 #include "fvSolution.H"
60 #include "volFieldsFwd.H"
61 #include "surfaceFieldsFwd.H"
62 #include "pointFieldsFwd.H"
63 #include "SlicedDimensionedField.H"
64 #include "slicedVolFieldsFwd.H"
65 #include "slicedSurfaceFieldsFwd.H"
66 #include "className.H"
67 #include "SolverPerformance.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 namespace Foam
72 {
73 
74 // Forward Declarations
75 class fvMeshLduAddressing;
76 class volMesh;
77 template<class Type> class fvMatrix;
78 
79 /*---------------------------------------------------------------------------*\
80  Class fvMesh Declaration
81 \*---------------------------------------------------------------------------*/
82 
83 class fvMesh
84 :
85  public polyMesh,
86  public lduMesh,
87  public fvSchemes,
88  public surfaceInterpolation, // needs input from fvSchemes
89  public fvSolution
90 {
91 protected:
92 
93  // Private data
94 
95  //- Boundary mesh
97 
98 
99  // Demand-driven data
100 
101  mutable fvMeshLduAddressing* lduPtr_;
102 
103  //- Current time index for cell volumes
104  // Note. The whole mechanism will be replaced once the
105  // dimensionedField is created and the dimensionedField
106  // will take care of the old-time levels.
107  mutable label curTimeIndex_;
108 
109  //- Cell volumes
111 
112  //- Cell volumes old time level
114 
115  //- Cell volumes old-old time level
118  //- Face area vectors
120 
121  //- Mag face area vectors
123 
124  //- Cell centres
125  mutable slicedVolVectorField* CPtr_;
126 
127  //- Face centres
129 
130  //- Face motion fluxes
131  mutable surfaceScalarField* phiPtr_;
133 
134  // Private Member Functions
135 
136  // Storage management
138  //- Clear geometry but not the old-time cell volumes
139  void clearGeomNotOldVol();
140 
141  //- Clear geometry like clearGeomNotOldVol but recreate any
142  // geometric demand-driven data that was set
143  void updateGeomNotOldVol();
144 
145  //- Clear local geometry
146  void clearGeom();
148  //- Clear local addressing
149  void clearAddressing(const bool isMeshUpdate = false);
150 
151  //- Clear local-only storage (geometry, addressing etc)
152  void clearOutLocal();
153 
154  //- Preserve old volume(s)
155  void storeOldVol(const scalarField&);
156 
157 
158  // Make geometric data
159 
160  void makeSf() const;
161  void makeMagSf() const;
162 
163  void makeC() const;
164  void makeCf() const;
165 
166 
167  //- No copy construct
168  fvMesh(const fvMesh&) = delete;
169 
170  //- No copy assignment
171  void operator=(const fvMesh&) = delete;
172 
173 
174 public:
175 
176  // Public Typedefs
177 
178  //- The mesh type
179  typedef fvMesh Mesh;
180 
181  //- The boundary type associated with the mesh
183 
184 
185  // Declare name of the class and its debug switch
186  ClassName("fvMesh");
187 
188 
189  // Constructors
190 
191  //- Construct from IOobject
192  explicit fvMesh(const IOobject& io, const bool doInit=true);
193 
194  //- Construct from IOobject or as zero-sized mesh
195  // Boundary is added using addFvPatches() member function
196  fvMesh(const IOobject& io, const Foam::zero, bool syncPar=true);
197 
198  //- Construct as copy (for dictionaries) and zero-sized components.
199  // Boundary is added using addFvPatches() member function
200  fvMesh
201  (
202  const IOobject& io,
203  const fvMesh& baseMesh,
204  const Foam::zero,
205  const bool syncPar = true
206  );
207 
208  //- Construct from components without boundary.
209  // Boundary is added using addFvPatches() member function
210  fvMesh
211  (
212  const IOobject& io,
213  pointField&& points,
215  labelList&& allOwner,
216  labelList&& allNeighbour,
217  const bool syncPar = true
218  );
220  //- Construct without boundary from cells rather than owner/neighbour.
221  // Boundary is added using addFvPatches() member function
222  fvMesh
223  (
224  const IOobject& io,
225  pointField&& points,
226  faceList&& faces,
227  cellList&& cells,
228  const bool syncPar = true
229  );
230 
231  //- Copy construct (for dictionaries) with components, without boundary.
232  // Boundary is added using addFvPatches() member function
233  fvMesh
234  (
235  const IOobject& io,
236  const fvMesh& baseMesh,
237  pointField&& points,
238  faceList&& faces,
239  labelList&& allOwner,
240  labelList&& allNeighbour,
241  const bool syncPar = true
242  );
243 
244  //- Copy construct (for dictionaries) with cells, without boundary.
245  // Boundary is added using addFvPatches() member function
246  fvMesh
247  (
248  const IOobject& io,
249  const fvMesh& baseMesh,
250  pointField&& points,
251  faceList&& faces,
252  cellList&& cells,
253  const bool syncPar = true
254  );
255 
256 
257  //- Destructor
258  virtual ~fvMesh();
259 
260 
261  // Member Functions
262 
263  // Helpers
264 
265  //- Initialise all non-demand-driven data
266  virtual bool init(const bool doInit);
267 
268  //- Add boundary patches. Constructor helper
269  void addFvPatches
270  (
271  polyPatchList& plist,
272  const bool validBoundary = true
273  );
274 
275  //- Add boundary patches. Constructor helper
276  void addFvPatches
277  (
278  const List<polyPatch*>& p,
279  const bool validBoundary = true
280  );
281 
282  //- Update the mesh based on the mesh files saved in time
283  // directories
284  virtual readUpdateState readUpdate();
285 
286 
287  // Access
288 
289  //- Return the top-level database
290  const Time& time() const
291  {
292  return polyMesh::time();
293  }
294 
295  //- Return true if thisDb() is a valid DB
296  virtual bool hasDb() const
297  {
298  return true;
299  }
300 
301  //- Return the object registry - resolve conflict polyMesh/lduMesh
302  virtual const objectRegistry& thisDb() const
303  {
304  return polyMesh::thisDb();
305  }
306 
307  //- Return reference to name
308  // Note: name() is currently ambiguous due to derivation from
309  // surfaceInterpolation
310  const word& name() const
311  {
312  return polyMesh::name();
313  }
314 
315  //- Return reference to boundary mesh
316  const fvBoundaryMesh& boundary() const noexcept
317  {
318  return boundary_;
319  }
320 
321  //- Return ldu addressing
322  virtual const lduAddressing& lduAddr() const;
323 
324  //- Return a list of pointers for each patch
325  // with only those pointing to interfaces being set
326  virtual lduInterfacePtrsList interfaces() const;
327 
328  //- Return communicator used for parallel communication
329  virtual label comm() const
330  {
331  return polyMesh::comm();
332  }
333 
334 
335  // Solution Control
336 
337  //- Non-null if fvSchemes exists (can test as bool).
338  const fvSchemes* hasSchemes() const;
339 
340  //- Non-null if fvSolution exists (can test as bool).
341  const fvSolution* hasSolution() const;
342 
343  //- Read-access to the fvSchemes controls
344  const fvSchemes& schemes() const;
345 
346  //- Read/write-access to the fvSchemes controls
347  fvSchemes& schemes();
348 
349  //- Read-access to the fvSolution controls
350  const fvSolution& solution() const;
351 
352  //- Read/write-access to the fvSolution controls
353  fvSolution& solution();
354 
355 
356  // Overlap
357 
358  //- Interpolate interpolationCells only
359  virtual void interpolate(volScalarField&) const
360  {}
361 
362  //- Interpolate interpolationCells only
363  virtual void interpolate(volVectorField&) const
364  {}
365 
366  //- Interpolate interpolationCells only
367  virtual void interpolate(volSphericalTensorField&) const
368  {}
369 
370  //- Interpolate interpolationCells only
371  virtual void interpolate(volSymmTensorField&) const
372  {}
373 
374  //- Interpolate interpolationCells only
375  virtual void interpolate(volTensorField&) const
376  {}
377 
378  //- Interpolate interpolationCells only. No bcs.
379  virtual void interpolate(scalarField&) const
380  {}
381 
382  //- Interpolate interpolationCells only. No bcs.
383  virtual void interpolate(vectorField&) const
384  {}
385 
386  //- Interpolate interpolationCells only. No bcs.
387  virtual void interpolate(sphericalTensorField&) const
388  {}
389 
390  //- Interpolate interpolationCells only. No bcs.
391  virtual void interpolate(symmTensorField&) const
392  {}
393 
394  //- Interpolate interpolationCells only. No bcs.
395  virtual void interpolate(tensorField&) const
396  {}
397 
398  //- Solve returning the solution statistics given convergence
399  //- tolerance. Use the given solver controls
401  (
403  const dictionary&
404  ) const;
405 
406  //- Solve returning the solution statistics given convergence
407  //- tolerance. Use the given solver controls
409  (
411  const dictionary&
412  ) const;
413 
414  //- Solve returning the solution statistics given convergence
415  //- tolerance. Use the given solver controls
417  (
419  const dictionary&
420  ) const;
421 
422  //- Solve returning the solution statistics given convergence
423  //- tolerance. Use the given solver controls
425  (
427  const dictionary&
428  ) const;
429 
430  //- Solve returning the solution statistics given convergence
431  //- tolerance. Use the given solver controls
433  (
435  const dictionary&
436  ) const;
437 
438 
439  //- Internal face owner. Note bypassing virtual mechanism so
440  // e.g. relaxation always gets done using original addressing
441  const labelUList& owner() const
442  {
443  return fvMesh::lduAddr().lowerAddr();
444  }
445 
446  //- Internal face neighbour
447  const labelUList& neighbour() const
448  {
449  return fvMesh::lduAddr().upperAddr();
450  }
451 
452  //- Return cell volumes
453  const DimensionedField<scalar, volMesh>& V() const;
454 
455  //- Return old-time cell volumes
456  const DimensionedField<scalar, volMesh>& V0() const;
457 
458  //- Return old-old-time cell volumes
460 
461  //- Return sub-cycle cell volumes
463 
464  //- Return sub-cycle old-time cell volumes
466 
467  //- Return cell face area vectors
468  const surfaceVectorField& Sf() const;
469 
470  //- Return cell face area magnitudes
471  const surfaceScalarField& magSf() const;
472 
473  //- Return cell face unit normals
475 
476  //- Return cell face motion fluxes
477  const surfaceScalarField& phi() const;
478 
479  //- Return cell centres as volVectorField
480  const volVectorField& C() const;
481 
482  //- Return face centres as surfaceVectorField
483  const surfaceVectorField& Cf() const;
484 
485  //- Return face deltas as surfaceVectorField
487 
488  //- Return a labelType of valid component indicators
489  // 1 : valid (solved)
490  // -1 : invalid (not solved)
491  template<class Type>
492  typename pTraits<Type>::labelType validComponents() const;
493 
494 
495  // Edit
496 
497  //- Clear all geometry and addressing
498  void clearOut();
499 
500  //- Update mesh corresponding to the given map
501  virtual void updateMesh(const mapPolyMesh& mpm);
502 
503  //- Avoid masking surfaceInterpolation method
505 
506  //- Move points, returns volumes swept by faces in motion
507  virtual void movePoints(const pointField&);
508 
509  //- Update all geometric data. This gets redirected up from
510  //- primitiveMesh level
511  virtual void updateGeom();
512 
513  //- Map all fields in time using given map.
514  virtual void mapFields(const mapPolyMesh& mpm);
515 
516  //- Remove boundary patches. Warning: fvPatchFields hold ref to
517  //- these fvPatches.
518  void removeFvBoundary();
519 
520  //- Return cell face motion fluxes (or null)
522 
523  //- Return old-time cell volumes
525 
526 
527  // Write
528 
529  //- Write the underlying polyMesh and other data
530  virtual bool writeObject
531  (
532  IOstreamOption streamOpt,
533  const bool writeOnProc
534  ) const;
535 
536  //- Write mesh using IO settings from time
537  virtual bool write(const bool writeOnProc = true) const;
538 
539 
540  // Member Operators
541 
542  //- Compares addresses
543  bool operator!=(const fvMesh& rhs) const;
544 
545  //- Compares addresses
546  bool operator==(const fvMesh& rhs) const;
547 };
548 
549 
550 template<>
552 fvMesh::validComponents<sphericalTensor>() const;
553 
554 
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
556 
557 } // End namespace Foam
558 
559 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
560 
561 #ifdef NoRepository
562  #include "fvMeshTemplates.C"
563  #include "fvPatchFvMeshTemplates.C"
564 #endif
565 
566 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
567 
568 #endif
569 
570 // ************************************************************************* //
label comm() const noexcept
The communicator used for parallel communication.
Definition: polyMesh.H:700
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:137
void operator=(const fvMesh &)=delete
No copy assignment.
const surfaceVectorField & Sf() const
Return cell face area vectors.
void clearAddressing()
Clear topological data.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition: fvMesh.C:692
surfaceScalarField * phiPtr_
Face motion fluxes.
Definition: fvMesh.H:147
void makeC() const
const fvSchemes & schemes() const
Read-access to the fvSchemes controls.
Definition: fvMesh.C:582
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1079
Specialisation of DimensionedField that holds a slice of a given field so that it acts as a Dimension...
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
DimensionedField< scalar, volMesh > * V0Ptr_
Cell volumes old time level.
Definition: fvMesh.H:117
void clearOutLocal()
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:215
const surfaceScalarField & phi() const
Return cell face motion fluxes.
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
virtual bool movePoints()
Do what is necessary if the mesh has moved.
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:157
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
const fvSolution & solution() const
Read-access to the fvSolution controls.
Definition: fvMesh.C:594
fvBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: fvMesh.H:219
const cellList & cells() const
Forwards and collection of common point field types.
A simple container for options an IOstream can normally have.
Cell to surface interpolation scheme. Included in fvMesh.
void makeCf() const
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:459
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
fvSchemes(const fvSchemes &)=delete
No copy construct.
label curTimeIndex_
Current time index for cell volumes.
Definition: fvMesh.H:107
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:69
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:84
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:580
fvMesh Mesh
The mesh type.
Definition: fvMesh.H:214
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:47
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:562
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:415
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:994
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:704
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
const Time & time() const noexcept
Return time registry.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1118
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:741
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:109
fvSolution(const fvSolution &)=delete
No copy construct.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
refPtr< surfaceScalarField > setPhi()
Return cell face motion fluxes (or null)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Solve returning the solution statistics given convergence tolerance. Use the given solver controls...
Definition: fvMesh.C:607
fvMesh(const fvMesh &)=delete
No copy construct.
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
DimensionedField< scalar, volMesh > * V00Ptr_
Cell volumes old-old time level.
Definition: fvMesh.H:122
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1124
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
const direction noexcept
Definition: Scalar.H:258
Foam::fvMeshLduAddressing.
void makeMagSf() const
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: fvMesh.H:368
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:662
fvMeshLduAddressing * lduPtr_
Definition: fvMesh.H:98
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:572
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:132
fvBoundaryMesh boundary_
Boundary mesh.
Definition: fvMesh.H:93
const fvSchemes * hasSchemes() const
Non-null if fvSchemes exists (can test as bool).
Definition: fvMesh.C:570
const word & name() const
Return reference to name.
Definition: fvMesh.H:387
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:758
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:51
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Foam::fvBoundaryMesh.
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const fvSolution * hasSolution() const
Non-null if fvSolution exists (can test as bool).
Definition: fvMesh.C:576
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:263
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:127
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:71
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition: polyMesh.H:689
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void updateGeom()
Update all geometric data. This gets redirected up from primitiveMesh level.
Definition: fvMesh.C:984
SlicedDimensionedField< scalar, volMesh > * VPtr_
Cell volumes.
Definition: fvMesh.H:112
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:90
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< surfaceVectorField > unitSf() const
Return cell face unit normals.
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:142
void makeSf() const
Namespace for OpenFOAM.
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:764
ClassName("fvMesh")