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. More...
#include <PopeIndex.H>


Public Member Functions | |
| TypeName ("PopeIndex") | |
| Runtime type information. | |
| PopeIndex (const word &name, const fvMesh &mesh, const dictionary &dict) | |
| Construct from components. | |
| PopeIndex (const PopeIndex &)=delete | |
| No copy construct. | |
| void | operator= (const PopeIndex &)=delete |
| No copy assignment. | |
| virtual | ~PopeIndex ()=default |
| virtual bool | read (const dictionary &dict) |
| Read the function-object dictionary. | |
| virtual bool | execute () |
| Execute the function-object operations. | |
| virtual bool | write () |
| Write the function-object results. | |
| Public Member Functions inherited from resolutionIndexModel | |
| TypeName ("resolutionIndexModel") | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, resolutionIndexModel, dictionary,(const word &name, const fvMesh &mesh, const dictionary &dict),(name, mesh, dict)) | |
| resolutionIndexModel (const word &name, const fvMesh &mesh, const dictionary &dict) | |
| Construct from components. | |
| resolutionIndexModel (const resolutionIndexModel &)=delete | |
| No copy construct. | |
| void | operator= (const resolutionIndexModel &)=delete |
| No copy assignment. | |
| virtual | ~resolutionIndexModel ()=default |
| Destructor. | |
| const fvMesh & | mesh () const noexcept |
| Return const reference to the mesh. | |
| const word & | resultName () const noexcept |
| Return const reference to the result name. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from resolutionIndexModel | |
| static autoPtr< resolutionIndexModel > | New (const word &name, const fvMesh &mesh, const dictionary &dict) |
| Return a reference to the selected resolutionIndex model. | |
| Protected Member Functions inherited from resolutionIndexModel | |
| template<class GeoFieldType> | |
| GeoFieldType & | getOrReadField (const word &fieldName) const |
| Return requested field from the object registry or read+register the field to the object registry. | |
| tmp< volScalarField > | V () const |
| Return cell volume field. | |
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.
![\[ \Gamma_{Pope}(\mathbf{x}, t) = \frac{k_{res}}{k_{tot}}
\]](form_361.png)
with
![\[ k_{tot} = k_{res} + k_{sgs} + |k_{num}|
\]](form_362.png)
where
![]() | = | Pope's index [-] |
![]() | = | Total turbulent kinetic energy [m^2/s^2] |
![]() | = | Resolved turbulent kinetic energy [m^2/s^2] |
![]() | = | Subgrid-scale turbulent kinetic energy [m^2/s^2] |
![]() | = | Numerical turbulent kinetic energy [m^2/s^2] |
Inclusion of 
true by default:
![\[ k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
\]](form_369.png)
where
![]() | = | Empirical constant [-] |
![]() | = | Characteristic length scale with ![]() |
![]() | = | Cell volume [m^3] |
![]() | = | Filter length scale [m] |
Typical criterion for acceptable-quality resolution:
![\[ \Gamma_{Pope}(\mathbf{x}) \geq 0.8
\]](form_373.png)
Required fields:
U | Velocity [m/s]
UMean | Mean velocity [m/s]
k | Subgrid-scale turbulent kinetic energy [m^2/s^2]
delta | Filter length [m]
References:
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
system/controlDict.functions: resolutionIndexFO
{
// Inherited entries
...
model PopeIndex;
// Optional entries
U <word>;
UMean <word>;
k <word>;
delta <word>;
includeKnum <bool>;
// Conditional entries
// when includeKnum = true
Cn <scalar>;
}
where the entries mean:
| 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 |


Definition at line 238 of file PopeIndex.H.
| PopeIndex | ( | const word & | name, |
| const fvMesh & | mesh, | ||
| const dictionary & | dict ) |
Construct from components.
Definition at line 54 of file PopeIndex.C.
References dict, resolutionIndexModel::mesh(), Foam::name(), read(), and resolutionIndexModel::resolutionIndexModel().
Referenced by operator=(), and PopeIndex().


|
delete |
| TypeName | ( | "PopeIndex" | ) |
Runtime type information.
References dict, resolutionIndexModel::mesh(), and Foam::name().

|
delete |
|
virtual |
Read the function-object dictionary.
Reimplemented from resolutionIndexModel.
Definition at line 75 of file PopeIndex.C.
References dict, and resolutionIndexModel::read().
Referenced by PopeIndex().


|
virtual |
Execute the function-object operations.
Implements resolutionIndexModel.
Definition at line 96 of file PopeIndex.C.
References resolutionIndexModel::getOrReadField(), Foam::mag(), Foam::magSqr(), Foam::max(), tmp< T >::ref(), resolutionIndexModel::resultName(), U, and UMean().

|
virtual |
Write the function-object results.
Implements resolutionIndexModel.
Definition at line 127 of file PopeIndex.C.
References Foam::endl(), resolutionIndexModel::getOrReadField(), Foam::Info, resolutionIndexModel::resultName(), and Foam::tab.
