43template<
class Type,
class WeightType = oneField>
48 const WeightType& weights =
oneField()
56 auto& scaling = tscaling.ref();
58 const auto& faceNeighbour =
mesh.faceNeighbour();
59 const auto& faceOwner =
mesh.faceOwner();
63 const label own = faceOwner[i];
64 const label nbr = faceNeighbour[i];
66 result[own] +=
field[nbr];
67 result[nbr] +=
field[own];
69 scaling[own] = scaling[own] + weights[nbr];
70 scaling[nbr] = scaling[nbr] + weights[own];
73 const auto&
pbm =
mesh.boundaryMesh();
77 const auto& pf =
field.boundaryField()[
pp.index()];
82 const scalarField nbrField(pf.patchNeighbourField());
87 result[celli] += nbrField[facei];
88 scaling[celli] = scaling[celli] + weights[celli];
113void Foam::LESModels::SLADelta::calcDelta()
138 *
sqrt(6.0)*
mag((S & vorticity)^vorticity)
141 vtm.correctBoundaryConditions();
157 vtmAve.primitiveFieldRef() /= tweights + 1;
171 (
nut +
nu)/(
max(tmagGradU, magGradUeps)*
sqr(kappa_*
y) + nuEps),
188 const point& cc = cellCentres[celli];
189 const vector& nv = nVecVort[celli];
191 scalar deltaMaxTmp = 0;
193 for (
const label facei : cFaces)
195 const point& fc = faceCentres[facei];
196 const scalar
tmp = 2.0*
mag(nv ^ (fc - cc));
198 if (
tmp > deltaMaxTmp)
212 + (FKHmax_ - FKHmin_)*(vtmAve[celli] - a1_)/(a2_ - a1_)
216 if ((shielding_ ==
true) && (fd[celli] < (1.0 - epsilon_)))
221 delta_[celli] = deltaCoeff_*deltaMaxTmp*FKH;
224 const label nD =
mesh.nGeometricD();
229 <<
"Case is 2D, LES is not strictly applicable" <<
nl
239 delta_.correctBoundaryConditions();
245Foam::LESModels::SLADelta::SLADelta
256 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
264 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<bool>
272 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<bool>
280 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
288 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
296 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
304 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
312 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
320 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
328 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
336 dict.optionalSubDict(
type() +
"Coeffs").getOrDefault<scalar>
351 dict.optionalSubDict(
"hmaxCoeffs"),
357 Info<<
"Employing " << maxDeltaxyz::typeName <<
" for hmax" <<
endl;
365 dict.optionalSubDict(
"hmaxCoeffs")
370 if (
mag(a2_ - a1_) < SMALL)
373 <<
"Model coefficients a1 = " << a1_
374 <<
", and a2 = " << a2_ <<
" cannot be equal."
386 coeffsDict.readIfPresent<scalar>(
"deltaCoeff", deltaCoeff_);
387 coeffsDict.readIfPresent<
bool>(
"requireUpdate", requireUpdate_);
388 coeffsDict.readIfPresent<scalar>(
"FKHmin", FKHmin_);
389 coeffsDict.readIfPresent<scalar>(
"FKHmax", FKHmax_);
390 coeffsDict.readIfPresent<scalar>(
"a1", a1_);
391 coeffsDict.readIfPresent<scalar>(
"a2", a2_);
392 coeffsDict.readIfPresent<scalar>(
"epsilon", epsilon_);
393 coeffsDict.readIfPresent<scalar>(
"kappa", kappa_);
394 coeffsDict.readIfPresent<scalar>(
"Cd1", Cd1_);
395 coeffsDict.readIfPresent<scalar>(
"Cd2", Cd2_);
403 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.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Generic GeometricField class.
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 and a flow-sensitised fun...
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_
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
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...
Smooth ATC in cells next to a set of patches supplied by type.
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 representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
A patch is a list of labels that address the faces in the global face list.
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const fvMesh & mesh() const
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#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 ...
Calculate the gradient of the given field.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
tmp< scalarField > sumNeighbours(const GeometricField< Type, fvPatchField, volMesh > &field, GeometricField< Type, fvPatchField, volMesh > &result, const WeightType &weights=oneField())
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
List< cell > cellList
List of cell.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
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)
UList< label > labelUList
A UList of labels.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.