mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: resolutionIndex: new function object to evaluate LES/DES resolution
Grid independency studies and grid adaptation for implicit LES/DES are nontrivial and intractable due to the inherent coupling between spatial resolution and subgrid-scale modelling. To enable assessments for LES/DES resolution, a function object of single-mesh resolution index with three submodels is introduced.
This commit is contained in:
@ -143,6 +143,13 @@ DMD/DMDModels/DMDModel/DMDModel.C
|
|||||||
DMD/DMDModels/DMDModel/DMDModelNew.C
|
DMD/DMDModels/DMDModel/DMDModelNew.C
|
||||||
DMD/DMDModels/derived/STDMD/STDMD.C
|
DMD/DMDModels/derived/STDMD/STDMD.C
|
||||||
|
|
||||||
|
resolutionIndex/resolutionIndex.C
|
||||||
|
resolutionIndex/resolutionIndexModels/resolutionIndexModel/resolutionIndexModel.C
|
||||||
|
resolutionIndex/resolutionIndexModels/resolutionIndexModel/resolutionIndexModelNew.C
|
||||||
|
resolutionIndex/resolutionIndexModels/PopeIndex/PopeIndex.C
|
||||||
|
resolutionIndex/resolutionIndexModels/CelikNuIndex/CelikNuIndex.C
|
||||||
|
resolutionIndex/resolutionIndexModels/CelikEtaIndex/CelikEtaIndex.C
|
||||||
|
|
||||||
age/age.C
|
age/age.C
|
||||||
comfort/comfort.C
|
comfort/comfort.C
|
||||||
|
|
||||||
|
|||||||
123
src/functionObjects/field/resolutionIndex/resolutionIndex.C
Normal file
123
src/functionObjects/field/resolutionIndex/resolutionIndex.C
Normal file
@ -0,0 +1,123 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "resolutionIndex.H"
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
#include "turbulenceModel.H"
|
||||||
|
#include "RASModel.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace functionObjects
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(resolutionIndex, 0);
|
||||||
|
addToRunTimeSelectionTable(functionObject, resolutionIndex, dictionary);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::functionObjects::resolutionIndex::resolutionIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const Time& runTime,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
fvMeshFunctionObject(name, runTime, dict),
|
||||||
|
resolutionIndexModelPtr_()
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::functionObjects::resolutionIndex::~resolutionIndex()
|
||||||
|
{} // resolutionIndexModel was forward declared
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::functionObjects::resolutionIndex::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (mesh_.foundObject<RASModelBaseName>(turbulenceModel::propertiesName))
|
||||||
|
{
|
||||||
|
FatalIOErrorInFunction(dict)
|
||||||
|
<< type() << " " << name()
|
||||||
|
<< " is not available for RANS-based turbulence models."
|
||||||
|
<< exit(FatalIOError);
|
||||||
|
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
if (!fvMeshFunctionObject::read(dict))
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
Info<< type() << " " << name() << ":" << endl;
|
||||||
|
|
||||||
|
resolutionIndexModelPtr_.reset
|
||||||
|
(
|
||||||
|
resolutionIndexModel::New(name(), mesh_, dict)
|
||||||
|
);
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::functionObjects::resolutionIndex::execute()
|
||||||
|
{
|
||||||
|
if (!resolutionIndexModelPtr_->execute())
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::functionObjects::resolutionIndex::write()
|
||||||
|
{
|
||||||
|
Info<< type() << " " << name() << " write:" << endl;
|
||||||
|
|
||||||
|
if (!resolutionIndexModelPtr_->write())
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
Info<< endl;
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
191
src/functionObjects/field/resolutionIndex/resolutionIndex.H
Normal file
191
src/functionObjects/field/resolutionIndex/resolutionIndex.H
Normal file
@ -0,0 +1,191 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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::functionObjects::resolutionIndex
|
||||||
|
|
||||||
|
Group
|
||||||
|
grpFieldFunctionObjects
|
||||||
|
|
||||||
|
Description
|
||||||
|
Computes a single-mesh resolution index according to the specified model,
|
||||||
|
which is used as a LES/DES quality/post-verification metric that does
|
||||||
|
not require any experimental or DNS data.
|
||||||
|
|
||||||
|
Operands:
|
||||||
|
\table
|
||||||
|
Operand | Type | Location
|
||||||
|
input | vol\<Type\>Field(s) | \<time\>/\<inpField\>
|
||||||
|
output file | - | -
|
||||||
|
output field | volScalarField | \<time\>/\<outField\>
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
Governing equation (tag:P):
|
||||||
|
Pope, S. B. (2000).
|
||||||
|
Turbulent flows.
|
||||||
|
Cambridge, UK: Cambridge Univ. Press
|
||||||
|
DOI:10.1017/CBO9780511840531
|
||||||
|
|
||||||
|
Governing equations (tag:CKJ):
|
||||||
|
Celik, I., Klein, M., & Janicka, J. (2009).
|
||||||
|
Assessment measures for engineering LES applications.
|
||||||
|
Journal of fluids engineering, 131(3).
|
||||||
|
DOI:10.1115/1.3059703
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
Usage
|
||||||
|
Minimal example by using \c system/controlDict.functions:
|
||||||
|
\verbatim
|
||||||
|
resolutionIndex1
|
||||||
|
{
|
||||||
|
// Mandatory entries
|
||||||
|
type resolutionIndex;
|
||||||
|
libs (fieldFunctionObjects);
|
||||||
|
model <word>;
|
||||||
|
|
||||||
|
// Conditional entries
|
||||||
|
|
||||||
|
// Option-1: when model == PopeIndex
|
||||||
|
|
||||||
|
// Option-2: when model == CelikNuIndex
|
||||||
|
|
||||||
|
// Option-3: when model == CelikEtaIndex
|
||||||
|
|
||||||
|
// Inherited entries
|
||||||
|
...
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
where the entries mean:
|
||||||
|
\table
|
||||||
|
Property | Description | Type | Reqd | Deflt
|
||||||
|
type | Type name: resolutionIndex | word | yes | -
|
||||||
|
libs | Library name: fieldFunctionObjects | word | yes | -
|
||||||
|
model | Name of the resolutionIndex model | word | yes | -
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
Options for the \c model entry:
|
||||||
|
\verbatim
|
||||||
|
PopeIndex | Resolution method proposed by Pope (2000)
|
||||||
|
CelikNuIndex | Resolution method proposed by Celik et al. (2009)
|
||||||
|
CelikEtaIndex | Resolution method proposed by Celik et al. (2009)
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
The inherited entries are elaborated in:
|
||||||
|
- \link functionObject.H \endlink
|
||||||
|
|
||||||
|
Note
|
||||||
|
- The \c resolutionIndex throws runtime error
|
||||||
|
when the turbulence model is RANS based.
|
||||||
|
- Resolution-estimator models strictly assume that
|
||||||
|
the flow is in a fully turbulent regime.
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
resolutionIndex.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_functionObjects_resolutionIndex_H
|
||||||
|
#define Foam_functionObjects_resolutionIndex_H
|
||||||
|
|
||||||
|
#include "fvMeshFunctionObject.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
// Forward Declarations
|
||||||
|
class resolutionIndexModel;
|
||||||
|
|
||||||
|
namespace functionObjects
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class resolutionIndex Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class resolutionIndex
|
||||||
|
:
|
||||||
|
public fvMeshFunctionObject
|
||||||
|
{
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Resolution index model
|
||||||
|
autoPtr<resolutionIndexModel> resolutionIndexModelPtr_;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("resolutionIndex");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from Time and dictionary
|
||||||
|
resolutionIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const Time& runTime,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
//- No copy construct
|
||||||
|
resolutionIndex(const resolutionIndex&) = delete;
|
||||||
|
|
||||||
|
//- No copy assignment
|
||||||
|
void operator=(const resolutionIndex&) = delete;
|
||||||
|
|
||||||
|
|
||||||
|
// Destructor
|
||||||
|
virtual ~resolutionIndex();
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Read the resolutionIndex data
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Calculate the result field
|
||||||
|
virtual bool execute();
|
||||||
|
|
||||||
|
//- Write the result field
|
||||||
|
virtual bool write();
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace functionObjects
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,185 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "CelikEtaIndex.H"
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(CelikEtaIndex, 0);
|
||||||
|
addToRunTimeSelectionTable(resolutionIndexModel, CelikEtaIndex, dictionary);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::eta() const
|
||||||
|
{
|
||||||
|
const auto& nu = getOrReadField<volScalarField>(nuName_);
|
||||||
|
tmp<volScalarField> tepsilon = epsilon();
|
||||||
|
|
||||||
|
const dimensionedScalar epsilonMin(tepsilon.cref().dimensions(), SMALL);
|
||||||
|
|
||||||
|
// (CKJ:Eq. 23)
|
||||||
|
return pow025(pow3(nu)/max(epsilonMin, tepsilon));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::epsilon() const
|
||||||
|
{
|
||||||
|
const auto& kSgs = getOrReadField<volScalarField>(kName_);
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
tmp<volScalarField> tnuEff = nuEff();
|
||||||
|
|
||||||
|
// (Derived based on CKJ:Eq. 25-26, p.031102-5)
|
||||||
|
return tnuEff*kSgs/(Ck_*sqr(Delta));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::nuEff() const
|
||||||
|
{
|
||||||
|
const auto& nu = getOrReadField<volScalarField>(nuName_);
|
||||||
|
const auto& nuSgs = getOrReadField<volScalarField>(nutName_);
|
||||||
|
tmp<volScalarField> tnuNum = nuNum();
|
||||||
|
|
||||||
|
// (CKJ:p. 031102-3)
|
||||||
|
return tnuNum + nuSgs + nu;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::nuNum() const
|
||||||
|
{
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
|
||||||
|
tmp<volScalarField> tkNum = kNum();
|
||||||
|
|
||||||
|
// (CKJ:Eq. 35)
|
||||||
|
return sign(tkNum.cref())*Cnu_*Delta*sqrt(mag(tkNum.cref()));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::kNum() const
|
||||||
|
{
|
||||||
|
const auto& kSgs = getOrReadField<volScalarField>(kName_);
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
|
||||||
|
tmp<volScalarField> th = cbrt(V());
|
||||||
|
|
||||||
|
// (CKJ:Eq. 28)
|
||||||
|
return Cn_*sqr(th/Delta)*kSgs;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::resolutionIndexModels::CelikEtaIndex::CelikEtaIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
resolutionIndexModel(name, mesh, dict),
|
||||||
|
alphaEta_(),
|
||||||
|
m_(),
|
||||||
|
Cnu_(),
|
||||||
|
Cn_(),
|
||||||
|
Ck_(),
|
||||||
|
kName_(),
|
||||||
|
deltaName_(),
|
||||||
|
nuName_(),
|
||||||
|
nutName_()
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikEtaIndex::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (!resolutionIndexModel::read(dict))
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
// (Default values from CKJ:p. 031102-3, 031102-5)
|
||||||
|
alphaEta_ = dict.getOrDefault<scalar>("alphaEta", 0.05);
|
||||||
|
m_ = dict.getOrDefault<scalar>("m", 0.5);
|
||||||
|
Cnu_ = dict.getOrDefault<scalar>("Cnu", 0.1);
|
||||||
|
Cn_ = dict.getOrDefault<scalar>("Cn", 1.0);
|
||||||
|
Ck_ = dict.getOrDefault<scalar>("Ck", 0.0376);
|
||||||
|
kName_ = dict.getOrDefault<word>("k", "k");
|
||||||
|
deltaName_ = dict.getOrDefault<word>("delta", "delta");
|
||||||
|
nuName_ = dict.getOrDefault<word>("nu", "nu");
|
||||||
|
nutName_ = dict.getOrDefault<word>("nut", "nut");
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikEtaIndex::execute()
|
||||||
|
{
|
||||||
|
// Calculate Kolmogorov and mesh length scales
|
||||||
|
tmp<volScalarField> teta = eta();
|
||||||
|
tmp<volScalarField> th = cbrt(V());
|
||||||
|
|
||||||
|
// Calculate index field
|
||||||
|
auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
// (CKJ:Eq. 9)
|
||||||
|
index = 1.0/(1.0 + alphaEta_*pow(th/teta, m_));
|
||||||
|
index.correctBoundaryConditions();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikEtaIndex::write()
|
||||||
|
{
|
||||||
|
const auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
Info<< tab << "writing field:" << index.name() << endl;
|
||||||
|
|
||||||
|
index.write();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,257 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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::resolutionIndexModels::CelikEtaIndex
|
||||||
|
|
||||||
|
Description
|
||||||
|
Computes a single-mesh resolution index according to Celik et al.'s index
|
||||||
|
using Kolmogorov length scale, which is used as a LES/DES quality/post
|
||||||
|
verification metric that does not require any experimental or DNS data.
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\Gamma_{Celik,\eta}(\mathbf{x}, t) =
|
||||||
|
\frac{1}{1 + \alpha_\eta \left(\frac{h}{\eta_{eff}}\right)^m}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
with
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\eta_{eff} = \left(\frac{\nu^3}{\epsilon}\right)^{1/4}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\epsilon = \nu_{eff} \frac{k_{sgs}}{C_k \Delta^2}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\nu_{eff} = \nu_{num} + \nu_{sgs} + \nu
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\nu_{num} = {sgn}(k_{num}) C_\nu \Delta \sqrt{k_{num}}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
where
|
||||||
|
\vartable
|
||||||
|
\Gamma_{Celik,\eta}(\mathbf{x}, t) | Celik et al.'s index [-]
|
||||||
|
\alpha_\eta | Empirical constant [-]
|
||||||
|
h | Characteristic length scale with \f$h = V^{1/3} \f$ [m]
|
||||||
|
V | Cell volume [m^3]
|
||||||
|
\eta_{eff} | Kolmogorov length scale [m]
|
||||||
|
m | Empirical exponent [-]
|
||||||
|
\nu | Kinematic viscosity of fluid [m^2/s]
|
||||||
|
\epsilon | Kinetic energy dissipation rate [m^2/s^3]
|
||||||
|
\nu_{eff} | Effective eddy viscosity [m^2/s]
|
||||||
|
\nu_{num} | Numerical eddy viscosity [m^2/s]
|
||||||
|
\nu_{sgs} | Subgrid-scale eddy viscosity [m^2/s]
|
||||||
|
k_{num} | Numerical turbulent kinetic energy [m^2/s^2]
|
||||||
|
C_\nu | Empirical constant [-]
|
||||||
|
\Delta | Filter length scale [m]
|
||||||
|
k_{sgs} | Subgrid-scale turbulent kinetic energy [m^2/s^2]
|
||||||
|
C_n | Empirical constant [-]
|
||||||
|
C_k | Empirical constant [-]
|
||||||
|
\endvartable
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
Governing equations (tag:CCY):
|
||||||
|
Celik, I. B., Cehreli Z. N., Yavuz I. (2005).
|
||||||
|
Index of resolution quality for large eddy simulations.
|
||||||
|
Journal of Fluids Engineering. 127:949–958.
|
||||||
|
DOI:10.1115/1.1990201
|
||||||
|
|
||||||
|
Governing equations (tag:CKJ):
|
||||||
|
Celik, I., Klein, M., & Janicka, J. (2009).
|
||||||
|
Assessment measures for engineering LES applications.
|
||||||
|
Journal of fluids engineering, 131(3).
|
||||||
|
DOI:10.1115/1.3059703
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
Usage
|
||||||
|
Minimal example by using \c system/controlDict.functions:
|
||||||
|
\verbatim
|
||||||
|
resolutionIndex1
|
||||||
|
{
|
||||||
|
// Inherited entries
|
||||||
|
...
|
||||||
|
model CelikEtaIndex;
|
||||||
|
|
||||||
|
// Optional entries
|
||||||
|
alphaEta <scalar>;
|
||||||
|
m <scalar>;
|
||||||
|
Cnu <scalar>;
|
||||||
|
Cn <scalar>;
|
||||||
|
Ck <scalar>;
|
||||||
|
k <word>;
|
||||||
|
delta <word>;
|
||||||
|
nu <word>;
|
||||||
|
nut <word>;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
where the entries mean:
|
||||||
|
\table
|
||||||
|
Property | Description | Type | Reqd | Deflt
|
||||||
|
model | Model name: CelikEtaIndex | word | yes | -
|
||||||
|
alphaEta | Empirical constant | scalar | no | 0.05
|
||||||
|
m | Empirical exponent | scalar | no | 0.5
|
||||||
|
Cnu | Empirical constant | scalar | no | 0.1
|
||||||
|
Cn | Empirical constant | scalar | no | 1.0
|
||||||
|
Ck | Empirical constant | scalar | no | 0.0376
|
||||||
|
k | Name of subgrid-scale turbulent kinetic energy field <!--
|
||||||
|
--> | word | no | k
|
||||||
|
delta | Name of filter field | word | no | delta
|
||||||
|
nu | Name of kinematic viscosity field | word | no | nu
|
||||||
|
nut | Name of subgrid-scale eddy viscosity field | word | no | nut
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
CelikEtaIndex.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_resolutionIndexModels_CelikEtaIndex_H
|
||||||
|
#define Foam_resolutionIndexModels_CelikEtaIndex_H
|
||||||
|
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class CelikEtaIndex Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class CelikEtaIndex
|
||||||
|
:
|
||||||
|
public resolutionIndexModel
|
||||||
|
{
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar alphaEta_;
|
||||||
|
|
||||||
|
//- Empirical exponent
|
||||||
|
scalar m_;
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar Cnu_;
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar Cn_;
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar Ck_;
|
||||||
|
|
||||||
|
//- Name of subgrid-scale turbulent kinetic energy field
|
||||||
|
word kName_;
|
||||||
|
|
||||||
|
//- Name of filter field
|
||||||
|
word deltaName_;
|
||||||
|
|
||||||
|
//- Name of kinematic viscosity field
|
||||||
|
word nuName_;
|
||||||
|
|
||||||
|
//- Name of subgrid-scale eddy viscosity field
|
||||||
|
word nutName_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Return Kolmogorov length scale field
|
||||||
|
tmp<volScalarField> eta() const;
|
||||||
|
|
||||||
|
//- Return kinetic energy dissipation rate field
|
||||||
|
tmp<volScalarField> epsilon() const;
|
||||||
|
|
||||||
|
//- Return effective eddy viscosity field
|
||||||
|
tmp<volScalarField> nuEff() const;
|
||||||
|
|
||||||
|
//- Return numerical eddy viscosity field
|
||||||
|
tmp<volScalarField> nuNum() const;
|
||||||
|
|
||||||
|
//- Return numerical turbulent kinetic energy field
|
||||||
|
tmp<volScalarField> kNum() const;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("CelikEtaIndex");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
CelikEtaIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
//- No copy construct
|
||||||
|
CelikEtaIndex(const CelikEtaIndex&) = delete;
|
||||||
|
|
||||||
|
//- No copy assignment
|
||||||
|
void operator=(const CelikEtaIndex&) = delete;
|
||||||
|
|
||||||
|
|
||||||
|
// Destructor
|
||||||
|
virtual ~CelikEtaIndex() = default;
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Read top-level dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Calculate the result field
|
||||||
|
virtual bool execute();
|
||||||
|
|
||||||
|
//- Write the result field
|
||||||
|
virtual bool write();
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace resolutionIndexModels
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,148 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "CelikNuIndex.H"
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(CelikNuIndex, 0);
|
||||||
|
addToRunTimeSelectionTable(resolutionIndexModel, CelikNuIndex, dictionary);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikNuIndex::nuNum() const
|
||||||
|
{
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
|
||||||
|
tmp<volScalarField> tkNum = kNum();
|
||||||
|
|
||||||
|
// (CKJ:Eq. 35)
|
||||||
|
return sign(tkNum.cref())*Cnu_*Delta*sqrt(tkNum.cref());
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::CelikNuIndex::kNum() const
|
||||||
|
{
|
||||||
|
const auto& kSgs = getOrReadField<volScalarField>(kName_);
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
|
||||||
|
tmp<volScalarField> th = cbrt(V());
|
||||||
|
|
||||||
|
// (CKJ:Eq. 28)
|
||||||
|
return Cn_*sqr(th/Delta)*kSgs;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::resolutionIndexModels::CelikNuIndex::CelikNuIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
resolutionIndexModel(name, mesh, dict),
|
||||||
|
alphaNu_(),
|
||||||
|
n_(),
|
||||||
|
Cnu_(),
|
||||||
|
Cn_(),
|
||||||
|
kName_(),
|
||||||
|
deltaName_(),
|
||||||
|
nuName_(),
|
||||||
|
nutName_()
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikNuIndex::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (!resolutionIndexModel::read(dict))
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
// (Default values from CKJ:p. 031102-3, 031102-5)
|
||||||
|
alphaNu_ = dict.getOrDefault<scalar>("alphaNu", 0.05);
|
||||||
|
n_ = dict.getOrDefault<scalar>("n", 0.53);
|
||||||
|
Cnu_ = dict.getOrDefault<scalar>("Cnu", 0.1);
|
||||||
|
Cn_ = dict.getOrDefault<scalar>("Cn", 1.0);
|
||||||
|
kName_ = dict.getOrDefault<word>("k", "k");
|
||||||
|
deltaName_ = dict.getOrDefault<word>("delta", "delta");
|
||||||
|
nuName_ = dict.getOrDefault<word>("nu", "nu");
|
||||||
|
nutName_ = dict.getOrDefault<word>("nut", "nut");
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikNuIndex::execute()
|
||||||
|
{
|
||||||
|
// Calculate effective eddy viscosity field
|
||||||
|
const auto& nu = getOrReadField<volScalarField>(nuName_);
|
||||||
|
const auto& nuSgs = getOrReadField<volScalarField>(nutName_);
|
||||||
|
tmp<volScalarField> tnuNum = nuNum();
|
||||||
|
tmp<volScalarField> tnuEff = tnuNum + nuSgs + nu;
|
||||||
|
|
||||||
|
// Calculate index field
|
||||||
|
auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
// (CKJ:Eq. 10)
|
||||||
|
index = 1.0/(1.0 + alphaNu_*pow(tnuEff/nu, n_));
|
||||||
|
index.correctBoundaryConditions();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::CelikNuIndex::write()
|
||||||
|
{
|
||||||
|
const auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
Info<< tab << "writing field:" << index.name() << endl;
|
||||||
|
|
||||||
|
index.write();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,246 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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::resolutionIndexModels::CelikNuIndex
|
||||||
|
|
||||||
|
Description
|
||||||
|
Computes a single-mesh resolution index according to Celik et al.'s index
|
||||||
|
using effective viscosity, which is used as a LES/DES quality/post
|
||||||
|
verification metric that does not require any experimental or DNS data.
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\Gamma_{Celik,\nu}(\mathbf{x}, t) =
|
||||||
|
\frac{1}{1 + \alpha_\nu \left(\frac{\nu_{eff}}{\nu}\right)^n}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
with
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\nu_{eff} = \nu_{num} + \nu_{sgs} + \nu
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\nu_{num} = {sgn}(k_{num}) C_\nu \Delta \sqrt{k_{num}}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
\f[
|
||||||
|
k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
where
|
||||||
|
\vartable
|
||||||
|
\Gamma_{Celik,\nu}(\mathbf{x}, t) | Celik et al.'s index [-]
|
||||||
|
\alpha_\nu | Empirical constant [-]
|
||||||
|
\nu_{eff} | Effective eddy viscosity [m^2/s]
|
||||||
|
\nu_{num} | Numerical eddy viscosity [m^2/s]
|
||||||
|
\nu_{sgs} | Subgrid-scale eddy viscosity [m^2/s]
|
||||||
|
\nu | Kinematic viscosity of fluid [m^2/s]
|
||||||
|
n | Empirical exponent [-]
|
||||||
|
k_{num} | Numerical turbulent kinetic energy [m^2/s^2]
|
||||||
|
C_\nu | Empirical constant [-]
|
||||||
|
\Delta | Filter length scale [m]
|
||||||
|
C_n | Empirical constant [-]
|
||||||
|
h | Characteristic length scale with \f$h = V^{1/3} \f$ [m]
|
||||||
|
V | Cell volume [m^3]
|
||||||
|
k_{sgs} | Subgrid-scale turbulent kinetic energy [m^2/s^2]
|
||||||
|
\endvartable
|
||||||
|
|
||||||
|
Criterion for acceptable-quality resolution:
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\Gamma_{Celik,\nu}(\mathbf{x}) \geq 0.8
|
||||||
|
\f]
|
||||||
|
|
||||||
|
Required fields:
|
||||||
|
\verbatim
|
||||||
|
k | Subgrid scale turbulent kinetic energy [m^2/s^2]
|
||||||
|
delta | Filter length [m]
|
||||||
|
nu | Kinematic viscosity [m^2/s]
|
||||||
|
nut | Subgrid-scale eddy viscosity [m^2/s]
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
Governing equations (tag:CCY):
|
||||||
|
Celik, I. B., Cehreli Z. N., Yavuz I. (2005).
|
||||||
|
Index of resolution quality for large eddy simulations.
|
||||||
|
Journal of Fluids Engineering. 127:949–958.
|
||||||
|
DOI:10.1115/1.1990201
|
||||||
|
|
||||||
|
Governing equations (tag:CKJ):
|
||||||
|
Celik, I., Klein, M., & Janicka, J. (2009).
|
||||||
|
Assessment measures for engineering LES applications.
|
||||||
|
Journal of fluids engineering, 131(3).
|
||||||
|
DOI:10.1115/1.3059703
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
Usage
|
||||||
|
Minimal example by using \c system/controlDict.functions:
|
||||||
|
\verbatim
|
||||||
|
resolutionIndex1
|
||||||
|
{
|
||||||
|
// Inherited entries
|
||||||
|
...
|
||||||
|
model CelikNuIndex;
|
||||||
|
|
||||||
|
// Optional entries
|
||||||
|
alphaNu <scalar>;
|
||||||
|
n <scalar>;
|
||||||
|
Cnu <scalar>;
|
||||||
|
Cn <scalar>;
|
||||||
|
k <word>;
|
||||||
|
delta <word>;
|
||||||
|
nu <word>;
|
||||||
|
nut <word>;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
where the entries mean:
|
||||||
|
\table
|
||||||
|
Property | Description | Type | Reqd | Deflt
|
||||||
|
model | Model name: CelikNuIndex | word | yes | -
|
||||||
|
alphaNu | Empirical constant | scalar | no | 0.05
|
||||||
|
n | Empirical exponent | scalar | no | 0.53
|
||||||
|
Cnu | Empirical constant | scalar | no | 0.1
|
||||||
|
Cn | Empirical constant | scalar | no | 1.0
|
||||||
|
k | Name of subgrid-scale turbulent kinetic energy field <!--
|
||||||
|
--> | word | no | k
|
||||||
|
delta | Name of filter field | word | no | delta
|
||||||
|
nu | Name of kinematic viscosity field | word | no | nu
|
||||||
|
nut | Name of subgrid-scale eddy viscosity field | word | no | nut
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
CelikNuIndex.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_resolutionIndexModels_CelikNuIndex_H
|
||||||
|
#define Foam_resolutionIndexModels_CelikNuIndex_H
|
||||||
|
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class CelikNuIndex Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class CelikNuIndex
|
||||||
|
:
|
||||||
|
public resolutionIndexModel
|
||||||
|
{
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar alphaNu_;
|
||||||
|
|
||||||
|
//- Empirical exponent
|
||||||
|
scalar n_;
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar Cnu_;
|
||||||
|
|
||||||
|
//- Empirical constant
|
||||||
|
scalar Cn_;
|
||||||
|
|
||||||
|
//- Name of subgrid-scale turbulent kinetic energy field
|
||||||
|
word kName_;
|
||||||
|
|
||||||
|
//- Name of filter field
|
||||||
|
word deltaName_;
|
||||||
|
|
||||||
|
//- Name of kinematic viscosity field
|
||||||
|
word nuName_;
|
||||||
|
|
||||||
|
//- Name of subgrid-scale eddy viscosity field
|
||||||
|
word nutName_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Return numerical eddy viscosity field
|
||||||
|
tmp<volScalarField> nuNum() const;
|
||||||
|
|
||||||
|
//- Return numerical turbulent kinetic energy field
|
||||||
|
tmp<volScalarField> kNum() const;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("CelikNuIndex");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
CelikNuIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
//- No copy construct
|
||||||
|
CelikNuIndex(const CelikNuIndex&) = delete;
|
||||||
|
|
||||||
|
//- No copy assignment
|
||||||
|
void operator=(const CelikNuIndex&) = delete;
|
||||||
|
|
||||||
|
|
||||||
|
// Destructor
|
||||||
|
virtual ~CelikNuIndex() = default;
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Read top-level dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Calculate the result field
|
||||||
|
virtual bool execute();
|
||||||
|
|
||||||
|
//- Write the result field
|
||||||
|
virtual bool write();
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace resolutionIndexModels
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,146 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "PopeIndex.H"
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
#include "addToRunTimeSelectionTable.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(PopeIndex, 0);
|
||||||
|
addToRunTimeSelectionTable(resolutionIndexModel, PopeIndex, dictionary);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField>
|
||||||
|
Foam::resolutionIndexModels::PopeIndex::kNum() const
|
||||||
|
{
|
||||||
|
const auto& kSgs = getOrReadField<volScalarField>(kName_);
|
||||||
|
const auto& Delta = getOrReadField<volScalarField>(deltaName_);
|
||||||
|
|
||||||
|
tmp<volScalarField> th = cbrt(V());
|
||||||
|
|
||||||
|
// (CKJ:Eq. 28)
|
||||||
|
return Cn_*sqr(th/Delta)*kSgs;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::resolutionIndexModels::PopeIndex::PopeIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
resolutionIndexModel(name, mesh, dict),
|
||||||
|
includeKnum_(),
|
||||||
|
Cn_(),
|
||||||
|
UName_(),
|
||||||
|
UMeanName_(),
|
||||||
|
kName_(),
|
||||||
|
deltaName_()
|
||||||
|
{
|
||||||
|
read(dict);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::PopeIndex::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
if (!resolutionIndexModel::read(dict))
|
||||||
|
{
|
||||||
|
return false;
|
||||||
|
}
|
||||||
|
|
||||||
|
includeKnum_ = dict.getOrDefault<bool>("includeKnum", true);
|
||||||
|
if (includeKnum_)
|
||||||
|
{
|
||||||
|
Cn_ = dict.getOrDefault<scalar>("Cnu", 1.0);
|
||||||
|
}
|
||||||
|
UName_ = dict.getOrDefault<word>("U", "U");
|
||||||
|
UMeanName_ = dict.getOrDefault<word>("UMean", "UMean");
|
||||||
|
kName_ = dict.getOrDefault<word>("k", "k");
|
||||||
|
deltaName_ = dict.getOrDefault<word>("delta", "delta");
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::PopeIndex::execute()
|
||||||
|
{
|
||||||
|
// Calculate resolved k field
|
||||||
|
const auto& U = getOrReadField<volVectorField>(UName_);
|
||||||
|
const auto& UMean = getOrReadField<volVectorField>(UMeanName_);
|
||||||
|
const volScalarField kRes(0.5*magSqr(U - UMean));
|
||||||
|
|
||||||
|
// Calculate subgrid-scale k field
|
||||||
|
const auto& kSgs = getOrReadField<volScalarField>(kName_);
|
||||||
|
|
||||||
|
// Calculate total k field
|
||||||
|
tmp<volScalarField> tkTot = kRes + kSgs;
|
||||||
|
if (includeKnum_)
|
||||||
|
{
|
||||||
|
tkTot.ref() += mag(kNum());
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// Calculate index field
|
||||||
|
auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
const dimensionedScalar kMin(kSgs.dimensions(), SMALL);
|
||||||
|
|
||||||
|
// (P:p. 560;CKJ:Eq. 11)
|
||||||
|
index = kRes/max(kMin, tkTot);
|
||||||
|
index.correctBoundaryConditions();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModels::PopeIndex::write()
|
||||||
|
{
|
||||||
|
const auto& index = getOrReadField<volScalarField>(resultName());
|
||||||
|
|
||||||
|
Info<< tab << "writing field:" << index.name() << endl;
|
||||||
|
|
||||||
|
index.write();
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,237 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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::resolutionIndexModels::PopeIndex
|
||||||
|
|
||||||
|
Description
|
||||||
|
Computes a single-mesh resolution index according to Pope's index,
|
||||||
|
which is used as a LES/DES quality/post-verification metric that does
|
||||||
|
not require any experimental or DNS data.
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\Gamma_{Pope}(\mathbf{x}, t) = \frac{k_{res}}{k_{tot}}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
with
|
||||||
|
|
||||||
|
\f[
|
||||||
|
k_{tot} = k_{res} + k_{sgs} + |k_{num}|
|
||||||
|
\f]
|
||||||
|
|
||||||
|
where
|
||||||
|
\vartable
|
||||||
|
\Gamma_{Pope}(\mathbf{x}, t) | Pope's index [-]
|
||||||
|
k_{tot} | Total turbulent kinetic energy [m^2/s^2]
|
||||||
|
k_{res} | Resolved turbulent kinetic energy [m^2/s^2]
|
||||||
|
k_{sgs} | Subgrid-scale turbulent kinetic energy [m^2/s^2]
|
||||||
|
k_{num} | Numerical turbulent kinetic energy [m^2/s^2]
|
||||||
|
\endvartable
|
||||||
|
|
||||||
|
Inclusion of \f$k_{num}\f$ is optional, and set as \c true by default:
|
||||||
|
|
||||||
|
\f[
|
||||||
|
k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
|
||||||
|
\f]
|
||||||
|
|
||||||
|
where
|
||||||
|
\vartable
|
||||||
|
C_n | Empirical constant [-]
|
||||||
|
h | Characteristic length scale with \f$h = V^{1/3} \f$ [m]
|
||||||
|
V | Cell volume [m^3]
|
||||||
|
\Delta | Filter length scale [m]
|
||||||
|
\endvartable
|
||||||
|
|
||||||
|
Typical criterion for acceptable-quality resolution:
|
||||||
|
|
||||||
|
\f[
|
||||||
|
\Gamma_{Pope}(\mathbf{x}) \geq 0.8
|
||||||
|
\f]
|
||||||
|
|
||||||
|
Required fields:
|
||||||
|
\verbatim
|
||||||
|
U | Velocity [m/s]
|
||||||
|
UMean | Mean velocity [m/s]
|
||||||
|
k | Subgrid-scale turbulent kinetic energy [m^2/s^2]
|
||||||
|
delta | Filter length [m]
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
References:
|
||||||
|
\verbatim
|
||||||
|
Governing equation (tag:P):
|
||||||
|
Pope, S. B. (2000).
|
||||||
|
Turbulent flows.
|
||||||
|
Cambridge, UK: Cambridge Univ. Press
|
||||||
|
DOI:10.1017/CBO9780511840531
|
||||||
|
|
||||||
|
Governing equations for the denominator kNum term (tag:CKJ):
|
||||||
|
Celik, I., Klein, M., & Janicka, J. (2009).
|
||||||
|
Assessment measures for engineering LES applications.
|
||||||
|
Journal of fluids engineering, 131(3).
|
||||||
|
DOI:10.1115/1.3059703
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
Usage
|
||||||
|
Minimal example by using \c system/controlDict.functions:
|
||||||
|
\verbatim
|
||||||
|
resolutionIndex1
|
||||||
|
{
|
||||||
|
// Inherited entries
|
||||||
|
...
|
||||||
|
model PopeIndex;
|
||||||
|
|
||||||
|
// Optional entries
|
||||||
|
U <word>;
|
||||||
|
UMean <word>;
|
||||||
|
k <word>;
|
||||||
|
delta <word>;
|
||||||
|
includeKnum <bool>;
|
||||||
|
|
||||||
|
// Conditional entries
|
||||||
|
// when includeKnum = true
|
||||||
|
Cn <scalar>;
|
||||||
|
}
|
||||||
|
\endverbatim
|
||||||
|
|
||||||
|
where the entries mean:
|
||||||
|
\table
|
||||||
|
Property | Description | Type | Reqd | Deflt
|
||||||
|
model | Model name: PopeIndex | word | yes | -
|
||||||
|
U | Name of velocity field | word | no | U
|
||||||
|
UMean | Name of mean velocity field | word | no | UMean
|
||||||
|
k | Name of subgrid-scale turbulent kinetic energy field <!--
|
||||||
|
--> | word | no | k
|
||||||
|
delta | Name of filter field | word | no | delta
|
||||||
|
includeKnum | Flag to include numerical k field | bool | no | true
|
||||||
|
Cn | Empirical constant | choice | no | 1.0
|
||||||
|
\endtable
|
||||||
|
|
||||||
|
Note
|
||||||
|
- Some studies measured \f$\Gamma_{Pope} > 1\f$ with \f$k_{res}\f$ comparisons
|
||||||
|
between a LES and a corresponding filtered DNS. Nonphysical results may
|
||||||
|
occur.
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
PopeIndex.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_resolutionIndexModels_PopeIndex_H
|
||||||
|
#define Foam_resolutionIndexModels_PopeIndex_H
|
||||||
|
|
||||||
|
#include "resolutionIndexModel.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
namespace resolutionIndexModels
|
||||||
|
{
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class PopeIndex Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class PopeIndex
|
||||||
|
:
|
||||||
|
public resolutionIndexModel
|
||||||
|
{
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Flag to include numerical turbulent kinetic energy field
|
||||||
|
bool includeKnum_;
|
||||||
|
|
||||||
|
//- Empirical constant in numerical turbulent kinetic energy estimation
|
||||||
|
scalar Cn_;
|
||||||
|
|
||||||
|
//- Name of velocity field
|
||||||
|
word UName_;
|
||||||
|
|
||||||
|
//- Name of mean velocity field
|
||||||
|
word UMeanName_;
|
||||||
|
|
||||||
|
//- Name of subgrid-scale turbulent kinetic energy field
|
||||||
|
word kName_;
|
||||||
|
|
||||||
|
//- Name of filter field
|
||||||
|
word deltaName_;
|
||||||
|
|
||||||
|
|
||||||
|
// Private Member Functions
|
||||||
|
|
||||||
|
//- Return numerical turbulent kinetic energy field
|
||||||
|
tmp<volScalarField> kNum() const;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("PopeIndex");
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
PopeIndex
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
//- No copy construct
|
||||||
|
PopeIndex(const PopeIndex&) = delete;
|
||||||
|
|
||||||
|
//- No copy assignment
|
||||||
|
void operator=(const PopeIndex&) = delete;
|
||||||
|
|
||||||
|
|
||||||
|
// Destructor
|
||||||
|
virtual ~PopeIndex() = default;
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
//- Read top-level dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Calculate the result field
|
||||||
|
virtual bool execute();
|
||||||
|
|
||||||
|
//- Write the result field
|
||||||
|
virtual bool write();
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace resolutionIndexModels
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,117 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "resolutionIndexModel.H"
|
||||||
|
#include "fvMesh.H"
|
||||||
|
#include "ListOps.H"
|
||||||
|
#include "turbulenceFields.H"
|
||||||
|
#include "turbulenceModel.H"
|
||||||
|
#include "zeroGradientFvPatchFields.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
defineTypeNameAndDebug(resolutionIndexModel, 0);
|
||||||
|
defineRunTimeSelectionTable(resolutionIndexModel, dictionary);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::tmp<Foam::volScalarField> Foam::resolutionIndexModel::V() const
|
||||||
|
{
|
||||||
|
auto tV = tmp<volScalarField>::New
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
"V",
|
||||||
|
mesh_.time().timeName(),
|
||||||
|
mesh_,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::NO_WRITE,
|
||||||
|
false
|
||||||
|
),
|
||||||
|
mesh_,
|
||||||
|
dimVolume,
|
||||||
|
zeroGradientFvPatchScalarField::typeName
|
||||||
|
);
|
||||||
|
|
||||||
|
tV.ref().primitiveFieldRef() = mesh_.V();
|
||||||
|
tV.ref().correctBoundaryConditions();
|
||||||
|
|
||||||
|
return tV;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::resolutionIndexModel::resolutionIndexModel
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
:
|
||||||
|
mesh_(mesh),
|
||||||
|
resultName_(word::null)
|
||||||
|
{}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
bool Foam::resolutionIndexModel::read(const dictionary& dict)
|
||||||
|
{
|
||||||
|
resultName_ = dict.getOrDefault("result", type());
|
||||||
|
|
||||||
|
auto* indexPtr = mesh_.getObjectPtr<volScalarField>(resultName_);
|
||||||
|
|
||||||
|
if (!indexPtr)
|
||||||
|
{
|
||||||
|
indexPtr = new volScalarField
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
resultName_,
|
||||||
|
mesh_.time().timeName(),
|
||||||
|
mesh_,
|
||||||
|
IOobject::READ_IF_PRESENT,
|
||||||
|
IOobject::NO_WRITE
|
||||||
|
),
|
||||||
|
mesh_,
|
||||||
|
dimensionedScalar(dimless, Zero),
|
||||||
|
zeroGradientFvPatchScalarField::typeName
|
||||||
|
);
|
||||||
|
|
||||||
|
mesh_.objectRegistry::store(indexPtr);
|
||||||
|
}
|
||||||
|
|
||||||
|
return true;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,183 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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/>.
|
||||||
|
|
||||||
|
Namespace
|
||||||
|
Foam::resolutionIndexModels
|
||||||
|
|
||||||
|
Description
|
||||||
|
A namespace for various \c resolutionIndex model implementations.
|
||||||
|
|
||||||
|
Class
|
||||||
|
Foam::resolutionIndexModel
|
||||||
|
|
||||||
|
Description
|
||||||
|
A base class for \c resolutionIndex models.
|
||||||
|
|
||||||
|
SourceFiles
|
||||||
|
resolutionIndexModel.C
|
||||||
|
resolutionIndexModelNew.C
|
||||||
|
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
#ifndef Foam_resolutionIndexModel_H
|
||||||
|
#define Foam_resolutionIndexModel_H
|
||||||
|
|
||||||
|
#include "volFields.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
namespace Foam
|
||||||
|
{
|
||||||
|
|
||||||
|
// Forward Declarations
|
||||||
|
class fvMesh;
|
||||||
|
|
||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
Class resolutionIndexModel Declaration
|
||||||
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
class resolutionIndexModel
|
||||||
|
{
|
||||||
|
// Private Data
|
||||||
|
|
||||||
|
//- Const reference to the mesh
|
||||||
|
const fvMesh& mesh_;
|
||||||
|
|
||||||
|
//- Name of result field
|
||||||
|
word resultName_;
|
||||||
|
|
||||||
|
|
||||||
|
protected:
|
||||||
|
|
||||||
|
// Protected Member Functions
|
||||||
|
|
||||||
|
//- Return requested field from the object registry
|
||||||
|
//- or read+register the field to the object registry
|
||||||
|
template<class GeoFieldType>
|
||||||
|
GeoFieldType& getOrReadField(const word& fieldName) const;
|
||||||
|
|
||||||
|
//- Return cell volume field
|
||||||
|
tmp<volScalarField> V() const;
|
||||||
|
|
||||||
|
|
||||||
|
public:
|
||||||
|
|
||||||
|
//- Runtime type information
|
||||||
|
TypeName("resolutionIndexModel");
|
||||||
|
|
||||||
|
|
||||||
|
// Declare runtime constructor selection table
|
||||||
|
|
||||||
|
declareRunTimeSelectionTable
|
||||||
|
(
|
||||||
|
autoPtr,
|
||||||
|
resolutionIndexModel,
|
||||||
|
dictionary,
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
),
|
||||||
|
(name, mesh, dict)
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Selectors
|
||||||
|
|
||||||
|
//- Return a reference to the selected resolutionIndex model
|
||||||
|
static autoPtr<resolutionIndexModel> New
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
|
// Constructors
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
resolutionIndexModel
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
);
|
||||||
|
|
||||||
|
//- No copy construct
|
||||||
|
resolutionIndexModel(const resolutionIndexModel&) = delete;
|
||||||
|
|
||||||
|
//- No copy assignment
|
||||||
|
void operator=(const resolutionIndexModel&) = delete;
|
||||||
|
|
||||||
|
|
||||||
|
//- Destructor
|
||||||
|
virtual ~resolutionIndexModel() = default;
|
||||||
|
|
||||||
|
|
||||||
|
// Member Functions
|
||||||
|
|
||||||
|
// Access
|
||||||
|
|
||||||
|
//- Return const reference to the mesh
|
||||||
|
const fvMesh& mesh() const noexcept
|
||||||
|
{
|
||||||
|
return mesh_;
|
||||||
|
}
|
||||||
|
|
||||||
|
//- Return const reference to the result name
|
||||||
|
const word& resultName() const noexcept
|
||||||
|
{
|
||||||
|
return resultName_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// I-O
|
||||||
|
|
||||||
|
//- Read top-level dictionary
|
||||||
|
virtual bool read(const dictionary& dict);
|
||||||
|
|
||||||
|
//- Calculate the result field
|
||||||
|
virtual bool execute() = 0;
|
||||||
|
|
||||||
|
//- Write the result field
|
||||||
|
virtual bool write() = 0;
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
} // End namespace Foam
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#ifdef NoRepository
|
||||||
|
#include "resolutionIndexModelTemplates.C"
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
#endif
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,61 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "resolutionIndexModel.H"
|
||||||
|
#include "fvMesh.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::autoPtr<Foam::resolutionIndexModel> Foam::resolutionIndexModel::New
|
||||||
|
(
|
||||||
|
const word& name,
|
||||||
|
const fvMesh& mesh,
|
||||||
|
const dictionary& dict
|
||||||
|
)
|
||||||
|
{
|
||||||
|
const word modelType(dict.get<word>("model"));
|
||||||
|
|
||||||
|
Info<< " Selecting model: " << modelType << nl << endl;
|
||||||
|
|
||||||
|
auto* ctorPtr = dictionaryConstructorTable(modelType);
|
||||||
|
|
||||||
|
if (!ctorPtr)
|
||||||
|
{
|
||||||
|
FatalIOErrorInLookup
|
||||||
|
(
|
||||||
|
dict,
|
||||||
|
"resolutionIndexModel",
|
||||||
|
modelType,
|
||||||
|
*dictionaryConstructorTablePtr_
|
||||||
|
) << exit(FatalIOError);
|
||||||
|
}
|
||||||
|
|
||||||
|
return autoPtr<resolutionIndexModel>(ctorPtr(name, mesh, dict));
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
@ -0,0 +1,62 @@
|
|||||||
|
/*---------------------------------------------------------------------------*\
|
||||||
|
========= |
|
||||||
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
|
\\ / O peration |
|
||||||
|
\\ / A nd | www.openfoam.com
|
||||||
|
\\/ M anipulation |
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
Copyright (C) 2022 OpenCFD Ltd.
|
||||||
|
-------------------------------------------------------------------------------
|
||||||
|
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 "resolutionIndexModel.H"
|
||||||
|
#include "volFields.H"
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class GeoFieldType>
|
||||||
|
GeoFieldType& Foam::resolutionIndexModel::getOrReadField
|
||||||
|
(
|
||||||
|
const word& fieldName
|
||||||
|
) const
|
||||||
|
{
|
||||||
|
auto* ptr = mesh_.getObjectPtr<GeoFieldType>(fieldName);
|
||||||
|
|
||||||
|
if (!ptr)
|
||||||
|
{
|
||||||
|
ptr = new GeoFieldType
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
fieldName,
|
||||||
|
mesh_.time().timeName(),
|
||||||
|
mesh_,
|
||||||
|
IOobject::MUST_READ,
|
||||||
|
IOobject::AUTO_WRITE
|
||||||
|
),
|
||||||
|
mesh_
|
||||||
|
);
|
||||||
|
mesh_.objectRegistry::store(ptr);
|
||||||
|
}
|
||||||
|
|
||||||
|
return *ptr;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// ************************************************************************* //
|
||||||
Reference in New Issue
Block a user