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-2024 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, template<class> class PatchField, class GeoMesh>
80 class GeometricBoundaryField;
81 
82 template<class Type>
83 class UCompactListList;
84 
85 namespace fvMeshTopoChangers
86 {
87  class list;
88 }
89 
90 
91 /*---------------------------------------------------------------------------*\
92  Class fvMesh Declaration
93 \*---------------------------------------------------------------------------*/
94 
95 class fvMesh
96 :
97  public polyMesh,
98  public lduMesh,
100 {
101 public:
102 
103  // Public data types
104 
105  //- Enumeration defining the state of the mesh after a read update.
106  // Used for post-processing applications, where the mesh
107  // needs to update based on the files written in time directories.
108  enum readUpdateState
109  {
113  TOPO_CHANGE,
115  };
116 
117  //- Extent to which to stitch on read and readUpdate. By default, a
118  // full geometric stitch is performed. A non-geometric stitch can be
119  // done as an optimisation in situations when finite volume geometry
120  // is not needed (e.g., decomposition). Stitching can also be
121  // prevented altogether if that is appropriate (e.g., if the mesh is
122  // loaded for mapping in an un-stitched state).
123  enum class stitchType
124  {
125  none,
126  nonGeometric,
127  geometric
128  };
129 
130 
131 private:
132 
133  // Private Data
134 
135  //- Boundary mesh
136  fvBoundaryMesh boundary_;
137 
138  //- The stitcher function class
139  autoPtr<fvMeshStitcher> stitcher_;
140 
141  //- The topo-changer function class
142  autoPtr<fvMeshTopoChanger> topoChanger_;
143 
144  //- The distributor function class
145  autoPtr<fvMeshDistributor> distributor_;
146 
147  //- The mover function class
148  autoPtr<fvMeshMover> mover_;
149 
150 
151  // Demand-driven data
152 
153  //- Matrix addressing
154  mutable fvMeshLduAddressing* lduPtr_;
155 
156  //- Face-poly-face IO
157  mutable IOobject* polyFacesBfIOPtr_;
158 
159  //- Face-poly-face addressing
161  polyFacesBfPtr_;
162 
163  //- Offsets for poly-bFace-patch and patch-face addressing
164  mutable labelList* polyBFaceOffsetsPtr_;
165 
166  //- Poly-bFace-patch addressing
167  mutable labelList* polyBFaceOffsetPatchesPtr_;
168 
169  //- Poly-bFace-patch-face addressing
170  mutable labelList* polyBFaceOffsetPatchFacesPtr_;
171 
172  //- Poly-bFace-patch addressing
173  mutable UCompactListList<label>* polyBFacePatchesPtr_;
174 
175  //- Poly-bFace-patch-face addressing
176  mutable UCompactListList<label>* polyBFacePatchFacesPtr_;
177 
178  //- Face-owner addressing
180  ownerBfPtr_;
181 
182  //- Current time index for cell volumes
183  // Note. The whole mechanism will be replaced once the
184  // dimensionedField is created and the dimensionedField
185  // will take care of the old-time levels.
186  mutable label curTimeIndex_;
187 
188  //- Cell volumes
190 
191  //- Cell volumes old time level
192  mutable DimensionedField<scalar, volMesh>* V0Ptr_;
193 
194  //- Cell volumes old-old time level
195  mutable DimensionedField<scalar, volMesh>* V00Ptr_;
196 
197  //- Face area vectors
198  mutable slicedSurfaceVectorField* SfSlicePtr_;
199 
200  //- Face area vectors
201  mutable surfaceVectorField* SfPtr_;
202 
203  //- Mag face area vectors
204  mutable slicedSurfaceScalarField* magSfSlicePtr_;
205 
206  //- Mag face area vectors
207  mutable surfaceScalarField* magSfPtr_;
208 
209  //- Cell centres
210  mutable slicedVolVectorField* CSlicePtr_;
211 
212  //- Cell centres
213  mutable volVectorField* CPtr_;
214 
215  //- Face centres
216  mutable slicedSurfaceVectorField* CfSlicePtr_;
217 
218  //- Face centres
219  mutable surfaceVectorField* CfPtr_;
220 
221  //- Face motion fluxes
222  mutable surfaceScalarField* phiPtr_;
223 
224  //- fvSchemes created on demand
225  mutable autoPtr<fvSchemes> fvSchemes_;
226 
227  //- fvSolution created on demand
228  mutable autoPtr<fvSolution> fvSolution_;
229 
230 
231  // Private Member Functions
232 
233  // Storage management
234 
235  //- Clear geometry but not the old-time cell volumes
236  void clearGeomNotOldVol();
237 
238  //- Clear geometry like clearGeomNotOldVol but recreate any
239  // geometric demand-driven data that was set
240  void updateGeomNotOldVol();
241 
242  //- Clear geometry
243  void clearGeom();
244 
245  //- Clear addressing
246  void clearAddressing(const bool isMeshUpdate = false);
247 
248  //- Store old times fields for all fields of a given type and class
249  template<class Type, template<class> class GeoField>
250  void storeOldTimeFields();
251 
252  //- Store old times fields for all fields of a given class
253  template<template<class> class GeoField>
254  void storeOldTimeFields();
255 
256  //- Store old times fields for all fields
257  // Used prior to field mapping following mesh-motion
258  // and NCC stitching
259  void storeOldTimeFields();
260 
261  //- Null oldest time field for all fields of a given type and class
262  template<class Type, template<class> class GeoField>
263  void nullOldestTimeFields();
264 
265  //- Null oldest time field for all fields of a given class
266  template<template<class> class GeoField>
267  void nullOldestTimeFields();
268 
269  //- Null oldest time field for all fields
270  // Used prior to field mapping following mesh topology change
271  void nullOldestTimeFields();
272 
273 
274  // Make topological data
275 
276  //- Get the patch types for the poly faces field
277  wordList polyFacesPatchTypes() const;
278 
279  //- Modify face-poly-face addressing
281  polyFacesBfRef();
282 
283 
284  // Make geometric data
285 
286  //- Make the face area vectors
287  void makeSf() const;
288 
289  //- Make the mag face area vectors
290  void makeMagSf() const;
291 
292  //- Make the cell centres
293  void makeC() const;
294 
295  //- Make the face centres
296  void makeCf() const;
297 
298  //- Modify cell face area vectors
299  surfaceVectorField& SfRef();
300 
301  //- Modify cell face area magnitudes
302  surfaceScalarField& magSfRef();
303 
304  //- Modify cell centres
305  volVectorField& CRef();
306 
307  //- Modify face centres
308  surfaceVectorField& CfRef();
309 
310  //- Modify cell face motion fluxes
311  surfaceScalarField& phiRef();
312 
313 
314 public:
315 
316  // Public Typedefs
317 
318  typedef fvMesh Mesh;
320 
322 
323 
324  // Declare name of the class and its debug switch
325  ClassName("fvMesh");
326 
327 
328  // Static data
329 
330  //- Set of names of registered geometric fields
331  // V, Sf, magSf, C, Cf
332  const static HashSet<word> geometryFields;
333 
334 
335  // Constructors
336 
337  //- Construct from IOobject. Optionally prevent post-construction, so
338  // that postConstruct may be called manually for finer control of
339  // construction of mesh changes and the level of stitching.
340  explicit fvMesh(const IOobject& io, const bool doPost = true);
341 
342  //- Construct from cellShapes with boundary.
343  fvMesh
344  (
345  const IOobject& io,
346  pointField&& points,
347  const cellShapeList& shapes,
348  const faceListList& boundaryFaces,
349  const wordList& boundaryPatchNames,
350  const PtrList<dictionary>& boundaryDicts,
351  const word& defaultBoundaryPatchName,
352  const word& defaultBoundaryPatchType,
353  const bool syncPar = true
354  );
355 
356  //- Construct from components without boundary.
357  // Boundary is added using addFvPatches() member function
358  fvMesh
359  (
360  const IOobject& io,
361  pointField&& points,
362  faceList&& faces,
363  labelList&& allOwner,
364  labelList&& allNeighbour,
365  const bool syncPar = true
366  );
367 
368  //- Construct without boundary from cells rather than owner/neighbour.
369  // Boundary is added using addFvPatches() member function
370  fvMesh
371  (
372  const IOobject& io,
373  pointField&& points,
374  faceList&& faces,
375  cellList&& cells,
376  const bool syncPar = true
377  );
378 
379  //- Move constructor
380  fvMesh(fvMesh&&);
381 
382  //- Disallow default bitwise copy construction
383  fvMesh(const fvMesh&) = delete;
384 
385 
386  //- Destructor
387  virtual ~fvMesh();
388 
389 
390  // Member Functions
391 
392  // Helpers
393 
394  //- Complete construction of the mesh
395  void postConstruct
396  (
397  const bool changers,
398  const stitchType stitch
399  );
400 
401  //- Add boundary patches. Constructor helper
402  void addFvPatches
403  (
404  const List<polyPatch*>&,
405  const bool validBoundary = true
406  );
407 
408  //- Update the mesh based on the mesh files saved in time
409  // directories
411  (
412  const stitchType stitch = stitchType::geometric
413  );
414 
415 
416  // Access
417 
418  //- Return the top-level database
419  const Time& time() const
420  {
421  return polyMesh::time();
422  }
423 
424  //- Return the object registry - resolve conflict polyMesh/lduMesh
425  virtual const objectRegistry& thisDb() const
426  {
427  return polyMesh::thisDb();
428  }
429 
430  //- Return reference to name
431  // Note: name() is currently ambiguous due to derivation from
432  // surfaceInterpolation
433  const word& name() const
434  {
435  return polyMesh::name();
436  }
437 
438  //- Return reference to boundary mesh
439  const fvBoundaryMesh& boundary() const;
440 
441  //- Return reference to polyMesh
442  const polyMesh& mesh() const
443  {
444  return *this;
445  }
446 
447  //- Return reference to polyMesh
448  const polyMesh& operator()() const
449  {
450  return *this;
451  }
452 
453  //- Return ldu addressing
454  virtual const lduAddressing& lduAddr() const;
455 
456  //- Return a list of pointers for each patch
457  // with only those pointing to interfaces being set
458  virtual lduInterfacePtrsList interfaces() const
459  {
460  return boundary().interfaces();
461  }
462 
463  //- Return communicator used for parallel communication
464  virtual label comm() const
465  {
466  return polyMesh::comm();
467  }
468 
469  //- Internal face owner
470  const labelUList& owner() const
471  {
472  return lduAddr().lowerAddr();
473  }
474 
475  //- Internal face neighbour
476  const labelUList& neighbour() const
477  {
478  return lduAddr().upperAddr();
479  }
480 
481  //- Return whether the fvMesh is conformal with the polyMesh
482  bool conformal() const;
483 
484  //- Return face-poly-face addressing
486  polyFacesBf() const;
487 
488  //- Return poly-bFace-patch addressing
490 
491  //- Return poly-bFace-patch-face addressing
493 
494  //- Return face-owner addressing
496  ownerBf() const;
497 
498  //- Return the stitcher function class
499  const fvMeshStitcher& stitcher() const;
500 
501  //- Return the stitcher function class
503 
504  //- Return the topo-changer function class
505  const fvMeshTopoChanger& topoChanger() const;
506 
507  //- Return the distributor function class
508  const fvMeshDistributor& distributor() const;
509 
510  //- Return the mover function class
511  const fvMeshMover& mover() const;
512 
513  //- Return cell volumes
514  const DimensionedField<scalar, volMesh>& V() const;
515 
516  //- Return old-time cell volumes
517  const DimensionedField<scalar, volMesh>& V0() const;
518 
519  //- Return old-old-time cell volumes
520  const DimensionedField<scalar, volMesh>& V00() const;
521 
522  //- Return sub-cycle cell volumes
524 
525  //- Return sub-cycle old-time cell volumes
527 
528  //- Return cell face area vectors
529  const surfaceVectorField& Sf() const;
530 
531  //- Return cell face area magnitudes
532  const surfaceScalarField& magSf() const;
533 
534  //- Return cell centres
535  const volVectorField& C() const;
536 
537  //- Return face centres
538  const surfaceVectorField& Cf() const;
539 
540  //- Return face deltas as surfaceVectorField
542 
543  //- Return cell face motion fluxes
544  const surfaceScalarField& phi() const;
545 
546  //- Return the list of fields of type GeoField
547  // i.e. fields that are not geometry
548  template<class GeoField>
549  UPtrList<GeoField> fields(const bool strict = false) const;
550 
551  //- Return the list of current fields of type GeoField
552  // i.e. fields that are not old-times or geometry
553  template<class GeoField>
555 
556  //- Return a labelType of valid component indicators
557  // 1 : valid (solved)
558  // -1 : invalid (not solved)
559  template<class Type>
560  typename pTraits<Type>::labelType validComponents() const;
561 
562  //- Return the fvSchemes
563  const fvSchemes& schemes() const;
564 
565  //- Return the fvSolution
566  const fvSolution& solution() const;
567 
568 
569  // Edit
570 
571  //- Does the mesh topology change?
572  bool topoChanging() const;
573 
574  //- Does the mesh get redistributed?
575  bool distributing() const;
576 
577  //- Is this mesh dynamic?
578  bool dynamic() const;
579 
580  //- Update the mesh for topology change, mesh to mesh mapping
581  // or redistribution
582  bool update();
583 
584  //- Move the mesh
585  bool move();
586 
587  //- Clear all geometry and addressing
588  void clearOut();
589 
590  //- Update mesh corresponding to the given map
591  virtual void topoChange(const polyTopoChangeMap& map);
592 
593  //- Update from another mesh using the given map
594  virtual void mapMesh(const polyMeshMap&);
595 
596  //- Redistribute or update using the given distribution map
597  virtual void distribute(const polyDistributionMap& map);
598 
599  //- Reset the points
600  // without storing old points or returning swept volumes
601  virtual void setPoints(const pointField&);
602 
603  //- Move points, returns volumes swept by faces in motion
604  virtual tmp<scalarField> movePoints(const pointField&);
605 
606  //- Get the current instance for the poly faces boundary field
607  const fileName& polyFacesBfInstance() const;
608 
609  //- Set the instance for the poly faces boundary field
610  void setPolyFacesBfInstance(const fileName&);
611 
612  //- Conform the fvMesh to the polyMesh
613  void conform
614  (
615  const surfaceScalarField& phi =
616  NullObjectRef<surfaceScalarField>()
617  );
618 
619  //- Unconform the fvMesh from the polyMesh
620  void unconform
621  (
623  polyFacesBf,
624  const surfaceVectorField& Sf,
625  const surfaceVectorField& Cf,
626  const surfaceScalarField& phi =
627  NullObjectRef<surfaceScalarField>(),
628  const bool sync = true
629  );
630 
631  //- Map all fields in time using given map.
632  void mapFields(const polyTopoChangeMap& map);
633 
634  //- Add/insert single patch. If validBoundary the new situation
635  // is consistent across processors.
636  virtual void addPatch
637  (
638  const label insertPatchi,
639  const polyPatch& patch,
640  const dictionary& patchFieldDict,
641  const word& defaultPatchFieldType,
642  const bool validBoundary
643  );
644 
645  //- Reorder and trim existing patches. If validBoundary the new
646  // situation is consistent across processors
647  virtual void reorderPatches
648  (
649  const labelUList& newToOld,
650  const bool validBoundary
651  );
652 
653  //- Remove boundary patches. Warning: fvPatchFields hold ref to
654  // these fvPatches.
655  void removeFvBoundary();
656 
657  //- Swap mesh
658  // For run-time mesh replacement and mesh to mesh mapping
659  void swap(fvMesh&);
660 
661 
662  // Write
663 
664  //- Write the underlying polyMesh and other data
665  virtual bool writeObject
666  (
670  const bool write = true
671  ) const;
672 
673  //- Write mesh using IO settings from time
674  virtual bool write(const bool write = true) const;
675 
676 
677  // Member Operators
678 
679  //- Disallow default bitwise assignment
680  void operator=(const fvMesh&) = delete;
681 
682  bool operator!=(const fvMesh&) const;
683 
684  bool operator==(const fvMesh&) const;
685 };
686 
687 
688 template<>
690 fvMesh::validComponents<sphericalTensor>() const;
691 
692 
693 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
694 
695 } // End namespace Foam
696 
697 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
698 
699 #ifdef NoRepository
700  #include "fvMeshTemplates.C"
701 #endif
702 
703 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
704 
705 #endif
706 
707 // ************************************************************************* //
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: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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const volVectorField & C() const
Return cell centres.
const fileName & polyFacesBfInstance() const
Get the current instance for the poly faces boundary field.
Definition: fvMesh.C:1381
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:469
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
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:123
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:424
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:863
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:1519
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:911
const fvMeshDistributor & distributor() const
Return the distributor function class.
Definition: fvMesh.C:1044
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
fvMesh Mesh
Definition: fvMesh.H:317
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:108
@ POINTS_MOVED
Definition: fvMesh.H:110
@ TOPO_CHANGE
Definition: fvMesh.H:112
@ TOPO_PATCH_CHANGE
Definition: fvMesh.H:113
virtual ~fvMesh()
Destructor.
Definition: fvMesh.C:530
pTraits< Type >::labelType validComponents() const
Return a labelType of valid component indicators.
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & ownerBf() const
Return face-owner addressing.
Definition: fvMesh.C:998
static const HashSet< word > geometryFields
Set of names of registered geometric fields.
Definition: fvMesh.H:331
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:1415
const fvMeshStitcher & stitcher() const
Return the stitcher function class.
Definition: fvMesh.C:1026
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1762
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:680
const surfaceVectorField & Cf() const
Return face centres.
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:714
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:880
bool operator==(const fvMesh &) const
Definition: fvMesh.C:1781
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:463
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:987
bool operator!=(const fvMesh &) const
Definition: fvMesh.C:1775
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:1676
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1317
void mapFields(const polyTopoChangeMap &map)
Map all fields in time using given map.
Definition: fvMesh.C:1056
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:1737
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:457
const surfaceVectorField & Sf() const
Return cell face area vectors.
bool update()
Update the mesh for topology change, mesh to mesh mapping.
Definition: fvMesh.C:618
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:475
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:1369
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
Definition: fvMesh.C:1343
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:698
bool distributing() const
Does the mesh get redistributed?
Definition: fvMesh.C:604
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
fvBoundaryMesh BoundaryMesh
Definition: fvMesh.H:318
const word & name() const
Return reference to name.
Definition: fvMesh.H:432
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: fvMesh.C:1651
ClassName("fvMesh")
bool dynamic() const
Is this mesh dynamic?
Definition: fvMesh.C:610
fvMesh(const IOobject &io, const bool doPost=true)
Construct from IOobject. Optionally prevent post-construction, so.
Definition: fvMesh.C:273
const fvMeshTopoChanger & topoChanger() const
Return the topo-changer function class.
Definition: fvMesh.C:1038
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1264
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
tmp< surfaceVectorField > delta() const
Return face deltas as surfaceVectorField.
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:598
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:223
bool move()
Move the mesh.
Definition: fvMesh.C:653
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1125
void postConstruct(const bool changers, const stitchType stitch)
Complete construction of the mesh.
Definition: fvMesh.C:538
const fvMeshMover & mover() const
Return the mover function class.
Definition: fvMesh.C:1050
const polyMesh & operator()() const
Return reference to polyMesh.
Definition: fvMesh.H:447
void conform(const surfaceScalarField &phi=NullObjectRef< surfaceScalarField >())
Conform the fvMesh to the polyMesh.
Definition: fvMesh.C:1392
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:874
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:1326
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:467
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1533
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:1313
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