84 class faMeshBoundaryHalo;
85 class faMeshLduAddressing;
95 public MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>,
134 for (
auto& tuples : list)
156 return (procNo() >= 0 && patchEdgei() >= 0);
160 label procNo()
const {
return (*
this)[0]; }
161 void procNo(label val) { (*this)[0] = val; }
164 label patchi()
const {
return (*
this)[1]; }
165 void patchi(label val) { (*this)[1] = val; }
169 label patchEdgei()
const {
return (*
this)[2]; }
170 void patchEdgei(label val) { (*this)[2] = val; }
173 label meshFacei()
const {
return (*
this)[3]; }
174 void meshFacei(label val) { (*this)[3] = val; }
177 label realPatchi()
const 179 const label
id = patchi();
180 return (
id < 0 ? -(
id + 1) :
id);
184 void faPatchi(label val)
192 return (patchi() < 0);
227 mutable label nPoints_;
230 mutable label nEdges_;
233 mutable label nInternalEdges_;
236 mutable label nFaces_;
248 mutable label curTimeIndex_;
254 mutable std::unique_ptr<uindirectPrimitivePatch> patchPtr_;
257 mutable std::unique_ptr<List<labelPair>> polyPatchFacesPtr_;
260 mutable std::unique_ptr<labelList> polyPatchIdsPtr_;
263 mutable std::unique_ptr<List<labelPair>> bndConnectPtr_;
302 mutable std::unique_ptr<vectorField> pointAreaNormalsPtr_;
311 mutable boolList* correctPatchPointNormalsPtr_;
320 mutable std::unique_ptr<faMeshBoundaryHalo> haloMapPtr_;
323 mutable std::unique_ptr<pointField> haloFaceCentresPtr_;
326 mutable std::unique_ptr<vectorField> haloFaceNormalsPtr_;
332 static const int quadricsFit_;
341 void operator=(
const faMesh&) =
delete;
345 void initPatch()
const;
349 void setPrimitiveMeshData();
356 void calcBoundaryConnections()
const;
360 void setBoundaryConnections
369 void calcLduAddressing()
const;
372 void calcPatchStarts()
const;
376 void calcWhichPatchFaces()
const;
389 void calcMagLe()
const;
393 void calcFaceCentres()
const;
397 void calcEdgeCentres()
const;
405 void calcFaceAreaNormals()
const;
409 void calcEdgeAreaNormals()
const;
413 void calcPointAreaNormals(
vectorField& result)
const;
416 void calcPointAreaNormalsByQuadricsFit(
vectorField& result)
const;
420 void calcFaceCurvatures()
const;
423 void calcEdgeTransformTensors()
const;
426 void clearGeomNotAreas()
const;
429 void clearHalo()
const;
432 void clearGeom()
const;
435 void clearAddressing()
const;
438 void clearOut()
const;
444 void calcHaloFaceGeometry()
const;
452 const word& patchName,
453 const word& patchType =
"" 460 const word& emptyPatchName =
"",
461 const dictionary* defaultPatchDefinition =
nullptr 466 void checkBoundaryEdgeLabelRange(
const labelUList& edgeLabels)
const;
478 checkBoundaryEdgeLabelRange(edgeLabels);
484 result[i] = bndField[edgeLabels[i] - nInternalEdges_];
493 static bool hasSystemFiles(
const polyMesh& pMesh);
496 static bool hasFiles(
const polyMesh& pMesh);
569 const bool doInit =
true 607 const bool validBoundary =
true 614 const bool validBoundary =
true 618 bool init(
const bool doInit);
699 virtual
bool hasDb() const;
724 inline label
whichFace(const label meshFacei) const;
760 return (edgeIndex < nInternalEdges_);
824 return bool(faceAreaNormalsPtr_);
830 return bool(edgeAreaNormalsPtr_);
836 return bool(pointAreaNormalsPtr_);
842 return bool(faceCurvaturesPtr_);
864 virtual void mapFields(
const faMeshMapper& mapper)
const;
867 virtual void mapOldAreas(
const faMeshMapper& mapper)
const;
894 const DimensionedField<scalar, areaMesh>&
S()
const;
897 const DimensionedField<scalar, areaMesh>&
S0()
const;
900 const DimensionedField<scalar, areaMesh>&
S00()
const;
937 virtual bool write(
const bool valid =
true)
const;
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
labelList internalPoints() const
Return internal point labels.
void size(const label n)
Older name for setAddressableSize.
virtual void updateMesh(const mapPolyMesh &)
Update after topo change.
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
A class for handling file names.
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
A 1D vector of objects of type <T> with a fixed length <N>.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch with only those pointing to interfaces being set...
virtual bool write(const bool valid=true) const
Write mesh.
const Time & time() const
Return reference to time.
const word & name() const noexcept
Return the object name.
virtual void mapFields(const faMeshMapper &mapper) const
Map all fields in time using given map.
bool hasAdgeCentres() const noexcept
Has edge centres: edgeCentres()
const edgeList::subList internalEdges() const
Sub-list of local internal edges.
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Face to edge interpolation scheme. Included in faMesh.
bool operator!=(const faMesh &m) const
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
label whichFace(const label meshFacei) const
The area-face corresponding to the mesh-face, -1 if not found.
Various mesh related information for a parallel run.
void stableSort(UList< T > &list)
Stable sort the list.
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Generic GeometricField class.
const faGlobalMeshData & globalData() const
Return parallel info.
const labelList & patchStarts() const
Return patch starts.
bool isInternalEdge(const label edgeIndex) const noexcept
True if given edge label is internal to the mesh.
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
const List< labelPair > & boundaryConnections() const
List of proc/face for the boundary edge neighbours using primitive patch edge numbering.
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
labelList boundaryProcs() const
Boundary edge neighbour processors (does not include own proc)
bool hasMagLe() const noexcept
Has edge length magnitudes: magLe()
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
const edgeScalarField & magLe() const
Return edge length magnitudes.
label nFaces() const noexcept
Number of patch faces.
UList< label > labelUList
A UList of labels.
label nPoints() const noexcept
Number of local mesh points.
A field of fields is a PtrList of fields with reference counting.
#define forAll(list, i)
Loop across all elements in list.
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
static int geometryOrder() noexcept
Return the current geometry treatment.
const polyMesh & operator()() const
Interface to referenced polyMesh (similar to GeoMesh)
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
const labelUList & neighbour() const
Internal face neighbour.
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
const faceList & faces() const
Return local faces.
List< labelPair > boundaryProcSizes() const
List of proc/size for the boundary edge neighbour processors (does not include own proc) ...
bool hasLe() const noexcept
Has edge length vectors: Le()
vectorField pointField
pointField is a vectorField.
Forwards for edge field types.
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
An ordered pair of two objects of type <T> with first() and second() elements.
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
virtual ~faMesh()
Destructor.
bool hasEdgeAreaNormals() const noexcept
Has edge area normals: edgeAreaNormals()
A class for handling words, derived from Foam::string.
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
bool operator==(const faMesh &m) const
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
void sort(UList< T > &list)
Sort the list.
const edgeVectorField & Le() const
Return edge length vectors.
bool hasFaceCurvatures() const noexcept
Has face curvatures: faceCurvatures()
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
const faMeshBoundaryHalo & boundaryHaloMap() const
Mapping/swapping for boundary halo neighbours.
bool hasPointAreaNormals() const noexcept
Has point area normals: pointAreaNormals()
bool hasAreaCentres() const noexcept
Has face centres: areaCentres()
static const word prefix
The prefix to local: finite-area.
bool hasInternalEdgeLabels() const noexcept
True if the internalEdges use an ordering that does not correspond 1-to-1 with the patch internalEdge...
const pointField & points() const
Return local points.
const vectorField & haloFaceNormals() const
Face unit-normals of boundary halo neighbours.
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
virtual void mapOldAreas(const faMeshMapper &mapper) const
Map face areas in time using given map.
virtual bool movePoints()
Update after mesh motion.
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
faBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Database for solution data, solver performance and other reduced data.
const polyMesh & mesh() const
Return access to polyMesh.
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
const areaScalarField & faceCurvatures() const
Return face curvatures.
label nEdges() const noexcept
Number of local mesh edges.
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
const labelUList & owner() const
Internal face owner.
label comm() const noexcept
Return communicator used for parallel communication.
bool moving() const noexcept
Is mesh moving.
labelList boundaryPoints() const
Return boundary point labels.
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
const fileName & facesInstance() const
Return the current instance directory for faces.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
label nInternalEdges() const noexcept
Number of internal faces.
Selector class for finite area differencing schemes. faMesh is derived from faSchemes so that all fie...
const edgeVectorField & edgeCentres() const
Return edge centres as edgeVectorField.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
bool moving() const
Is mesh moving.
static autoPtr< faMesh > TryNew(const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
bool hasFaceAreas() const noexcept
Has face areas: S()
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
const word & name() const
Name function is needed to disambiguate those inherited from base classes.
const fileName & pointsInstance() const
Return the current instance directory for points.
Finite area boundary mesh.
Selector class for finite area solution. faMesh is derived from faSolution so that all fields have ac...
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
bool init(const bool doInit)
Initialise non-demand-driven data etc.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
labelList faceCells() const
The volume (owner) cells associated with the area-mesh.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static int geometryOrder_
Geometry treatment.
faMesh Mesh
The mesh type.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Mesh consisting of general polyhedral cells.
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
List< label > labelList
A List of labels.
A class for managing temporary objects.
TypeName("faMesh")
Runtime type information.
void removeFiles() const
Remove all files from mesh instance()
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
label nBoundaryEdges() const noexcept
Number of boundary edges (== nEdges - nInternalEdges)
bool hasFaceAreaNormals() const noexcept
Has face area normals: faceAreaNormals()
const vectorField & pointAreaNormals() const
Return point area normals.
Defines the attributes of an object for which implicit objectRegistry management is supported...
const areaVectorField & faceAreaNormals() const
Return face area normals.
List< bool > boolList
A List of bools.
GeometricField< vector, faPatchField, areaMesh > areaVectorField
const labelList & edgeOwner() const noexcept
Edge owner addressing.
lduAddressing wrapper for faMesh
Forwards and collection of common area field types.
const DimensionedField< scalar, areaMesh > & S0() const
Return old-time face areas.
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.