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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
36 namespace Foam
37 {
38 namespace heatTransferCoeffModels
39 {
40  defineTypeNameAndDebug(ReynoldsAnalogy, 0);
42  (
43  heatTransferCoeffModel,
44  ReynoldsAnalogy,
45  dictionary
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  {
63  const volScalarField& rho =
65  return rho.boundaryField()[patchi];
66  }
67 
69  << "Unable to set rho for patch " << patchi
70  << exit(FatalError);
71 
72  return nullptr;
73 }
74 
75 
78 {
79  if (CpName_ == "CpInf")
80  {
81  const label n = mesh_.boundary()[patchi].size();
82  return tmp<Field<scalar>>::New(n, CpRef_);
83  }
84  else if (mesh_.foundObject<fluidThermo>(fluidThermo::typeName))
85  {
86  const fluidThermo& thermo =
87  mesh_.lookupObject<fluidThermo>(fluidThermo::typeName);
88 
89  const scalarField& pp = thermo.p().boundaryField()[patchi];
90  const scalarField& Tp = thermo.T().boundaryField()[patchi];
91 
92  return thermo.Cp(pp, Tp, patchi);
93  }
94 
96  << "Unable to set Cp for patch " << patchi
97  << exit(FatalError);
98 
99  return nullptr;
100 }
101 
102 
105 {
106  typedef compressible::turbulenceModel cmpTurbModel;
107  typedef incompressible::turbulenceModel icoTurbModel;
108 
109  if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
110  {
111  const cmpTurbModel& turb =
112  mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
113 
114  return turb.devRhoReff()/turb.rho();
115  }
116  else if (mesh_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
117  {
119  mesh_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
120 
121  return turb.devReff();
122  }
123  else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
124  {
125  const fluidThermo& thermo =
126  mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
127 
128  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
129 
130  return -thermo.nu()*dev(twoSymm(fvc::grad(U)));
131  }
132  else if (mesh_.foundObject<transportModel>("transportProperties"))
133  {
134  const transportModel& laminarT =
135  mesh_.lookupObject<transportModel>("transportProperties");
136 
137  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
138 
139  return -laminarT.nu()*dev(twoSymm(fvc::grad(U)));
140  }
141  else if (mesh_.foundObject<dictionary>("transportProperties"))
142  {
144  mesh_.lookupObject<dictionary>("transportProperties");
145 
147 
148  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
149 
150  return -nu*dev(twoSymm(fvc::grad(U)));
151  }
152  else
153  {
155  << "No valid model for viscous stress calculation"
156  << exit(FatalError);
157 
158  return volSymmTensorField::null();
159  }
160 }
161 
162 
165 {
166  const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
167  const volVectorField::Boundary& Ubf = U.boundaryField();
168 
169  auto tCf = tmp<FieldField<Field, scalar>>::New(Ubf.size());
170  auto& Cf = tCf.ref();
171 
172  forAll(Cf, patchi)
173  {
174  Cf.set(patchi, new Field<scalar>(Ubf[patchi].size(), Zero));
175  }
176 
177  const volSymmTensorField R(devReff());
178  const volSymmTensorField::Boundary& Rbf = R.boundaryField();
179 
180  for (const label patchi : patchSet_)
181  {
182  const fvPatchVectorField& Up = Ubf[patchi];
183 
184  const symmTensorField& Rp = Rbf[patchi];
185 
186  const vectorField nHat(Up.patch().nf());
187 
188  const scalarField tauByRhop(mag(nHat & Rp));
189 
190  Cf[patchi] = 2*tauByRhop/magSqr(URef_);
191  }
192 
193  return tCf;
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
199 Foam::heatTransferCoeffModels::ReynoldsAnalogy::ReynoldsAnalogy
200 (
201  const dictionary& dict,
202  const fvMesh& mesh,
203  const word& TName
204 )
205 :
207  UName_("U"),
208  URef_(Zero),
209  rhoName_("rho"),
210  rhoRef_(0.0),
211  CpName_("Cp"),
212  CpRef_(0.0)
213 {
214  read(dict);
215 }
216 
217 
218 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
219 
221 (
222  const dictionary& dict
223 )
224 {
226  {
227  dict.readEntry("UInf", URef_);
228 
229  dict.readIfPresent("Cp", CpName_);
230  if (CpName_ == "CpInf")
231  {
232  dict.readEntry("CpInf", CpRef_);
233  }
234 
235  dict.readIfPresent("rho", rhoName_);
236  if (rhoName_ == "rhoInf")
237  {
238  dict.readEntry("rhoInf", rhoRef_);
239  }
240 
241  return true;
242  }
243 
244  return false;
245 }
246 
247 
249 {
250  const FieldField<Field, scalar> CfBf(Cf());
251  const scalar magU = mag(URef_);
252 
253  volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
254 
255  for (const label patchi : patchSet_)
256  {
257  const scalarField rhop(rho(patchi));
258  const scalarField Cpp(Cp(patchi));
259 
260  htcBf[patchi] = 0.5*rhop*Cpp*magU*CfBf[patchi];
261  }
262 }
263 
264 
265 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::heatTransferCoeffModels::ReynoldsAnalogy::devReff
virtual tmp< volSymmTensorField > devReff() const
Definition: ReynoldsAnalogy.C:104
Foam::heatTransferCoeffModels::ReynoldsAnalogy::rhoRef_
scalar rhoRef_
Reference density.
Definition: ReynoldsAnalogy.H:152
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::heatTransferCoeffModel
An abstract base class for heat transfer coeffcient models.
Definition: heatTransferCoeffModel.H:63
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::transportModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::heatTransferCoeffModel::mesh_
const fvMesh & mesh_
Definition: heatTransferCoeffModel.H:79
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::heatTransferCoeffModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(heatTransferCoeffModel, fixedReferenceTemperature, dictionary)
fluidThermo.H
Foam::heatTransferCoeffModels::ReynoldsAnalogy::htc
virtual void htc(volScalarField &htc)
Set the heat transfer coefficient.
Definition: ReynoldsAnalogy.C:248
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
transportProperties
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
turbulentTransportModel.H
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
rho
rho
Definition: readInitialConditions.H:96
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cf
tmp< FieldField< Field, scalar > > Cf() const
Definition: ReynoldsAnalogy.C:164
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:120
R
#define R(A, B, C, D, E, F, K, M)
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
Foam::heatTransferCoeffModel::New
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
Definition: heatTransferCoeffModelNew.C:35
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cp
virtual tmp< Field< scalar > > Cp(const label patchi) const
Definition: ReynoldsAnalogy.C:77
Foam::heatTransferCoeffModels::ReynoldsAnalogy::rho
virtual tmp< Field< scalar > > rho(const label patchi) const
Definition: ReynoldsAnalogy.C:54
ReynoldsAnalogy.H
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::heatTransferCoeffModels::ReynoldsAnalogy::rhoName_
word rhoName_
Name of density field.
Definition: ReynoldsAnalogy.H:149
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
U
U
Definition: pEqn.H:72
Foam::heatTransferCoeffModels::ReynoldsAnalogy::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: ReynoldsAnalogy.C:221
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::heatTransferCoeffModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: heatTransferCoeffModel.C:134
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::GeometricField< symmTensor, fvPatchField, volMesh >::null
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:32
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
turbulentFluidThermoModel.H
Foam::heatTransferCoeffModels::defineTypeNameAndDebug
defineTypeNameAndDebug(fixedReferenceTemperature, 0)
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106