ENH: Added new nutUBlendedWallFunction

This forms part of what is termed 'automatic wall treatment' in the
reference:

    Menter, F., Carregal Ferreira, J., Esch, T., Konno, B. (2003).
    The SST Turbulence Model with Improved Wall Treatment
    for Heat Transfer Predictions in Gas Turbines.
    Proceedings of the International Gas Turbine Congress 2003 Tokyo

Note
    The full 'automatic wall treatment' description also requires use of
    the  Foam::omegaWallFunction with the \c blended flag set to 'on'
This commit is contained in:
Andrew Heather
2017-12-12 08:51:47 +00:00
parent 3b70a82bb7
commit 0af97856f1
3 changed files with 432 additions and 0 deletions

View File

@ -40,6 +40,7 @@ $(nutWallFunctions)/nutkWallFunction/nutkWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutkRoughWallFunction/nutkRoughWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutkAtmRoughWallFunction/nutkAtmRoughWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutUBlendedWallFunction/nutUBlendedWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutUWallFunction/nutUWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutUSpaldingWallFunction/nutUSpaldingWallFunctionFvPatchScalarField.C
$(nutWallFunctions)/nutUTabulatedWallFunction/nutUTabulatedWallFunctionFvPatchScalarField.C

View File

@ -0,0 +1,225 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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 "nutUBlendedWallFunctionFvPatchScalarField.H"
#include "turbulenceModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::scalarField>
Foam::nutUBlendedWallFunctionFvPatchScalarField::calcNut() const
{
const label patchi = patch().index();
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
const scalarField magGradU(mag(Uw.snGrad()));
const tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
return max
(
scalar(0),
sqr(calcUTau(magGradU))/(magGradU + ROOTVSMALL) - nuw
);
}
Foam::tmp<Foam::scalarField>
Foam::nutUBlendedWallFunctionFvPatchScalarField::calcUTau
(
const scalarField& magGradU
) const
{
const label patchi = patch().index();
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
const scalarField& y = turbModel.y()[patchi];
const tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
const vectorField n(patch().nf());
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
vectorField Up(Uw.patchInternalField() - Uw);
Up -= n*(n & Up);
const scalarField magUp(mag(Up));
tmp<scalarField> tuTaup(new scalarField(patch().size(), 0.0));
scalarField& uTaup = tuTaup.ref();
const scalarField& nutw = *this;
forAll(uTaup, facei)
{
scalar ut = sqrt((nutw[facei] + nuw[facei])*magGradU[facei]);
scalar error = GREAT;
label iter = 0;
while (iter++ < 10 && error > 0.001)
{
scalar yPlus = y[facei]*ut/nuw[facei];
scalar uTauVis = magUp[facei]/yPlus;
scalar uTauLog = kappa_*magUp[facei]/log(E_*yPlus);
scalar utNew = pow(pow(uTauVis, n_) + pow(uTauLog, n_), 1.0/n_);
error = mag(ut - utNew)/(ut + ROOTVSMALL);
ut = 0.5*(ut + utNew);
}
uTaup[facei] = ut;
}
return tuTaup;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::nutUBlendedWallFunctionFvPatchScalarField::
nutUBlendedWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
nutWallFunctionFvPatchScalarField(p, iF),
n_(4)
{}
Foam::nutUBlendedWallFunctionFvPatchScalarField::
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
n_(ptf.n_)
{}
Foam::nutUBlendedWallFunctionFvPatchScalarField::
nutUBlendedWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
nutWallFunctionFvPatchScalarField(p, iF, dict),
n_(dict.lookupOrDefault<scalar>("n", 4))
{}
Foam::nutUBlendedWallFunctionFvPatchScalarField::
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField& wfpsf
)
:
nutWallFunctionFvPatchScalarField(wfpsf),
n_(wfpsf.n_)
{}
Foam::nutUBlendedWallFunctionFvPatchScalarField::
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField& wfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
nutWallFunctionFvPatchScalarField(wfpsf, iF),
n_(wfpsf.n_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::scalarField>
Foam::nutUBlendedWallFunctionFvPatchScalarField::yPlus() const
{
const label patchi = patch().index();
const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
(
IOobject::groupName
(
turbulenceModel::propertiesName,
internalField().group()
)
);
const scalarField& y = turbModel.y()[patchi];
const tmp<scalarField> tnuw = turbModel.nu(patchi);
const scalarField& nuw = tnuw();
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
const scalarField magGradU(mag(Uw.snGrad()));
return y*calcUTau(magGradU)/nuw;
}
void Foam::nutUBlendedWallFunctionFvPatchScalarField::write(Ostream& os) const
{
fvPatchField<scalar>::write(os);
writeLocalEntries(os);
writeEntry("value", os);
os.writeEntry("n", n_);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
nutUBlendedWallFunctionFvPatchScalarField
);
}
// ************************************************************************* //

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016-2017 OpenCFD Ltd.
\\/ M anipulation |
-------------------------------------------------------------------------------
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/>.
Class
Foam::nutUBlendedWallFunctionFvPatchScalarField
Group
grpWallFunctions
Description
This boundary condition provides a turbulent kinematic viscosity condition
when using wall functions, based on a blending of laminar sub-layer and
log region contributions.
\f[
uTau = (uTau_v^n + uTau_l^n)^(1/n)
\f]
where
\vartable
uTau | friction velocity
uTau_v | friction velocity in the viscous region
uTau_l | friction velocity in the log region
\endvartable
Usage
Example of the boundary condition specification:
\verbatim
myPatch
{
type nutUBlendedWallFunction;
}
\endverbatim
Reference:
See the section that describes 'automatic wall treatment'
\verbatim
Menter, F., Carregal Ferreira, J., Esch, T., Konno, B. (2003).
The SST Turbulence Model with Improved Wall Treatment
for Heat Transfer Predictions in Gas Turbines.
Proceedings of the International Gas Turbine Congress 2003 Tokyo
\endverbatim
Note
The full 'automatic wall treatment' description also requires use of the
Foam::omegaWallFunction with the \c blended flag set to 'on'
SeeAlso
Foam::nutWallFunctionFvPatchScalarField
Foam::omegaWallFunctionFvPatchScalarField
SourceFiles
nutUBlendedWallFunctionFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef nutUBlendedWallFunctionFvPatchScalarField_H
#define nutUBlendedWallFunctionFvPatchScalarField_H
#include "nutWallFunctionFvPatchScalarField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class nutUBlendedWallFunctionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class nutUBlendedWallFunctionFvPatchScalarField
:
public nutWallFunctionFvPatchScalarField
{
protected:
// Protected data
//- Model coefficient; default = 4
scalar n_;
// Protected Member Functions
//- Calculate the turbulence viscosity
virtual tmp<scalarField> calcNut() const;
//- Calculate the friction velocity
virtual tmp<scalarField> calcUTau(const scalarField& magGradU) const;
public:
//- Runtime type information
TypeName("nutUBlendedWallFunction");
// Constructors
//- Construct from patch and internal field
nutUBlendedWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
nutUBlendedWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// nutUBlendedWallFunctionFvPatchScalarField
// onto a new patch
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField>
(
new nutUBlendedWallFunctionFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
nutUBlendedWallFunctionFvPatchScalarField
(
const nutUBlendedWallFunctionFvPatchScalarField&,
const DimensionedField<scalar, volMesh>&
);
//- Construct and return a clone setting internal field reference
virtual tmp<fvPatchScalarField> clone
(
const DimensionedField<scalar, volMesh>& iF
) const
{
return tmp<fvPatchScalarField>
(
new nutUBlendedWallFunctionFvPatchScalarField(*this, iF)
);
}
// Member functions
// Evaluation functions
//- Calculate and return the yPlus at the boundary
virtual tmp<scalarField> yPlus() const;
// I-O
//- Write
virtual void write(Ostream& os) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //