ENH: Added new incompressible epsilon low-Re wall function

This commit is contained in:
andy
2012-11-02 15:57:08 +00:00
parent 0763bc2312
commit 291bf72175
2 changed files with 397 additions and 0 deletions

View File

@ -0,0 +1,206 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "epsilonLowReWallFunctionFvPatchScalarField.H"
#include "RASModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
scalar 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;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
epsilonWallFunctionFvPatchScalarField(p, iF),
yPlusLam_(yPlusLam(kappa_, E_))
{}
epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
epsilonWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
yPlusLam_(ptf.yPlusLam_)
{}
epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
epsilonWallFunctionFvPatchScalarField(p, iF, dict),
yPlusLam_(yPlusLam(kappa_, E_))
{}
epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf
)
:
epsilonWallFunctionFvPatchScalarField(ewfpsf),
yPlusLam_(ewfpsf.yPlusLam_)
{}
epsilonLowReWallFunctionFvPatchScalarField::
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
epsilonWallFunctionFvPatchScalarField(ewfpsf, iF),
yPlusLam_(ewfpsf.yPlusLam_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void epsilonLowReWallFunctionFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
const label patchI = patch().index();
const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
const scalarField& y = rasModel.y()[patchI];
volScalarField& G =
const_cast<volScalarField&>(db().lookupObject<volScalarField>(GName_));
DimensionedField<scalar, volMesh>& epsilon =
const_cast<DimensionedField<scalar, volMesh>&>
(
dimensionedInternalField()
);
const tmp<volScalarField> tk = rasModel.k();
const volScalarField& k = tk();
const tmp<volScalarField> tnu = rasModel.nu();
const scalarField& nuw = tnu().boundaryField()[patchI];
const tmp<volScalarField> tnut = rasModel.nut();
const volScalarField& nut = tnut();
const scalarField& nutw = nut.boundaryField()[patchI];
const fvPatchVectorField& Uw = rasModel.U().boundaryField()[patchI];
const scalarField magGradUw(mag(Uw.snGrad()));
const scalar Cmu25 = pow025(Cmu_);
const scalar Cmu75 = pow(Cmu_, 0.75);
// Set epsilon and G
forAll(nutw, faceI)
{
label faceCellI = patch().faceCells()[faceI];
scalar yPlus = Cmu25*sqrt(k[faceCellI])*y[faceI]/nuw[faceI];
if (yPlus > yPlusLam_)
{
epsilon[faceCellI] = Cmu75*pow(k[faceCellI], 1.5)/(kappa_*y[faceI]);
}
else
{
epsilon[faceCellI] = 2.0*Cmu25*pow(k[faceCellI], 1.5)/y[faceI];
}
G[faceCellI] =
(nutw[faceI] + nuw[faceI])
*magGradUw[faceI]
*Cmu25*sqrt(k[faceCellI])
/(kappa_*y[faceI]);
}
fixedInternalValueFvPatchField<scalar>::updateCoeffs();
// TODO: perform averaging for cells sharing more than one boundary face
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField
(
fvPatchScalarField,
epsilonLowReWallFunctionFvPatchScalarField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,191 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::epsilonLowReWallFunctionFvPatchScalarField
Group
grpIcoWallFunctions
Description
This boundary condition provides a turbulence dissipation wall function
condition for low- and high-Reynolds number turbulent flow cases.
The condition can be applied to wall boundaries, whereby it inserts near
wall epsilon values directly into the epsilon equation to act as a
constraint.
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
\endtable
Example of the boundary condition specification:
\verbatim
myPatch
{
type epsilonLowReWallFunction;
}
\endverbatim
SeeAlso
Foam::epsilonWallFunctionFvPatchScalarField
SourceFiles
epsilonLowReWallFunctionFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef epsilonLowReWallFunctionFvPatchScalarField_H
#define epsilonLowReWallFunctionFvPatchScalarField_H
#include "epsilonWallFunctionFvPatchScalarField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class epsilonLowReWallFunctionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class epsilonLowReWallFunctionFvPatchScalarField
:
public epsilonWallFunctionFvPatchScalarField
{
protected:
// Protected data
//- Y+ at the edge of the laminar sublayer
scalar yPlusLam_;
// Protected Member Functions
//- Calculate the Y+ at the edge of the laminar sublayer
scalar yPlusLam(const scalar kappa, const scalar E);
public:
//- Runtime type information
TypeName("epsilonLowReWallFunction");
// Constructors
//- Construct from patch and internal field
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
epsilonLowReWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// epsilonLowReWallFunctionFvPatchScalarField
// onto a new patch
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField>
(
new epsilonLowReWallFunctionFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
epsilonLowReWallFunctionFvPatchScalarField
(
const epsilonLowReWallFunctionFvPatchScalarField&,
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 epsilonLowReWallFunctionFvPatchScalarField(*this, iF)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //