Loading...
Searching...
No Matches
ReynoldsAnalogy.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "ReynoldsAnalogy.H"
29#include "fluidThermo.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39{
42 (
46 );
48}
49
50
51// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52
55{
56 if (rhoName_ == "rhoInf")
57 {
58 const label n = mesh_.boundary()[patchi].size();
60 }
61 else if (mesh_.foundObject<volScalarField>(rhoName_, false))
62 {
64 return rho.boundaryField()[patchi];
65 }
66
68 << "Unable to set rho for patch " << patchi
70
71 return nullptr;
72}
73
74
77{
78 if (CpName_ == "CpInf")
79 {
80 const label n = mesh_.boundary()[patchi].size();
81 return tmp<scalarField>::New(n, CpRef_);
82 }
83 else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
84 {
85 const auto& thermo =
86 mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
87
88 const scalarField& pp = thermo.p().boundaryField()[patchi];
89 const scalarField& Tp = thermo.T().boundaryField()[patchi];
90
91 return thermo.Cp(pp, Tp, patchi);
92 }
93
95 << "Unable to set Cp for patch " << patchi
97
98 return nullptr;
99}
100
101
104{
105 typedef compressible::turbulenceModel cmpTurbModel;
106 typedef incompressible::turbulenceModel icoTurbModel;
107
108 if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
109 {
110 const auto& turb =
111 mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
112
113 return turb.devRhoReff()/turb.rho();
114 }
115 else if (mesh_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
116 {
117 const auto& turb =
118 mesh_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
119
120 return turb.devReff();
121 }
122 else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
123 {
124 const auto& thermo =
125 mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
126
127 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
128
129 return -thermo.nu()*devTwoSymm(fvc::grad(U));
130 }
131 else if (mesh_.foundObject<transportModel>("transportProperties"))
132 {
133 const auto& laminarT =
134 mesh_.lookupObject<transportModel>("transportProperties");
135
136 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
137
138 return -laminarT.nu()*devTwoSymm(fvc::grad(U));
139 }
140 else if (mesh_.foundObject<dictionary>("transportProperties"))
141 {
142 const auto& transportProperties =
143 mesh_.lookupObject<dictionary>("transportProperties");
144
146
147 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
148
149 return -nu*devTwoSymm(fvc::grad(U));
150 }
151
153 << "No valid model for viscous stress calculation"
155
156 return nullptr;
157}
158
159
162{
163 const auto& U = mesh_.lookupObject<volVectorField>(UName_);
164 const volVectorField::Boundary& Ubf = U.boundaryField();
165
166 auto tCf = tmp<FieldField<Field, scalar>>::New(Ubf.size());
167 auto& Cf = tCf.ref();
168
169 forAll(Cf, patchi)
170 {
171 Cf.set(patchi, new Field<scalar>(Ubf[patchi].size(), Zero));
172 }
173
174 const volSymmTensorField R(devReff());
175 const volSymmTensorField::Boundary& Rbf = R.boundaryField();
176
177 for (const label patchi : patchIDs_)
178 {
179 const fvPatchVectorField& Up = Ubf[patchi];
180
181 const symmTensorField& Rp = Rbf[patchi];
182
183 tmp<vectorField> tnHat = Up.patch().nf();
184
185 tmp<scalarField> ttauByRhop = mag(tnHat & Rp);
186
187 Cf[patchi] = 2*ttauByRhop/magSqr(URef_);
188 }
189
190 return tCf;
191}
192
193
195(
196 volScalarField& htc,
198)
199{
200 const FieldField<Field, scalar> CfBf(Cf());
201 const scalar magU = mag(URef_);
202
204
205 for (const label patchi : patchIDs_)
206 {
207 tmp<scalarField> trhop = rho(patchi);
208 tmp<scalarField> tCpp = Cp(patchi);
209
210 htcBf[patchi] = 0.5*trhop*tCpp*magU*CfBf[patchi];
211 }
212}
213
214
215// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
216
218(
219 const dictionary& dict,
220 const fvMesh& mesh,
221 const word& TName
222)
223:
225 UName_("U"),
226 URef_(Zero),
227 rhoName_("rho"),
228 rhoRef_(0),
229 CpName_("Cp"),
230 CpRef_(0)
232 read(dict);
233}
234
235
236// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
237
239(
240 const dictionary& dict
241)
242{
244 {
245 return false;
246 }
247
248 dict.readIfPresent("U", UName_);
249 dict.readEntry("UInf", URef_);
250
251 dict.readIfPresent("Cp", CpName_);
252 if (CpName_ == "CpInf")
253 {
254 dict.readEntry("CpInf", CpRef_);
255 }
256
257 dict.readIfPresent("rho", rhoName_);
258 if (rhoName_ == "rhoInf")
259 {
260 dict.readEntry("rhoInf", rhoRef_);
261 }
262
263 return true;
264}
265
266
267// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
label n
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())
compressible::turbulenceModel & turb
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const word dictName
The dictionary name ("thermophysicalProperties").
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
const fvPatch & patch() const noexcept
Return the patch.
tmp< vectorField > nf() const
Return face unit normals, like the fvMesh::unitSf() method Same as unitSf().
Definition fvPatch.C:143
A base class for heat transfer coefficient models.
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
const fvMesh & mesh_
Const reference to the mesh.
heatTransferCoeffModel(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
tmp< FieldField< Field, scalar > > q() const
Return boundary fields of heat-flux field.
const word & TName() const noexcept
Return const reference to name of temperature field.
labelList patchIDs_
List of (wall) patches to process (selected by name).
Heat transfer coefficient calculation based on Reynolds Analogy, which is used to relate turbulent mo...
scalar CpRef_
Reference specific heat capacity.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual tmp< scalarField > Cp(const label patchi) const
Return heat capacity at constant pressure [J/kg/K].
ReynoldsAnalogy(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
word rhoName_
Name of fluid density field.
virtual tmp< scalarField > rho(const label patchi) const
Return fluid density field [kg/m^3].
word CpName_
Name of specific heat capacity field.
tmp< FieldField< Field, scalar > > Cf() const
Return skin friction coefficient field [-].
virtual void htc(volScalarField &htc, const FieldField< Field, scalar > &q)
Set the heat transfer coefficient.
bool foundObject(const word &name, const bool recursive=false) const
Contains the named Type?
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const volScalarField & Cp
Definition EEqn.H:7
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
A namespace for various heat transfer coefficient model implementations.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
const dimensionSet dimViscosity
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
volScalarField & nu
dictionary dict
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299