fvPatch: Added polyFaceFraction calculation

This returns the fraction of the poly-face that each fv-face in the
patch covers. This is useful for mapping some calculations from the
polyMesh to the fvMesh in cases with non-conformal couplings.
This commit is contained in:
Will Bainbridge
2023-12-08 16:34:18 +00:00
parent b33fcacc79
commit 42ff732d58
8 changed files with 40 additions and 52 deletions

View File

@ -31,7 +31,7 @@ License
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
Foam::scalar Foam::epsilonWallFunctionFvPatchScalarField::tolerance_ = 1e-1;
Foam::scalar Foam::epsilonWallFunctionFvPatchScalarField::tol_ = 1e-1;
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
@ -361,12 +361,10 @@ void Foam::epsilonWallFunctionFvPatchScalarField::updateCoeffs()
internalField()
);
scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
forAll(weights, facei)
{
scalar& w = weights[facei];
w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
}
const scalarField weights
(
max((patch().polyFaceFraction() - tol_)/(1 - tol_), scalar(0))
);
forAll(weights, facei)
{
@ -393,20 +391,11 @@ void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
return;
}
const scalarField& epsilon0 = this->epsilon();
scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
forAll(weights, facei)
{
scalar& w = weights[facei];
w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
}
matrix.setValues
(
patch().faceCells(),
UIndirectList<scalar>(epsilon0, patch().faceCells()),
weights
UIndirectList<scalar>(this->epsilon(), patch().faceCells()),
max((patch().polyFaceFraction() - tol_)/(1 - tol_), scalar(0))
);
fvPatchField<scalar>::manipulateMatrix(matrix);

View File

@ -86,7 +86,7 @@ protected:
// Protected data
//- Tolerance used in weighted calculations
static scalar tolerance_;
static scalar tol_;
//- Local copy of turbulence G field
scalarField G_;

View File

@ -36,7 +36,7 @@ namespace Foam
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
scalar omegaWallFunctionFvPatchScalarField::tolerance_ = 1e-1;
scalar omegaWallFunctionFvPatchScalarField::tol_ = 1e-1;
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
@ -402,12 +402,10 @@ void omegaWallFunctionFvPatchScalarField::updateCoeffs()
internalField()
);
scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
forAll(weights, facei)
{
scalar& w = weights[facei];
w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
}
const scalarField weights
(
max((patch().polyFaceFraction() - tol_)/(1 - tol_), scalar(0))
);
forAll(weights, facei)
{
@ -434,20 +432,11 @@ void omegaWallFunctionFvPatchScalarField::manipulateMatrix
return;
}
const scalarField& omega0 = this->omega();
scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
forAll(weights, facei)
{
scalar& w = weights[facei];
w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
}
matrix.setValues
(
patch().faceCells(),
UIndirectList<scalar>(omega0, patch().faceCells()),
weights
UIndirectList<scalar>(this->omega(), patch().faceCells()),
max((patch().polyFaceFraction() - tol_)/(1 - tol_), scalar(0))
);
fvPatchField<scalar>::manipulateMatrix(matrix);

View File

@ -86,7 +86,7 @@ protected:
// Protected data
//- Tolerance used in weighted calculations
static scalar tolerance_;
static scalar tol_;
//- beta1 coefficient
scalar beta1_;

View File

@ -29,7 +29,7 @@ License
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
Foam::scalar Foam::epsilonmWallFunctionFvPatchScalarField::tolerance_ = 1e-1;
Foam::scalar Foam::epsilonmWallFunctionFvPatchScalarField::tol_ = 1e-1;
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -85,20 +85,11 @@ void Foam::epsilonmWallFunctionFvPatchScalarField::manipulateMatrix
return;
}
const DimensionedField<scalar, volMesh>& epsilon = internalField();
scalarField weights(patch().magSf()/patch().patch().magFaceAreas());
forAll(weights, facei)
{
scalar& w = weights[facei];
w = w <= tolerance_ ? 0 : (w - tolerance_)/(1 - tolerance_);
}
matrix.setValues
(
patch().faceCells(),
UIndirectList<scalar>(epsilon, patch().faceCells()),
weights
UIndirectList<scalar>(internalField(), patch().faceCells()),
max((patch().polyFaceFraction() - tol_)/(1 - tol_), scalar(0))
);
fvPatchField<scalar>::manipulateMatrix(matrix);

View File

@ -54,7 +54,7 @@ class epsilonmWallFunctionFvPatchScalarField
// Private member data
//- Tolerance used in weighted calculations
static scalar tolerance_;
static scalar tol_;
public:

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -153,6 +153,20 @@ Foam::tmp<Foam::vectorField> Foam::fvPatch::delta() const
}
Foam::tmp<Foam::scalarField> Foam::fvPatch::polyFaceFraction() const
{
return
boundaryMesh().mesh().conformal()
? tmp<scalarField>(new scalarField(size(), scalar(1)))
: magSf()
/scalarField
(
boundaryMesh().mesh().magFaceAreas(),
boundaryMesh().mesh().polyFacesBf()[patch().index()]
);
}
void Foam::fvPatch::makeWeights(scalarField& w) const
{
w = 1.0;

View File

@ -221,6 +221,11 @@ public:
// to coupled-cell-centre vector is returned
virtual tmp<vectorField> delta() const;
//- Return the fraction of the poly-face that each fv-face in this
// patch covers. Will be equal to one, unless the mesh is
// non-conformal.
tmp<scalarField> polyFaceFraction() const;
// Access functions for demand driven data