The turbulentDFSEMInlet
is a synthesised-eddy based velocity inlet boundary condition to generate synthetic turbulence-alike time-series from a given set of turbulence statistics for LES and hybrid RANS-LES computations.
More...
Public Member Functions | |
TypeName ("turbulentDFSEMInlet") | |
Runtime type information. More... | |
turbulentDFSEMInletFvPatchVectorField (const fvPatch &, const DimensionedField< vector, volMesh > &) | |
Construct from patch and internal field. More... | |
turbulentDFSEMInletFvPatchVectorField (const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &) | |
Construct from patch, internal field and dictionary. More... | |
turbulentDFSEMInletFvPatchVectorField (const turbulentDFSEMInletFvPatchVectorField &, const fvPatch &, const DimensionedField< vector, volMesh > &, const fvPatchFieldMapper &) | |
turbulentDFSEMInletFvPatchVectorField (const turbulentDFSEMInletFvPatchVectorField &) | |
Construct as copy. More... | |
virtual tmp< fvPatchVectorField > | clone () const |
Construct and return a clone. More... | |
turbulentDFSEMInletFvPatchVectorField (const turbulentDFSEMInletFvPatchVectorField &, const DimensionedField< vector, volMesh > &) | |
Construct as copy setting internal field reference. More... | |
virtual tmp< fvPatchVectorField > | clone (const DimensionedField< vector, volMesh > &iF) const |
Construct and return a clone setting internal field reference. More... | |
virtual | ~turbulentDFSEMInletFvPatchVectorField ()=default |
Destructor. More... | |
virtual void | autoMap (const fvPatchFieldMapper &m) |
Map (and resize as needed) from self given a mapping object. More... | |
virtual void | rmap (const fvPatchVectorField &ptf, const labelList &addr) |
Reverse map the given fvPatchField onto this fvPatchField. More... | |
virtual void | updateCoeffs () |
Update the coefficients associated with the patch field. More... | |
virtual void | write (Ostream &) const |
Write. More... | |
Static Public Member Functions | |
static void | checkStresses (const symmTensorField &R) |
Check if input Reynolds stresses are valid. More... | |
static void | checkStresses (const scalarField &R) |
Check if input Reynolds stresses are valid. More... | |
The turbulentDFSEMInlet
is a synthesised-eddy based velocity inlet boundary condition to generate synthetic turbulence-alike time-series from a given set of turbulence statistics for LES and hybrid RANS-LES computations.
Reference:
Standard model (tag:PCR): Poletto, R., Craft, T., & Revell, A. (2013). A new divergence free synthetic eddy method for the reproduction of inlet flow conditions for LES. Flow, turbulence and combustion, 91(3), 519-539. DOI:10.1007/s10494-013-9488-2 Expression for the average length scale (tag:SST): Shur, M., Strelets, M., Travin, A., Probst, A., Probst, S., Schwamborn, D., ... & Revell, A. (2018). Improved embedded approaches. In: Mockett C., Haase W., Schwamborn D. (eds) Go4Hybrid: Grey area mitigation for hybrid RANS-LES methods. Notes on Numerical Fluid Mechanics and Multidisciplinary Design. p. 51-87. Springer, Cham. DOI:10.1007/978-3-319-52995-0_3
<patchName> { // Mandatory entries type turbulentDFSEMInlet; delta <scalar>; R <PatchFunction1>; U <PatchFunction1>; L <PatchFunction1>; // e.g. // R uniform (<Rxx> <Rxy> <Rxz> <Ryy> <Ryz> <Rzz>); // U uniform (<Ux> <Uy> <Uz>); // L uniform <L>; // Optional entries d <scalar>; nCellPerEddy <label>; kappa <scalar>; Uref <scalar>; Lref <scalar>; scale <scalar>; m <scalar>; writeEddies <bool>; // Inherited entries ... }
where the entries mean:
Property | Description | Type | Reqd | Deflt |
---|---|---|---|---|
type | Type name: turbulentDFSEMInlet | word | yes | - |
delta | Characteristic length scale | scalar | yes | - |
R | Reynolds-stress tensor field | PatchFunction1<symmTensor> | yes | - |
U | Mean velocity field | PatchFunction1<vector> | yes | - |
L | Integral-length scale field | PatchFunction1<scalar> | yes | - |
d | Ratio of sum of eddy volumes to eddy box volume i.e. eddy density (fill fraction) | scalar | no | 1.0 |
nCellPerEddy | Minimum eddy length in units of number of cells | label | no | 5 |
kappa | von Karman constant | scalar | no | 0.41 |
Uref | Normalisation factor for Reynolds-stress tensor and mean velocity | scalar | no | 1.0 |
Lref | Normalisation factor for integral-length scale | scalar | no | 1.0 |
scale | Heuristic scaling factor being applied on the normalisation factor C1 | scalar | no | 1.0 |
m | The power of V defined in C1 | scalar | no | 0.5 |
writeEddies | Flag to write eddies as OBJ file | bool | no | false |
The inherited entries are elaborated in:
delta
value typically represents the characteristic scale of flow or flow domain, e.g. a channel half height for plane channel flows.nCellPerEddy
and scale
entries are heuristic entries which do not exist in the standard method, yet are necessary to reproduce the results published in the original journal paper.C1
in Eq. 11 is not dimensionless. It is not clear whether this dimensionality issue was intentional. To alleviate this matter, users can alter the input entry m
, which is the power of the eddy-box volume defined in the C1
, to obtain a dimensionless C1
coefficient. The default value of m
is 0.5, which is the value stated in the original journal paper, and m=1/3
leads to a dimensionless C1
.Definition at line 258 of file turbulentDFSEMInletFvPatchVectorField.H.
turbulentDFSEMInletFvPatchVectorField | ( | const fvPatch & | p, |
const DimensionedField< vector, volMesh > & | iF | ||
) |
Construct from patch and internal field.
Definition at line 603 of file turbulentDFSEMInletFvPatchVectorField.C.
turbulentDFSEMInletFvPatchVectorField | ( | const fvPatch & | p, |
const DimensionedField< vector, volMesh > & | iF, | ||
const dictionary & | dict | ||
) |
Construct from patch, internal field and dictionary.
Definition at line 685 of file turbulentDFSEMInletFvPatchVectorField.C.
References turbulentDFSEMInletFvPatchVectorField::checkStresses(), eddy::debug, R, and Foam::sqr().
turbulentDFSEMInletFvPatchVectorField | ( | const turbulentDFSEMInletFvPatchVectorField & | ptf, |
const fvPatch & | p, | ||
const DimensionedField< vector, volMesh > & | iF, | ||
const fvPatchFieldMapper & | mapper | ||
) |
Construct by mapping given turbulentDFSEMInletFvPatchVectorField onto a new patch
Definition at line 643 of file turbulentDFSEMInletFvPatchVectorField.C.
Construct as copy.
Definition at line 733 of file turbulentDFSEMInletFvPatchVectorField.C.
turbulentDFSEMInletFvPatchVectorField | ( | const turbulentDFSEMInletFvPatchVectorField & | ptf, |
const DimensionedField< vector, volMesh > & | iF | ||
) |
Construct as copy setting internal field reference.
Definition at line 772 of file turbulentDFSEMInletFvPatchVectorField.C.
|
virtualdefault |
Destructor.
TypeName | ( | "turbulentDFSEMInlet" | ) |
Runtime type information.
|
inlinevirtual |
Construct and return a clone.
Definition at line 429 of file turbulentDFSEMInletFvPatchVectorField.H.
|
inlinevirtual |
Construct and return a clone setting internal field reference.
Definition at line 445 of file turbulentDFSEMInletFvPatchVectorField.H.
|
static |
Check if input Reynolds stresses are valid.
Realizability conditions (tag:S): Schumann, U. (1977). Realizability of Reynoldsâstress turbulence models. The Physics of Fluids, 20(5), 721-725. DOI:10.1063/1.861942
Definition at line 814 of file turbulentDFSEMInletFvPatchVectorField.C.
References Foam::det(), Foam::diff(), Foam::endl(), forAll, Foam::Info, R, Foam::sqr(), WarningInFunction, SymmTensor< Cmpt >::xx(), SymmTensor< Cmpt >::xy(), and SymmTensor< Cmpt >::yy().
Referenced by turbulentDFSEMInletFvPatchVectorField::turbulentDFSEMInletFvPatchVectorField(), and turbulentDigitalFilterInletFvPatchField< Type >::turbulentDigitalFilterInletFvPatchField().
|
static |
Check if input Reynolds stresses are valid.
Definition at line 872 of file turbulentDFSEMInletFvPatchVectorField.C.
References Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::min(), and R.
|
virtual |
Map (and resize as needed) from self given a mapping object.
Definition at line 887 of file turbulentDFSEMInletFvPatchVectorField.C.
References atmBoundaryLayer::autoMap().
|
virtual |
Reverse map the given fvPatchField onto this fvPatchField.
Definition at line 911 of file turbulentDFSEMInletFvPatchVectorField.C.
References atmBoundaryLayer::rmap().
|
virtual |
Update the coefficients associated with the patch field.
Definition at line 939 of file turbulentDFSEMInletFvPatchVectorField.C.
References Foam::endl(), forAll, Foam::gMax(), Foam::gMin(), Foam::gSum(), Foam::Info, n, PstreamBuffers::nProcs(), UPstream::parRun(), Foam::pow(), Foam::returnReduce(), UList< T >::size(), Foam::sqrt(), timeIndex, U, and UMean().
|
virtual |
Write.
Definition at line 1055 of file turbulentDFSEMInletFvPatchVectorField.C.
References os(), ObukhovLength::write(), Ostream::writeEntry(), and Ostream::writeEntryIfDifferent().