ENH: polyPatch cached areaFraction setter with uniform fraction

STYLE: polyPatch cached areaFraction as std::unique_ptr

- more consistent with other demand-driven data.
  Getter now returns tmp field instead of const reference.
This commit is contained in:
Mark Olesen
2024-05-31 12:42:52 +01:00
parent 323149ad24
commit c61466b6e5
10 changed files with 63 additions and 83 deletions

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2018-2023 OpenCFD Ltd.
Copyright (C) 2018-2024 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -95,9 +95,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(start),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{
if (constraintType(patchType))
{
@ -124,9 +122,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(start),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{}
@ -151,9 +147,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(dict.get<label>("startFace")),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{
if (constraintType(patchType))
{
@ -180,9 +174,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(pp.start()),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{}
@ -207,9 +199,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(newStart),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{}
@ -234,9 +224,7 @@ Foam::polyPatch::polyPatch
bm.mesh().points()
),
start_(newStart),
boundaryMesh_(bm),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(bm)
{}
@ -245,9 +233,7 @@ Foam::polyPatch::polyPatch(const polyPatch& p)
patchIdentifier(p),
primitivePatch(p),
start_(p.start_),
boundaryMesh_(p.boundaryMesh_),
faceCellsPtr_(nullptr),
mePtr_(nullptr)
boundaryMesh_(p.boundaryMesh_)
{}
@ -362,24 +348,41 @@ Foam::tmp<Foam::scalarField> Foam::polyPatch::areaFraction
forAll(*this, facei)
{
const face& curFace = this->operator[](facei);
fraction[facei] =
mag(faceAreas[facei])/(curFace.mag(points) + ROOTVSMALL);
const face& f = this->operator[](facei);
fraction[facei] = faceAreas[facei].mag()/(f.mag(points) + ROOTVSMALL);
}
return tfraction;
}
const Foam::tmp<Foam::scalarField>& Foam::polyPatch::areaFraction() const
Foam::tmp<Foam::scalarField> Foam::polyPatch::areaFraction() const
{
return areaFraction_;
if (areaFractionPtr_)
{
return tmp<scalarField>(*areaFractionPtr_);
}
return nullptr;
}
void Foam::polyPatch::areaFraction(const scalar fraction)
{
areaFractionPtr_ = std::make_unique<scalarField>(size(), fraction);
}
void Foam::polyPatch::areaFraction(const tmp<scalarField>& fraction)
{
areaFraction_ = fraction;
if (fraction)
{
// Steal or clone
areaFractionPtr_.reset(fraction.ptr());
}
else
{
areaFractionPtr_.reset(nullptr);
}
}
@ -427,6 +430,7 @@ void Foam::polyPatch::clearAddressing()
primitivePatch::clearPatchMeshAddr();
faceCellsPtr_.reset(nullptr);
mePtr_.reset(nullptr);
areaFractionPtr_.reset(nullptr);
}

View File

@ -6,7 +6,7 @@
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2015 OpenFOAM Foundation
Copyright (C) 2015-2023 OpenCFD Ltd.
Copyright (C) 2015-2024 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -89,7 +89,7 @@ class polyPatch
mutable std::unique_ptr<labelList> mePtr_;
//- Cached area fraction
tmp<scalarField> areaFraction_;
std::unique_ptr<scalarField> areaFractionPtr_;
protected:
@ -445,15 +445,18 @@ public:
tmp<vectorField> faceCellCentres() const;
//- Calculate the area fraction as the ratio of the stored face
// area and the area given by the face points.
//- area and the area given by the face points.
tmp<scalarField> areaFraction(const pointField& points) const;
//- Return the cached area fraction. Usually only set for the
//- non-overlap patches on ACMI.
const tmp<scalarField>& areaFraction() const;
//- Return the cached area fraction.
//- Usually only set for the non-overlap patches on ACMI.
tmp<scalarField> areaFraction() const;
//- Set cached area fraction
void areaFraction(const tmp<scalarField>&);
//- Set uniform cached area fraction
void areaFraction(const scalar fraction);
//- Set cached area fraction (non-uniform)
void areaFraction(const tmp<scalarField>& fraction);
// Addressing into mesh

View File

@ -317,16 +317,10 @@ void Foam::activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs()
Info<< "Open fraction = " << openFraction_ << endl;
}
scalar areaFraction = 0.0;
if (opening_)
{
areaFraction = openFraction_;
}
else
{
areaFraction = 1 - openFraction_;
}
const scalar areaFraction =
(
opening_ ? openFraction_ : (1 - openFraction_)
);
if (patch().boundaryMesh().mesh().moving())
{
@ -352,40 +346,19 @@ void Foam::activePressureForceBaffleVelocityFvPatchVectorField::updateCoeffs()
}
const_cast<scalarField&>(patch().magSf()) = mag(patch().Sf());
// Cache fraction
const_cast<polyPatch&>(patch().patch()).areaFraction
(
tmp<scalarField>::New
(
patch().patch().size(),
(1-areaFraction)
)
);
const_cast<polyPatch&>(patch().patch()).areaFraction(1-areaFraction);
// Update owner side of cyclic
const_cast<vectorField&>(cyclicPatch.Sf()) = areaFraction*initCyclicSf_;
const_cast<scalarField&>(cyclicPatch.magSf()) = mag(cyclicPatch.Sf());
const_cast<polyPatch&>(cyclicPatch.patch()).areaFraction
(
tmp<scalarField>::New
(
cyclicPatch.patch().size(),
areaFraction
)
);
const_cast<polyPatch&>(cyclicPatch.patch()).areaFraction(areaFraction);
// Update neighbour side of cyclic
const_cast<vectorField&>(nbrPatch.Sf()) = areaFraction*nbrCyclicSf_;
const_cast<scalarField&>(nbrPatch.magSf()) = mag(nbrPatch.Sf());
const_cast<polyPatch&>(cyclicPatch.patch()).areaFraction
(
tmp<scalarField>::New
(
nbrPatch.patch().size(),
areaFraction
)
);
const_cast<polyPatch&>(nbrPatch.patch()).areaFraction(areaFraction);
curTimeIndex_ = this->db().time().timeIndex();
}

