ENH: Added new incompressible k low-Re wall function

This commit is contained in:
andy
2012-11-02 15:57:19 +00:00
parent 291bf72175
commit 44ca03d30f
2 changed files with 465 additions and 0 deletions

View File

@ -0,0 +1,255 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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 "kLowReWallFunctionFvPatchScalarField.H"
#include "RASModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
#include "wallFvPatch.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void kLowReWallFunctionFvPatchScalarField::checkType()
{
if (!isA<wallFvPatch>(patch()))
{
FatalErrorIn("kLowReWallFunctionFvPatchScalarField::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);
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
scalar kLowReWallFunctionFvPatchScalarField::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;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchField<scalar>(p, iF),
Cmu_(0.09),
kappa_(0.41),
E_(9.8),
Ceps2_(1.9),
yPlusLam_(yPlusLam(kappa_, E_))
{
checkType();
}
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
Cmu_(ptf.Cmu_),
kappa_(ptf.kappa_),
E_(ptf.E_),
Ceps2_(ptf.Ceps2_),
yPlusLam_(ptf.yPlusLam_)
{
checkType();
}
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<scalar>(p, iF, dict),
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
Ceps2_(dict.lookupOrDefault<scalar>("Ceps2", 1.9)),
yPlusLam_(yPlusLam(kappa_, E_))
{
checkType();
}
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField& kwfpsf
)
:
fixedValueFvPatchField<scalar>(kwfpsf),
Cmu_(kwfpsf.Cmu_),
kappa_(kwfpsf.kappa_),
E_(kwfpsf.E_),
Ceps2_(kwfpsf.Ceps2_),
yPlusLam_(kwfpsf.yPlusLam_)
{
checkType();
}
kLowReWallFunctionFvPatchScalarField::kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField& kwfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchField<scalar>(kwfpsf, iF),
Cmu_(kwfpsf.Cmu_),
kappa_(kwfpsf.kappa_),
E_(kwfpsf.E_),
Ceps2_(kwfpsf.Ceps2_),
yPlusLam_(kwfpsf.yPlusLam_)
{
checkType();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void kLowReWallFunctionFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
const label patchI = patch().index();
const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
const scalarField& y = rasModel.y()[patchI];
const tmp<volScalarField> tk = rasModel.k();
const volScalarField& k = tk();
const tmp<volScalarField> tnu = rasModel.nu();
const scalarField& nuw = tnu().boundaryField()[patchI];
const scalar Cmu25 = pow025(Cmu_);
scalarField& kw = *this;
// Set k wall values
forAll(kw, faceI)
{
label faceCellI = patch().faceCells()[faceI];
scalar uTau = Cmu25*sqrt(k[faceCellI]);
scalar yPlus = uTau*y[faceI]/nuw[faceI];
if (yPlus > yPlusLam_)
{
scalar Ck = -0.416;
scalar Bk = 8.366;
kw[faceI] = Ck/kappa_*log(yPlus) + Bk;
}
else
{
scalar C = 11.0;
scalar Cf = (1.0/sqr(yPlus + C) + 2.0*yPlus/pow3(C) - 1.0/sqr(C));
kw[faceI] = 2400.0/sqr(Ceps2_)*Cf;
}
kw[faceI] *= sqr(uTau);
}
fixedValueFvPatchField<scalar>::updateCoeffs();
// TODO: perform averaging for cells sharing more than one boundary face
}
void kLowReWallFunctionFvPatchScalarField::evaluate
(
const Pstream::commsTypes commsType
)
{
fixedValueFvPatchField<scalar>::evaluate(commsType);
}
void kLowReWallFunctionFvPatchScalarField::write(Ostream& os) const
{
fixedValueFvPatchField<scalar>::write(os);
os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
os.writeKeyword("Ceps2") << Ceps2_ << token::END_STATEMENT << nl;
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField
(
fvPatchScalarField,
kLowReWallFunctionFvPatchScalarField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,210 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ 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::incompressible::RASModels::kLowReWallFunctionFvPatchScalarField
Group
grpIcoWallFunctions
Description
This boundary condition provides a turbulence kinetic energy wall function
condition for low- and high-Reynolds number turbulent flow cases.
The model operates in two modes, based on the computed laminar-to-turbulent
switch-over y+ value derived from kappa and E.
\heading Patch usage
\table
Property | Description | Required | Default value
Cmu | model coefficient | no | 0.09
kappa | Von Karman constant | no | 0.41
E | model coefficient | no | 9.8
Ceps2 | model coefficient | no | 1.9
\endtable
Example of the boundary condition specification:
\verbatim
myPatch
{
type kLowReWallFunction;
}
\endverbatim
SeeAlso
Foam::fixedValueFvPatchField
SourceFiles
kLowReWallFunctionFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef kLowReWallFunctionFvPatchScalarField_H
#define kLowReWallFunctionFvPatchScalarField_H
#include "fixedValueFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class kLowReWallFunctionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class kLowReWallFunctionFvPatchScalarField
:
public fixedValueFvPatchField<scalar>
{
protected:
// Protected data
//- Cmu coefficient
scalar Cmu_;
//- Von Karman constant
scalar kappa_;
//- E coefficient
scalar E_;
//- Ceps2 coefficient
scalar Ceps2_;
//- Y+ at the edge of the laminar sublayer
scalar yPlusLam_;
// Protected Member Functions
//- Check the type of the patch
virtual void checkType();
//- Calculate the Y+ at the edge of the laminar sublayer
scalar yPlusLam(const scalar kappa, const scalar E);
public:
//- Runtime type information
TypeName("kLowReWallFunction");
// Constructors
//- Construct from patch and internal field
kLowReWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
kLowReWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given kLowReWallFunctionFvPatchScalarField
// onto a new patch
kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField>
(
new kLowReWallFunctionFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
kLowReWallFunctionFvPatchScalarField
(
const kLowReWallFunctionFvPatchScalarField&,
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 kLowReWallFunctionFvPatchScalarField(*this, iF)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Evaluate the patchField
virtual void evaluate(const Pstream::commsTypes);
// I-O
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //