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-2019 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 /*---------------------------------------------------------------------------*\
82  Class fvMesh Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class fvMesh
86 :
87  public polyMesh,
88  public lduMesh,
89  public surfaceInterpolation,
90  public fvSchemes,
91  public fvSolution,
92  public data
93 {
94  // Private data
95 
96  //- Boundary mesh
97  fvBoundaryMesh boundary_;
98 
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
114  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
115 
116  //- Cell volumes old-old time level
117  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
118 
119  //- Face area vectors
120  mutable slicedSurfaceVectorField* SfPtr_;
121 
122  //- Mag face area vectors
123  mutable surfaceScalarField* magSfPtr_;
124 
125  //- Cell centres
126  mutable slicedVolVectorField* CPtr_;
127 
128  //- Face centres
129  mutable slicedSurfaceVectorField* CfPtr_;
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);
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 
217  //- Destructor
218  virtual ~fvMesh();
219 
220 
221  // Member Functions
222 
223  // Helpers
224 
225  //- Add boundary patches. Constructor helper
226  void addFvPatches
227  (
228  PtrList<polyPatch>& plist,
229  const bool validBoundary = true
230  );
231 
232  //- Add boundary patches. Constructor helper
233  void addFvPatches
234  (
235  const List<polyPatch*>& p,
236  const bool validBoundary = true
237  );
238 
239  //- Update the mesh based on the mesh files saved in time
240  // directories
241  virtual readUpdateState readUpdate();
242 
243 
244  // Access
245 
246  //- Return the top-level database
247  const Time& time() const
248  {
249  return polyMesh::time();
250  }
251 
252  //- Return true if thisDb() is a valid DB
253  virtual bool hasDb() const
254  {
255  return true;
256  }
257 
258  //- Return the object registry - resolve conflict polyMesh/lduMesh
259  virtual const objectRegistry& thisDb() const
260  {
261  return polyMesh::thisDb();
262  }
263 
264  //- Return reference to name
265  // Note: name() is currently ambiguous due to derivation from
266  // surfaceInterpolation
267  const word& name() const
268  {
269  return polyMesh::name();
270  }
271 
272  //- Return reference to boundary mesh
273  const fvBoundaryMesh& boundary() const;
274 
275  //- Return ldu addressing
276  virtual const lduAddressing& lduAddr() const;
277 
278  //- Return a list of pointers for each patch
279  // with only those pointing to interfaces being set
280  virtual lduInterfacePtrsList interfaces() const
281  {
282  return boundary().interfaces();
283  }
284 
285  //- Return communicator used for parallel communication
286  virtual label comm() const
287  {
288  return polyMesh::comm();
289  }
290 
291 
292  // Overlap
293 
294  //- Interpolate interpolationCells only
295  virtual void interpolate(volScalarField&) const
296  {}
297 
298  //- Interpolate interpolationCells only
299  virtual void interpolate(volVectorField&) const
300  {}
301 
302  //- Interpolate interpolationCells only
303  virtual void interpolate(volSphericalTensorField&) const
304  {}
305 
306  //- Interpolate interpolationCells only
307  virtual void interpolate(volSymmTensorField&) const
308  {}
309 
310  //- Interpolate interpolationCells only
311  virtual void interpolate(volTensorField&) const
312  {}
313 
314  //- Interpolate interpolationCells only. No bcs.
315  virtual void interpolate(scalarField&) const
316  {}
317 
318  //- Interpolate interpolationCells only. No bcs.
319  virtual void interpolate(vectorField&) const
320  {}
321 
322  //- Interpolate interpolationCells only. No bcs.
323  virtual void interpolate(sphericalTensorField&) const
324  {}
325 
326  //- Interpolate interpolationCells only. No bcs.
327  virtual void interpolate(symmTensorField&) const
328  {}
329 
330  //- Interpolate interpolationCells only. No bcs.
331  virtual void interpolate(tensorField&) const
332  {}
333 
334  //- Solve returning the solution statistics given convergence
335  //- tolerance. Use the given solver controls
337  (
339  const dictionary&
340  ) const;
341 
342  //- Solve returning the solution statistics given convergence
343  //- tolerance. Use the given solver controls
345  (
347  const dictionary&
348  ) const;
349 
350  //- Solve returning the solution statistics given convergence
351  //- tolerance. Use the given solver controls
353  (
355  const dictionary&
356  ) const;
357 
358  //- Solve returning the solution statistics given convergence
359  //- tolerance. Use the given solver controls
361  (
363  const dictionary&
364  ) const;
365 
366  //- Solve returning the solution statistics given convergence
367  //- tolerance. Use the given solver controls
369  (
371  const dictionary&
372  ) const;
373 
374 
375  //- Internal face owner. Note bypassing virtual mechanism so
376  // e.g. relaxation always gets done using original addressing
377  const labelUList& owner() const
378  {
379  return fvMesh::lduAddr().lowerAddr();
380  }
381 
382  //- Internal face neighbour
383  const labelUList& neighbour() const
384  {
385  return fvMesh::lduAddr().upperAddr();
386  }
387 
388  //- Return cell volumes
389  const DimensionedField<scalar, volMesh>& V() const;
390 
391  //- Return old-time cell volumes
392  const DimensionedField<scalar, volMesh>& V0() const;
393 
394  //- Return old-old-time cell volumes
395  const DimensionedField<scalar, volMesh>& V00() const;
396 
397  //- Return sub-cycle cell volumes
399 
400  //- Return sub-cycle old-time cell volumes
402 
403  //- Return cell face area vectors
404  const surfaceVectorField& Sf() const;
405 
406  //- Return cell face area magnitudes
407  const surfaceScalarField& magSf() const;
408 
409  //- Return cell face motion fluxes
410  const surfaceScalarField& phi() const;
411 
412  //- Return cell centres as volVectorField
413  const volVectorField& C() const;
414 
415  //- Return face centres as surfaceVectorField
416  const surfaceVectorField& Cf() const;
417 
418  //- Return face deltas as surfaceVectorField
420 
421  //- Return a labelType of valid component indicators
422  // 1 : valid (solved)
423  // -1 : invalid (not solved)
424  template<class Type>
425  typename pTraits<Type>::labelType validComponents() const;
426 
427 
428  // Edit
429 
430  //- Clear all geometry and addressing
431  void clearOut();
432 
433  //- Update mesh corresponding to the given map
434  virtual void updateMesh(const mapPolyMesh& mpm);
435 
436  //- Move points, returns volumes swept by faces in motion
437  virtual tmp<scalarField> movePoints(const pointField&);
438 
439  //- Map all fields in time using given map.
440  virtual void mapFields(const mapPolyMesh& mpm);
441 
442  //- Remove boundary patches. Warning: fvPatchFields hold ref to
443  // these fvPatches.
444  void removeFvBoundary();
445 
446  //- Return cell face motion fluxes
448 
449  //- Return old-time cell volumes
451 
452 
453  // Write
454 
455  //- Write the underlying polyMesh and other data
456  virtual bool writeObject
457  (
461  const bool valid
462  ) const;
463 
464  //- Write mesh using IO settings from time
465  virtual bool write(const bool valid = true) const;
466 
467 
468  // Member Operators
469 
470  //- Compares addresses
471  bool operator!=(const fvMesh& rhs) const;
472 
473  //- Compares addresses
474  bool operator==(const fvMesh& rhs) const;
475 };
476 
477 
478 template<>
479 typename pTraits<sphericalTensor>::labelType
480 fvMesh::validComponents<sphericalTensor>() const;
481 
482 
483 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
484 
485 } // End namespace Foam
486 
487 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
488 
489 #ifdef NoRepository
490  #include "fvMeshTemplates.C"
491  #include "fvPatchFvMeshTemplates.C"
492 #endif
493 
494 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
495 
496 #endif
497 
498 // ************************************************************************* //
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
volFieldsFwd.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::fvMesh::~fvMesh
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:428
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:322
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:46
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:939
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:143
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::fvMesh::Vsc0
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
Definition: fvMeshGeometry.C:307
Foam::fvMesh::thisDb
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:258
Foam::fvMesh::writeObject
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:914
DimensionedField.H
Foam::fvSchemes
Selector class for finite volume differencing schemes. fvMesh is derived from fvShemes so that all fi...
Definition: fvSchemes.H:52
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:613
Foam::fvMesh::setPhi
surfaceScalarField & setPhi()
Return cell face motion fluxes.
Definition: fvMeshGeometry.C:450
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:326
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:522
primitiveMesh.H
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:955
fvMeshTemplates.C
Foam::fvMesh::V00
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
Definition: fvMeshGeometry.C:247
polyMesh.H
Foam::fvMesh::operator==
bool operator==(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:961
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:306
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:51
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
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:538
Foam::fvMesh::V0
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
Definition: fvMeshGeometry.C:221
Foam::fvMesh::solve
virtual SolverPerformance< scalar > solve(fvMatrix< scalar > &, const dictionary &) const
Definition: fvMesh.C:437
Foam::fvMesh::comm
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:285
Foam::fvMesh::interpolate
virtual void interpolate(scalarField &) const
Interpolate interpolationCells only. No bcs.
Definition: fvMesh.H:314
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:298
Foam::label
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:62
Foam::IOstreamOption::versionNumber
Representation of a major/minor version number.
Definition: IOstreamOption.H:79
Foam::Field< scalar >
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:347
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:595
Foam::fvMesh::interpolate
virtual void interpolate(volTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:310
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:279
Foam::fvSolution
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:49
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:851
fvPatchFvMeshTemplates.C
Foam::fvBoundaryMesh::interfaces
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvBoundaryMesh.C:119
Foam::UPtrList< const lduInterface >
Foam::polyMesh::thisDb
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:498
SolverPerformance.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
fvSchemes.H
Foam::polyMesh::comm
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1259
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:358
fvBoundaryMesh.H
Foam::IOstreamOption::streamFormat
streamFormat
Data format (ascii | binary)
Definition: IOstreamOption.H:64
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:382
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:492
Foam::fvMesh::Mesh
fvMesh Mesh
Definition: fvMesh.H:175
fvSolution.H
Foam::surfaceInterpolation::movePoints
bool movePoints()
Do what is necessary if the mesh has moved.
Definition: surfaceInterpolation.C:125
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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:318
Foam::fvMesh::owner
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:376
Foam::fvMesh::phi
const surfaceScalarField & phi() const
Return cell face motion fluxes.
Definition: fvMeshGeometry.C:428
Foam::fvMesh::hasDb
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: fvMesh.H:252
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
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:330
Foam::List< face >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:52
pointFieldsFwd.H
Forwards and collection of common point field types.
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:279
Foam::UList< label >
surfaceFieldsFwd.H
Foam::surfaceInterpolation
Cell to surface interpolation scheme. Included in fvMesh.
Definition: surfaceInterpolation.H:54
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:160
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:234
Foam::IOstreamOption::compressionType
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
Definition: IOstreamOption.H:71
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
data.H
Foam::fvMesh::interpolate
virtual void interpolate(volSphericalTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:302
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:369
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:54
Foam::fvMesh::delta
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
Definition: fvMeshGeometry.C:380
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:234
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:266
Foam::fvMesh::interpolate
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:294
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:190
Foam::zero
A class representing the concept of 0 (zero), which can be used to avoid manipulating objects that ar...
Definition: zero.H:61
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:336