mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: IDDESDelta - updated the hmax delta to be run-time selectable
The delta used for hmax can be set using the optional 'hmax' entry - if not supplied, the maxDeltaxyz delta is used (backwards compatibility)
This commit is contained in:
@ -3,7 +3,7 @@
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
|
||||
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
|
||||
\\/ M anipulation | Copyright (C) 2016-2017 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -26,6 +26,7 @@ License
|
||||
#include "IDDESDelta.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "wallDist.H"
|
||||
#include "maxDeltaxyz.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
@ -43,7 +44,7 @@ namespace LESModels
|
||||
|
||||
void Foam::LESModels::IDDESDelta::calcDelta()
|
||||
{
|
||||
const volScalarField& hmax = hmax_;
|
||||
const volScalarField& hmax = hmaxPtr_();
|
||||
const fvMesh& mesh = turbulenceModel_.mesh();
|
||||
|
||||
// Wall-normal vectors
|
||||
@ -143,12 +144,7 @@ Foam::LESModels::IDDESDelta::IDDESDelta
|
||||
)
|
||||
:
|
||||
LESdelta(name, turbulence),
|
||||
hmax_
|
||||
(
|
||||
IOobject::groupName("hmax", turbulence.U().group()),
|
||||
turbulence,
|
||||
dict
|
||||
),
|
||||
hmaxPtr_(nullptr),
|
||||
Cw_
|
||||
(
|
||||
dict.optionalSubDict(type() + "Coeffs").lookupOrDefault<scalar>
|
||||
@ -158,6 +154,33 @@ Foam::LESModels::IDDESDelta::IDDESDelta
|
||||
)
|
||||
)
|
||||
{
|
||||
if (dict.optionalSubDict(type() + "Coeffs").found("hmax"))
|
||||
{
|
||||
// User-defined hmax
|
||||
hmaxPtr_ =
|
||||
LESdelta::New
|
||||
(
|
||||
IOobject::groupName("hmax", turbulence.U().group()),
|
||||
turbulence,
|
||||
dict.optionalSubDict(type() + "Coeffs"),
|
||||
"hmax"
|
||||
);
|
||||
}
|
||||
else
|
||||
{
|
||||
Info<< "Employing " << maxDeltaxyz::typeName << " for hmax" << endl;
|
||||
|
||||
hmaxPtr_.reset
|
||||
(
|
||||
new maxDeltaxyz
|
||||
(
|
||||
IOobject::groupName("hmax", turbulence.U().group()),
|
||||
turbulence,
|
||||
dict.optionalSubDict(type() + "Coeffs")
|
||||
)
|
||||
);
|
||||
}
|
||||
|
||||
calcDelta();
|
||||
}
|
||||
|
||||
@ -178,7 +201,7 @@ void Foam::LESModels::IDDESDelta::correct()
|
||||
{
|
||||
if (turbulenceModel_.mesh().changing())
|
||||
{
|
||||
hmax_.correct();
|
||||
hmaxPtr_->correct();
|
||||
calcDelta();
|
||||
}
|
||||
}
|
||||
|
||||
@ -3,7 +3,7 @@
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
|
||||
\\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd.
|
||||
\\/ M anipulation | Copyright (C) 2016-2017 OpenCFD Ltd.
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of OpenFOAM.
|
||||
@ -37,7 +37,7 @@ SourceFiles
|
||||
#ifndef LESModels_IDDESDelta_H
|
||||
#define LESModels_IDDESDelta_H
|
||||
|
||||
#include "maxDeltaxyz.H"
|
||||
#include "LESdelta.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -56,7 +56,9 @@ class IDDESDelta
|
||||
{
|
||||
// Private data
|
||||
|
||||
maxDeltaxyz hmax_;
|
||||
//- Run-time selectable delta for hmax
|
||||
// Defaults to the maxDeltaXYZ model if not supplied
|
||||
autoPtr<LESdelta> hmaxPtr_;
|
||||
|
||||
scalar Cw_;
|
||||
|
||||
@ -102,7 +104,7 @@ public:
|
||||
//- Return the hmax delta field
|
||||
const volScalarField& hmax() const
|
||||
{
|
||||
return hmax_;
|
||||
return hmaxPtr_();
|
||||
}
|
||||
|
||||
// Correct values
|
||||
|
||||
@ -66,12 +66,13 @@ Foam::autoPtr<Foam::LESdelta> Foam::LESdelta::New
|
||||
(
|
||||
const word& name,
|
||||
const turbulenceModel& turbulence,
|
||||
const dictionary& dict
|
||||
const dictionary& dict,
|
||||
const word& lookupName
|
||||
)
|
||||
{
|
||||
const word deltaType(dict.lookup("delta"));
|
||||
const word deltaType(dict.lookup(lookupName));
|
||||
|
||||
Info<< "Selecting LES delta type " << deltaType << endl;
|
||||
Info<< "Selecting LES " << lookupName << " type " << deltaType << endl;
|
||||
|
||||
auto cstrIter = dictionaryConstructorTablePtr_->cfind(deltaType);
|
||||
|
||||
@ -94,12 +95,13 @@ Foam::autoPtr<Foam::LESdelta> Foam::LESdelta::New
|
||||
const word& name,
|
||||
const turbulenceModel& turbulence,
|
||||
const dictionary& dict,
|
||||
const dictionaryConstructorTable& additionalConstructors
|
||||
const dictionaryConstructorTable& additionalConstructors,
|
||||
const word& lookupName
|
||||
)
|
||||
{
|
||||
const word deltaType(dict.lookup("delta"));
|
||||
const word deltaType(dict.lookup(lookupName));
|
||||
|
||||
Info<< "Selecting LES delta type " << deltaType << endl;
|
||||
Info<< "Selecting LES " << lookupName << " type " << deltaType << endl;
|
||||
|
||||
// First any additional ones
|
||||
{
|
||||
|
||||
@ -108,7 +108,8 @@ public:
|
||||
(
|
||||
const word& name,
|
||||
const turbulenceModel& turbulence,
|
||||
const dictionary& dict
|
||||
const dictionary& dict,
|
||||
const word& lookupName = "delta"
|
||||
);
|
||||
|
||||
//- Return a reference to the selected LES delta
|
||||
@ -117,7 +118,8 @@ public:
|
||||
const word& name,
|
||||
const turbulenceModel& turbulence,
|
||||
const dictionary& dict,
|
||||
const dictionaryConstructorTable&
|
||||
const dictionaryConstructorTable& additionalConstructors,
|
||||
const word& lookupName = "delta"
|
||||
);
|
||||
|
||||
|
||||
|
||||
Reference in New Issue
Block a user