Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\/ 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 "epsilonLowReWallFunctionFvPatchScalarField.H"
#include "turbulenceModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Henry
committed
Foam::scalar Foam::epsilonLowReWallFunctionFvPatchScalarField::yPlusLam
(
const scalar kappa,
const scalar E
)
{
scalar ypl = 11.0;
for (int i=0; i<10; i++)
{
ypl = log(max(E*ypl, 1))/kappa;
}
return ypl;
}
Henry
committed
void Foam::epsilonLowReWallFunctionFvPatchScalarField::calculate
Henry
committed
const turbulenceModel& turbModel,
const List<scalar>& cornerWeights,
const fvPatch& patch,
Henry
committed
scalarField& G0,
scalarField& epsilon0
Henry
committed
const scalarField& y = turbModel.y()[patchi];
const scalar Cmu25 = pow025(Cmu_);
const scalar Cmu75 = pow(Cmu_, 0.75);
Henry
committed
const tmp<volScalarField> tk = turbModel.k();
Henry
committed
const tmp<scalarField> tnuw = turbModel.nu(patchi);
Henry
committed
const tmp<scalarField> tnutw = turbModel.nut(patchi);
Henry
committed
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
const scalarField magGradUw(mag(Uw.snGrad()));
Henry
committed
typedef DimensionedField<scalar, volMesh> FieldType;
const DimensionedField<scalar, volMesh>& G =
db().lookupObject<DimensionedField<scalar, volMesh> >
(
turbModel.GName()
);
scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
Henry
committed
epsilon0[celli] += w*Cmu75*pow(k[celli], 1.5)/(kappa_*y[facei]);
G0[celli] +=
w
*(nutw[facei] + nuw[facei])
*magGradUw[facei]
*Cmu25*sqrt(k[celli])
/(kappa_*y[facei]);
Henry
committed
epsilon0[celli] += w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
G0[celli] += G[celli];
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Henry
committed
Foam::epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
epsilonWallFunctionFvPatchScalarField(p, iF),
yPlusLam_(yPlusLam(kappa_, E_))
{}
Henry
committed
Foam::epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
epsilonWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
yPlusLam_(ptf.yPlusLam_)
{}
Henry
committed
Foam::epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
epsilonWallFunctionFvPatchScalarField(p, iF, dict),
yPlusLam_(yPlusLam(kappa_, E_))
{}
Henry
committed
Foam::epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf
)
:
epsilonWallFunctionFvPatchScalarField(ewfpsf),
yPlusLam_(ewfpsf.yPlusLam_)
{}
Henry
committed
Foam::epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
epsilonWallFunctionFvPatchScalarField(ewfpsf, iF),
yPlusLam_(ewfpsf.yPlusLam_)
{}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Henry
committed
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
epsilonLowReWallFunctionFvPatchScalarField
);
}
// ************************************************************************* //