/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2016, 2019 OpenFOAM Foundation
Copyright (C) 2017-2024 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 .
\*---------------------------------------------------------------------------*/
#include "omegaWallFunctionFvPatchScalarField.H"
#include "nutWallFunctionFvPatchScalarField.H"
#include "turbulenceModel.H"
#include "fvMatrix.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
Foam::scalar Foam::omegaWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
void Foam::omegaWallFunctionFvPatchScalarField::setMaster()
{
if (master_ != -1)
{
return;
}
const auto& omega =
static_cast(this->internalField());
const volScalarField::Boundary& bf = omega.boundaryField();
label master = -1;
forAll(bf, patchi)
{
if (isA(bf[patchi]))
{
omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
if (master == -1)
{
master = patchi;
}
opf.master() = master;
}
}
}
void Foam::omegaWallFunctionFvPatchScalarField::createAveragingWeights()
{
const auto& omega =
static_cast(this->internalField());
const volScalarField::Boundary& bf = omega.boundaryField();
const fvMesh& mesh = omega.mesh();
if (initialised_ && !mesh.changing())
{
return;
}
volScalarField weights
(
IOobject
(
"weights",
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
IOobject::NO_REGISTER
),
mesh,
dimensionedScalar(dimless, Zero)
);
DynamicList omegaPatches(bf.size());
forAll(bf, patchi)
{
if (isA(bf[patchi]))
{
omegaPatches.append(patchi);
const labelUList& faceCells = bf[patchi].patch().faceCells();
for (const auto& celli : faceCells)
{
++weights[celli];
}
}
}
cornerWeights_.setSize(bf.size());
for (const auto& patchi : omegaPatches)
{
const fvPatchScalarField& wf = weights.boundaryField()[patchi];
cornerWeights_[patchi] = 1.0/wf.patchInternalField();
}
G_.setSize(internalField().size(), 0.0);
omega_.setSize(internalField().size(), 0.0);
initialised_ = true;
}
Foam::omegaWallFunctionFvPatchScalarField&
Foam::omegaWallFunctionFvPatchScalarField::omegaPatch
(
const label patchi
)
{
const auto& omega =
static_cast(this->internalField());
const volScalarField::Boundary& bf = omega.boundaryField();
const auto& opf =
refCast(bf[patchi]);
return const_cast(opf);
}
void Foam::omegaWallFunctionFvPatchScalarField::calculateTurbulenceFields
(
const turbulenceModel& turbModel,
scalarField& G0,
scalarField& omega0
)
{
// accumulate all of the G and omega contributions
forAll(cornerWeights_, patchi)
{
if (!cornerWeights_[patchi].empty())
{
omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
const List& w = cornerWeights_[patchi];
opf.calculate(turbModel, w, opf.patch(), G0, omega0);
}
}
// apply zero-gradient condition for omega
forAll(cornerWeights_, patchi)
{
if (!cornerWeights_[patchi].empty())
{
omegaWallFunctionFvPatchScalarField& opf = omegaPatch(patchi);
opf == scalarField(omega0, opf.patch().faceCells());
}
}
}
void Foam::omegaWallFunctionFvPatchScalarField::calculate
(
const turbulenceModel& turbModel,
const List& cornerWeights,
const fvPatch& patch,
scalarField& G0,
scalarField& omega0
)
{
const label patchi = patch.index();
const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
const scalar kappa = wallCoeffs_.kappa();
const scalar yPlusLam = wallCoeffs_.yPlusLam();
const scalarField& y = turbModel.y()[patchi];
const labelUList& faceCells = patch.faceCells();
const tmp tnutw = turbModel.nut(patchi);
const scalarField& nutw = tnutw();
const tmp tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
const tmp tk = turbModel.k();
const volScalarField& k = tk();
// Calculate y-plus
const auto yPlus = [&](const label facei) -> scalar
{
return
(
Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nuw[facei]
);
};
// Contribution from the viscous sublayer
const auto omegaVis = [&](const label facei) -> scalar
{
return
(
6.0*nuw[facei]/(beta1_*sqr(y[facei]))
);
};
// Contribution from the inertial sublayer
const auto omegaLog = [&](const label facei) -> scalar
{
return
(
sqrt(k[faceCells[facei]])
/ (Cmu25*kappa*y[facei])
);
};
switch (blender_)
{
case blenderType::STEPWISE:
{
forAll(faceCells, facei)
{
if (yPlus(facei) > yPlusLam)
{
omega0[faceCells[facei]] +=
cornerWeights[facei]
* omegaLog(facei);
}
else
{
omega0[faceCells[facei]] +=
cornerWeights[facei]
* omegaVis(facei);
}
}
break;
}
case blenderType::BINOMIAL:
{
forAll(faceCells, facei)
{
omega0[faceCells[facei]] +=
cornerWeights[facei]
* pow
(
pow(omegaVis(facei), n_) + pow(omegaLog(facei), n_),
scalar(1)/n_
);
}
break;
}
case blenderType::MAX:
{
forAll(faceCells, facei)
{
// (PH:Eq. 27)
omega0[faceCells[facei]] +=
cornerWeights[facei]
* max(omegaVis(facei), omegaLog(facei));
}
break;
}
case blenderType::EXPONENTIAL:
{
forAll(faceCells, facei)
{
// (PH:Eq. 31)
const scalar yPlusFace = yPlus(facei);
const scalar Gamma = 0.01*pow4(yPlusFace)/(1 + 5*yPlusFace);
const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
omega0[faceCells[facei]] +=
cornerWeights[facei]
* (
omegaVis(facei)*exp(-Gamma)
+ omegaLog(facei)*exp(-invGamma)
);
}
break;
}
case blenderType::TANH:
{
forAll(faceCells, facei)
{
// (KAS:Eqs. 33-34)
const scalar omegaVisFace = omegaVis(facei);
const scalar omegaLogFace = omegaLog(facei);
const scalar b1 = omegaVisFace + omegaLogFace;
const scalar b2 =
pow
(
pow(omegaVisFace, 1.2) + pow(omegaLogFace, 1.2),
1.0/1.2
);
const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
omega0[faceCells[facei]] +=
cornerWeights[facei]
*(phiTanh*b1 + (1 - phiTanh)*b2);
}
break;
}
}
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
const scalarField magGradUw(mag(Uw.snGrad()));
forAll(faceCells, facei)
{
if (!(blender_ == blenderType::STEPWISE) || yPlus(facei) > yPlusLam)
{
G0[faceCells[facei]] +=
cornerWeights[facei]
*(nutw[facei] + nuw[facei])
*magGradUw[facei]
*Cmu25*sqrt(k[faceCells[facei]])
/(kappa*y[facei]);
}
}
}
void Foam::omegaWallFunctionFvPatchScalarField::writeLocalEntries
(
Ostream& os
) const
{
wallFunctionBlenders::writeEntries(os);
os.writeEntryIfDifferent("beta1", 0.075, beta1_);
wallCoeffs_.writeEntries(os);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField& iF
)
:
fixedValueFvPatchField(p, iF),
wallFunctionBlenders(),
initialised_(false),
master_(-1),
beta1_(0.075),
wallCoeffs_(),
G_(),
omega_(),
cornerWeights_()
{}
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
(
const omegaWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField(ptf, p, iF, mapper),
wallFunctionBlenders(ptf),
initialised_(false),
master_(-1),
beta1_(ptf.beta1_),
wallCoeffs_(ptf.wallCoeffs_),
G_(),
omega_(),
cornerWeights_()
{}
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField& iF,
const dictionary& dict
)
:
fixedValueFvPatchField(p, iF, dict),
wallFunctionBlenders(dict, blenderType::BINOMIAL, scalar(2)),
initialised_(false),
master_(-1),
beta1_(dict.getOrDefault("beta1", 0.075)),
wallCoeffs_(dict),
G_(),
omega_(),
cornerWeights_()
{
// Apply zero-gradient condition on start-up
this->extrapolateInternal();
}
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
(
const omegaWallFunctionFvPatchScalarField& owfpsf
)
:
fixedValueFvPatchField(owfpsf),
wallFunctionBlenders(owfpsf),
initialised_(false),
master_(-1),
beta1_(owfpsf.beta1_),
wallCoeffs_(owfpsf.wallCoeffs_),
G_(),
omega_(),
cornerWeights_()
{}
Foam::omegaWallFunctionFvPatchScalarField::omegaWallFunctionFvPatchScalarField
(
const omegaWallFunctionFvPatchScalarField& owfpsf,
const DimensionedField& iF
)
:
fixedValueFvPatchField(owfpsf, iF),
wallFunctionBlenders(owfpsf),
initialised_(false),
master_(-1),
beta1_(owfpsf.beta1_),
wallCoeffs_(owfpsf.wallCoeffs_),
G_(),
omega_(),
cornerWeights_()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::scalarField& Foam::omegaWallFunctionFvPatchScalarField::G
(
bool init
)
{
if (patch().index() == master_)
{
if (init)
{
G_ = 0.0;
}
return G_;
}
return omegaPatch(master_).G();
}
Foam::scalarField& Foam::omegaWallFunctionFvPatchScalarField::omega
(
bool init
)
{
if (patch().index() == master_)
{
if (init)
{
omega_ = 0.0;
}
return omega_;
}
return omegaPatch(master_).omega(init);
}
void Foam::omegaWallFunctionFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
const auto& turbModel = db().lookupObject
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
setMaster();
if (patch().index() == master_)
{
createAveragingWeights();
calculateTurbulenceFields(turbModel, G(true), omega(true));
}
const scalarField& G0 = this->G();
const scalarField& omega0 = this->omega();
typedef DimensionedField FieldType;
FieldType& G = db().lookupObjectRef(turbModel.GName());
FieldType& omega = const_cast(internalField());
forAll(*this, facei)
{
const label celli = patch().faceCells()[facei];
G[celli] = G0[celli];
omega[celli] = omega0[celli];
}
fvPatchField::updateCoeffs();
}
void Foam::omegaWallFunctionFvPatchScalarField::updateWeightedCoeffs
(
const scalarField& weights
)
{
if (updated())
{
return;
}
const auto& turbModel = db().lookupObject
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
setMaster();
if (patch().index() == master_)
{
createAveragingWeights();
calculateTurbulenceFields(turbModel, G(true), omega(true));
}
const scalarField& G0 = this->G();
const scalarField& omega0 = this->omega();
typedef DimensionedField FieldType;
FieldType& G = db().lookupObjectRef(turbModel.GName());
FieldType& omega = const_cast(internalField());
scalarField& omegaf = *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];
omega[celli] = (1.0 - w)*omega[celli] + w*omega0[celli];
omegaf[facei] = omega[celli];
}
}
fvPatchField::updateCoeffs();
}
void Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix
(
fvMatrix& matrix
)
{
if (manipulatedMatrix())
{
return;
}
matrix.setValues(patch().faceCells(), patchInternalField());
fvPatchField::manipulateMatrix(matrix);
}
void Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix
(
fvMatrix& matrix,
const Field& weights
)
{
if (manipulatedMatrix())
{
return;
}
DynamicList constraintCells(weights.size());
DynamicList constraintValues(weights.size());
const labelUList& faceCells = patch().faceCells();
const DimensionedField& fld = internalField();
forAll(weights, facei)
{
// only set the values if the weights are > tolerance
if (tolerance_ < weights[facei])
{
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::manipulateMatrix(matrix);
}
void Foam::omegaWallFunctionFvPatchScalarField::write
(
Ostream& os
) const
{
fvPatchField::write(os);
writeLocalEntries(os);
fvPatchField::writeValueEntry(os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
omegaWallFunctionFvPatchScalarField
);
}
// ************************************************************************* //