mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Added compressible form of vanDriest wall damping.
This commit is contained in:
@ -9,6 +9,8 @@ dynOneEqEddy/dynOneEqEddy.C
|
||||
DeardorffDiffStress/DeardorffDiffStress.C
|
||||
SpalartAllmaras/SpalartAllmaras.C
|
||||
|
||||
vanDriestDelta/vanDriestDelta.C
|
||||
|
||||
/* Wall functions */
|
||||
wallFunctions=derivedFvPatchFields/wallFunctions
|
||||
|
||||
|
||||
@ -0,0 +1,160 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-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 "vanDriestDelta.H"
|
||||
#include "LESModel.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "wallDistData.H"
|
||||
#include "wallPointYPlus.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace compressible
|
||||
{
|
||||
namespace LESModels
|
||||
{
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
defineTypeNameAndDebug(vanDriestDelta, 0);
|
||||
addToRunTimeSelectionTable(LESdelta, vanDriestDelta, dictionary);
|
||||
|
||||
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||
|
||||
void vanDriestDelta::calcDelta()
|
||||
{
|
||||
const LESModel& lesModel = mesh_.lookupObject<LESModel>("LESProperties");
|
||||
|
||||
const volVectorField& U = lesModel.U();
|
||||
const volScalarField& rho = lesModel.rho();
|
||||
const volScalarField& mu = lesModel.mu();
|
||||
tmp<volScalarField> muSgs = lesModel.muSgs();
|
||||
|
||||
volScalarField ystar
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"ystar",
|
||||
mesh_.time().constant(),
|
||||
mesh_
|
||||
),
|
||||
mesh_,
|
||||
dimensionedScalar("ystar", dimLength, GREAT)
|
||||
);
|
||||
|
||||
const fvPatchList& patches = mesh_.boundary();
|
||||
forAll(patches, patchi)
|
||||
{
|
||||
if (isType<wallFvPatch>(patches[patchi]))
|
||||
{
|
||||
const fvPatchVectorField& Uw = U.boundaryField()[patchi];
|
||||
const scalarField& rhow = rho.boundaryField()[patchi];
|
||||
const scalarField& muw = mu.boundaryField()[patchi];
|
||||
const scalarField& muSgsw = muSgs().boundaryField()[patchi];
|
||||
|
||||
ystar.boundaryField()[patchi] =
|
||||
muw/(rhow*sqrt((muw + muSgsw)*mag(Uw.snGrad())/rhow + VSMALL));
|
||||
}
|
||||
}
|
||||
|
||||
wallPointYPlus::yPlusCutOff = 500;
|
||||
wallDistData<wallPointYPlus> y(mesh_, ystar);
|
||||
|
||||
delta_ = min
|
||||
(
|
||||
static_cast<const volScalarField&>(geometricDelta_()),
|
||||
(kappa_/Cdelta_)*((scalar(1) + SMALL) - exp(-y/ystar/Aplus_))*y
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
vanDriestDelta::vanDriestDelta
|
||||
(
|
||||
const word& name,
|
||||
const fvMesh& mesh,
|
||||
const dictionary& dd
|
||||
)
|
||||
:
|
||||
LESdelta(name, mesh),
|
||||
geometricDelta_
|
||||
(
|
||||
LESdelta::New("geometricDelta", mesh, dd.subDict(type() + "Coeffs"))
|
||||
),
|
||||
kappa_(dd.lookupOrDefault<scalar>("kappa", 0.4187)),
|
||||
Aplus_
|
||||
(
|
||||
dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Aplus", 26.0)
|
||||
),
|
||||
Cdelta_
|
||||
(
|
||||
dd.subDict(type() + "Coeffs").lookupOrDefault<scalar>("Cdelta", 0.158)
|
||||
),
|
||||
calcInterval_
|
||||
(
|
||||
dd.subDict(type() + "Coeffs").lookupOrDefault<label>("calcInterval", 1)
|
||||
)
|
||||
{
|
||||
delta_ = geometricDelta_();
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
void vanDriestDelta::read(const dictionary& d)
|
||||
{
|
||||
const dictionary& dd(d.subDict(type() + "Coeffs"));
|
||||
|
||||
geometricDelta_().read(dd);
|
||||
d.readIfPresent<scalar>("kappa", kappa_);
|
||||
dd.readIfPresent<scalar>("Aplus", Aplus_);
|
||||
dd.readIfPresent<scalar>("Cdelta", Cdelta_);
|
||||
dd.readIfPresent<label>("calcInterval", calcInterval_);
|
||||
calcDelta();
|
||||
}
|
||||
|
||||
|
||||
void vanDriestDelta::correct()
|
||||
{
|
||||
if (mesh().time().timeIndex() % calcInterval_ == 0)
|
||||
{
|
||||
geometricDelta_().correct();
|
||||
calcDelta();
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace LESModels
|
||||
} // End namespace compressible
|
||||
} // End namespace Foam
|
||||
|
||||
// ************************************************************************* //
|
||||
@ -0,0 +1,114 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-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
|
||||
|
||||
Class
|
||||
Foam::compressible::LESModels::vanDriestDelta
|
||||
|
||||
Description
|
||||
Simple cube-root of cell volume delta used in compressible LES models.
|
||||
|
||||
SourceFiles
|
||||
vanDriestDelta.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef vanDriestDelta_H
|
||||
#define vanDriestDelta_H
|
||||
|
||||
#include "LESdelta.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace compressible
|
||||
{
|
||||
namespace LESModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class vanDriestDelta Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class vanDriestDelta
|
||||
:
|
||||
public LESdelta
|
||||
{
|
||||
// Private data
|
||||
|
||||
autoPtr<LESdelta> geometricDelta_;
|
||||
scalar kappa_;
|
||||
scalar Aplus_;
|
||||
scalar Cdelta_;
|
||||
label calcInterval_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Disallow default bitwise copy construct and assignment
|
||||
vanDriestDelta(const vanDriestDelta&);
|
||||
void operator=(const vanDriestDelta&);
|
||||
|
||||
// Calculate the delta values
|
||||
void calcDelta();
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("vanDriest");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from name, mesh and IOdictionary
|
||||
vanDriestDelta(const word& name, const fvMesh& mesh, const dictionary&);
|
||||
|
||||
|
||||
//- Destructor
|
||||
virtual ~vanDriestDelta()
|
||||
{}
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Read the LESdelta dictionary
|
||||
virtual void read(const dictionary&);
|
||||
|
||||
// Correct values
|
||||
virtual void correct();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace LESModels
|
||||
} // End namespace compressible
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
Reference in New Issue
Block a user