117#ifndef Foam_DEShybrid_H
118#define Foam_DEShybrid_H
121#include "surfaceInterpolate.H"
141 typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
142 typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
183 DEShybrid(
const DEShybrid&) =
delete;
186 void operator=(
const DEShybrid&) =
delete;
194 if (U0_.value() <= 0)
197 <<
"U0 coefficient must be > 0. "
200 if (L0_.value() <= 0)
203 <<
"L0 coefficient must be > 0. "
209 <<
"sigmaMin coefficient must be >= 0. "
215 <<
"sigmaMax coefficient must be >= 0. "
221 <<
"sigmaMin coefficient must be <= 1. "
227 <<
"sigmaMax coefficient must be <= 1. "
234 <<
" delta : " << deltaName_ <<
nl
235 <<
" CDES : " << CDES_ <<
nl
236 <<
" U0 : " << U0_.value() <<
nl
237 <<
" L0 : " << L0_.value() <<
nl
238 <<
" sigmaMin : " << sigmaMin_ <<
nl
239 <<
" sigmaMax : " << sigmaMax_ <<
nl
240 <<
" OmegaLim : " << OmegaLim_ <<
nl
241 <<
" nutLim : " << nutLim_ <<
nl
242 <<
" CH1 : " << CH1_ <<
nl
243 <<
" CH2 : " << CH2_ <<
nl
244 <<
" CH3 : " << CH3_ <<
nl
245 <<
" Cs : " << Cs_ <<
nl
254 const VolFieldType& vf,
268 CH3_*Omega*
max(S, Omega)
282 /(
pow(0.09, 1.5)*tK),
292 CDES_*
delta/
max(tlTurb*tg, SMALL*L0_) - 0.5
303 mesh.time().timeName(),
310 if (blendedSchemeBaseName::debug)
327 auto& factor = *factorPtr;
329 factor =
max(sigmaMax_*
tanh(
pow(
A, CH1_)), sigmaMin_);
333 vf.name() +
"BlendingFactor",
349 vf.name() +
"BlendingFactor",
367 DEShybrid(
const fvMesh&
mesh, Istream& is)
369 surfaceInterpolationScheme<Type>(
mesh),
370 tScheme1_(surfaceInterpolationScheme<Type>::
New(
mesh, is)),
371 tScheme2_(surfaceInterpolationScheme<Type>::
New(
mesh, is)),
373 CDES_(readScalar(is)),
376 sigmaMin_(readScalar(is)),
377 sigmaMax_(readScalar(is)),
378 OmegaLim_(readScalar(is)),
379 nutLim_(readScalarOrDefault(is, scalar(1))),
396 surfaceInterpolationScheme<Type>(
mesh),
399 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
403 surfaceInterpolationScheme<Type>::
New(
mesh, faceFlux, is)
406 CDES_(readScalar(is)),
409 sigmaMin_(readScalar(is)),
410 sigmaMax_(readScalar(is)),
411 OmegaLim_(readScalar(is)),
412 nutLim_(readScalarOrDefault(is, scalar(1))),
437 const auto* modelPtr =
445 const auto& model = *modelPtr;
447 return calcBlendingFactor
449 vf, model.nut(), model.nu(), model.U(),
delta
454 <<
"Scheme requires a turbulence model to be present. "
455 <<
"Unable to retrieve turbulence model from the mesh "
463 tmp<surfaceScalarField>
weights(
const VolFieldType& vf)
const
468 (scalar(1) - bf)*tScheme1_().weights(vf)
469 + bf*tScheme2_().weights(vf);
480 (scalar(1) - bf)*tScheme1_().interpolate(vf)
481 + bf*tScheme2_().interpolate(vf);
488 return tScheme1_().corrected() || tScheme2_().corrected();
494 virtual tmp<SurfaceFieldType>
correction(
const VolFieldType& vf)
const
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Improved hybrid convection scheme of Travin et al. for hybrid RAS/LES calculations with enhanced Grey...
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
DEShybrid(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< surfaceScalarField > weights(const VolFieldType &vf) const
Return the interpolation weighting factors.
virtual tmp< SurfaceFieldType > correction(const VolFieldType &vf) const
Return the explicit correction to the face-interpolate for the given field.
DEShybrid(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
tmp< SurfaceFieldType > interpolate(const VolFieldType &vf) const
Return the face-interpolate of the given cell field with explicit correction.
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
TypeName("DEShybrid")
Runtime type information.
Generic GeometricField class.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
bool registerObject() const noexcept
Should objects created with this IOobject be registered?
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Base class for blended schemes to provide access to the blending factor surface field.
blendedSchemeBase()=default
Constructor.
Mesh data needed to do the Finite Volume discretisation.
bool store()
Register object with its registry and transfer ownership to the registry.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the gradient of the given field.
Namespace for handling debugging switches.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
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).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
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.
dimensionedScalar pow4(const dimensionedScalar &ds)
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...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedTensor skew(const dimensionedTensor &dt)
constexpr char nl
The newline '\n' character (0x0a).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.