Go to the documentation of this file.
38 namespace heatTransferCoeffModels
43 heatTransferCoeffModel,
64 return rho.boundaryField()[patchi];
68 <<
"Unable to set rho for patch " << patchi
78 if (CpName_ ==
"CpInf")
80 const label
n = mesh_.boundary()[patchi].size();
91 return thermo.Cp(pp, Tp, patchi);
95 <<
"Unable to set Cp for patch " << patchi
108 if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
111 mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
113 return turb.devRhoReff()/
turb.rho();
115 else if (mesh_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
118 mesh_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
120 return turb.devReff();
131 else if (mesh_.foundObject<
transportModel>(
"transportProperties"))
133 const auto& laminarT =
140 else if (mesh_.foundObject<
dictionary>(
"transportProperties"))
143 mesh_.lookupObject<
dictionary>(
"transportProperties");
153 <<
"No valid model for viscous stress calculation"
164 const volVectorField::Boundary& Ubf =
U.boundaryField();
167 auto& Cf = tCf.ref();
175 const volSymmTensorField::Boundary& Rbf =
R.boundaryField();
177 for (
const label patchi : patchSet_)
187 Cf[patchi] = 2*tauByRhop/
magSqr(URef_);
224 dict.readEntry(
"UInf", URef_);
226 dict.readIfPresent(
"Cp", CpName_);
227 if (CpName_ ==
"CpInf")
229 dict.readEntry(
"CpInf", CpRef_);
232 dict.readIfPresent(
"rho", rhoName_);
233 if (rhoName_ ==
"rhoInf")
235 dict.readEntry(
"rhoInf", rhoRef_);
252 const scalar magU =
mag(URef_);
256 for (
const label patchi : patchSet_)
261 htcBf[patchi] = 0.5*rhop*Cpp*magU*CfBf[patchi];
virtual tmp< volSymmTensorField > devReff() const
scalar rhoRef_
Reference density.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
An abstract base class for heat transfer coeffcient models.
A class for handling words, derived from Foam::string.
A field of fields is a PtrList of fields with reference counting.
A class for managing temporary objects.
const fvMesh & mesh_
Mesh reference.
static constexpr const zero Zero
Global zero (0)
addToRunTimeSelectionTable(heatTransferCoeffModel, fixedReferenceTemperature, dictionary)
ReynoldsAnalogy(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
bool read(const char *buf, int32_t &val)
Same as readInt32.
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Fundamental fluid thermodynamic properties.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
#define forAll(list, i)
Loop across all elements in list.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
tmp< FieldField< Field, scalar > > Cf() const
tmp< vectorField > nf() const
Return face normals.
#define R(A, B, C, D, E, F, K, M)
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
virtual tmp< Field< scalar > > Cp(const label patchi) const
virtual tmp< Field< scalar > > rho(const label patchi) const
const Type & lookupObject(const word &name, const bool recursive=false) const
const dimensionSet dimViscosity
word rhoName_
Name of density field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Base-class for all transport models used by the incompressible turbulence models.
virtual void htc(volScalarField &htc, const FieldField< Field, scalar > &q)
Set the heat transfer coefficient.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
virtual bool read(const dictionary &dict)
Read from dictionary.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual bool read(const dictionary &dict)
Read from dictionary.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volScalarField & Cp
Templated abstract base class for single-phase incompressible turbulence models.
const fvPatch & patch() const
Return patch.
compressible::turbulenceModel & turb
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
defineTypeNameAndDebug(fixedReferenceTemperature, 0)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
word dictName() const
The local dictionary name (final part of scoped name)