openfoam/src/turbulenceModels/compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C

273 lines
8.0 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2008-2009 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "epsilonWallFunctionFvPatchScalarField.H"
#include "RASModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
#include "wallFvPatch.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace compressible
{
namespace RASModels
{
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void epsilonWallFunctionFvPatchScalarField::checkType()
{
if (!isA<wallFvPatch>(patch()))
{
FatalErrorIn("epsilonWallFunctionFvPatchScalarField::checkType()")
<< "Invalid wall function specification" << nl
<< " Patch type for patch " << patch().name()
<< " must be wall" << nl
<< " Current patch type is " << patch().type() << nl << endl
<< abort(FatalError);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedInternalValueFvPatchField<scalar>(p, iF),
UName_("U"),
kName_("k"),
GName_("RASModel::G"),
rhoName_("rho"),
muName_("mu"),
mutName_("mut"),
Cmu_(turbulenceModel::standardCmu),
kappa_(turbulenceModel::standardKappa),
E_(turbulenceModel::standardE)
{
checkType();
}
epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
(
const epsilonWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedInternalValueFvPatchField<scalar>(ptf, p, iF, mapper),
UName_(ptf.UName_),
kName_(ptf.kName_),
GName_(ptf.GName_),
rhoName_(ptf.rhoName_),
muName_(ptf.muName_),
mutName_(ptf.mutName_),
Cmu_(ptf.Cmu_),
kappa_(ptf.kappa_),
E_(ptf.E_)
{
checkType();
}
epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedInternalValueFvPatchField<scalar>(p, iF, dict),
UName_(dict.lookupOrDefault<word>("U", "U")),
kName_(dict.lookupOrDefault<word>("k", "k")),
GName_(dict.lookupOrDefault<word>("G", "RASModel::G")),
rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
muName_(dict.lookupOrDefault<word>("mu", "mu")),
mutName_(dict.lookupOrDefault<word>("mut", "mut")),
Cmu_(dict.lookupOrDefault<scalar>("Cmu", turbulenceModel::standardCmu)),
kappa_(dict.lookupOrDefault<scalar>("kappa", turbulenceModel::standardKappa)),
E_(dict.lookupOrDefault<scalar>("E", turbulenceModel::standardE))
{
checkType();
}
epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
(
const epsilonWallFunctionFvPatchScalarField& ewfpsf
)
:
fixedInternalValueFvPatchField<scalar>(ewfpsf),
UName_(ewfpsf.UName_),
kName_(ewfpsf.kName_),
GName_(ewfpsf.GName_),
rhoName_(ewfpsf.rhoName_),
muName_(ewfpsf.muName_),
mutName_(ewfpsf.mutName_),
Cmu_(ewfpsf.Cmu_),
kappa_(ewfpsf.kappa_),
E_(ewfpsf.E_)
{
checkType();
}
epsilonWallFunctionFvPatchScalarField::epsilonWallFunctionFvPatchScalarField
(
const epsilonWallFunctionFvPatchScalarField& ewfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedInternalValueFvPatchField<scalar>(ewfpsf, iF),
UName_(ewfpsf.UName_),
kName_(ewfpsf.kName_),
GName_(ewfpsf.GName_),
rhoName_(ewfpsf.rhoName_),
muName_(ewfpsf.muName_),
mutName_(ewfpsf.mutName_),
Cmu_(ewfpsf.Cmu_),
kappa_(ewfpsf.kappa_),
E_(ewfpsf.E_)
{
checkType();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void epsilonWallFunctionFvPatchScalarField::updateCoeffs()
{
const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
const scalar Cmu25 = pow(Cmu_, 0.25);
const scalar Cmu75 = pow(Cmu_, 0.75);
const scalar yPlusLam = rasModel.yPlusLam(kappa_, E_);
const scalarField& y = rasModel.y()[patch().index()];
volScalarField& G = const_cast<volScalarField&>
(db().lookupObject<volScalarField>(GName_));
volScalarField& epsilon = const_cast<volScalarField&>
(db().lookupObject<volScalarField>(dimensionedInternalField().name()));
const volScalarField& k = db().lookupObject<volScalarField>(kName_);
const scalarField& rhow =
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
const scalarField& muw =
patch().lookupPatchField<volScalarField, scalar>(muName_);
const scalarField& mutw =
patch().lookupPatchField<volScalarField, scalar>(mutName_);
const fvPatchVectorField& Uw =
patch().lookupPatchField<volVectorField, vector>(UName_);
const scalarField magGradUw = mag(Uw.snGrad());
// Set epsilon and G
forAll(mutw, faceI)
{
label faceCellI = patch().faceCells()[faceI];
scalar yPlus =
Cmu25*y[faceI]*sqrt(k[faceCellI])
/(muw[faceI]/rhow[faceI]);
epsilon[faceCellI] = Cmu75*pow(k[faceCellI], 1.5)/(kappa_*y[faceI]);
if (yPlus > yPlusLam)
{
G[faceCellI] =
(mutw[faceI] + muw[faceI])
*magGradUw[faceI]
*Cmu25*sqrt(k[faceCellI])
/(kappa_*y[faceI]);
}
else
{
G[faceCellI] = 0.0;
}
}
// TODO: perform averaging for cells sharing more than one boundary face
}
void epsilonWallFunctionFvPatchScalarField::evaluate
(
const Pstream::commsTypes commsType
)
{
fixedInternalValueFvPatchField<scalar>::evaluate(commsType);
}
void epsilonWallFunctionFvPatchScalarField::write(Ostream& os) const
{
fixedInternalValueFvPatchField<scalar>::write(os);
writeEntryIfDifferent<word>(os, "U", "U", UName_);
writeEntryIfDifferent<word>(os, "k", "k", kName_);
writeEntryIfDifferent<word>(os, "G", "RASModel::G", GName_);
writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
writeEntryIfDifferent<word>(os, "mu", "mu", muName_);
writeEntryIfDifferent<word>(os, "mut", "mut", mutName_);
writeEntryIfDifferent<scalar>(os, "Cmu", turbulenceModel::standardCmu, Cmu_);
writeEntryIfDifferent<scalar>(os, "kappa", turbulenceModel::standardKappa, kappa_);
writeEntryIfDifferent<scalar>(os, "E", turbulenceModel::standardE, E_);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField
(
fvPatchScalarField,
epsilonWallFunctionFvPatchScalarField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace compressible
} // End namespace Foam
// ************************************************************************* //