Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows, with a limit imposed on the turbulent viscosity given by Davidson et al. More...


Public Types | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Types inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
| typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
| typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
Public Types inherited from linearViscousStress< RASModel< BasicTurbulenceModel > > | |
| typedef RASModel< BasicTurbulenceModel > ::alphaField | alphaField |
| typedef RASModel< BasicTurbulenceModel > ::rhoField | rhoField |
| typedef RASModel< BasicTurbulenceModel > ::transportModel | transportModel |
Public Types inherited from RASModel< BasicTurbulenceModel > | |
| typedef BasicTurbulenceModel::alphaField | alphaField |
| typedef BasicTurbulenceModel::rhoField | rhoField |
| typedef BasicTurbulenceModel::transportModel | transportModel |
Public Member Functions | |
| TypeName ("v2f") | |
| Runtime type information. More... | |
| v2f (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName) | |
| Construct from components. More... | |
| virtual | ~v2f () |
| Destructor. More... | |
| virtual bool | read () |
| Read RASProperties dictionary. More... | |
| tmp< volScalarField > | DkEff () const |
| Return the effective diffusivity for k. More... | |
| tmp< volScalarField > | DepsilonEff () const |
| Return the effective diffusivity for epsilon. More... | |
| virtual tmp< volScalarField > | k () const |
| Return the turbulence kinetic energy. More... | |
| virtual tmp< volScalarField > | epsilon () const |
| Return the turbulence kinetic energy dissipation rate. More... | |
| virtual tmp< volScalarField > | v2 () const |
| Return turbulence stress normal to streamlines. More... | |
| virtual tmp< volScalarField > | f () const |
| Return the damping function. More... | |
| virtual void | correct () |
| Solve the turbulence equations and correct the turbulence viscosity. More... | |
Public Member Functions inherited from eddyViscosity< RASModel< BasicTurbulenceModel > > | |
| eddyViscosity (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
| Construct from components. More... | |
| virtual | ~eddyViscosity () |
| Destructor. More... | |
| virtual tmp< volScalarField > | nut () const |
| Return the turbulence viscosity. More... | |
| virtual tmp< scalarField > | nut (const label patchi) const |
| Return the turbulence viscosity on patch. More... | |
| virtual tmp< volSymmTensorField > | R () const |
| Return the Reynolds stress tensor. More... | |
| virtual void | validate () |
| Validate the turbulence fields after construction. More... | |
Public Member Functions inherited from linearViscousStress< RASModel< BasicTurbulenceModel > > | |
| linearViscousStress (const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
| Construct from components. More... | |
| virtual | ~linearViscousStress () |
| Destructor. More... | |
| virtual tmp< volSymmTensorField > | devRhoReff () const |
| Return the effective stress tensor. More... | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
| virtual tmp< fvVectorMatrix > | divDevRhoReff (const volScalarField &rho, volVectorField &U) const |
| Return the source term for the momentum equation. More... | |
Public Member Functions inherited from RASModel< BasicTurbulenceModel > | |
| TypeName ("RAS") | |
| Runtime type information. More... | |
| declareRunTimeSelectionTable (autoPtr, RASModel, dictionary,(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName),(alpha, rho, U, alphaRhoPhi, phi, transport, propertiesName)) | |
| RASModel (const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName) | |
| Construct from components. More... | |
| virtual | ~RASModel () |
| Destructor. More... | |
| const dimensionedScalar & | kMin () const |
| Return the lower allowable limit for k (default: SMALL) More... | |
| const dimensionedScalar & | epsilonMin () const |
| Return the lower allowable limit for epsilon (default: SMALL) More... | |
| const dimensionedScalar & | omegaMin () const |
| Return the lower allowable limit for omega (default: SMALL) More... | |
| dimensionedScalar & | kMin () |
| Allow kMin to be changed. More... | |
| dimensionedScalar & | epsilonMin () |
| Allow epsilonMin to be changed. More... | |
| dimensionedScalar & | omegaMin () |
| Allow omegaMin to be changed. More... | |
| virtual const dictionary & | coeffDict () const |
| Const access to the coefficients dictionary. More... | |
| virtual tmp< volScalarField > | nuEff () const |
| Return the effective viscosity. More... | |
| virtual tmp< scalarField > | nuEff (const label patchi) const |
| Return the effective viscosity on patch. More... | |
Public Member Functions inherited from v2fBase | |
| TypeName ("v2fBase") | |
| Runtime type information. More... | |
| v2fBase () | |
| virtual | ~v2fBase () |
| Destructor. More... | |
Protected Member Functions | |
| virtual void | correctNut () |
| tmp< volScalarField > | Ts () const |
| Return time scale, Ts. More... | |
| tmp< volScalarField > | Ls () const |
| Return length scale, Ls. More... | |
Protected Member Functions inherited from RASModel< BasicTurbulenceModel > | |
| virtual void | printCoeffs (const word &type) |
| Print model coefficients. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from RASModel< BasicTurbulenceModel > | |
| static autoPtr< RASModel > | New (const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName) |
| Return a reference to the selected RAS model. More... | |
Lien and Kalitzin's v2-f turbulence model for incompressible and compressible flows, with a limit imposed on the turbulent viscosity given by Davidson et al.
The model solves for turbulence kinetic energy k and turbulence dissipation rate epsilon, with additional equations for the turbulence stress normal to streamlines, v2, and elliptic damping function, f.
The variant implemented employs N=6, such that f=0 on walls.
Wall boundary conditions are:
k = kLowReWallFunction epsilon = epsilonWallFunction v2 = v2WallFunction f = fWallFunction
These are applicable to both low- and high-Reynolds number flows.
Inlet values can be approximated by:
v2 = 2/3 k f = zero-gradient
References:
Lien, F. S., & Kalitzin, G. (2001).
Computations of transonic flow with the v2f turbulence model.
International Journal of Heat and Fluid Flow, 22(1), 53-61.
Davidson, L., Nielsen, P., & Sveningsson, A. (2003).
Modifications of the v2-f model for computing the flow in a
3D wall jet.
Turbulence, Heat and Mass Transfer, 4, 577-584
The default model coefficients are
v2fCoeffs
{
Cmu 0.22;
CmuKEps 0.09;
C1 1.4;
C2 0.3;
CL 0.23;
Ceta 70;
Ceps2 1.9;
Ceps3 -0.33;
sigmaEps 1.3;
sigmaK 1;
}
| typedef BasicTurbulenceModel::alphaField alphaField |
| typedef BasicTurbulenceModel::transportModel transportModel |
| v2f | ( | const alphaField & | alpha, |
| const rhoField & | rho, | ||
| const volVectorField & | U, | ||
| const surfaceScalarField & | alphaRhoPhi, | ||
| const surfaceScalarField & | phi, | ||
| const transportModel & | transport, | ||
| const word & | propertiesName = turbulenceModel::propertiesName, |
||
| const word & | type = typeName |
||
| ) |
Construct from components.
Definition at line 70 of file v2f.C.
References Foam::bound().
Referenced by v2f< BasicTurbulenceModel >::correctNut().


|
inlinevirtual |
Destructor.
Definition at line 204 of file v2f.H.
References v2f< BasicTurbulenceModel >::read().

|
protectedvirtual |
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 56 of file v2f.C.
References optionList::correct(), Foam::min(), options::New(), Foam::sqr(), and v2f< BasicTurbulenceModel >::v2f().

|
protected |
Return time scale, Ts.
Definition at line 40 of file v2f.C.
References Foam::max(), nu, and Foam::sqrt().

|
protected |
Return length scale, Ls.
Definition at line 47 of file v2f.C.
References Foam::max(), nu, Foam::pow(), Foam::pow025(), and Foam::pow3().

| TypeName | ( | "v2f< BasicTurbulenceModel >" | ) |
Runtime type information.
|
virtual |
Read RASProperties dictionary.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 251 of file v2f.C.
References Foam::read().
Referenced by v2f< BasicTurbulenceModel >::~v2f().


|
inline |
Return the effective diffusivity for k.
Definition at line 214 of file v2f.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inline |
Return the effective diffusivity for epsilon.
Definition at line 227 of file v2f.H.
References nu, and eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_.
|
inlinevirtual |
Return the turbulence kinetic energy.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 240 of file v2f.H.
References v2f< BasicTurbulenceModel >::k_.
|
inlinevirtual |
Return the turbulence kinetic energy dissipation rate.
Definition at line 246 of file v2f.H.
References v2f< BasicTurbulenceModel >::epsilon_.
|
inlinevirtual |
Return turbulence stress normal to streamlines.
Implements v2fBase.
Definition at line 252 of file v2f.H.
References v2f< BasicTurbulenceModel >::v2_.
|
inlinevirtual |
Return the damping function.
Implements v2fBase.
Definition at line 258 of file v2f.H.
References v2f< BasicTurbulenceModel >::correct(), and v2f< BasicTurbulenceModel >::f_.

|
virtual |
Solve the turbulence equations and correct the turbulence viscosity.
Implements eddyViscosity< RASModel< BasicTurbulenceModel > >.
Definition at line 276 of file v2f.C.
References Foam::fvc::absolute(), Foam::bound(), optionList::constrain(), correct, optionList::correct(), Foam::fvm::ddt(), Foam::dev(), Foam::dimless, Foam::fvm::div(), Foam::fvc::div(), divU, fvOptions, Foam::constant::universal::G, Foam::fvc::grad(), Foam::fvm::laplacian(), Foam::magSqr(), Foam::min(), N, options::New(), nut, phi, tmp< T >::ref(), Foam::solve(), Foam::fvm::Sp(), Foam::sqr(), Foam::sqrt(), Foam::fvm::SuSp(), Foam::symm(), and Foam::type().
Referenced by v2f< BasicTurbulenceModel >::f().


|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Turbulence kinetic energy.
Definition at line 147 of file v2f.H.
Referenced by v2f< BasicTurbulenceModel >::k().
|
protected |
Turbulence dissipation.
Definition at line 150 of file v2f.H.
Referenced by v2f< BasicTurbulenceModel >::epsilon().
|
protected |
Turbulence stress normal to streamlines.
Definition at line 153 of file v2f.H.
Referenced by v2f< BasicTurbulenceModel >::v2().
|
protected |
Damping function.
Definition at line 156 of file v2f.H.
Referenced by v2f< BasicTurbulenceModel >::f().
|
protected |
|
protected |
1.8.13