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  //- Store old times fields for all fields of a given type and class
235  template<class Type, template<class> class GeoField>
236  void storeOldTimeFields();
237 
238  //- Store old times fields for all fields of a given class
239  template<template<class> class GeoField>
240  void storeOldTimeFields();
241 
242  //- Store old times fields for all fields
243  // Used prior to field mapping following mesh-motion
244  // and NCC stitching
245  void storeOldTimeFields();
246 
247  //- Null oldest time field for all fields of a given type and class
248  template<class Type, template<class> class GeoField>
249  void nullOldestTimeFields();
250 
251  //- Null oldest time field for all fields of a given class
252  template<template<class> class GeoField>
253  void nullOldestTimeFields();
254 
255  //- Null oldest time field for all fields
256  // Used prior to field mapping following mesh topology change
257  void nullOldestTimeFields();
258 
259 
260  // Make topological data
261 
262  //- Get the patch types for the poly faces field
263  wordList polyFacesPatchTypes() const;
264 
265  //- Modify face-poly-face addressing
267  polyFacesBfRef();
268 
269 
270  // Make geometric data
271 
272  //- Make the face area vectors
273  void makeSf() const;
274 
275  //- Make the mag face area vectors
276  void makeMagSf() const;
277 
278  //- Make the cell centres
279  void makeC() const;
280 
281  //- Make the face centres
282  void makeCf() const;
283 
284  //- Modify cell face area vectors
285  surfaceVectorField& SfRef();
286 
287  //- Modify cell face area magnitudes
288  surfaceScalarField& magSfRef();
289 
290  //- Modify cell centres
291  volVectorField& CRef();
292 
293  //- Modify face centres
294  surfaceVectorField& CfRef();
295 
296  //- Modify cell face motion fluxes
297  surfaceScalarField& phiRef();
298 
299 
300 public:
301 
302  // Public Typedefs
303 
304  typedef fvMesh Mesh;
306 
308 
309 
310  // Declare name of the class and its debug switch
311  ClassName("fvMesh");
312 
313 
314  // Static data
315 
316  //- Set of names of registered geometric fields
317  // V, Sf, magSf, C, Cf
318  const static HashSet<word> geometryFields;
319 
320 
321  // Constructors
322 
323  //- Construct from IOobject
324  // with the option to not instantiate the mesh changers and the option
325  // to prevent some or all of the stitching
326  explicit fvMesh
327  (
328  const IOobject& io,
329  const bool changers = true,
330  const stitchType stitch = stitchType::geometric
331  );
332 
333  //- Construct from cellShapes with boundary.
334  fvMesh
335  (
336  const IOobject& io,
337  pointField&& points,
338  const cellShapeList& shapes,
339  const faceListList& boundaryFaces,
340  const wordList& boundaryPatchNames,
341  const PtrList<dictionary>& boundaryDicts,
342  const word& defaultBoundaryPatchName,
343  const word& defaultBoundaryPatchType,
344  const bool syncPar = true
345  );
346 
347  //- Construct from components without boundary.
348  // Boundary is added using addFvPatches() member function
349  fvMesh
350  (
351  const IOobject& io,
352  pointField&& points,
353  faceList&& faces,
354  labelList&& allOwner,
355  labelList&& allNeighbour,
356  const bool syncPar = true
357  );
358 
359  //- Construct without boundary from cells rather than owner/neighbour.
360  // Boundary is added using addPatches() member function
361  fvMesh
362  (
363  const IOobject& io,
364  pointField&& points,
365  faceList&& faces,
366  cellList&& cells,
367  const bool syncPar = true
368  );
369 
370  //- Move constructor
371  fvMesh(fvMesh&&);
372 
373  //- Disallow default bitwise copy construction
374  fvMesh(const fvMesh&) = delete;
375 
376 
377  //- Destructor
378  virtual ~fvMesh();
379 
380 
381  // Member Functions
382 
383  // Helpers
384 
385  //- Add boundary patches. Constructor helper
386  void addFvPatches
387  (
388  const List<polyPatch*>&,
389  const bool validBoundary = true
390  );
391 
392  //- Update the mesh based on the mesh files saved in time
393  // directories
395  (
396  const stitchType stitch = stitchType::geometric
397  );
398 
399 
400  // Access
401 
402  //- Return the top-level database
403  const Time& time() const
404  {
405  return polyMesh::time();
406  }
407 
408  //- Return the object registry - resolve conflict polyMesh/lduMesh
409  virtual const objectRegistry& thisDb() const
410  {
411  return polyMesh::thisDb();
412  }
413 
414  //- Return reference to name
415  // Note: name() is currently ambiguous due to derivation from
416  // surfaceInterpolation
417  const word& name() const
418  {
419  return polyMesh::name();
420  }
421 
422  //- Return reference to boundary mesh
423  const fvBoundaryMesh& boundary() const;
424 
425  //- Return reference to polyMesh
426  const polyMesh& operator()() const
427  {
428  return *this;
429  }
430 
431  //- Return ldu addressing
432  virtual const lduAddressing& lduAddr() const;
433 
434  //- Return a list of pointers for each patch
435  // with only those pointing to interfaces being set
436  virtual lduInterfacePtrsList interfaces() const
437  {
438  return boundary().interfaces();
439  }
440 
441  //- Return communicator used for parallel communication
442  virtual label comm() const
443  {
444  return polyMesh::comm();
445  }
446 
447  //- Internal face owner
448  const labelUList& owner() const
449  {
450  return lduAddr().lowerAddr();
451  }
452 
453  //- Internal face neighbour
454  const labelUList& neighbour() const
455  {
456  return lduAddr().upperAddr();
457  }
458 
459  //- Return whether the fvMesh is conformal with the polyMesh
460  bool conformal() const;
461 
462  //- Get the poly faces IO object
464 
465  //- Return face-poly-face addressing
467  polyFacesBf() const;
468 
469  //- Return poly-bFace-patch addressing
471 
472  //- Return poly-bFace-patch-face addressing
474 
475  //- Return the stitcher function class
476  const fvMeshStitcher& stitcher() const;
477 
478  //- Return the topo-changer function class
479  const fvMeshTopoChanger& topoChanger() const;
480 
481  //- Return the distributor function class
482  const fvMeshDistributor& distributor() const;
483 
484  //- Return the mover function class
485  const fvMeshMover& mover() const;
486 
487  //- Return cell volumes
488  const DimensionedField<scalar, volMesh>& V() const;
489 
490  //- Return old-time cell volumes
491  const DimensionedField<scalar, volMesh>& V0() const;
492 
493  //- Return old-old-time cell volumes
494  const DimensionedField<scalar, volMesh>& V00() const;
495 
496  //- Return sub-cycle cell volumes
498 
499  //- Return sub-cycle old-time cell volumes
501 
502  //- Return cell face area vectors
503  const surfaceVectorField& Sf() const;
504 
505  //- Return cell face area magnitudes
506  const surfaceScalarField& magSf() const;
507 
508  //- Return cell centres
509  const volVectorField& C() const;
510 
511  //- Return face centres
512  const surfaceVectorField& Cf() const;
513 
514  //- Return face deltas as surfaceVectorField
516 
517  //- Return cell face motion fluxes
518  const surfaceScalarField& phi() const;
519 
520  //- Return the list of fields of type GeoField
521  // i.e. fields that are not geometry
522  template<class GeoField>
523  UPtrList<GeoField> fields(const bool strict = false) const;
524 
525  //- Return the list of current fields of type GeoField
526  // i.e. fields that are not old-times or geometry
527  template<class GeoField>
529 
530  //- Return a labelType of valid component indicators
531  // 1 : valid (solved)
532  // -1 : invalid (not solved)
533  template<class Type>
534  typename pTraits<Type>::labelType validComponents() const;
535 
536  //- Return the fvSchemes
537  const fvSchemes& schemes() const;
538 
539  //- Return the fvSchemes
540  const fvSolution& solution() const;
541 
542 
543  // Edit
544 
545  //- Does the mesh topology change?
546  bool topoChanging() const;
547 
548  //- Is this mesh dynamic?
549  bool dynamic() const;
550 
551  //- Update the mesh for topology change, mesh to mesh mapping
552  // or redistribution
553  bool update();
554 
555  //- Move the mesh
556  bool move();
557 
558  //- Clear all geometry and addressing
559  void clearOut();
560 
561  //- Update mesh corresponding to the given map
562  virtual void topoChange(const polyTopoChangeMap& map);
563 
564  //- Update from another mesh using the given map
565  virtual void mapMesh(const polyMeshMap&);
566 
567  //- Redistribute or update using the given distribution map
568  virtual void distribute(const polyDistributionMap& map);
569 
570  //- Reset the points
571  // without storing old points or returning swept volumes
572  virtual void setPoints(const pointField&);
573 
574  //- Move points, returns volumes swept by faces in motion
575  virtual tmp<scalarField> movePoints(const pointField&);
576 
577  //- Conform the fvMesh to the polyMesh
578  void conform
579  (
580  const surfaceScalarField& phi =
581  NullObjectRef<surfaceScalarField>()
582  );
583 
584  //- Unconform the fvMesh from the polyMesh
585  void unconform
586  (
588  polyFacesBf,
589  const surfaceVectorField& Sf,
590  const surfaceVectorField& Cf,
591  const surfaceScalarField& phi =
592  NullObjectRef<surfaceScalarField>(),
593  const bool sync = true
594  );
595 
596  //- Map all fields in time using given map.
597  void mapFields(const polyTopoChangeMap& map);
598 
599  //- Add/insert single patch. If validBoundary the new situation
600  // is consistent across processors.
601  virtual void addPatch
602  (
603  const label insertPatchi,
604  const polyPatch& patch,
605  const dictionary& patchFieldDict,
606  const word& defaultPatchFieldType,
607  const bool validBoundary
608  );
609 
610  //- Reorder and trim existing patches. If validBoundary the new
611  // situation is consistent across processors
612  virtual void reorderPatches
613  (
614  const labelUList& newToOld,
615  const bool validBoundary
616  );
617 
618  //- Remove boundary patches. Warning: fvPatchFields hold ref to
619  // these fvPatches.
620  void removeFvBoundary();
621 
622  //- Swap mesh
623  // For run-time mesh replacement and mesh to mesh mapping
624  void swap(fvMesh&);
625 
626 
627  // Write
628 
629  //- Write the underlying polyMesh and other data
630  virtual bool writeObject
631  (
635  const bool write = true
636  ) const;
637 
638  //- Write mesh using IO settings from time
639  virtual bool write(const bool write = true) const;
640 
641 
642  // Member Operators
643 
644  //- Disallow default bitwise assignment
645  void operator=(const fvMesh&) = delete;
646 
647  bool operator!=(const fvMesh&) const;
648 
649  bool operator==(const fvMesh&) const;
650 };
651 
652 
653 template<>
655 fvMesh::validComponents<sphericalTensor>() const;
656 
657 
658 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
659 
660 } // End namespace Foam
661 
662 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
663 
664 #ifdef NoRepository
665  #include "fvMeshTemplates.C"
666 #endif
667 
668 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
669 
670 #endif
671 
672 // ************************************************************************* //
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
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
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:402
const volVectorField & C() const
Return cell centres.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:447
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:895
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: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:408
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:901
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:1523
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:953
const fvMeshDistributor & distributor() const
Return the distributor function class.
Definition: fvMesh.C:1052
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1750
fvMesh Mesh
Definition: fvMesh.H:303
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:639
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
static const HashSet< word > geometryFields
Set of names of registered geometric fields.
Definition: fvMesh.H:317
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:1422
const fvMeshStitcher & stitcher() const
Return the stitcher function class.
Definition: fvMesh.C:1040
const fvSolution & solution() const
Return the fvSchemes.
Definition: fvMesh.C:1761
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:727
const surfaceVectorField & Cf() const
Return face centres.
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:761
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:934
bool operator==(const fvMesh &) const
Definition: fvMesh.C:1780
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:441
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1029
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:1774
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:1680
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1349
void mapFields(const polyTopoChangeMap &map)
Map all fields in time using given map.
Definition: fvMesh.C:1064
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:1736
IOobject polyFacesBfIO(IOobject::readOption r) const
Get the poly faces IO object.
Definition: fvMesh.C:918
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:435
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool update()
Update the mesh for topology change, mesh to mesh mapping.
Definition: fvMesh.C:661
UPtrList< GeoField > fields(const bool strict=false) const
Return the list of fields of type GeoField.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:453
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:1374
UPtrList< GeoField > curFields() const
Return the list of current fields of type GeoField.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:745
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:304
const word & name() const
Return reference to name.
Definition: fvMesh.H:416
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1655
fvMesh(const IOobject &io, const bool changers=true, const stitchType stitch=stitchType::geometric)
Construct from IOobject.
Definition: fvMesh.C:335
ClassName("fvMesh")
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:653
const fvMeshTopoChanger & topoChanger() const
Return the topo-changer function class.
Definition: fvMesh.C:1046
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1236
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:647
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:286
bool move()
Move the mesh.
Definition: fvMesh.C:700
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1133
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1058
const polyMesh & operator()() const
Return reference to polyMesh.
Definition: fvMesh.H:425
void conform(const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >())
Conform the fvMesh to the polyMesh.
Definition: fvMesh.C:1399
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:912
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:1373
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:468
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1580
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:1360
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