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