mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: extend polyBoundaryMesh patch/face query
- whichPatchFace() returns the (patchi, patchFacei) tuple, whichPatch() simply wraps whichPatchFace() - groupNames() : similar to zones ENH: simplify calls to faPatch/fvPatch patchField, lookupPatchField - make second (ununsed) template parameter optional. Was previously needed for old compilers (2008 and earlier).
This commit is contained in:
@ -586,6 +586,12 @@ Foam::List<Foam::labelRange> Foam::polyBoundaryMesh::patchRanges() const
|
||||
}
|
||||
|
||||
|
||||
Foam::wordList Foam::polyBoundaryMesh::groupNames() const
|
||||
{
|
||||
return this->groupPatchIDs().sortedToc();
|
||||
}
|
||||
|
||||
|
||||
Foam::label Foam::polyBoundaryMesh::start() const noexcept
|
||||
{
|
||||
return mesh_.nInternalFaces();
|
||||
@ -787,26 +793,29 @@ Foam::label Foam::polyBoundaryMesh::findPatchID
|
||||
}
|
||||
|
||||
|
||||
Foam::label Foam::polyBoundaryMesh::whichPatch(const label faceIndex) const
|
||||
Foam::labelPair
|
||||
Foam::polyBoundaryMesh::whichPatchFace(const label faceIndex) const
|
||||
{
|
||||
// Find out which patch the current face belongs to by comparing label
|
||||
// with patch start labels.
|
||||
// If the face is internal, return -1;
|
||||
// if it is off the end of the list, abort
|
||||
if (faceIndex < mesh().nInternalFaces())
|
||||
{
|
||||
return -1;
|
||||
// Internal face: return (-1, meshFace)
|
||||
return labelPair(-1, faceIndex);
|
||||
}
|
||||
else if (faceIndex >= mesh().nFaces())
|
||||
{
|
||||
// Bounds error: abort
|
||||
FatalErrorInFunction
|
||||
<< "Face " << faceIndex
|
||||
<< " out of bounds. Number of geometric faces " << mesh().nFaces()
|
||||
<< abort(FatalError);
|
||||
|
||||
return labelPair(-1, faceIndex);
|
||||
}
|
||||
|
||||
|
||||
// Patches are ordered, use binary search
|
||||
// Find out which patch index and local patch face the specified
|
||||
// mesh face belongs to by comparing label with patch start labels.
|
||||
|
||||
const polyPatchList& patches = *this;
|
||||
|
||||
@ -831,10 +840,23 @@ Foam::label Foam::polyBoundaryMesh::whichPatch(const label faceIndex) const
|
||||
<< " total number of faces:" << mesh().nFaces()
|
||||
<< abort(FatalError);
|
||||
|
||||
return -1;
|
||||
return labelPair(-1, faceIndex);
|
||||
}
|
||||
|
||||
return patchi;
|
||||
// (patch, local face index)
|
||||
return labelPair(patchi, faceIndex - patches[patchi].start());
|
||||
}
|
||||
|
||||
|
||||
Foam::labelPairList
|
||||
Foam::polyBoundaryMesh::whichPatchFace(const labelUList& faceIndices) const
|
||||
{
|
||||
labelPairList output(faceIndices.size());
|
||||
forAll(faceIndices, i)
|
||||
{
|
||||
output[i] = whichPatchFace(faceIndices[i]);
|
||||
}
|
||||
return output;
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -190,6 +190,9 @@ public:
|
||||
//- Return a list of patch ranges
|
||||
List<labelRange> patchRanges() const;
|
||||
|
||||
//- A list of the group names (if any)
|
||||
wordList groupNames() const;
|
||||
|
||||
//- The start label of boundary faces in the polyMesh face list
|
||||
// Same as polyMesh::nInternalFaces()
|
||||
label start() const noexcept;
|
||||
@ -242,8 +245,18 @@ public:
|
||||
template<class Type>
|
||||
labelHashSet findPatchIDs() const;
|
||||
|
||||
//- Return patch index for a given face label
|
||||
label whichPatch(const label faceIndex) const;
|
||||
//- Lookup mesh face index and return (patchi, patchFacei) tuple
|
||||
//- or (-1, meshFacei) for internal faces
|
||||
labelPair whichPatchFace(const label faceIndex) const;
|
||||
|
||||
//- Lookup mesh face indices and return (patchi, patchFacei) tuples
|
||||
labelPairList whichPatchFace(const labelUList& faceIndices) const;
|
||||
|
||||
//- Return patch index for a given mesh face index
|
||||
label whichPatch(const label faceIndex) const
|
||||
{
|
||||
return whichPatchFace(faceIndex).first();
|
||||
}
|
||||
|
||||
//- Per boundary face label the patch index
|
||||
const labelList& patchID() const;
|
||||
|
||||
@ -421,6 +421,12 @@ Foam::List<Foam::labelRange> Foam::faBoundaryMesh::patchRanges() const
|
||||
}
|
||||
|
||||
|
||||
Foam::wordList Foam::faBoundaryMesh::groupNames() const
|
||||
{
|
||||
return this->groupPatchIDs().sortedToc();
|
||||
}
|
||||
|
||||
|
||||
Foam::label Foam::faBoundaryMesh::start() const
|
||||
{
|
||||
return mesh_.nInternalEdges();
|
||||
|
||||
@ -160,6 +160,9 @@ public:
|
||||
//- Return a list of patch ranges
|
||||
List<labelRange> patchRanges() const;
|
||||
|
||||
//- A list of the group names (if any)
|
||||
wordList groupNames() const;
|
||||
|
||||
//- The start label of the edges in the faMesh edges list
|
||||
// Same as mesh.nInternalEdges()
|
||||
label start() const;
|
||||
|
||||
@ -414,12 +414,21 @@ public:
|
||||
|
||||
// Evaluation
|
||||
|
||||
//- Extract internal field next to patch using edgeFaces mapping
|
||||
template<class Type>
|
||||
inline void patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
const labelUList& edgeFaces,
|
||||
Field<Type>& pfld
|
||||
) const;
|
||||
|
||||
//- Return given internal field next to patch as patch field
|
||||
template<class Type>
|
||||
tmp<Field<Type>> patchInternalField(const UList<Type>&) const;
|
||||
|
||||
//- Return given internal field next to patch as patch field
|
||||
// providing addressing
|
||||
//- using provided addressing
|
||||
template<class Type>
|
||||
tmp<Foam::Field<Type>> patchInternalField
|
||||
(
|
||||
@ -427,24 +436,28 @@ public:
|
||||
const labelUList& edgeFaces
|
||||
) const;
|
||||
|
||||
//- Return the corresponding patchField of the named field
|
||||
template<class GeometricField, class Type>
|
||||
//- Extract internal field next to patch as patch field
|
||||
template<class Type>
|
||||
void patchInternalField(const UList<Type>&, Field<Type>&) const;
|
||||
|
||||
//- Return the patch field of the GeometricField
|
||||
//- corresponding to this patch.
|
||||
template<class GeometricField, class AnyType = bool>
|
||||
const typename GeometricField::Patch& patchField
|
||||
(
|
||||
const GeometricField&
|
||||
const GeometricField& gf
|
||||
) const;
|
||||
|
||||
//- Lookup and return the patchField of the named field from the
|
||||
//- local objectRegistry.
|
||||
//- Lookup the named field from the local registry and
|
||||
//- return the patch field corresponding to this patch.
|
||||
// N.B. The dummy pointer arguments are used if this function is
|
||||
// instantiated within a templated function to avoid a bug in gcc.
|
||||
// See inletOutletFvPatchField.C and outletInletFvPatchField.C
|
||||
template<class GeometricField, class Type>
|
||||
template<class GeometricField, class AnyType = bool>
|
||||
const typename GeometricField::Patch& lookupPatchField
|
||||
(
|
||||
const word& name,
|
||||
const GeometricField* = nullptr,
|
||||
const Type* = nullptr
|
||||
const AnyType* = nullptr
|
||||
) const;
|
||||
|
||||
|
||||
|
||||
@ -6,7 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2016-2017 Wikki Ltd
|
||||
Copyright (C) 2021 OpenCFD Ltd.
|
||||
Copyright (C) 2021-2022 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -30,19 +30,18 @@ License
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class GeometricField, class Type>
|
||||
template<class GeometricField, class AnyType>
|
||||
const typename GeometricField::Patch& Foam::faPatch::lookupPatchField
|
||||
(
|
||||
const word& name,
|
||||
const GeometricField*,
|
||||
const Type*
|
||||
const AnyType*
|
||||
) const
|
||||
{
|
||||
return patchField<GeometricField, Type>
|
||||
(
|
||||
boundaryMesh().mesh().mesh().objectRegistry::template
|
||||
return
|
||||
boundaryMesh().mesh().thisDb().template
|
||||
lookupObject<GeometricField>(name)
|
||||
);
|
||||
.boundaryField()[this->index()];
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -6,7 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2016-2017 Wikki Ltd
|
||||
Copyright (C) 2019 OpenCFD Ltd.
|
||||
Copyright (C) 2019-2022 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -30,13 +30,32 @@ License
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class Type>
|
||||
void Foam::faPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
const labelUList& edgeFaces,
|
||||
Field<Type>& pfld
|
||||
) const
|
||||
{
|
||||
pfld.resize(size());
|
||||
|
||||
forAll(pfld, i)
|
||||
{
|
||||
pfld[i] = f[edgeFaces[i]];
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
Foam::tmp<Foam::Field<Type>> Foam::faPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f
|
||||
) const
|
||||
{
|
||||
return patchInternalField(f, this->edgeFaces());
|
||||
auto tpfld = tmp<Field<Type>>::New(size());
|
||||
patchInternalField(f, this->edgeFaces(), tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
@ -47,25 +66,30 @@ Foam::tmp<Foam::Field<Type>> Foam::faPatch::patchInternalField
|
||||
const labelUList& edgeFaces
|
||||
) const
|
||||
{
|
||||
auto tpif = tmp<Field<Type>>::New(size());
|
||||
auto& pif = tpif.ref();
|
||||
|
||||
forAll(pif, facei)
|
||||
{
|
||||
pif[facei] = f[edgeFaces[facei]];
|
||||
}
|
||||
|
||||
return tpif;
|
||||
auto tpfld = tmp<Field<Type>>::New(size());
|
||||
patchInternalField(f, edgeFaces, tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
template<class GeometricField, class Type>
|
||||
template<class Type>
|
||||
void Foam::faPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
Field<Type>& pfld
|
||||
) const
|
||||
{
|
||||
patchInternalField(f, this->edgeFaces(), pfld);
|
||||
}
|
||||
|
||||
|
||||
template<class GeometricField, class AnyType>
|
||||
const typename GeometricField::Patch& Foam::faPatch::patchField
|
||||
(
|
||||
const GeometricField& gf
|
||||
) const
|
||||
{
|
||||
return gf.boundaryField()[index()];
|
||||
return gf.boundaryField()[this->index()];
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -28,6 +28,7 @@ License
|
||||
|
||||
#include "fvBoundaryMesh.H"
|
||||
#include "fvMesh.H"
|
||||
#include "PtrListOps.H"
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
@ -99,6 +100,8 @@ Foam::label Foam::fvBoundaryMesh::findPatchID(const word& patchName) const
|
||||
return -1;
|
||||
}
|
||||
|
||||
// OR: return PtrListOps::firstMatching(*this, patchName);
|
||||
|
||||
const fvPatchList& patches = *this;
|
||||
|
||||
forAll(patches, patchi)
|
||||
|
||||
@ -135,10 +135,10 @@ public:
|
||||
using fvPatchList::operator[];
|
||||
|
||||
//- Return const reference to fvPatch by name.
|
||||
const fvPatch& operator[](const word&) const;
|
||||
const fvPatch& operator[](const word& patchName) const;
|
||||
|
||||
//- Return reference to fvPatch by name.
|
||||
fvPatch& operator[](const word&);
|
||||
fvPatch& operator[](const word& patchName);
|
||||
|
||||
|
||||
// Housekeeping
|
||||
|
||||
@ -6,7 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2011-2018 OpenFOAM Foundation
|
||||
Copyright (C) 2020 OpenCFD Ltd.
|
||||
Copyright (C) 2020-2022 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -257,18 +257,27 @@ public:
|
||||
const scalarField& deltaCoeffs() const;
|
||||
|
||||
|
||||
// Evaluation functions
|
||||
// Evaluation Functions
|
||||
|
||||
//- Extract internal field next to patch using faceCells mapping
|
||||
template<class Type>
|
||||
inline void patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
const labelUList& faceCells,
|
||||
Field<Type>& pfld
|
||||
) const;
|
||||
|
||||
//- Return given internal field next to patch as patch field
|
||||
template<class Type>
|
||||
tmp<Field<Type>> patchInternalField(const UList<Type>&) const;
|
||||
|
||||
//- Return given internal field next to patch as patch field
|
||||
//- using faceCells mapping
|
||||
//- using provided addressing
|
||||
template<class Type>
|
||||
tmp<Field<Type>> patchInternalField
|
||||
(
|
||||
const UList<Type>&,
|
||||
const UList<Type>& f,
|
||||
const labelUList& faceCells
|
||||
) const;
|
||||
|
||||
@ -276,24 +285,24 @@ public:
|
||||
template<class Type>
|
||||
void patchInternalField(const UList<Type>&, Field<Type>&) const;
|
||||
|
||||
//- Return the corresponding patchField of the named field
|
||||
template<class GeometricField, class Type>
|
||||
//- Return the patch field of the GeometricField
|
||||
//- corresponding to this patch.
|
||||
template<class GeometricField, class AnyType = bool>
|
||||
const typename GeometricField::Patch& patchField
|
||||
(
|
||||
const GeometricField&
|
||||
const GeometricField& gf
|
||||
) const;
|
||||
|
||||
//- Lookup and return the patchField of the named field from the
|
||||
//- local objectRegistry.
|
||||
//- Lookup the named field from the local registry and
|
||||
//- return the patch field corresponding to this patch.
|
||||
// N.B. The dummy pointer arguments are used if this function is
|
||||
// instantiated within a templated function to avoid a bug in gcc.
|
||||
// See inletOutletFvPatchField.C and outletInletFvPatchField.C
|
||||
template<class GeometricField, class Type>
|
||||
template<class GeometricField, class AnyType = bool>
|
||||
const typename GeometricField::Patch& lookupPatchField
|
||||
(
|
||||
const word& name,
|
||||
const GeometricField* = nullptr,
|
||||
const Type* = nullptr
|
||||
const AnyType* = nullptr
|
||||
) const;
|
||||
};
|
||||
|
||||
|
||||
@ -6,6 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2011-2016 OpenFOAM Foundation
|
||||
Copyright (C) 2022 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -29,19 +30,18 @@ License
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class GeometricField, class Type>
|
||||
template<class GeometricField, class AnyType>
|
||||
const typename GeometricField::Patch& Foam::fvPatch::lookupPatchField
|
||||
(
|
||||
const word& name,
|
||||
const GeometricField*,
|
||||
const Type*
|
||||
const AnyType*
|
||||
) const
|
||||
{
|
||||
return patchField<GeometricField, Type>
|
||||
(
|
||||
boundaryMesh().mesh().objectRegistry::template
|
||||
return
|
||||
boundaryMesh().mesh().thisDb().template
|
||||
lookupObject<GeometricField>(name)
|
||||
);
|
||||
.boundaryField()[this->index()];
|
||||
}
|
||||
|
||||
|
||||
|
||||
@ -6,7 +6,7 @@
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
Copyright (C) 2011-2016 OpenFOAM Foundation
|
||||
Copyright (C) 2019 OpenCFD Ltd.
|
||||
Copyright (C) 2019-2022 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -30,13 +30,32 @@ License
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
template<class Type>
|
||||
void Foam::fvPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
const labelUList& faceCells,
|
||||
Field<Type>& pfld
|
||||
) const
|
||||
{
|
||||
pfld.resize(size());
|
||||
|
||||
forAll(pfld, i)
|
||||
{
|
||||
pfld[i] = f[faceCells[i]];
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
template<class Type>
|
||||
Foam::tmp<Foam::Field<Type>> Foam::fvPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f
|
||||
) const
|
||||
{
|
||||
return patchInternalField(f, this->faceCells());
|
||||
auto tpfld = tmp<Field<Type>>::New(size());
|
||||
patchInternalField(f, this->faceCells(), tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
@ -47,15 +66,9 @@ Foam::tmp<Foam::Field<Type>> Foam::fvPatch::patchInternalField
|
||||
const labelUList& faceCells
|
||||
) const
|
||||
{
|
||||
auto tpif = tmp<Field<Type>>::New(size());
|
||||
auto& pif = tpif.ref();
|
||||
|
||||
forAll(pif, facei)
|
||||
{
|
||||
pif[facei] = f[faceCells[facei]];
|
||||
}
|
||||
|
||||
return tpif;
|
||||
auto tpfld = tmp<Field<Type>>::New(size());
|
||||
patchInternalField(f, faceCells, tpfld.ref());
|
||||
return tpfld;
|
||||
}
|
||||
|
||||
|
||||
@ -63,27 +76,20 @@ template<class Type>
|
||||
void Foam::fvPatch::patchInternalField
|
||||
(
|
||||
const UList<Type>& f,
|
||||
Field<Type>& pif
|
||||
Field<Type>& pfld
|
||||
) const
|
||||
{
|
||||
pif.resize(size());
|
||||
|
||||
const labelUList& faceCells = this->faceCells();
|
||||
|
||||
forAll(pif, facei)
|
||||
{
|
||||
pif[facei] = f[faceCells[facei]];
|
||||
}
|
||||
patchInternalField(f, this->faceCells(), pfld);
|
||||
}
|
||||
|
||||
|
||||
template<class GeometricField, class Type>
|
||||
template<class GeometricField, class AnyType>
|
||||
const typename GeometricField::Patch& Foam::fvPatch::patchField
|
||||
(
|
||||
const GeometricField& gf
|
||||
) const
|
||||
{
|
||||
return gf.boundaryField()[index()];
|
||||
return gf.boundaryField()[this->index()];
|
||||
}
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user