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-2020 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 namespace Foam
73 {
74 
75 class fvMeshLduAddressing;
76 class volMesh;
77 template<class Type>
78 class fvMatrix;
79 
80 /*---------------------------------------------------------------------------*\
81  Class fvMesh Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class fvMesh
85 :
86  public polyMesh,
87  public lduMesh,
88  public fvSchemes,
89  public surfaceInterpolation, // needs input from fvSchemes
90  public fvSolution,
91  public data
92 {
93 protected:
94 
95  // Private data
96 
97  //- Boundary mesh
99 
100  // Demand-driven data
101 
102  mutable fvMeshLduAddressing* lduPtr_;
103 
104  //- Current time index for cell volumes
105  // Note. The whole mechanism will be replaced once the
106  // dimensionedField is created and the dimensionedField
107  // will take care of the old-time levels.
108  mutable label curTimeIndex_;
109 
110  //- Cell volumes old time level
111  mutable void* VPtr_;
112 
113  //- Cell volumes old time level
115 
116  //- Cell volumes old-old time level
118 
119  //- Face area vectors
121 
122  //- Mag face area vectors
123  mutable surfaceScalarField* magSfPtr_;
124 
125  //- Cell centres
126  mutable slicedVolVectorField* CPtr_;
127 
128  //- Face centres
130 
131  //- Face motion fluxes
132  mutable surfaceScalarField* phiPtr_;
133 
134 
135  // Private Member Functions
136 
137  // Storage management
138 
139  //- Clear geometry but not the old-time cell volumes
140  void clearGeomNotOldVol();
141 
142  //- Clear geometry like clearGeomNotOldVol but recreate any
143  // geometric demand-driven data that was set
144  void updateGeomNotOldVol();
145 
146  //- Clear geometry
147  void clearGeom();
148 
149  //- Clear addressing
150  void clearAddressing(const bool isMeshUpdate = false);
151 
152  //- Preserve old volume(s)
153  void storeOldVol(const scalarField&);
154 
155 
156  // Make geometric data
157 
158  void makeSf() const;
159  void makeMagSf() const;
160 
161  void makeC() const;
162  void makeCf() const;
163 
164 
165  //- No copy construct
166  fvMesh(const fvMesh&) = delete;
167 
168  //- No copy assignment
169  void operator=(const fvMesh&) = delete;
170 
171 
172 public:
173 
174  // Public typedefs
175 
176  typedef fvMesh Mesh;
178 
179 
180  // Declare name of the class and its debug switch
181  ClassName("fvMesh");
182 
183 
184  // Constructors
185 
186  //- Construct from IOobject
187  explicit fvMesh(const IOobject& io, const bool doInit=true);
188 
189  //- Construct from IOobject or as zero-sized mesh
190  // Boundary is added using addFvPatches() member function
191  fvMesh(const IOobject& io, const zero, bool syncPar=true);
192 
193  //- Construct from components without boundary.
194  // Boundary is added using addFvPatches() member function
195  fvMesh
196  (
197  const IOobject& io,
198  pointField&& points,
199  faceList&& faces,
200  labelList&& allOwner,
201  labelList&& allNeighbour,
202  const bool syncPar = true
203  );
204 
205  //- Construct without boundary from cells rather than owner/neighbour.
206  // Boundary is added using addFvPatches() member function
207  fvMesh
208  (
209  const IOobject& io,
210  pointField&& points,
211  faceList&& faces,
212  cellList&& cells,
213  const bool syncPar = true
214  );
215 
216  //- Construct as copy (for dictionaries) and components without
217  // boundary. Boundary is added using addFvPatches() member function
218  fvMesh
219  (
220  const IOobject& io,
221  const fvMesh& baseMesh,
222  pointField&& points,
223  faceList&& faces,
224  labelList&& allOwner,
225  labelList&& allNeighbour,
226  const bool syncPar = true
227  );
228 
229  //- Construct as copy (for dictionaries) without boundary from cells
230  // rather than owner/neighbour. Boundary is added using addFvPatches()
231  // member function
232  fvMesh
233  (
234  const IOobject& io,
235  const fvMesh& baseMesh,
236  pointField&& points,
237  faceList&& faces,
238  cellList&& cells,
239  const bool syncPar = true
240  );
241 
242 
243  //- Destructor
244  virtual ~fvMesh();
245 
246 
247  // Member Functions
248 
249  // Helpers
250 
251  //- Initialise all non-demand-driven data
252  virtual bool init(const bool doInit);
253 
254  //- Add boundary patches. Constructor helper
255  void addFvPatches
256  (
257  PtrList<polyPatch>& plist,
258  const bool validBoundary = true
259  );
260 
261  //- Add boundary patches. Constructor helper
262  void addFvPatches
263  (
264  const List<polyPatch*>& p,
265  const bool validBoundary = true
266  );
267 
268  //- Update the mesh based on the mesh files saved in time
269  // directories
270  virtual readUpdateState readUpdate();
271 
272 
273  // Access
274 
275  //- Return the top-level database
276  const Time& time() const
277  {
278  return polyMesh::time();
279  }
280 
281  //- Return true if thisDb() is a valid DB
282  virtual bool hasDb() const
283  {
284  return true;
285  }
286 
287  //- Return the object registry - resolve conflict polyMesh/lduMesh
288  virtual const objectRegistry& thisDb() const
289  {
290  return polyMesh::thisDb();
291  }
292 
293  //- Return reference to name
294  // Note: name() is currently ambiguous due to derivation from
295  // surfaceInterpolation
296  const word& name() const
297  {
298  return polyMesh::name();
299  }
300 
301  //- Return reference to boundary mesh
302  const fvBoundaryMesh& boundary() const;
303 
304  //- Return ldu addressing
305  virtual const lduAddressing& lduAddr() const;
306 
307  //- Return a list of pointers for each patch
308  // with only those pointing to interfaces being set
309  virtual lduInterfacePtrsList interfaces() const
310  {
311  return boundary().interfaces();
312  }
313 
314  //- Return communicator used for parallel communication
315  virtual label comm() const
316  {
317  return polyMesh::comm();
318  }
319 
320 
321  // Overlap
322 
323  //- Interpolate interpolationCells only
324  virtual void interpolate(volScalarField&) const
325  {}
326 
327  //- Interpolate interpolationCells only
328  virtual void interpolate(volVectorField&) const
329  {}
330 
331  //- Interpolate interpolationCells only
332  virtual void interpolate(volSphericalTensorField&) const
333  {}
334 
335  //- Interpolate interpolationCells only
336  virtual void interpolate(volSymmTensorField&) const
337  {}
338 
339  //- Interpolate interpolationCells only
340  virtual void interpolate(volTensorField&) const
341  {}
342 
343  //- Interpolate interpolationCells only. No bcs.
344  virtual void interpolate(scalarField&) const
345  {}
346 
347  //- Interpolate interpolationCells only. No bcs.
348  virtual void interpolate(vectorField&) const
349  {}
350 
351  //- Interpolate interpolationCells only. No bcs.
352  virtual void interpolate(sphericalTensorField&) const
353  {}
354 
355  //- Interpolate interpolationCells only. No bcs.
356  virtual void interpolate(symmTensorField&) const
357  {}
358 
359  //- Interpolate interpolationCells only. No bcs.
360  virtual void interpolate(tensorField&) const
361  {}
362 
363  //- Solve returning the solution statistics given convergence
364  //- tolerance. Use the given solver controls
366  (
368  const dictionary&
369  ) const;
370 
371  //- Solve returning the solution statistics given convergence
372  //- tolerance. Use the given solver controls
374  (
376  const dictionary&
377  ) const;
378 
379  //- Solve returning the solution statistics given convergence
380  //- tolerance. Use the given solver controls
382  (
384  const dictionary&
385  ) const;
386 
387  //- Solve returning the solution statistics given convergence
388  //- tolerance. Use the given solver controls
390  (
392  const dictionary&
393  ) const;
394 
395  //- Solve returning the solution statistics given convergence
396  //- tolerance. Use the given solver controls
398  (
400  const dictionary&
401  ) const;
402 
403 
404  //- Internal face owner. Note bypassing virtual mechanism so
405  // e.g. relaxation always gets done using original addressing
406  const labelUList& owner() const
407  {
408  return fvMesh::lduAddr().lowerAddr();
409  }
410 
411  //- Internal face neighbour
412  const labelUList& neighbour() const
413  {
414  return fvMesh::lduAddr().upperAddr();
415  }
416 
417  //- Return cell volumes
418  const DimensionedField<scalar, volMesh>& V() const;
419 
420  //- Return old-time cell volumes
421  const DimensionedField<scalar, volMesh>& V0() const;
422 
423  //- Return old-old-time cell volumes
424  const DimensionedField<scalar, volMesh>& V00() const;
425 
426  //- Return sub-cycle cell volumes
428 
429  //- Return sub-cycle old-time cell volumes
431 
432  //- Return cell face area vectors
433  const surfaceVectorField& Sf() const;
434 
435  //- Return cell face area magnitudes
436  const surfaceScalarField& magSf() const;
437 
438  //- Return cell face motion fluxes
439  const surfaceScalarField& phi() const;
440 
441  //- Return cell centres as volVectorField
442  const volVectorField& C() const;
443 
444  //- Return face centres as surfaceVectorField
445  const surfaceVectorField& Cf() const;
446 
447  //- Return face deltas as surfaceVectorField
449 
450  //- Return a labelType of valid component indicators
451  // 1 : valid (solved)
452  // -1 : invalid (not solved)
453  template<class Type>
454  typename pTraits<Type>::labelType validComponents() const;
455 
456 
457  // Edit
458 
459  //- Clear all geometry and addressing
460  void clearOut();
461 
462  //- Update mesh corresponding to the given map
463  virtual void updateMesh(const mapPolyMesh& mpm);
464 
465  //- Avoid masking surfaceInterpolation method
467 
468  //- Move points, returns volumes swept by faces in motion
469  virtual tmp<scalarField> movePoints(const pointField&);
470 
471  //- Update all geometric data. This gets redirected up from
472  //- primitiveMesh level
473  virtual void updateGeom();
474 
475  //- Map all fields in time using given map.
476  virtual void mapFields(const mapPolyMesh& mpm);
477 
478  //- Remove boundary patches. Warning: fvPatchFields hold ref to
479  //- these fvPatches.
480  void removeFvBoundary();
481 
482  //- Return cell face motion fluxes
484 
485  //- Return old-time cell volumes
487 
488 
489  // Write
490 
491  //- Write the underlying polyMesh and other data
492  virtual bool writeObject
493  (
494  IOstreamOption streamOpt,
495  const bool valid
496  ) const;
497 
498  //- Write mesh using IO settings from time
499  virtual bool write(const bool valid = true) const;
500 
501 
502  // Member Operators
503 
504  //- Compares addresses
505  bool operator!=(const fvMesh& rhs) const;
506 
507  //- Compares addresses
508  bool operator==(const fvMesh& rhs) const;
509 };
510 
511 
512 template<>
513 typename pTraits<sphericalTensor>::labelType
514 fvMesh::validComponents<sphericalTensor>() const;
515 
516 
517 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
518 
519 } // End namespace Foam
520 
521 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
522 
523 #ifdef NoRepository
524  #include "fvMeshTemplates.C"
525  #include "fvPatchFvMeshTemplates.C"
526 #endif
527 
528 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
529 
530 #endif
531 
532 // ************************************************************************* //
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:536
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::fvMesh::interpolate
virtual void interpolate(sphericalTensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:351
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::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:70
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fvMesh::BoundaryMesh
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:176
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1027
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:287
DimensionedField.H
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:53
Foam::fvMesh::VPtr_
void * VPtr_
Cell volumes old time level.
Definition: fvMesh.H:110
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
Foam::fvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:700
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:355
Foam::fvMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: fvMesh.C:116
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Definition: fvMesh.C:630
primitiveMesh.H
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:1043
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:1049
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::fvMesh::interpolate
virtual void interpolate(volSymmTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:335
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:51
Foam::fvMesh::magSfPtr_
surfaceScalarField * magSfPtr_
Mag face area vectors.
Definition: fvMesh.H:122
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:56
Foam::fvMesh::SfPtr_
slicedSurfaceVectorField * SfPtr_
Face area vectors.
Definition: fvMesh.H:119
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:643
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:545
Foam::fvMesh::comm
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:314
Foam::fvMesh::CfPtr_
slicedSurfaceVectorField * CfPtr_
Face centres.
Definition: fvMesh.H:128
Foam::fvMesh::V0Ptr_
DimensionedField< scalar, volMesh > * V0Ptr_
Cell volumes old time level.
Definition: fvMesh.H:113
Foam::fvMesh::interpolate
virtual void interpolate(scalarField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:343
Foam::fvMesh::updateGeom
virtual void updateGeom()
Definition: fvMesh.C:931
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:327
Foam::Field< scalar >
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:685
Foam::fvMesh::interpolate
virtual void interpolate(volTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:339
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:939
fvPatchFvMeshTemplates.C
Foam::fvBoundaryMesh::interfaces
lduInterfacePtrsList interfaces() const
Definition: fvBoundaryMesh.C:124
Foam::fvMesh::operator=
void operator=(const fvMesh &)=delete
No copy assignment.
Foam::UPtrList< const lduInterface >
Foam::polyMesh::thisDb
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:507
SolverPerformance.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
fvSchemes.H
Foam::polyMesh::comm
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1313
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:121
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:411
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:600
Foam::fvMesh::Mesh
fvMesh Mesh
Definition: fvMesh.H:175
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:83
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:347
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:405
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:281
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:679
Foam::fvMesh::CPtr_
slicedVolVectorField * CPtr_
Cell centres.
Definition: fvMesh.H:125
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:271
Foam::fvMesh::lduPtr_
fvMeshLduAddressing * lduPtr_
Definition: fvMesh.H:101
Foam::lduAddressing::lowerAddr
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Foam::fvMesh::interpolate
virtual void interpolate(tensorField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:359
Foam::fvMesh::boundary_
fvBoundaryMesh boundary_
Boundary mesh.
Definition: fvMesh.H:97
Foam::List< face >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:54
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:116
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:308
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:76
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:275
data.H
Foam::fvMesh::curTimeIndex_
label curTimeIndex_
Current time index for cell volumes.
Definition: fvMesh.H:107
Foam::fvMesh::interpolate
virtual void interpolate(volSphericalTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:331
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:51
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:55
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:227
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:295
Foam::fvMesh::makeMagSf
void makeMagSf() const
Definition: fvMeshGeometry.C:75
Foam::fvMesh::interpolate
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:323
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:131
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:1004