turbulentDFSEMInletFvPatchVectorField Class Reference

Velocity boundary condition including synthesised eddies for use with LES and DES turbulent flows. More...

Inheritance diagram for turbulentDFSEMInletFvPatchVectorField:
[legend]
Collaboration diagram for turbulentDFSEMInletFvPatchVectorField:
[legend]

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< fvPatchVectorFieldclone () const
 Construct and return a clone. More...
 
 turbulentDFSEMInletFvPatchVectorField (const turbulentDFSEMInletFvPatchVectorField &, const DimensionedField< vector, volMesh > &)
 Construct as copy setting internal field reference. More...
 
virtual tmp< fvPatchVectorFieldclone (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...
 
template<class Type >
Foam::tmp< Foam::Field< Type > > interpolateOrRead (const word &fieldName, const dictionary &dict, bool &interpolateField) const
 
template<class Type >
Foam::tmp< Foam::Field< Type > > interpolateBoundaryData (const word &fieldName) const
 

Static Public Member Functions

static bool checkStresses (const symmTensorField &Rf)
 Helper function to check that Reynold stresses are valid. More...
 

Detailed Description

Velocity boundary condition including synthesised eddies for use with LES and DES turbulent flows.

Reference:

    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

Reynolds stress, velocity and turbulence length scale values can either be specified directly, or mapped. If mapping, the values should be entered in the same form as the timeVaryingMappedFixedValue condition, except that no interpolation in time is supported. These should be located in the directory:

    \$FOAM_CASE/constant/boundaryData/\<patchName\>/points
    \$FOAM_CASE/constant/boundaryData/\<patchName\>/0/\{R|U|L\}
Usage
Property Description Required Default value
value Restart value yes
delta Local limiting length scale yes
R Reynolds stress field no
U Velocity field no
L Turbulence length scale field no
d Eddy density (fill fraction) no 1
kappa Von Karman constant no 0.41
mapMethod Method to map reference values no nearestCell
perturb Point perturbation for interpolation no 1e-5
interpolateR Flag to interpolate the R field no false
interpolateL Flag to interpolate the L field no false
interpolateU Flag to interpolate the U field no false
writeEddies Flag to write eddies as OBJ file no no
Note
  • The delta value typically represents the characteristic scale of flow or flow domain, e.g. a channel half-height
  • For R, U and L specification: if the entry is not user input, it is assumed that the data will be mapped
See also
timeVaryingMappedFixedValueFvPatchField turbulentDigitalFilterInlet
Source files

Definition at line 179 of file turbulentDFSEMInletFvPatchVectorField.H.

Constructor & Destructor Documentation

◆ turbulentDFSEMInletFvPatchVectorField() [1/5]

Construct from patch and internal field.

Definition at line 735 of file turbulentDFSEMInletFvPatchVectorField.C.

Referenced by turbulentDFSEMInletFvPatchVectorField::clone().

Here is the caller graph for this function:

◆ turbulentDFSEMInletFvPatchVectorField() [2/5]

turbulentDFSEMInletFvPatchVectorField ( const fvPatch p,
const DimensionedField< vector, volMesh > &  iF,
const dictionary dict 
)

Construct from patch, internal field and dictionary.

Definition at line 825 of file turbulentDFSEMInletFvPatchVectorField.C.

References Foam::expressions::patchExpr::debug, eddy::debug, Foam::gSum(), and Foam::foamVersion::patch.

Here is the call graph for this function:

◆ turbulentDFSEMInletFvPatchVectorField() [3/5]

Construct by mapping given turbulentDFSEMInletFvPatchVectorField onto a new patch

Definition at line 779 of file turbulentDFSEMInletFvPatchVectorField.C.

◆ turbulentDFSEMInletFvPatchVectorField() [4/5]

◆ turbulentDFSEMInletFvPatchVectorField() [5/5]

Construct as copy setting internal field reference.

Definition at line 920 of file turbulentDFSEMInletFvPatchVectorField.C.

◆ ~turbulentDFSEMInletFvPatchVectorField()

virtual ~turbulentDFSEMInletFvPatchVectorField ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "turbulentDFSEMInlet"  )

Runtime type information.

◆ clone() [1/2]

virtual tmp<fvPatchVectorField> clone ( ) const
inlinevirtual

Construct and return a clone.

Definition at line 383 of file turbulentDFSEMInletFvPatchVectorField.H.

References turbulentDFSEMInletFvPatchVectorField::turbulentDFSEMInletFvPatchVectorField().

Here is the call graph for this function:

◆ clone() [2/2]

virtual tmp<fvPatchVectorField> clone ( const DimensionedField< vector, volMesh > &  iF) const
inlinevirtual

Construct and return a clone setting internal field reference.

Definition at line 400 of file turbulentDFSEMInletFvPatchVectorField.H.

References turbulentDFSEMInletFvPatchVectorField::turbulentDFSEMInletFvPatchVectorField().

Here is the call graph for this function:

◆ checkStresses()

bool checkStresses ( const symmTensorField Rf)
static

Helper function to check that Reynold stresses are valid.

Definition at line 965 of file turbulentDFSEMInletFvPatchVectorField.C.

References Foam::expressions::patchExpr::debug, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::Pout, R, Foam::sqr(), and Foam::sqrt().

Here is the call graph for this function:

◆ autoMap()

void autoMap ( const fvPatchFieldMapper m)
virtual

Map (and resize as needed) from self given a mapping object.

Definition at line 1033 of file turbulentDFSEMInletFvPatchVectorField.C.

References fvPatchField< Type >::autoMap().

Here is the call graph for this function:

◆ rmap()

void rmap ( const fvPatchVectorField ptf,
const labelList addr 
)
virtual

Reverse map the given fvPatchField onto this fvPatchField.

Definition at line 1050 of file turbulentDFSEMInletFvPatchVectorField.C.

References fvPatchField< Type >::rmap().

Here is the call graph for this function:

◆ updateCoeffs()

void updateCoeffs ( )
virtual

Update the coefficients associated with the patch field.

Definition at line 1071 of file turbulentDFSEMInletFvPatchVectorField.C.

References Foam::constant::universal::c, Foam::expressions::patchExpr::debug, Foam::endl(), forAll, Foam::gMax(), Foam::gMin(), Foam::gSum(), Foam::Info, n, UPstream::nProcs(), UPstream::parRun(), Foam::foamVersion::patch, Foam::returnReduce(), Foam::sqrt(), timeIndex, and U.

Here is the call graph for this function:

◆ write()

void write ( Ostream os) const
virtual

Write.

Definition at line 1183 of file turbulentDFSEMInletFvPatchVectorField.C.

References Foam::constant::electromagnetic::e, fvPatchField< Type >::write(), Ostream::writeEntry(), and Ostream::writeEntryIfDifferent().

Here is the call graph for this function:

◆ interpolateOrRead()

Foam::tmp<Foam::Field<Type> > interpolateOrRead ( const word fieldName,
const dictionary dict,
bool interpolateField 
) const

◆ interpolateBoundaryData()

Foam::tmp<Foam::Field<Type> > interpolateBoundaryData ( const word fieldName) const

Definition at line 70 of file turbulentDFSEMInletFvPatchVectorFieldTemplates.C.

References Foam::endl(), Foam::Info, and Foam::foamVersion::patch.

Here is the call graph for this function:

The documentation for this class was generated from the following files: