48void Foam::LESModels::DeltaOmegaTildeDelta::calcDelta()
75 const point& cc = cellCentres[celli];
76 const vector& nv = nvecvort[celli];
78 scalar deltaMaxTmp = 0;
80 for (
const label facei : cFaces)
82 const point& fc = faceCentres[facei];
83 const scalar
tmp = 2.0*
mag(nv ^ (fc - cc));
85 if (
tmp > deltaMaxTmp)
91 delta_[celli] = deltaCoeff_*deltaMaxTmp;
94 const label nD =
mesh.nGeometricD();
99 <<
"Case is 2D, LES is not strictly applicable" <<
nl
109 delta_.correctBoundaryConditions();
115Foam::LESModels::DeltaOmegaTildeDelta::DeltaOmegaTildeDelta
126 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
134 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<bool>
149 dict.optionalSubDict(
"hmaxCoeffs"),
155 Info<<
"Employing " << maxDeltaxyz::typeName <<
" for hmax" <<
endl;
163 dict.optionalSubDict(
"hmaxCoeffs")
179 coeffsDict.
readIfPresent<
bool>(
"requireUpdate", requireUpdate_);
187 if (turbulenceModel_.mesh().changing() && requireUpdate_)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Delta formulation that accounts for the orientation of the vorticity vector. In "2D-regions" (i....
virtual void read(const dictionary &)
Read the LESdelta dictionary.
Abstract base class for LES deltas.
LESdelta(const LESdelta &)=delete
No copy construct.
static autoPtr< LESdelta > New(const word &name, const turbulenceModel &turbulence, const dictionary &dict, const word &lookupName="delta")
Return a reference to the selected LES delta.
const turbulenceModel & turbulence() const
Return turbulenceModel reference.
const turbulenceModel & turbulenceModel_
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Mesh data needed to do the Finite Volume discretisation.
Delta calculated by taking the maximum distance between the cell centre and any face centre....
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the curl of the given volField by constructing the Hodge-dual of the symmetric part of the ...
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
List< cell > cellList
List of cell.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.