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-2020 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"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace heatTransferCoeffModels
39{
42 (
46 );
47}
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
69 << exit(FatalError);
70
71 return nullptr;
72}
73
74
77{
78 if (CpName_ == "CpInf")
79 {
80 const label n = mesh_.boundary()[patchi].size();
81 return tmp<Field<scalar>>::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
96 << exit(FatalError);
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 {
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()*dev(twoSymm(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()*dev(twoSymm(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*dev(twoSymm(fvc::grad(U)));
150 }
151
153 << "No valid model for viscous stress calculation"
154 << exit(FatalError);
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 : patchSet_)
178 {
179 const fvPatchVectorField& Up = Ubf[patchi];
180
181 const symmTensorField& Rp = Rbf[patchi];
182
183 const vectorField nHat(Up.patch().nf());
184
185 const scalarField tauByRhop(mag(nHat & Rp));
186
187 Cf[patchi] = 2*tauByRhop/magSqr(URef_);
188 }
189
190 return tCf;
191}
192
193
194// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
195
197(
198 const dictionary& dict,
199 const fvMesh& mesh,
200 const word& TName
201)
202:
204 UName_("U"),
205 URef_(Zero),
206 rhoName_("rho"),
207 rhoRef_(0),
208 CpName_("Cp"),
209 CpRef_(0)
210{
211 read(dict);
212}
213
214
215// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
216
218(
219 const dictionary& dict
220)
221{
223 {
224 dict.readEntry("UInf", URef_);
225
226 dict.readIfPresent("Cp", CpName_);
227 if (CpName_ == "CpInf")
228 {
229 dict.readEntry("CpInf", CpRef_);
230 }
231
232 dict.readIfPresent("rho", rhoName_);
233 if (rhoName_ == "rhoInf")
234 {
235 dict.readEntry("rhoInf", rhoRef_);
236 }
237
238 return true;
239 }
240
241 return false;
242}
243
244
246(
247 volScalarField& htc,
249)
250{
251 const FieldField<Field, scalar> CfBf(Cf());
252 const scalar magU = mag(URef_);
253
255
256 for (const label patchi : patchSet_)
257 {
258 const scalarField rhop(rho(patchi));
259 const scalarField Cpp(Cp(patchi));
260
261 htcBf[patchi] = 0.5*rhop*Cpp*magU*CfBf[patchi];
262 }
263}
264
265
266// ************************************************************************* //
#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.
compressible::turbulenceModel & turb
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:144
An abstract base class for heat transfer coeffcient models.
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
const fvMesh & mesh_
Mesh reference.
Heat transfer coefficient calculation based on Reynolds Analogy, which is used to relate turbulent mo...
virtual tmp< volSymmTensorField > devReff() const
virtual bool read(const dictionary &dict)
Read from dictionary.
tmp< FieldField< Field, scalar > > Cf() const
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
Is the named Type found?
const Type & lookupObject(const word &name, const bool recursive=false) const
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:65
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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:453
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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:333