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-2023 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 "data.H"
58 #include "DimensionedField.H"
59 #include "SlicedDimensionedField.H"
60 #include "volFieldsFwd.H"
61 #include "surfaceFieldsFwd.H"
62 #include "pointFieldsFwd.H"
63 #include "slicedVolFieldsFwd.H"
64 #include "slicedSurfaceFieldsFwd.H"
65 #include "className.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 
72 class fvMeshLduAddressing;
73 class fvMeshStitcher;
74 class fvMeshTopoChanger;
75 class fvMeshDistributor;
76 class fvMeshMover;
77 class volMesh;
78 class polyDistributionMap;
79 
80 template<class Type, template<class> class PatchField, class GeoMesh>
81 class GeometricBoundaryField;
82 
83 template<class Type>
84 class UCompactListList;
85 
86 namespace fvMeshTopoChangers
87 {
88  class list;
89 }
90 
91 
92 /*---------------------------------------------------------------------------*\
93  Class fvMesh Declaration
94 \*---------------------------------------------------------------------------*/
95 
96 class fvMesh
97 :
98  public polyMesh,
99  public lduMesh,
100  public surfaceInterpolation,
101  public data
102 {
103 public:
104 
105  // Public data types
106 
107  //- Extent to which to stitch on read and readUpdate. By default, a
108  // full geometric stitch is performed. A non-geometric stitch can be
109  // done as an optimisation in situations when finite volume geometry
110  // is not needed (e.g., decomposition). Stitching can also be
111  // prevented altogether if that is appropriate (e.g., if the mesh is
112  // loaded for mapping in an un-stitched state).
113  enum class stitchType
114  {
115  none,
116  nonGeometric,
117  geometric
118  };
119 
120 
121 private:
122 
123  // Private Data
124 
125  //- Boundary mesh
126  fvBoundaryMesh boundary_;
127 
128  //- The stitcher function class
129  autoPtr<fvMeshStitcher> stitcher_;
130 
131  //- The topo-changer function class
132  autoPtr<fvMeshTopoChanger> topoChanger_;
133 
134  //- The distributor function class
135  autoPtr<fvMeshDistributor> distributor_;
136 
137  //- The mover function class
138  autoPtr<fvMeshMover> mover_;
139 
140 
141  // Demand-driven data
142 
143  //- Matrix addressing
144  mutable fvMeshLduAddressing* lduPtr_;
145 
146  //- Face-poly-face addressing
148  polyFacesBfPtr_;
149 
150  //- Offsets for poly-bFace-patch and patch-face addressing
151  mutable labelList* polyBFaceOffsetsPtr_;
152 
153  //- Poly-bFace-patch addressing
154  mutable labelList* polyBFaceOffsetPatchesPtr_;
155 
156  //- Poly-bFace-patch-face addressing
157  mutable labelList* polyBFaceOffsetPatchFacesPtr_;
158 
159  //- Poly-bFace-patch addressing
160  mutable UCompactListList<label>* polyBFacePatchesPtr_;
161 
162  //- Poly-bFace-patch-face addressing
163  mutable UCompactListList<label>* polyBFacePatchFacesPtr_;
164 
165  //- Current time index for cell volumes
166  // Note. The whole mechanism will be replaced once the
167  // dimensionedField is created and the dimensionedField
168  // will take care of the old-time levels.
169  mutable label curTimeIndex_;
170 
171  //- Cell volumes
173 
174  //- Cell volumes old time level
175  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
176 
177  //- Cell volumes old-old time level
178  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
179 
180  //- Face area vectors
181  mutable slicedSurfaceVectorField* SfSlicePtr_;
182 
183  //- Face area vectors
184  mutable surfaceVectorField* SfPtr_;
185 
186  //- Mag face area vectors
187  mutable slicedSurfaceScalarField* magSfSlicePtr_;
188 
189  //- Mag face area vectors
190  mutable surfaceScalarField* magSfPtr_;
191 
192  //- Cell centres
193  mutable slicedVolVectorField* CSlicePtr_;
194 
195  //- Cell centres
196  mutable volVectorField* CPtr_;
197 
198  //- Face centres
199  mutable slicedSurfaceVectorField* CfSlicePtr_;
200 
201  //- Face centres
202  mutable surfaceVectorField* CfPtr_;
203 
204  //- Face motion fluxes
205  mutable surfaceScalarField* phiPtr_;
206 
207  //- fvSchemes created on demand
208  mutable autoPtr<fvSchemes> fvSchemes_;
209 
210  //- fvSolution created on demand
211  mutable autoPtr<fvSolution> fvSolution_;
212 
213 
214  // Private Member Functions
215 
216  // Storage management
217 
218  //- Clear geometry but not the old-time cell volumes
219  void clearGeomNotOldVol();
220 
221  //- Clear geometry like clearGeomNotOldVol but recreate any
222  // geometric demand-driven data that was set
223  void updateGeomNotOldVol();
224 
225  //- Clear geometry
226  void clearGeom();
227 
228  //- Clear addressing
229  void clearAddressing(const bool isMeshUpdate = false);
230 
231  //- Preserve old volume(s)
232  void storeOldVol(const scalarField&);
233 
234  //- Return the list of current fields of type GeoField<Type>
235  // i.e. fields that are not old-time or old-old-time etc.
236  template<class Type, template<class> class GeoField>
237  UPtrList<GeoField<Type>> curFields();
238 
239  //- Store old times fields for all fields of a given type and class
240  template<class Type, template<class> class GeoField>
241  void storeOldTimeFields();
242 
243  //- Store old times fields for all fields of a given class
244  template<template<class> class GeoField>
245  void storeOldTimeFields();
246 
247  //- Store old times fields for all fields
248  // Used prior to field mapping following mesh-motion
249  // and NCC stitching
250  void storeOldTimeFields();
251 
252  //- Null oldest time field for all fields of a given type and class
253  template<class Type, template<class> class GeoField>
254  void nullOldestTimeFields();
255 
256  //- Null oldest time field for all fields of a given class
257  template<template<class> class GeoField>
258  void nullOldestTimeFields();
259 
260  //- Null oldest time field for all fields
261  // Used prior to field mapping following mesh topology change
262  void nullOldestTimeFields();
263 
264 
265  // Make topological data
266 
267  //- Get the patch types for the poly faces field
268  wordList polyFacesPatchTypes() const;
269 
270  //- Modify face-poly-face addressing
272  polyFacesBfRef();
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 old-time cell volumes
291 
292  //- Modify old-old-time cell volumes
294 
295  //- Modify cell face area vectors
296  surfaceVectorField& SfRef();
297 
298  //- Modify cell face area magnitudes
299  surfaceScalarField& magSfRef();
300 
301  //- Modify cell centres
302  volVectorField& CRef();
303 
304  //- Modify face centres
305  surfaceVectorField& CfRef();
306 
307  //- Modify cell face motion fluxes
308  surfaceScalarField& phiRef();
309 
310 
311 public:
312 
313  // Public Typedefs
314 
315  typedef fvMesh Mesh;
317 
319  friend fvMeshDistributor;
320 
321  // Declare name of the class and its debug switch
322  ClassName("fvMesh");
323 
324 
325  // Constructors
326 
327  //- Construct from IOobject
328  // with the option to not instantiate the mesh changers and the option
329  // to prevent some or all of the stitching
330  explicit fvMesh
331  (
332  const IOobject& io,
333  const bool changers = true,
334  const stitchType stitch = stitchType::geometric
335  );
336 
337  //- Construct from cellShapes with boundary.
338  fvMesh
339  (
340  const IOobject& io,
341  pointField&& points,
342  const cellShapeList& shapes,
343  const faceListList& boundaryFaces,
344  const wordList& boundaryPatchNames,
345  const PtrList<dictionary>& boundaryDicts,
346  const word& defaultBoundaryPatchName,
347  const word& defaultBoundaryPatchType,
348  const bool syncPar = true
349  );
350 
351  //- Construct from components without boundary.
352  // Boundary is added using addFvPatches() member function
353  fvMesh
354  (
355  const IOobject& io,
356  pointField&& points,
357  faceList&& faces,
358  labelList&& allOwner,
359  labelList&& allNeighbour,
360  const bool syncPar = true
361  );
362 
363  //- Construct without boundary from cells rather than owner/neighbour.
364  // Boundary is added using addPatches() member function
365  fvMesh
366  (
367  const IOobject& io,
368  pointField&& points,
369  faceList&& faces,
370  cellList&& cells,
371  const bool syncPar = true
372  );
373 
374  //- Move constructor
375  fvMesh(fvMesh&&);
376 
377  //- Disallow default bitwise copy construction
378  fvMesh(const fvMesh&) = delete;
379 
380 
381  //- Destructor
382  virtual ~fvMesh();
383 
384 
385  // Member Functions
386 
387  // Helpers
388 
389  //- Add boundary patches. Constructor helper
390  void addFvPatches
391  (
392  const List<polyPatch*>&,
393  const bool validBoundary = true
394  );
395 
396  //- Update the mesh based on the mesh files saved in time
397  // directories
399  (
400  const stitchType stitch = stitchType::geometric
401  );
402 
403 
404  // Access
405 
406  //- Return the top-level database
407  const Time& time() const
408  {
409  return polyMesh::time();
410  }
411 
412  //- Return the object registry - resolve conflict polyMesh/lduMesh
413  virtual const objectRegistry& thisDb() const
414  {
415  return polyMesh::thisDb();
416  }
417 
418  //- Return reference to name
419  // Note: name() is currently ambiguous due to derivation from
420  // surfaceInterpolation
421  const word& name() const
422  {
423  return polyMesh::name();
424  }
425 
426  //- Return reference to boundary mesh
427  const fvBoundaryMesh& boundary() const;
428 
429  //- Return reference to polyMesh
430  const polyMesh& operator()() const
431  {
432  return *this;
433  }
434 
435  //- Return ldu addressing
436  virtual const lduAddressing& lduAddr() const;
437 
438  //- Return a list of pointers for each patch
439  // with only those pointing to interfaces being set
440  virtual lduInterfacePtrsList interfaces() const
441  {
442  return boundary().interfaces();
443  }
444 
445  //- Return communicator used for parallel communication
446  virtual label comm() const
447  {
448  return polyMesh::comm();
449  }
450 
451  //- Internal face owner
452  const labelUList& owner() const
453  {
454  return lduAddr().lowerAddr();
455  }
456 
457  //- Internal face neighbour
458  const labelUList& neighbour() const
459  {
460  return lduAddr().upperAddr();
461  }
462 
463  //- Return whether the fvMesh is conformal with the polyMesh
464  bool conformal() const;
465 
466  //- Get the poly faces IO object
468 
469  //- Return face-poly-face addressing
471  polyFacesBf() const;
472 
473  //- Return poly-bFace-patch addressing
475 
476  //- Return poly-bFace-patch-face addressing
478 
479  //- Return the stitcher function class
480  const fvMeshStitcher& stitcher() const;
481 
482  //- Return the topo-changer function class
483  const fvMeshTopoChanger& topoChanger() const;
484 
485  //- Return the distributor function class
486  const fvMeshDistributor& distributor() const;
487 
488  //- Return the mover function class
489  const fvMeshMover& mover() const;
490 
491  //- Return cell volumes
492  const DimensionedField<scalar, volMesh>& V() const;
493 
494  //- Return old-time cell volumes
495  const DimensionedField<scalar, volMesh>& V0() const;
496 
497  //- Return old-old-time cell volumes
498  const DimensionedField<scalar, volMesh>& V00() const;
499 
500  //- Return sub-cycle cell volumes
502 
503  //- Return sub-cycle old-time cell volumes
505 
506  //- Return cell face area vectors
507  const surfaceVectorField& Sf() const;
508 
509  //- Return cell face area magnitudes
510  const surfaceScalarField& magSf() const;
511 
512  //- Return cell centres
513  const volVectorField& C() const;
514 
515  //- Return face centres
516  const surfaceVectorField& Cf() const;
517 
518  //- Return face deltas as surfaceVectorField
520 
521  //- Return cell face motion fluxes
522  const surfaceScalarField& phi() const;
523 
524  //- Return a labelType of valid component indicators
525  // 1 : valid (solved)
526  // -1 : invalid (not solved)
527  template<class Type>
528  typename pTraits<Type>::labelType validComponents() const;
529 
530  //- Return the fvSchemes
531  const fvSchemes& schemes() const;
532 
533  //- Return the fvSchemes
534  const fvSolution& solution() const;
535 
536 
537  // Edit
538 
539  //- Does the mesh topology change?
540  bool topoChanging() const;
541 
542  //- Is this mesh dynamic?
543  bool dynamic() const;
544 
545  //- Update the mesh for topology change, mesh to mesh mapping
546  // or redistribution
547  bool update();
548 
549  //- Move the mesh
550  bool move();
551 
552  //- Clear all geometry and addressing
553  void clearOut();
554 
555  //- Update mesh corresponding to the given map
556  virtual void topoChange(const polyTopoChangeMap& map);
557 
558  //- Update from another mesh using the given map
559  virtual void mapMesh(const polyMeshMap&);
560 
561  //- Redistribute or update using the given distribution map
562  virtual void distribute(const polyDistributionMap& map);
563 
564  //- Reset the points
565  // without storing old points or returning swept volumes
566  virtual void setPoints(const pointField&);
567 
568  //- Move points, returns volumes swept by faces in motion
569  virtual tmp<scalarField> movePoints(const pointField&);
570 
571  //- Conform the fvMesh to the polyMesh
572  void conform
573  (
574  const surfaceScalarField& phi =
575  NullObjectRef<surfaceScalarField>()
576  );
577 
578  //- Unconform the fvMesh from the polyMesh
579  void unconform
580  (
582  polyFacesBf,
583  const surfaceVectorField& Sf,
584  const surfaceVectorField& Cf,
585  const surfaceScalarField& phi =
586  NullObjectRef<surfaceScalarField>(),
587  const bool sync = true
588  );
589 
590  //- Map all fields in time using given map.
591  void mapFields(const polyTopoChangeMap& map);
592 
593  //- Add/insert single patch. If validBoundary the new situation
594  // is consistent across processors.
595  virtual void addPatch
596  (
597  const label insertPatchi,
598  const polyPatch& patch,
599  const dictionary& patchFieldDict,
600  const word& defaultPatchFieldType,
601  const bool validBoundary
602  );
603 
604  //- Reorder and trim existing patches. If validBoundary the new
605  // situation is consistent across processors
606  virtual void reorderPatches
607  (
608  const labelUList& newToOld,
609  const bool validBoundary
610  );
611 
612  //- Remove boundary patches. Warning: fvPatchFields hold ref to
613  // these fvPatches.
614  void removeFvBoundary();
615 
616  //- Swap mesh
617  // For run-time mesh replacement and mesh to mesh mapping
618  void swap(fvMesh&);
619 
620 
621  // Write
622 
623  //- Write the underlying polyMesh and other data
624  virtual bool writeObject
625  (
629  const bool write = true
630  ) const;
631 
632  //- Write mesh using IO settings from time
633  virtual bool write(const bool write = true) const;
634 
635 
636  // Member Operators
637 
638  //- Disallow default bitwise assignment
639  void operator=(const fvMesh&);
640 
641  bool operator!=(const fvMesh&) const;
642  bool operator==(const fvMesh&) const;
643 };
644 
645 
646 template<>
648 fvMesh::validComponents<sphericalTensor>() const;
649 
650 
651 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
652 
653 } // End namespace Foam
654 
655 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
656 
657 #ifdef NoRepository
658  #include "fvMeshTemplates.C"
659 #endif
660 
661 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
662 
663 #endif
664 
665 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
const word & name() const
Return name.
Definition: IOobject.H:310
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
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Database for solution and other reduced data.
Definition: data.H:54
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
const volVectorField & C() const
Return cell centres.
friend fvMeshDistributor
Definition: fvMesh.H:318
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:451
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const DimensionedField< scalar, volMesh > & V00() const
Return old-old-time cell volumes.
stitchType
Extent to which to stitch on read and readUpdate. By default, a.
Definition: fvMesh.H:113
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:412
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:899
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: fvMesh.C:1521
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:951
const fvMeshDistributor & distributor() const
Return the distributor function class.
Definition: fvMesh.C:1050
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
fvMesh Mesh
Definition: fvMesh.H:314
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:627
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
void operator=(const fvMesh &)
Disallow default bitwise assignment.
void unconform(const GeometricBoundaryField< label, fvsPatchField, 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:1420
const fvMeshStitcher & stitcher() const
Return the stitcher function class.
Definition: fvMesh.C:1038
const fvSolution & solution() const
Return the fvSchemes.
Definition: fvMesh.C:1759
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:725
const surfaceVectorField & Cf() const
Return face centres.
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:759
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:932
bool operator==(const fvMesh &) const
Definition: fvMesh.C:1778
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:445
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1027
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:1772
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:1678
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1347
void mapFields(const polyTopoChangeMap &map)
Map all fields in time using given map.
Definition: fvMesh.C:1062
const DimensionedField< scalar, volMesh > & V0() const
Return old-time cell volumes.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1734
IOobject polyFacesBfIO(IOobject::readOption r) const
Get the poly faces IO object.
Definition: fvMesh.C:916
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:439
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool update()
Update the mesh for topology change, mesh to mesh mapping.
Definition: fvMesh.C:647
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:457
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: fvMesh.C:1372
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:743
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:315
const word & name() const
Return reference to name.
Definition: fvMesh.H:420
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1653
fvMesh(const IOobject &io, const bool changers=true, const stitchType stitch=stitchType::geometric)
Construct from IOobject.
Definition: fvMesh.C:323
ClassName("fvMesh")
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:641
const fvMeshTopoChanger & topoChanger() const
Return the topo-changer function class.
Definition: fvMesh.C:1044
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1234
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:635
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:274
bool move()
Move the mesh.
Definition: fvMesh.C:698
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1131
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1056
const polyMesh & operator()() const
Return reference to polyMesh.
Definition: fvMesh.H:429
void conform(const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >())
Conform the fvMesh to the polyMesh.
Definition: fvMesh.C:1397
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:910
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:1374
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:466
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1581
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:99
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
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.
#define nullOldestTimeFields(Type, nullArg)
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