Please see the Extended Code Guide and header files for the details. ENH: add wall-function blending treatments to epsilonWallFunc adds `stepwise`, `max`, `binomial`, and `exponential` wall-function blending treatments: COMP: check backward compatibility for: `serial` `parallel` `serial restart` `parallel restart` computations in comparison to the following tutorials from v1906: `circuitBoardCooling condensatingVessel evaporationMultiComponent hotRoom movingBox multiRegionHeaterRadiation reverseBurner solidQuenching2D` STYLE: simplify `forAll`s ENH: add new wall blending approaches into omegaWallFunction adds stepwise, max, binomial, and exponential wall function blending apprs. ensures/forces the backward compatibility: - blended = true (blending:BINOMIAL2) (default) - blended = false (blending:STEPWISE) simplifies forAlls deprecates objects: - "blended" - use "blending::" options - "includeG" - as was hardcoded ENH: add wall-func blending treatments into nutWallFuncs - nutWallFunction - nutUWallFunction - nutkWallFunction COMP: check backward compatibility for: - serial - parallel - serial restart - parallel restart computations in comparison to the following tutorials from v1906: - heatTransfer/buoyantSimpleFoam/buoyantCavity - compressible/rhoSimpleFoam/gasMixing/injectorPipe DOC: modify header docs in wallFuncs - nutUTabulatedWallFunction - nutUSpaldingWallFunction - nutURoughWallFunction - nutUBlendedWallFunction - REVERT: change write order - nutLowReWallFunction - kLowReWallFunction: - ENH: protect against zero-division error through 'Ceps2' entry - STYLE: remove few redundant empty lines
645 lines
16 KiB
C
645 lines
16 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | www.openfoam.com
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
Copyright (C) 2011-2019 OpenFOAM Foundation
|
|
Copyright (C) 2017-2020 OpenCFD Ltd.
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "epsilonWallFunctionFvPatchScalarField.H"
|
|
#include "nutWallFunctionFvPatchScalarField.H"
|
|
#include "turbulenceModel.H"
|
|
#include "fvMatrix.H"
|
|
#include "addToRunTimeSelectionTable.H"
|
|
|
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
|
|
|
const Foam::Enum
|
|
<
|
|
Foam::epsilonWallFunctionFvPatchScalarField::blendingType
|
|
>
|
|
Foam::epsilonWallFunctionFvPatchScalarField::blendingTypeNames
|
|
({
|
|
{ blendingType::STEPWISE , "stepwise" },
|
|
{ blendingType::MAX , "max" },
|
|
{ blendingType::BINOMIAL , "binomial" },
|
|
{ blendingType::EXPONENTIAL, "exponential" }
|
|
});
|
|
|
|
Foam::scalar Foam::epsilonWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
|
|
|
|
// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::setMaster()
|
|
{
|
|
if (master_ != -1)
|
|
{
|
|
return;
|
|
}
|
|
|
|
const volScalarField& epsilon =
|
|
static_cast<const volScalarField&>(this->internalField());
|
|
|
|
const volScalarField::Boundary& bf = epsilon.boundaryField();
|
|
|
|
label master = -1;
|
|
forAll(bf, patchi)
|
|
{
|
|
if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
|
|
{
|
|
epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
|
|
|
|
if (master == -1)
|
|
{
|
|
master = patchi;
|
|
}
|
|
|
|
epf.master() = master;
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::createAveragingWeights()
|
|
{
|
|
const volScalarField& epsilon =
|
|
static_cast<const volScalarField&>(this->internalField());
|
|
|
|
const volScalarField::Boundary& bf = epsilon.boundaryField();
|
|
|
|
const fvMesh& mesh = epsilon.mesh();
|
|
|
|
if (initialised_ && !mesh.changing())
|
|
{
|
|
return;
|
|
}
|
|
|
|
volScalarField weights
|
|
(
|
|
IOobject
|
|
(
|
|
"weights",
|
|
mesh.time().timeName(),
|
|
mesh,
|
|
IOobject::NO_READ,
|
|
IOobject::NO_WRITE,
|
|
false // do not register
|
|
),
|
|
mesh,
|
|
dimensionedScalar(dimless, Zero)
|
|
);
|
|
|
|
DynamicList<label> epsilonPatches(bf.size());
|
|
forAll(bf, patchi)
|
|
{
|
|
if (isA<epsilonWallFunctionFvPatchScalarField>(bf[patchi]))
|
|
{
|
|
epsilonPatches.append(patchi);
|
|
|
|
const labelUList& faceCells = bf[patchi].patch().faceCells();
|
|
for (const auto& faceCell : faceCells)
|
|
{
|
|
++weights[faceCell];
|
|
}
|
|
}
|
|
}
|
|
|
|
cornerWeights_.setSize(bf.size());
|
|
|
|
for (const auto& patchi : epsilonPatches)
|
|
{
|
|
const fvPatchScalarField& wf = weights.boundaryField()[patchi];
|
|
cornerWeights_[patchi] = 1.0/wf.patchInternalField();
|
|
}
|
|
|
|
G_.setSize(internalField().size(), Zero);
|
|
epsilon_.setSize(internalField().size(), Zero);
|
|
|
|
initialised_ = true;
|
|
}
|
|
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField&
|
|
Foam::epsilonWallFunctionFvPatchScalarField::epsilonPatch
|
|
(
|
|
const label patchi
|
|
)
|
|
{
|
|
const volScalarField& epsilon =
|
|
static_cast<const volScalarField&>(this->internalField());
|
|
|
|
const volScalarField::Boundary& bf = epsilon.boundaryField();
|
|
|
|
const epsilonWallFunctionFvPatchScalarField& epf =
|
|
refCast<const epsilonWallFunctionFvPatchScalarField>(bf[patchi]);
|
|
|
|
return const_cast<epsilonWallFunctionFvPatchScalarField&>(epf);
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::calculateTurbulenceFields
|
|
(
|
|
const turbulenceModel& turbulence,
|
|
scalarField& G0,
|
|
scalarField& epsilon0
|
|
)
|
|
{
|
|
// Accumulate all of the G and epsilon contributions
|
|
forAll(cornerWeights_, patchi)
|
|
{
|
|
if (!cornerWeights_[patchi].empty())
|
|
{
|
|
epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
|
|
|
|
const List<scalar>& w = cornerWeights_[patchi];
|
|
|
|
epf.calculate(turbulence, w, epf.patch(), G0, epsilon0);
|
|
}
|
|
}
|
|
|
|
// Apply zero-gradient condition for epsilon
|
|
forAll(cornerWeights_, patchi)
|
|
{
|
|
if (!cornerWeights_[patchi].empty())
|
|
{
|
|
epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
|
|
|
|
epf == scalarField(epsilon0, epf.patch().faceCells());
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::calculate
|
|
(
|
|
const turbulenceModel& turbModel,
|
|
const List<scalar>& cornerWeights,
|
|
const fvPatch& patch,
|
|
scalarField& G0,
|
|
scalarField& epsilon0
|
|
)
|
|
{
|
|
const label patchi = patch.index();
|
|
|
|
const nutWallFunctionFvPatchScalarField& nutw =
|
|
nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
|
|
|
|
const scalarField& y = turbModel.y()[patchi];
|
|
|
|
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
|
const scalarField& nuw = tnuw();
|
|
|
|
const tmp<volScalarField> tk = turbModel.k();
|
|
const volScalarField& k = tk();
|
|
|
|
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
|
|
|
const scalarField magGradUw(mag(Uw.snGrad()));
|
|
|
|
const scalar Cmu25 = pow025(nutw.Cmu());
|
|
const scalar Cmu75 = pow(nutw.Cmu(), 0.75);
|
|
|
|
// Set epsilon and G
|
|
forAll(nutw, facei)
|
|
{
|
|
const label celli = patch.faceCells()[facei];
|
|
|
|
const scalar yPlus = Cmu25*y[facei]*sqrt(k[celli])/nuw[facei];
|
|
|
|
const scalar w = cornerWeights[facei];
|
|
|
|
scalar epsilonBlended = 0.0;
|
|
|
|
// Contribution from the viscous sublayer
|
|
const scalar epsilonVis = w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
|
|
|
|
// Contribution from the inertial sublayer
|
|
const scalar epsilonLog =
|
|
w*Cmu75*pow(k[celli], 1.5)/(nutw.kappa()*y[facei]);
|
|
|
|
switch (blending_)
|
|
{
|
|
case blendingType::STEPWISE:
|
|
{
|
|
if (lowReCorrection_ && yPlus < nutw.yPlusLam())
|
|
{
|
|
epsilonBlended = epsilonVis;
|
|
}
|
|
else
|
|
{
|
|
epsilonBlended = epsilonLog;
|
|
}
|
|
break;
|
|
}
|
|
|
|
case blendingType::MAX:
|
|
{
|
|
// (PH:Eq. 27)
|
|
epsilonBlended = max(epsilonVis, epsilonLog);
|
|
break;
|
|
}
|
|
|
|
case blendingType::BINOMIAL:
|
|
{
|
|
// (ME:Eqs. 15-16)
|
|
epsilonBlended =
|
|
pow
|
|
(
|
|
pow(epsilonVis, n_) + pow(epsilonLog, n_),
|
|
1.0/n_
|
|
);
|
|
break;
|
|
}
|
|
|
|
case blendingType::EXPONENTIAL:
|
|
{
|
|
// (PH:p. 193)
|
|
const scalar Gamma = 0.001*pow4(yPlus)/(1.0 + yPlus);
|
|
const scalar invGamma = 1.0/(Gamma + ROOTVSMALL);
|
|
epsilonBlended =
|
|
epsilonVis*exp(-Gamma) + epsilonLog*exp(-invGamma);
|
|
break;
|
|
}
|
|
}
|
|
|
|
epsilon0[celli] += epsilonBlended;
|
|
|
|
if (!(lowReCorrection_ && yPlus < nutw.yPlusLam()))
|
|
{
|
|
G0[celli] +=
|
|
w
|
|
*(nutw[facei] + nuw[facei])
|
|
*magGradUw[facei]
|
|
*Cmu25*sqrt(k[celli])
|
|
/(nutw.kappa()*y[facei]);
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField::
|
|
epsilonWallFunctionFvPatchScalarField
|
|
(
|
|
const fvPatch& p,
|
|
const DimensionedField<scalar, volMesh>& iF
|
|
)
|
|
:
|
|
fixedValueFvPatchField<scalar>(p, iF),
|
|
blending_(blendingType::STEPWISE),
|
|
n_(2.0),
|
|
lowReCorrection_(false),
|
|
initialised_(false),
|
|
master_(-1),
|
|
G_(),
|
|
epsilon_(),
|
|
cornerWeights_()
|
|
{}
|
|
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField::
|
|
epsilonWallFunctionFvPatchScalarField
|
|
(
|
|
const epsilonWallFunctionFvPatchScalarField& ptf,
|
|
const fvPatch& p,
|
|
const DimensionedField<scalar, volMesh>& iF,
|
|
const fvPatchFieldMapper& mapper
|
|
)
|
|
:
|
|
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
|
|
blending_(ptf.blending_),
|
|
n_(ptf.n_),
|
|
lowReCorrection_(ptf.lowReCorrection_),
|
|
initialised_(false),
|
|
master_(-1),
|
|
G_(),
|
|
epsilon_(),
|
|
cornerWeights_()
|
|
{}
|
|
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField::
|
|
epsilonWallFunctionFvPatchScalarField
|
|
(
|
|
const fvPatch& p,
|
|
const DimensionedField<scalar, volMesh>& iF,
|
|
const dictionary& dict
|
|
)
|
|
:
|
|
fixedValueFvPatchField<scalar>(p, iF, dict),
|
|
blending_
|
|
(
|
|
blendingTypeNames.getOrDefault
|
|
(
|
|
"blending",
|
|
dict,
|
|
blendingType::STEPWISE
|
|
)
|
|
),
|
|
n_
|
|
(
|
|
dict.getCheckOrDefault<scalar>
|
|
(
|
|
"n",
|
|
2.0,
|
|
scalarMinMax::ge(0)
|
|
)
|
|
),
|
|
lowReCorrection_(dict.getOrDefault("lowReCorrection", false)),
|
|
initialised_(false),
|
|
master_(-1),
|
|
G_(),
|
|
epsilon_(),
|
|
cornerWeights_()
|
|
{
|
|
// Apply zero-gradient condition on start-up
|
|
this->operator==(patchInternalField());
|
|
}
|
|
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField::
|
|
epsilonWallFunctionFvPatchScalarField
|
|
(
|
|
const epsilonWallFunctionFvPatchScalarField& ewfpsf
|
|
)
|
|
:
|
|
fixedValueFvPatchField<scalar>(ewfpsf),
|
|
blending_(ewfpsf.blending_),
|
|
n_(ewfpsf.n_),
|
|
lowReCorrection_(ewfpsf.lowReCorrection_),
|
|
initialised_(false),
|
|
master_(-1),
|
|
G_(),
|
|
epsilon_(),
|
|
cornerWeights_()
|
|
{}
|
|
|
|
|
|
Foam::epsilonWallFunctionFvPatchScalarField::
|
|
epsilonWallFunctionFvPatchScalarField
|
|
(
|
|
const epsilonWallFunctionFvPatchScalarField& ewfpsf,
|
|
const DimensionedField<scalar, volMesh>& iF
|
|
)
|
|
:
|
|
fixedValueFvPatchField<scalar>(ewfpsf, iF),
|
|
blending_(ewfpsf.blending_),
|
|
n_(ewfpsf.n_),
|
|
lowReCorrection_(ewfpsf.lowReCorrection_),
|
|
initialised_(false),
|
|
master_(-1),
|
|
G_(),
|
|
epsilon_(),
|
|
cornerWeights_()
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
Foam::scalarField& Foam::epsilonWallFunctionFvPatchScalarField::G
|
|
(
|
|
bool init
|
|
)
|
|
{
|
|
if (patch().index() == master_)
|
|
{
|
|
if (init)
|
|
{
|
|
G_ = 0.0;
|
|
}
|
|
|
|
return G_;
|
|
}
|
|
|
|
return epsilonPatch(master_).G();
|
|
}
|
|
|
|
|
|
Foam::scalarField& Foam::epsilonWallFunctionFvPatchScalarField::epsilon
|
|
(
|
|
bool init
|
|
)
|
|
{
|
|
if (patch().index() == master_)
|
|
{
|
|
if (init)
|
|
{
|
|
epsilon_ = 0.0;
|
|
}
|
|
|
|
return epsilon_;
|
|
}
|
|
|
|
return epsilonPatch(master_).epsilon(init);
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::updateCoeffs()
|
|
{
|
|
if (updated())
|
|
{
|
|
return;
|
|
}
|
|
|
|
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
|
|
(
|
|
IOobject::groupName
|
|
(
|
|
turbulenceModel::propertiesName,
|
|
internalField().group()
|
|
)
|
|
);
|
|
|
|
setMaster();
|
|
|
|
if (patch().index() == master_)
|
|
{
|
|
createAveragingWeights();
|
|
calculateTurbulenceFields(turbModel, G(true), epsilon(true));
|
|
}
|
|
|
|
const scalarField& G0 = this->G();
|
|
const scalarField& epsilon0 = this->epsilon();
|
|
|
|
typedef DimensionedField<scalar, volMesh> FieldType;
|
|
|
|
FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
|
|
|
|
FieldType& epsilon = const_cast<FieldType&>(internalField());
|
|
|
|
forAll(*this, facei)
|
|
{
|
|
const label celli = patch().faceCells()[facei];
|
|
|
|
G[celli] = G0[celli];
|
|
epsilon[celli] = epsilon0[celli];
|
|
}
|
|
|
|
fvPatchField<scalar>::updateCoeffs();
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::updateWeightedCoeffs
|
|
(
|
|
const scalarField& weights
|
|
)
|
|
{
|
|
if (updated())
|
|
{
|
|
return;
|
|
}
|
|
|
|
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
|
|
(
|
|
IOobject::groupName
|
|
(
|
|
turbulenceModel::propertiesName,
|
|
internalField().group()
|
|
)
|
|
);
|
|
|
|
setMaster();
|
|
|
|
if (patch().index() == master_)
|
|
{
|
|
createAveragingWeights();
|
|
calculateTurbulenceFields(turbModel, G(true), epsilon(true));
|
|
}
|
|
|
|
const scalarField& G0 = this->G();
|
|
const scalarField& epsilon0 = this->epsilon();
|
|
|
|
typedef DimensionedField<scalar, volMesh> FieldType;
|
|
|
|
FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
|
|
|
|
FieldType& epsilon = const_cast<FieldType&>(internalField());
|
|
|
|
scalarField& epsilonf = *this;
|
|
|
|
// Only set the values if the weights are > tolerance
|
|
forAll(weights, facei)
|
|
{
|
|
const scalar w = weights[facei];
|
|
|
|
if (w > tolerance_)
|
|
{
|
|
const label celli = patch().faceCells()[facei];
|
|
|
|
G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
|
|
epsilon[celli] = (1.0 - w)*epsilon[celli] + w*epsilon0[celli];
|
|
epsilonf[facei] = epsilon[celli];
|
|
}
|
|
}
|
|
|
|
fvPatchField<scalar>::updateCoeffs();
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
|
|
(
|
|
fvMatrix<scalar>& matrix
|
|
)
|
|
{
|
|
if (manipulatedMatrix())
|
|
{
|
|
return;
|
|
}
|
|
|
|
matrix.setValues(patch().faceCells(), patchInternalField());
|
|
|
|
fvPatchField<scalar>::manipulateMatrix(matrix);
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix
|
|
(
|
|
fvMatrix<scalar>& matrix,
|
|
const Field<scalar>& weights
|
|
)
|
|
{
|
|
if (manipulatedMatrix())
|
|
{
|
|
return;
|
|
}
|
|
|
|
DynamicList<label> constraintCells(weights.size());
|
|
DynamicList<scalar> constraintValues(weights.size());
|
|
const labelUList& faceCells = patch().faceCells();
|
|
|
|
const DimensionedField<scalar, volMesh>& fld = internalField();
|
|
|
|
forAll(weights, facei)
|
|
{
|
|
// Only set the values if the weights are > tolerance
|
|
if (weights[facei] > tolerance_)
|
|
{
|
|
const label celli = faceCells[facei];
|
|
|
|
constraintCells.append(celli);
|
|
constraintValues.append(fld[celli]);
|
|
}
|
|
}
|
|
|
|
if (debug)
|
|
{
|
|
Pout<< "Patch: " << patch().name()
|
|
<< ": number of constrained cells = " << constraintCells.size()
|
|
<< " out of " << patch().size()
|
|
<< endl;
|
|
}
|
|
|
|
matrix.setValues(constraintCells, constraintValues);
|
|
|
|
fvPatchField<scalar>::manipulateMatrix(matrix);
|
|
}
|
|
|
|
|
|
void Foam::epsilonWallFunctionFvPatchScalarField::write
|
|
(
|
|
Ostream& os
|
|
) const
|
|
{
|
|
os.writeEntry("lowReCorrection", lowReCorrection_);
|
|
os.writeKeyword("blending") << blendingTypeNames[blending_]
|
|
<< token::END_STATEMENT << nl;
|
|
os.writeEntry("n", n_);
|
|
fixedValueFvPatchField<scalar>::write(os);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
makePatchTypeField
|
|
(
|
|
fvPatchScalarField,
|
|
epsilonWallFunctionFvPatchScalarField
|
|
);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|