View File

@ -154,7 +154,7 @@ void Foam::wallDistAddressing::correct(volScalarField& y)
for (const label patchi : patchIDs_)
{
const auto& fc = C.boundaryField()[patchi];
const auto& areaFraction = fc.patch().patch().areaFraction();
const auto areaFraction(fc.patch().patch().areaFraction());
forAll(fc, patchFacei)
{

View File

@ -299,7 +299,7 @@ void Foam::InteractionLists<ParticleType>::buildInteractionLists()
{
if (isA<wallPolyPatch>(patch))
{
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
forAll(patch, facei)
{

View File

@ -558,12 +558,12 @@ void Foam::cyclicAMIPolyPatch::movePoints
polyPatch::movePoints
-> primitivePatch::movePoints
-> primitivePatch::clearGeom:
deleteDemandDrivenData(localPointsPtr_);
deleteDemandDrivenData(faceCentresPtr_);
deleteDemandDrivenData(faceAreasPtr_);
deleteDemandDrivenData(magFaceAreasPtr_);
deleteDemandDrivenData(faceNormalsPtr_);
deleteDemandDrivenData(pointNormalsPtr_);
localPointsPtr_.reset(nullptr);
faceCentresPtr_.reset(nullptr);
faceAreasPtr_.reset(nullptr);
magFaceAreasPtr_.reset(nullptr);
faceNormalsPtr_.reset(nullptr);
pointNormalsPtr_.reset(nullptr);
*/
}

View File

@ -822,7 +822,7 @@ void Foam::FaceCellWave<Type, TrackingData>::handleAMICyclicPatches()
// Merge into global storage
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
forAll(receiveInfo, i)
{

View File

@ -249,7 +249,7 @@ void Foam::cellDistFuncs::correctBoundaryFaceCells
if (patchIDs.found(patchi))
{
const polyPatch& patch = pbm[patchi];
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
// Check cells with face on wall
forAll(patch, patchFacei)
@ -306,7 +306,7 @@ void Foam::cellDistFuncs::correctBoundaryPointCells
bitSet isWallPoint(meshPoints.size(), true);
{
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
// Check cells with face on wall
forAll(patch, patchFacei)

View File

@ -53,7 +53,7 @@ void Foam::patchDataWave<TransferType, TrackingData>::setChangedFaces
{
const polyPatch& patch = mesh.boundaryMesh()[patchi];
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
const auto faceCentres(patch.faceCentres());

View File

@ -48,7 +48,7 @@ void Foam::patchWave::setChangedFaces
{
const polyPatch& patch = mesh.boundaryMesh()[patchi];
const auto& areaFraction = patch.areaFraction();
const auto areaFraction(patch.areaFraction());
const auto faceCentres(patch.faceCentres());