fvMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fvMesh
26 
27 Description
28  Mesh data needed to do the Finite Volume discretisation.
29 
30  NOTE ON USAGE:
31  fvMesh contains all the topological and geometric information
32  related to the mesh. It is also responsible for keeping the data
33  up-to-date. This is done by deleting the cell volume, face area,
34  cell/face centre, addressing and other derived information as
35  required and recalculating it as necessary. The fvMesh therefore
36  reserves the right to delete the derived information upon every
37  topological (mesh refinement/morphing) or geometric change (mesh
38  motion). It is therefore unsafe to keep local references to the
39  derived data outside of the time loop.
40 
41 SourceFiles
42  fvMesh.C
43  fvMeshGeometry.C
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #ifndef fvMesh_H
48 #define fvMesh_H
49 
50 #include "polyMesh.H"
51 #include "lduMesh.H"
52 #include "primitiveMesh.H"
53 #include "fvBoundaryMesh.H"
54 #include "surfaceInterpolation.H"
55 #include "fvSchemes.H"
56 #include "fvSolution.H"
57 #include "DimensionedField.H"
58 #include "SlicedDimensionedField.H"
59 #include "volFieldsFwd.H"
60 #include "surfaceFieldsFwd.H"
61 #include "pointFieldsFwd.H"
62 #include "slicedVolFieldsFwd.H"
63 #include "slicedSurfaceFieldsFwd.H"
64 #include "className.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 
71 class fvMeshLduAddressing;
72 class fvMeshStitcher;
73 class fvMeshTopoChanger;
74 class fvMeshDistributor;
75 class fvMeshMover;
76 class volMesh;
77 class polyDistributionMap;
78 
79 template<class Type>
80 class UCompactListList;
81 
82 namespace fvMeshTopoChangers
83 {
84  class list;
85 }
86 
87 
88 /*---------------------------------------------------------------------------*\
89  Class fvMesh Declaration
90 \*---------------------------------------------------------------------------*/
91 
92 class fvMesh
93 :
94  public polyMesh,
95  public lduMesh,
97 {
98 public:
99 
100  // Public data types
101 
102  //- Enumeration defining the state of the mesh after a read update.
103  // Used for post-processing applications, where the mesh
104  // needs to update based on the files written in time directories.
105  enum readUpdateState
106  {
110  TOPO_CHANGE,
112  };
113 
114  //- Extent to which to stitch on read and readUpdate. By default, a
115  // full geometric stitch is performed. A non-geometric stitch can be
116  // done as an optimisation in situations when finite volume geometry
117  // is not needed (e.g., decomposition). Stitching can also be
118  // prevented altogether if that is appropriate (e.g., if the mesh is
119  // loaded for mapping in an un-stitched state).
120  enum class stitchType
121  {
122  none,
123  nonGeometric,
124  geometric
125  };
126 
127 
128 private:
129 
130  // Private Data
131 
132  //- Boundary mesh
133  fvBoundaryMesh boundary_;
134 
135  //- The stitcher function class
136  autoPtr<fvMeshStitcher> stitcher_;
137 
138  //- The topo-changer function class
139  autoPtr<fvMeshTopoChanger> topoChanger_;
140 
141  //- The distributor function class
142  autoPtr<fvMeshDistributor> distributor_;
143 
144  //- The mover function class
145  autoPtr<fvMeshMover> mover_;
146 
147 
148  // Demand-driven data
149 
150  //- Matrix addressing
151  mutable fvMeshLduAddressing* lduPtr_;
152 
153  //- Face-poly-face IO
154  mutable IOobject* polyFacesBfIOPtr_;
155 
156  //- Face-poly-face addressing
157  mutable GeometricBoundaryField<label, surfaceMesh>* polyFacesBfPtr_;
158 
159  //- Offsets for poly-bFace-patch and patch-face addressing
160  mutable labelList* polyBFaceOffsetsPtr_;
161 
162  //- Poly-bFace-patch addressing
163  mutable labelList* polyBFaceOffsetPatchesPtr_;
164 
165  //- Poly-bFace-patch-face addressing
166  mutable labelList* polyBFaceOffsetPatchFacesPtr_;
167 
168  //- Poly-bFace-patch addressing
169  mutable UCompactListList<label>* polyBFacePatchesPtr_;
170 
171  //- Poly-bFace-patch-face addressing
172  mutable UCompactListList<label>* polyBFacePatchFacesPtr_;
173 
174  //- Face-owner addressing
175  mutable GeometricBoundaryField<label, surfaceMesh>* ownerBfPtr_;
176 
177  //- Current time index for cell volumes
178  // Note. The whole mechanism will be replaced once the
179  // dimensionedField is created and the dimensionedField
180  // will take care of the old-time levels.
181  mutable label curTimeIndex_;
182 
183  //- Cell volumes
185 
186  //- Cell volumes old time level
187  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
188 
189  //- Cell volumes old-old time level
190  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
191 
192  //- Face area vectors
193  mutable slicedSurfaceVectorField* SfSlicePtr_;
194 
195  //- Face area vectors
196  mutable surfaceVectorField* SfPtr_;
197 
198  //- Mag face area vectors
199  mutable slicedSurfaceScalarField* magSfSlicePtr_;
200 
201  //- Mag face area vectors
202  mutable surfaceScalarField* magSfPtr_;
203 
204  //- Cell centres
205  mutable slicedVolVectorField* CSlicePtr_;
206 
207  //- Cell centres
208  mutable volVectorField* CPtr_;
209 
210  //- Face centres
211  mutable slicedSurfaceVectorField* CfSlicePtr_;
212 
213  //- Face centres
214  mutable surfaceVectorField* CfPtr_;
215 
216  //- Face motion fluxes
217  mutable surfaceScalarField* phiPtr_;
218 
219  //- fvSchemes created on demand
220  mutable autoPtr<fvSchemes> fvSchemes_;
221 
222  //- fvSolution created on demand
223  mutable autoPtr<fvSolution> fvSolution_;
224 
225 
226  // Private Member Functions
227 
228  // Storage management
229 
230  //- Clear FV geometry but not the old-time cell volumes
231  void clearFvGeomNotOldVol();
232 
233  //- Clear FV geometry
234  void clearFvGeom();
235 
236  //- Clear geometry like clearGeomNotOldVol but recreate any
237  // geometric demand-driven data that was set
238  void updateGeomNotOldVol();
239 
240  //- Store old times fields for all fields of a given type and class
241  template<class Type, template<class> class GeoField>
242  void storeOldTimeFields();
243 
244  //- Store old times fields for all fields of a given class
245  template<template<class> class GeoField>
246  void storeOldTimeFields();
247 
248  //- Store old times fields for all fields
249  // Used prior to field mapping following mesh-motion
250  // and NCC stitching
251  void storeOldTimeFields();
252 
253  //- Null oldest time field for all fields of a given type and class
254  template<class Type, template<class> class GeoField>
255  void nullOldestTimeFields();
256 
257  //- Null oldest time field for all fields of a given class
258  template<template<class> class GeoField>
259  void nullOldestTimeFields();
260 
261  //- Null oldest time field for all fields
262  // Used prior to field mapping following mesh topology change
263  void nullOldestTimeFields();
264 
265 
266  // Make topological data
267 
268  //- Get the patch types for the poly faces field
269  wordList polyFacesPatchTypes() const;
270 
271  //- Modify face-poly-face addressing
273 
274 
275  // Make geometric data
276 
277  //- Make the face area vectors
278  void makeSf() const;
279 
280  //- Make the mag face area vectors
281  void makeMagSf() const;
282 
283  //- Make the cell centres
284  void makeC() const;
285 
286  //- Make the face centres
287  void makeCf() const;
288 
289  //- Modify cell face area vectors
290  surfaceVectorField& SfRef();
291 
292  //- Modify cell face area magnitudes
293  surfaceScalarField& magSfRef();
294 
295  //- Modify cell centres
296  volVectorField& CRef();
297 
298  //- Modify face centres
299  surfaceVectorField& CfRef();
300 
301  //- Modify cell face motion fluxes
302  surfaceScalarField& phiRef();
303 
304 
305 public:
306 
307  // Public Typedefs
308 
309  typedef fvMesh Mesh;
311 
313 
314 
315  // Declare name of the class and its debug switch
316  ClassName("fvMesh");
317 
318 
319  // Static data
320 
321  //- Set of names of registered geometric fields
322  // Vc, Vc0, Vc00, Sf, magSf, C, Cf, meshPhi, meshPhi_0
323  const static HashSet<word> geometryFields;
324 
325  //- Set of names of registered current-time geometric fields
326  // Vc, Sf, magSf, C, Cf
327  const static HashSet<word> curGeometryFields;
328 
329 
330  // Constructors
331 
332  //- Construct from IOobject.
333  // Optionally prevent post-construction, so that postConstruct may be
334  // called manually for finer control of construction of mesh changes
335  // and the level of stitching.
336  // Optionally prevent zone generation
337  explicit fvMesh
338  (
339  const IOobject& io,
340  const bool doPost = true
341  );
342 
343  //- Construct from cellShapes with boundary.
344  fvMesh
345  (
346  const IOobject& io,
347  pointField&& points,
348  const cellShapeList& shapes,
349  const faceListList& boundaryFaces,
350  const wordList& boundaryPatchNames,
351  const PtrList<dictionary>& boundaryDicts,
352  const word& defaultBoundaryPatchName,
353  const word& defaultBoundaryPatchType,
354  const bool syncPar = true
355  );
356 
357  //- Construct from components without boundary.
358  // Boundary is added using addFvPatches() member function
359  fvMesh
360  (
361  const IOobject& io,
362  pointField&& points,
363  faceList&& faces,
364  labelList&& allOwner,
365  labelList&& allNeighbour,
366  const bool syncPar = true
367  );
368 
369  //- Construct without boundary from cells rather than owner/neighbour.
370  // Boundary is added using addFvPatches() member function
371  fvMesh
372  (
373  const IOobject& io,
374  pointField&& points,
375  faceList&& faces,
376  cellList&& cells,
377  const bool syncPar = true
378  );
379 
380  //- Move constructor
381  fvMesh(fvMesh&&);
382 
383  //- Disallow default bitwise copy construction
384  fvMesh(const fvMesh&) = delete;
385 
386 
387  //- Destructor
388  virtual ~fvMesh();
389 
390 
391  // Member Functions
392 
393  // Helpers
394 
395  //- Complete construction of the mesh
396  void postConstruct
397  (
398  const bool changers,
399  const bool zones,
400  const stitchType stitch
401  );
402 
403  //- Add boundary patches. Constructor helper
404  void addFvPatches
405  (
406  const List<polyPatch*>&,
407  const bool validBoundary = true
408  );
409 
410  //- Update the mesh based on the mesh files saved in time
411  // directories
413  (
414  const stitchType stitch = stitchType::geometric
415  );
416 
417 
418  // Access
419 
420  //- Return the top-level database
421  const Time& time() const
422  {
423  return polyMesh::time();
424  }
425 
426  //- Return the object registry - resolve conflict polyMesh/lduMesh
427  virtual const objectRegistry& thisDb() const
428  {
429  return polyMesh::thisDb();
430  }
431 
432  //- Return reference to name
433  // Note: name() is currently ambiguous due to derivation from
434  // surfaceInterpolation
435  const word& name() const
436  {
437  return polyMesh::name();
438  }
439 
440  //- Return reference to boundary mesh
441  const fvBoundaryMesh& boundary() const;
442 
443  //- Return reference to polyMesh
444  const polyMesh& mesh() const
445  {
446  return *this;
447  }
448 
449  //- Return reference to polyMesh
450  const polyMesh& operator()() const
451  {
452  return *this;
453  }
454 
455  //- Return ldu addressing
456  virtual const lduAddressing& lduAddr() const;
457 
458  //- Return a list of pointers for each patch
459  // with only those pointing to interfaces being set
460  virtual lduInterfacePtrsList interfaces() const
461  {
462  return boundary().interfaces();
463  }
464 
465  //- Return communicator used for parallel communication
466  virtual label comm() const
467  {
468  return polyMesh::comm();
469  }
470 
471  //- Internal face owner
472  const labelUList& owner() const
473  {
474  return lduAddr().lowerAddr();
475  }
476 
477  //- Internal face neighbour
478  const labelUList& neighbour() const
479  {
480  return lduAddr().upperAddr();
481  }
482 
483  //- Return whether the fvMesh is conformal with the polyMesh
484  bool conformal() const;
485 
486  //- Return face-poly-face addressing
488  polyFacesBf() const;
489 
490  //- Return poly-bFace-patch addressing
492 
493  //- Return poly-bFace-patch-face addressing
495 
496  //- Return face-owner addressing
498 
499  //- Return the stitcher function class
500  const fvMeshStitcher& stitcher() const;
501 
502  //- Return the stitcher function class
504 
505  //- Return the topo-changer function class
506  const fvMeshTopoChanger& topoChanger() const;
507 
508  //- Return the distributor function class
509  const fvMeshDistributor& distributor() const;
510 
511  //- Return the mover function class
512  const fvMeshMover& mover() const;
513 
514  //- Return cell volumes
515  const DimensionedField<scalar, volMesh>& V() const;
516 
517  //- Return old-time cell volumes
518  const DimensionedField<scalar, volMesh>& V0() const;
519 
520  //- Return old-old-time cell volumes
521  const DimensionedField<scalar, volMesh>& V00() const;
522 
523  //- Return sub-cycle cell volumes
525 
526  //- Return sub-cycle old-time cell volumes
528 
529  //- Return cell face area vectors
530  const surfaceVectorField& Sf() const;
531 
532  //- Return cell face area magnitudes
533  const surfaceScalarField& magSf() const;
534 
535  //- Return cell face normal vectors
536  tmp<surfaceVectorField> nf() const;
537 
538  //- Return cell centres
539  const volVectorField& C() const;
540 
541  //- Return face centres
542  const surfaceVectorField& Cf() const;
543 
544  //- Return face deltas as surfaceVectorField
546 
547  //- Return cell face motion fluxes
548  const surfaceScalarField& phi() const;
549 
550  //- Return the list of fields of type GeoField
551  // i.e. fields that are not geometry
552  template<class GeoField>
554  (
555  bool strict = false,
557  ) const;
558 
559  //- Return the list of current fields of type GeoField
560  // i.e. fields that are not old-times or geometry
561  template<class GeoField>
563  (
564  const bool strict = false,
566  ) const;
567 
568  //- Return the list of fields and old-time-geometry of type GeoField
569  // i.e. fields that are not current-time geometry
570  template<class GeoField>
572  (
573  const bool strict = false
574  ) const;
575 
576  //- Return a labelType of valid component indicators
577  // 1 : valid (solved)
578  // -1 : invalid (not solved)
579  template<class Type>
580  typename pTraits<Type>::labelType validComponents() const;
581 
582  //- Return the fvSchemes
583  const fvSchemes& schemes() const;
584 
585  //- Return the fvSolution
586  const fvSolution& solution() const;
587 
588 
589  // Edit
590 
591  //- Does the mesh topology change?
592  bool topoChanging() const;
593 
594  //- Does the mesh get redistributed?
595  bool distributing() const;
596 
597  //- Is this mesh dynamic?
598  bool dynamic() const;
599 
600  //- Update the mesh for topology change, mesh to mesh mapping
601  // or redistribution
602  bool update();
603 
604  //- Move the mesh
605  bool move();
606 
607  //- Print a list of all the currently allocated mesh data
608  void printAllocated() const;
609 
610  //- Clear geometry
611  void clearGeom();
612 
613  //- Clear addressing
614  void clearAddressing(const bool isMeshUpdate = false);
615 
616  //- Clear all geometry and addressing
617  void clearOut();
618 
619  //- Prepare for a mesh change
620  void preChange();
621 
622  //- Update mesh corresponding to the given map
623  virtual void topoChange(const polyTopoChangeMap& map);
624 
625  //- Update from another mesh using the given map
626  virtual void mapMesh(const polyMeshMap&);
627 
628  //- Redistribute or update using the given distribution map
629  virtual void distribute(const polyDistributionMap& map);
630 
631  //- Reset the points
632  // without storing old points or returning swept volumes
633  virtual void setPoints(const pointField&);
634 
635  //- Move points, returns volumes swept by faces in motion
636  virtual tmp<scalarField> movePoints(const pointField&);
637 
638  //- Get the current instance for the poly faces boundary field
639  const fileName& polyFacesBfInstance() const;
640 
641  //- Set the instance for the poly faces boundary field
642  void setPolyFacesBfInstance(const fileName&);
643 
644  //- Conform the fvMesh to the polyMesh
645  void conform
646  (
647  const surfaceScalarField& phi =
648  NullObjectRef<surfaceScalarField>()
649  );
650 
651  //- Unconform the fvMesh from the polyMesh
652  void unconform
653  (
655  const surfaceVectorField& Sf,
656  const surfaceVectorField& Cf,
657  const surfaceScalarField& phi =
658  NullObjectRef<surfaceScalarField>(),
659  const bool sync = true
660  );
661 
662  //- Map all fields in time using given map.
663  void mapFields(const polyTopoChangeMap& map);
664 
665  //- Add/insert single patch
666  virtual void addPatch
667  (
668  const label insertPatchi,
669  const polyPatch& patch
670  );
671 
672  //- Reorder and trim existing patches. If validBoundary the new
673  // situation is consistent across processors
674  virtual void reorderPatches
675  (
676  const labelUList& newToOld,
677  const bool validBoundary
678  );
679 
680  //- Remove boundary patches. Warning: fvPatchFields hold ref to
681  // these fvPatches.
682  void removeFvBoundary();
683 
684  //- Swap mesh
685  // For run-time mesh replacement and mesh to mesh mapping
686  void swap(fvMesh&);
687 
688 
689  // Write
690 
691  //- Write the underlying polyMesh and other data
692  virtual bool writeObject
693  (
697  const bool write = true
698  ) const;
699 
700  //- Write mesh using IO settings from time
701  virtual bool write(const bool write = true) const;
702 
703 
704  // Member Operators
705 
706  //- Disallow default bitwise assignment
707  void operator=(const fvMesh&) = delete;
708 
709  bool operator!=(const fvMesh&) const;
710 
711  bool operator==(const fvMesh&) const;
712 };
713 
714 
715 template<>
717 fvMesh::validComponents<sphericalTensor>() const;
718 
719 
720 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
721 
722 } // End namespace Foam
723 
724 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
725 
726 #ifdef NoRepository
727  #include "fvMeshTemplates.C"
728 #endif
729 
730 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
731 
732 #endif
733 
734 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Specialisation of DimensionedField which holds a slice of a given complete field in such a form that ...
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A class for handling file names.
Definition: fileName.H:82
Foam::fvBoundaryMesh.
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Abstract base class for fvMesh movers.
Foam::fvMeshLduAddressing.
Abstract base class for fvMesh movers.
Definition: fvMeshMover.H:53
Mesh manipulator that uses the intersection provided by the cyclic non-conformal poly patches to crea...
Abstract base class for fvMesh topology changers.
Sequence of mesh topology changes applied in order.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const volVectorField & C() const
Return cell centres.
static const HashSet< word > curGeometryFields
Set of names of registered current-time geometric fields.
Definition: fvMesh.H:326
const fileName & polyFacesBfInstance() const
Get the current instance for the poly faces boundary field.
Definition: fvMesh.C:1507
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:471
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
void operator=(const fvMesh &)=delete
Disallow default bitwise assignment.
stitchType
Extent to which to stitch on read and readUpdate. By default, a.
Definition: fvMesh.H:120
void postConstruct(const bool changers, const bool zones, const stitchType stitch)
Complete construction of the mesh.
Definition: fvMesh.C:634
const surfaceScalarField & phi() const
Return cell face motion fluxes.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:426
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:968
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:1016
const fvMeshDistributor & distributor() const
Return the distributor function class.
Definition: fvMesh.C:1149
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
void clearGeom()
Clear geometry.
Definition: fvMesh.C:252
fvMesh Mesh
Definition: fvMesh.H:308
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:105
@ POINTS_MOVED
Definition: fvMesh.H:107
@ TOPO_CHANGE
Definition: fvMesh.H:109
@ TOPO_PATCH_CHANGE
Definition: fvMesh.H:110
UPtrList< GeoField > fieldsAndOldTimeGeometry(const bool strict=false) const
Return the list of fields and old-time-geometry of type GeoField.
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:625
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
UPtrList< GeoField > curFields(const bool strict=false, const HashSet< word > &geometryFields=fvMesh::geometryFields) const
Return the list of current fields of type GeoField.
const GeometricBoundaryField< label, surfaceMesh > & ownerBf() const
Return face-owner addressing.
Definition: fvMesh.C:1103
static const HashSet< word > geometryFields
Set of names of registered geometric fields.
Definition: fvMesh.H:322
const fvMeshStitcher & stitcher() const
Return the stitcher function class.
Definition: fvMesh.C:1131
void preChange()
Prepare for a mesh change.
Definition: fvMesh.C:1222
tmp< surfaceVectorField > nf() const
Return cell face normal vectors.
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:785
const surfaceVectorField & Cf() const
Return face centres.
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:819
const GeometricBoundaryField< label, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:985
bool operator==(const fvMesh &) const
Definition: fvMesh.C:1829
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:465
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1092
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:1823
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1724
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1429
void mapFields(const polyTopoChangeMap &map)
Map all fields in time using given map.
Definition: fvMesh.C:1161
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
UPtrList< GeoField > fields(bool strict=false, const HashSet< word > &geometryFields=fvMesh::geometryFields) const
Return the list of fields of type GeoField.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:459
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool update()
Update the mesh for topology change, mesh to mesh mapping.
Definition: fvMesh.C:725
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:477
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
void setPolyFacesBfInstance(const fileName &)
Set the instance for the poly faces boundary field.
Definition: fvMesh.C:1495
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: fvMesh.C:1462
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:803
bool distributing() const
Does the mesh get redistributed?
Definition: fvMesh.C:711
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
void printAllocated() const
Print a list of all the currently allocated mesh data.
Definition: fvMesh.C:186
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:309
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1702
ClassName("fvMesh")
void unconform(const GeometricBoundaryField< label, surfaceMesh > &polyFacesBf, const surfaceVectorField &Sf, const surfaceVectorField &Cf, const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >(), const bool sync=true)
Unconform the fvMesh from the polyMesh.
Definition: fvMesh.C:1541
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:717
fvMesh(const IOobject &io, const bool doPost=true)
Construct from IOobject.
Definition: fvMesh.C:364
const fvMeshTopoChanger & topoChanger() const
Return the topo-changer function class.
Definition: fvMesh.C:1143
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:705
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
virtual void addPatch(const label insertPatchi, const polyPatch &patch)
Add/insert single patch.
Definition: fvMesh.C:1645
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:313
bool move()
Move the mesh.
Definition: fvMesh.C:765
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1228
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1155
const polyMesh & operator()() const
Return reference to polyMesh.
Definition: fvMesh.H:449
void conform(const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >())
Conform the fvMesh to the polyMesh.
Definition: fvMesh.C:1518
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:979
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:53
Selector class for finite volume solution solution. fvMesh is derived from fvSolution so that all fie...
Definition: fvSolution.H:50
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:60
Registry of regIOobjects.
const Time & time() const
Return time.
Traits class for primitives.
Definition: pTraits.H:53
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:464
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1539
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:124
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void clearAddressing()
Clear topological data.
const cellList & cells() const
Cell to surface interpolation scheme. Included in fvMesh.
bool movePoints()
Do what is necessary if the mesh has moved.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Namespace for OpenFOAM.
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:59