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 /*---------------------------------------------------------------------------*\
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 protected:
95 
96  // Private data
97 
98  //- Boundary mesh
100 
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 geometry
149  void clearGeom();
150 
151  //- Clear addressing
152  void clearAddressing(const bool isMeshUpdate = false);
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  typedef fvMesh Mesh;
180 
181 
182  // Declare name of the class and its debug switch
183  ClassName("fvMesh");
184 
185 
186  // Constructors
187 
188  //- Construct from IOobject
189  explicit fvMesh(const IOobject& io);
190 
191  //- Construct from IOobject or as zero-sized mesh
192  // Boundary is added using addFvPatches() member function
193  fvMesh(const IOobject& io, const zero, bool syncPar=true);
194 
195  //- Construct from components without boundary.
196  // Boundary is added using addFvPatches() member function
197  fvMesh
198  (
199  const IOobject& io,
200  pointField&& points,
201  faceList&& faces,
202  labelList&& allOwner,
203  labelList&& allNeighbour,
204  const bool syncPar = true
205  );
206 
207  //- Construct without boundary from cells rather than owner/neighbour.
208  // Boundary is added using addFvPatches() member function
209  fvMesh
210  (
211  const IOobject& io,
212  pointField&& points,
213  faceList&& faces,
214  cellList&& cells,
215  const bool syncPar = true
216  );
217 
218 
219  //- Destructor
220  virtual ~fvMesh();
221 
222 
223  // Member Functions
224 
225  // Helpers
226 
227  //- Add boundary patches. Constructor helper
228  void addFvPatches
229  (
230  PtrList<polyPatch>& plist,
231  const bool validBoundary = true
232  );
233 
234  //- Add boundary patches. Constructor helper
235  void addFvPatches
236  (
237  const List<polyPatch*>& p,
238  const bool validBoundary = true
239  );
240 
241  //- Update the mesh based on the mesh files saved in time
242  // directories
243  virtual readUpdateState readUpdate();
244 
245 
246  // Access
247 
248  //- Return the top-level database
249  const Time& time() const
250  {
251  return polyMesh::time();
252  }
253 
254  //- Return true if thisDb() is a valid DB
255  virtual bool hasDb() const
256  {
257  return true;
258  }
259 
260  //- Return the object registry - resolve conflict polyMesh/lduMesh
261  virtual const objectRegistry& thisDb() const
262  {
263  return polyMesh::thisDb();
264  }
265 
266  //- Return reference to name
267  // Note: name() is currently ambiguous due to derivation from
268  // surfaceInterpolation
269  const word& name() const
270  {
271  return polyMesh::name();
272  }
273 
274  //- Return reference to boundary mesh
275  const fvBoundaryMesh& boundary() const;
276 
277  //- Return ldu addressing
278  virtual const lduAddressing& lduAddr() const;
279 
280  //- Return a list of pointers for each patch
281  // with only those pointing to interfaces being set
282  virtual lduInterfacePtrsList interfaces() const
283  {
284  return boundary().interfaces();
285  }
286 
287  //- Return communicator used for parallel communication
288  virtual label comm() const
289  {
290  return polyMesh::comm();
291  }
292 
293 
294  // Overlap
295 
296  //- Interpolate interpolationCells only
297  virtual void interpolate(volScalarField&) const
298  {}
299 
300  //- Interpolate interpolationCells only
301  virtual void interpolate(volVectorField&) const
302  {}
303 
304  //- Interpolate interpolationCells only
305  virtual void interpolate(volSphericalTensorField&) const
306  {}
307 
308  //- Interpolate interpolationCells only
309  virtual void interpolate(volSymmTensorField&) const
310  {}
311 
312  //- Interpolate interpolationCells only
313  virtual void interpolate(volTensorField&) const
314  {}
315 
316  //- Interpolate interpolationCells only. No bcs.
317  virtual void interpolate(scalarField&) const
318  {}
319 
320  //- Interpolate interpolationCells only. No bcs.
321  virtual void interpolate(vectorField&) const
322  {}
323 
324  //- Interpolate interpolationCells only. No bcs.
325  virtual void interpolate(sphericalTensorField&) const
326  {}
327 
328  //- Interpolate interpolationCells only. No bcs.
329  virtual void interpolate(symmTensorField&) const
330  {}
331 
332  //- Interpolate interpolationCells only. No bcs.
333  virtual void interpolate(tensorField&) const
334  {}
335 
336  //- Solve returning the solution statistics given convergence
337  //- tolerance. Use the given solver controls
339  (
341  const dictionary&
342  ) const;
343 
344  //- Solve returning the solution statistics given convergence
345  //- tolerance. Use the given solver controls
347  (
349  const dictionary&
350  ) const;
351 
352  //- Solve returning the solution statistics given convergence
353  //- tolerance. Use the given solver controls
355  (
357  const dictionary&
358  ) const;
359 
360  //- Solve returning the solution statistics given convergence
361  //- tolerance. Use the given solver controls
363  (
365  const dictionary&
366  ) const;
367 
368  //- Solve returning the solution statistics given convergence
369  //- tolerance. Use the given solver controls
371  (
373  const dictionary&
374  ) const;
375 
376 
377  //- Internal face owner. Note bypassing virtual mechanism so
378  // e.g. relaxation always gets done using original addressing
379  const labelUList& owner() const
380  {
381  return fvMesh::lduAddr().lowerAddr();
382  }
383 
384  //- Internal face neighbour
385  const labelUList& neighbour() const
386  {
387  return fvMesh::lduAddr().upperAddr();
388  }
389 
390  //- Return cell volumes
391  const DimensionedField<scalar, volMesh>& V() const;
392 
393  //- Return old-time cell volumes
394  const DimensionedField<scalar, volMesh>& V0() const;
395 
396  //- Return old-old-time cell volumes
397  const DimensionedField<scalar, volMesh>& V00() const;
398 
399  //- Return sub-cycle cell volumes
401 
402  //- Return sub-cycle old-time cell volumes
404 
405  //- Return cell face area vectors
406  const surfaceVectorField& Sf() const;
407 
408  //- Return cell face area magnitudes
409  const surfaceScalarField& magSf() const;
410 
411  //- Return cell face motion fluxes
412  const surfaceScalarField& phi() const;
413 
414  //- Return cell centres as volVectorField
415  const volVectorField& C() const;
416 
417  //- Return face centres as surfaceVectorField
418  const surfaceVectorField& Cf() const;
419 
420  //- Return face deltas as surfaceVectorField
422 
423  //- Return a labelType of valid component indicators
424  // 1 : valid (solved)
425  // -1 : invalid (not solved)
426  template<class Type>
427  typename pTraits<Type>::labelType validComponents() const;
428 
429 
430  // Edit
431 
432  //- Clear all geometry and addressing
433  void clearOut();
434 
435  //- Update mesh corresponding to the given map
436  virtual void updateMesh(const mapPolyMesh& mpm);
437 
438  //- Move points, returns volumes swept by faces in motion
439  virtual tmp<scalarField> movePoints(const pointField&);
440 
441  //- Map all fields in time using given map.
442  virtual void mapFields(const mapPolyMesh& mpm);
443 
444  //- Remove boundary patches. Warning: fvPatchFields hold ref to
445  // these fvPatches.
446  void removeFvBoundary();
447 
448  //- Return cell face motion fluxes
450 
451  //- Return old-time cell volumes
453 
454 
455  // Write
456 
457  //- Write the underlying polyMesh and other data
458  virtual bool writeObject
459  (
460  IOstreamOption streamOpt,
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
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:1038
Foam::fvMesh::~fvMesh
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:412
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:324
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:178
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:895
Foam::fvMesh::makeCf
void makeCf() const
Definition: fvMeshGeometry.C:145
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::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:260
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::fvMesh::VPtr_
void * VPtr_
Cell volumes old time level.
Definition: fvMesh.H:112
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:576
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:328
Foam::fvMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: fvMesh.C:116
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:506
primitiveMesh.H
Foam::fvMesh::operator!=
bool operator!=(const fvMesh &rhs) const
Compares addresses.
Definition: fvMesh.C:911
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:917
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:308
Foam::fvMeshLduAddressing
Foam::fvMeshLduAddressing.
Definition: fvMeshLduAddressing.H:51
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:519
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:421
Foam::fvMesh::comm
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:287
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:316
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:300
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:561
Foam::fvMesh::interpolate
virtual void interpolate(volTensorField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:312
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:49
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:807
fvPatchFvMeshTemplates.C
Foam::fvBoundaryMesh::interfaces
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvBoundaryMesh.C:119
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:498
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:1252
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:384
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:476
Foam::fvMesh::Mesh
fvMesh Mesh
Definition: fvMesh.H:177
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:320
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:378
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:254
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:555
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:1063
Foam::fvMesh::lduPtr_
fvMeshLduAddressing * lduPtr_
Definition: fvMesh.H:103
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:332
Foam::fvMesh::boundary_
fvBoundaryMesh boundary_
Boundary mesh.
Definition: fvMesh.H:98
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:118
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:281
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:220
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:248
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:304
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:54
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:268
Foam::fvMesh::makeMagSf
void makeMagSf() const
Definition: fvMeshGeometry.C:75
Foam::fvMesh::interpolate
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:296
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:872