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 OpenFOAM Foundation
9  Copyright (C) 2016-2021 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 fvMesh_H
51 #define 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 "data.H"
61 #include "DimensionedField.H"
62 #include "volFieldsFwd.H"
63 #include "surfaceFieldsFwd.H"
64 #include "pointFieldsFwd.H"
65 #include "slicedVolFieldsFwd.H"
66 #include "slicedSurfaceFieldsFwd.H"
67 #include "className.H"
68 #include "SolverPerformance.H"
69 
70 
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 namespace Foam
75 {
76 
77 class fvMeshLduAddressing;
78 class volMesh;
79 template<class Type>
80 class fvMatrix;
81 
82 /*---------------------------------------------------------------------------*\
83  Class fvMesh Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 class fvMesh
87 :
88  public polyMesh,
89  public lduMesh,
90  public fvSchemes,
91  public surfaceInterpolation, // needs input from fvSchemes
92  public fvSolution,
93  public data
94 {
95 protected:
96 
97  // Private data
98 
99  //- Boundary mesh
101 
102  // Demand-driven data
103 
104  mutable fvMeshLduAddressing* lduPtr_;
105 
106  //- Current time index for cell volumes
107  // Note. The whole mechanism will be replaced once the
108  // dimensionedField is created and the dimensionedField
109  // will take care of the old-time levels.
110  mutable label curTimeIndex_;
111 
112  //- Cell volumes old time level
113  mutable void* VPtr_;
114 
115  //- Cell volumes old time level
117 
118  //- Cell volumes old-old time level
120 
121  //- Face area vectors
123 
124  //- Mag face area vectors
125  mutable surfaceScalarField* magSfPtr_;
126 
127  //- Cell centres
128  mutable slicedVolVectorField* CPtr_;
129 
130  //- Face centres
132 
133  //- Face motion fluxes
134  mutable surfaceScalarField* phiPtr_;
135 
136 
137  // Private Member Functions
138 
139  // Storage management
140 
141  //- Clear geometry but not the old-time cell volumes
142  void clearGeomNotOldVol();
143 
144  //- Clear geometry like clearGeomNotOldVol but recreate any
145  // geometric demand-driven data that was set
146  void updateGeomNotOldVol();
147 
148  //- Clear local geometry
149  void clearGeom();
150 
151  //- Clear local addressing
152  void clearAddressing(const bool isMeshUpdate = false);
153 
154  //- Clear local-only storage (geometry, addressing etc)
155  void clearOutLocal();
156 
157  //- Preserve old volume(s)
158  void storeOldVol(const scalarField&);
159 
160 
161  // Make geometric data
162 
163  void makeSf() const;
164  void makeMagSf() const;
165 
166  void makeC() const;
167  void makeCf() const;
168 
169 
170  //- No copy construct
171  fvMesh(const fvMesh&) = delete;
172 
173  //- No copy assignment
174  void operator=(const fvMesh&) = delete;
175 
176 
177 public:
178 
179  // Public typedefs
180 
181  typedef fvMesh 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 zero, bool syncPar=true);
197 
198  //- Construct from components without boundary.
199  // Boundary is added using addFvPatches() member function
200  fvMesh
201  (
202  const IOobject& io,
203  pointField&& points,
204  faceList&& faces,
205  labelList&& allOwner,
206  labelList&& allNeighbour,
207  const bool syncPar = true
208  );
209 
210  //- Construct without boundary from cells rather than owner/neighbour.
211  // Boundary is added using addFvPatches() member function
212  fvMesh
213  (
214  const IOobject& io,
215  pointField&& points,
216  faceList&& faces,
217  cellList&& cells,
218  const bool syncPar = true
219  );
220 
221  //- Construct as copy (for dictionaries) and components without
222  // boundary. Boundary is added using addFvPatches() member function
223  fvMesh
224  (
225  const IOobject& io,
226  const fvMesh& baseMesh,
227  pointField&& points,
228  faceList&& faces,
229  labelList&& allOwner,
230  labelList&& allNeighbour,
231  const bool syncPar = true
232  );
233 
234  //- Construct as copy (for dictionaries) without boundary from cells
235  // rather than owner/neighbour. Boundary is added using addFvPatches()
236  // member function
237  fvMesh
238  (
239  const IOobject& io,
240  const fvMesh& baseMesh,
241  pointField&& points,
242  faceList&& faces,
243  cellList&& cells,
244  const bool syncPar = true
245  );
246 
247 
248  //- Destructor
249  virtual ~fvMesh();
250 
251 
252  // Member Functions
253 
254  // Helpers
255 
256  //- Initialise all non-demand-driven data
257  virtual bool init(const bool doInit);
258 
259  //- Add boundary patches. Constructor helper
260  void addFvPatches
261  (
262  PtrList<polyPatch>& plist,
263  const bool validBoundary = true
264  );
265 
266  //- Add boundary patches. Constructor helper
267  void addFvPatches
268  (
269  const List<polyPatch*>& p,
270  const bool validBoundary = true
271  );
272 
273  //- Update the mesh based on the mesh files saved in time
274  // directories
275  virtual readUpdateState readUpdate();
276 
277 
278  // Access
279 
280  //- Return the top-level database
281  const Time& time() const
282  {
283  return polyMesh::time();
284  }
285 
286  //- Return true if thisDb() is a valid DB
287  virtual bool hasDb() const
288  {
289  return true;
290  }
291 
292  //- Return the object registry - resolve conflict polyMesh/lduMesh
293  virtual const objectRegistry& thisDb() const
294  {
295  return polyMesh::thisDb();
296  }
297 
298  //- Return reference to name
299  // Note: name() is currently ambiguous due to derivation from
300  // surfaceInterpolation
301  const word& name() const
302  {
303  return polyMesh::name();
304  }
305 
306  //- Return reference to boundary mesh
307  const fvBoundaryMesh& boundary() const;
308 
309  //- Return ldu addressing
310  virtual const lduAddressing& lduAddr() const;
311 
312  //- Return a list of pointers for each patch
313  // with only those pointing to interfaces being set
314  virtual lduInterfacePtrsList interfaces() const;
315 
316  //- Return communicator used for parallel communication
317  virtual label comm() const
318  {
319  return polyMesh::comm();
320  }
321 
322 
323  // Overlap
324 
325  //- Interpolate interpolationCells only
326  virtual void interpolate(volScalarField&) const
327  {}
328 
329  //- Interpolate interpolationCells only
330  virtual void interpolate(volVectorField&) const
331  {}
332 
333  //- Interpolate interpolationCells only
334  virtual void interpolate(volSphericalTensorField&) const
335  {}
336 
337  //- Interpolate interpolationCells only
338  virtual void interpolate(volSymmTensorField&) const
339  {}
340 
341  //- Interpolate interpolationCells only
342  virtual void interpolate(volTensorField&) const
343  {}
344 
345  //- Interpolate interpolationCells only. No bcs.
346  virtual void interpolate(scalarField&) const
347  {}
348 
349  //- Interpolate interpolationCells only. No bcs.
350  virtual void interpolate(vectorField&) const
351  {}
352 
353  //- Interpolate interpolationCells only. No bcs.
354  virtual void interpolate(sphericalTensorField&) const
355  {}
356 
357  //- Interpolate interpolationCells only. No bcs.
358  virtual void interpolate(symmTensorField&) const
359  {}
360 
361  //- Interpolate interpolationCells only. No bcs.
362  virtual void interpolate(tensorField&) const
363  {}
364 
365  //- Solve returning the solution statistics given convergence
366  //- tolerance. Use the given solver controls
368  (
370  const dictionary&
371  ) const;
372 
373  //- Solve returning the solution statistics given convergence
374  //- tolerance. Use the given solver controls
376  (
378  const dictionary&
379  ) const;
380 
381  //- Solve returning the solution statistics given convergence
382  //- tolerance. Use the given solver controls
384  (
386  const dictionary&
387  ) const;
388 
389  //- Solve returning the solution statistics given convergence
390  //- tolerance. Use the given solver controls
392  (
394  const dictionary&
395  ) const;
396 
397  //- Solve returning the solution statistics given convergence
398  //- tolerance. Use the given solver controls
400  (
402  const dictionary&
403  ) const;
404 
405 
406  //- Internal face owner. Note bypassing virtual mechanism so
407  // e.g. relaxation always gets done using original addressing
408  const labelUList& owner() const
409  {
410  return fvMesh::lduAddr().lowerAddr();
411  }
412 
413  //- Internal face neighbour
414  const labelUList& neighbour() const
415  {
416  return fvMesh::lduAddr().upperAddr();
417  }
418 
419  //- Return cell volumes
420  const DimensionedField<scalar, volMesh>& V() const;
421 
422  //- Return old-time cell volumes
423  const DimensionedField<scalar, volMesh>& V0() const;
424 
425  //- Return old-old-time cell volumes
426  const DimensionedField<scalar, volMesh>& V00() const;
427 
428  //- Return sub-cycle cell volumes
430 
431  //- Return sub-cycle old-time cell volumes
433 
434  //- Return cell face area vectors
435  const surfaceVectorField& Sf() const;
436 
437  //- Return cell face area magnitudes
438  const surfaceScalarField& magSf() const;
439 
440  //- Return cell face motion fluxes
441  const surfaceScalarField& phi() const;
442 
443  //- Return cell centres as volVectorField
444  const volVectorField& C() const;
445 
446  //- Return face centres as surfaceVectorField
447  const surfaceVectorField& Cf() const;
448 
449  //- Return face deltas as surfaceVectorField
451 
452  //- Return a labelType of valid component indicators
453  // 1 : valid (solved)
454  // -1 : invalid (not solved)
455  template<class Type>
456  typename pTraits<Type>::labelType validComponents() const;
457 
458 
459  // Edit
460 
461  //- Clear all geometry and addressing
462  void clearOut();
463 
464  //- Update mesh corresponding to the given map
465  virtual void updateMesh(const mapPolyMesh& mpm);
466 
467  //- Avoid masking surfaceInterpolation method
469 
470  //- Move points, returns volumes swept by faces in motion
471  virtual tmp<scalarField> movePoints(const pointField&);
472 
473  //- Update all geometric data. This gets redirected up from
474  //- primitiveMesh level
475  virtual void updateGeom();
476 
477  //- Map all fields in time using given map.
478  virtual void mapFields(const mapPolyMesh& mpm);
479 
480  //- Remove boundary patches. Warning: fvPatchFields hold ref to
481  //- these fvPatches.
482  void removeFvBoundary();
483 
484  //- Return cell face motion fluxes
486 
487  //- Return old-time cell volumes
489 
490 
491  // Write
492 
493  //- Write the underlying polyMesh and other data
494  virtual bool writeObject
495  (
496  IOstreamOption streamOpt,
497  const bool valid
498  ) const;
499 
500  //- Write mesh using IO settings from time
501  virtual bool write(const bool valid = true) const;
502 
503 
504  // Member Operators
505 
506  //- Compares addresses
507  bool operator!=(const fvMesh& rhs) const;
508 
509  //- Compares addresses
510  bool operator==(const fvMesh& rhs) const;
511 };
512 
513 
514 template<>
515 typename pTraits<sphericalTensor>::labelType
516 fvMesh::validComponents<sphericalTensor>() const;
517 
518 
519 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
520 
521 } // End namespace Foam
522 
523 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
524 
525 #ifdef NoRepository
526  #include "fvMeshTemplates.C"
527  #include "fvPatchFvMeshTemplates.C"
528 #endif
529 
530 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
531 
532 #endif
533 
534 // ************************************************************************* //
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::fvMesh::storeOldVol
void storeOldVol(const scalarField &)
Preserve old volume(s)
Definition: fvMesh.C:164
volFieldsFwd.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::fvMesh::~fvMesh
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:541
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::fvMesh::interpolate
virtual void interpolate(sphericalTensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:353
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fvMesh::BoundaryMesh
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:181
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1041
Foam::fvMesh::makeCf
void makeCf() const
Definition: fvMeshGeometry.C:145
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:157
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fvMesh::updateGeomNotOldVol
void updateGeomNotOldVol()
Clear geometry like clearGeomNotOldVol but recreate any.
Definition: fvMesh.C:82
Foam::fvMesh::Vsc0
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
Definition: fvMeshGeometry.C:290
Foam::fvMesh::thisDb
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:292
Foam::fvMesh::clearOutLocal
void clearOutLocal()
Clear local-only storage (geometry, addressing etc)
Definition: fvMesh.C:227
DimensionedField.H
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:52
Foam::fvMesh::VPtr_
void * VPtr_
Cell volumes old time level.
Definition: fvMesh.H:112
Foam::fvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:714
Foam::fvMesh::clearGeomNotOldVol
void clearGeomNotOldVol()
Clear geometry but not the old-time cell volumes.
Definition: fvMesh.C:54
Foam::fvMesh::setPhi
surfaceScalarField & setPhi()
Return cell face motion fluxes.
Definition: fvMeshGeometry.C:430
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::fvMesh::interpolate
virtual void interpolate(symmTensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:357
Foam::fvMesh::clearGeom
void clearGeom()
Clear local geometry.
Definition: fvMesh.C:116
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Definition: fvMesh.C:635
primitiveMesh.H
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1057
fvMeshTemplates.C
Foam::fvMesh::V00
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Definition: fvMeshGeometry.C:233
Foam::fvMesh::makeC
void makeC() const
Definition: fvMeshGeometry.C:108
polyMesh.H
Foam::fvMesh::operator==
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1063
Foam::lduAddressing::upperAddr
virtual const labelUList & upperAddr() const =0
Return upper addressing.
slicedSurfaceFieldsFwd.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
slicedVolFieldsFwd.H
Foam::polyMesh::comm
label comm() const noexcept
Return communicator used for parallel communication.
Definition: polyMesh.C:1313
Foam::fvMesh::interpolate
virtual void interpolate(volSymmTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:337
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:52
Foam::fvMesh::magSfPtr_
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:124
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::fvMesh::SfPtr_
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:121
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:648
Foam::fvMesh::V0
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
Definition: fvMeshGeometry.C:207
Foam::fvMesh::solve
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Definition: fvMesh.C:550
Foam::fvMesh::comm
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:316
Foam::fvMesh::CfPtr_
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:130
Foam::fvMesh::V0Ptr_
DimensionedField< scalar, volMesh > * V0Ptr_
Cell volumes old time level.
Definition: fvMesh.H:115
Foam::fvMesh::interpolate
virtual void interpolate(scalarField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:345
Foam::fvMesh::updateGeom
virtual void updateGeom()
Definition: fvMesh.C:945
Foam::fvMesh::makeSf
void makeSf() const
Definition: fvMeshGeometry.C:41
Foam::fvMesh::validComponents
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
Foam::fvMesh::interpolate
virtual void interpolate(volVectorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:329
Foam::Field< scalar >
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:708
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:691
Foam::fvMesh::interpolate
virtual void interpolate(volTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:341
className.H
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Foam::fvMesh::Vsc
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
Definition: fvMeshGeometry.C:262
Foam::fvSolution
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:50
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:953
fvPatchFvMeshTemplates.C
Foam::fvMesh::operator=
void operator=(const fvMesh &)=delete
No copy assignment.
Foam::UPtrList< const lduInterface >
SolverPerformance.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
fvSchemes.H
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
fvBoundaryMesh.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:413
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:605
Foam::fvMesh::Mesh
fvMesh Mesh
Definition: fvMesh.H:180
fvSolution.H
Foam::surfaceInterpolation::movePoints
virtual bool movePoints()
Do what is necessary if the mesh has moved.
Definition: surfaceInterpolation.C:151
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
surfaceInterpolation.H
Foam::fvMesh::interpolate
virtual void interpolate(vectorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:349
Foam::fvMesh::fvMesh
fvMesh(const fvMesh &)=delete
No copy construct.
Foam::fvMesh::owner
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:407
Foam::fvMesh::phi
const surfaceScalarField & phi() const
Return cell face motion fluxes.
Definition: fvMeshGeometry.C:408
Foam::fvMesh::hasDb
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: fvMesh.H:286
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::fvMesh::CPtr_
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:127
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::fvMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: fvMesh.C:276
Foam::fvMesh::lduPtr_
fvMeshLduAddressing * lduPtr_
Definition: fvMesh.H:103
Foam::lduAddressing::lowerAddr
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Foam::polyMesh::thisDb
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition: polyMesh.H:507
Foam::fvMesh::interpolate
virtual void interpolate(tensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:361
Foam::fvMesh::boundary_
fvBoundaryMesh boundary_
Boundary mesh.
Definition: fvMesh.H:99
Foam::List< face >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
pointFieldsFwd.H
Forwards and collection of common point field types.
Foam::fvMesh::V00Ptr_
DimensionedField< scalar, volMesh > * V00Ptr_
Cell volumes old-old time level.
Definition: fvMesh.H:118
Foam::UList< label >
surfaceFieldsFwd.H
Foam::surfaceInterpolation
Cell to surface interpolation scheme. Included in fvMesh.
Definition: surfaceInterpolation.H:58
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::SlicedGeometricField
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
Definition: slicedSurfaceFieldsFwd.H:58
Foam::fvMesh::setV0
DimensionedField< scalar, volMesh > & setV0()
Return old-time cell volumes.
Definition: fvMeshGeometry.C:220
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
data.H
Foam::fvMesh::curTimeIndex_
label curTimeIndex_
Current time index for cell volumes.
Definition: fvMesh.H:109
Foam::fvMesh::interpolate
virtual void interpolate(volSphericalTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:333
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
lduMesh.H
Foam::fvMesh::ClassName
ClassName("fvMesh")
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:352
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:52
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::fvMesh::delta
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
Definition: fvMeshGeometry.C:363
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62
Foam::lduInterfacePtrsList
UPtrList< const lduInterface > lduInterfacePtrsList
List of coupled interface fields to be used in coupling.
Definition: lduInterfacePtrsList.H:44
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:239
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:300
Foam::fvMesh::makeMagSf
void makeMagSf() const
Definition: fvMeshGeometry.C:75
Foam::fvMesh::interpolate
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:325
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179
Foam::fvMesh::phiPtr_
surfaceScalarField * phiPtr_
Face motion fluxes.
Definition: fvMesh.H:133
Foam::zero
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:62
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:319
Foam::fvMesh::writeObject
virtual bool writeObject(IOstreamOption streamOpt, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1018