ReynoldsAnalogy.H
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-2022 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 Class
27  Foam::heatTransferCoeffModels::ReynoldsAnalogy
28 
29 Description
30  Heat transfer coefficient calculation based on Reynolds Analogy,
31  which is used to relate turbulent momentum and heat transfer.
32 
33  The heat transfer coefficient is derived
34  from the skin friction coefficient:
35 
36  \f[
37  C_f = \frac{\tau_w}{0.5 \rho_{ref} |U|^2}
38  \f]
39 
40  as:
41 
42  \f[
43  h = 0.5 \rho_{ref} c_{p,ref} |U_{ref}| C_f
44  \f]
45 
46  where
47  \vartable
48  h | Heat transfer coefficient [W/m^2/K]
49  \rho_{ref} | Reference fluid density [kg/m^3]
50  c_{p,ref} | Reference specific heat capacity at constant pressure [J/kg/K]
51  U_{ref} | Reference velocity [m/s]
52  C_f | Skin friction coefficient [-]
53  \tau_w | Wall shear stress [m^2/s^2]
54  \endvartable
55 
56 Usage
57  Minimal example by using \c system/controlDict.functions:
58  \verbatim
59  heatTransferCoeff1
60  {
61  // Inherited entries
62  ...
63 
64  // Mandatory entries
65  htcModel ReynoldsAnalogy;
66  UInf <vector>;
67 
68  // Optional entries
69  U <word>;
70  Cp <word>;
71  rho <word>;
72 
73  // Conditional mandatory entries
74 
75  // when Cp == CpInf
76  CpInf <scalar>;
77 
78  // when rho == rhoInf
79  rhoInf <scalar>;
80  }
81  \endverbatim
82 
83  where the entries mean:
84  \table
85  Property | Description | Type | Reqd | Deflt
86  type | Model name: ReynoldsAnalogy | word | yes | -
87  UInf | Reference velocity | vector | yes | -
88  U | Name of velocity field | word | no | U
89  Cp | Name of reference specific heat capacity | word | no | Cp
90  CpInf | Reference specific heat capacity value | scalar | choice | -
91  rho | Name of fluid density field | word | no | rho
92  rhoInf | Reference fluid density value | scalar | choice | -
93  \endtable
94 
95 Note
96  - In order to use a reference \c Cp, set \c Cp to \c CpInf.
97  - In order to use a reference \c rho, set \c rho to \c rhoInf.
98 
99 SourceFiles
100  ReynoldsAnalogy.C
101 
102 \*---------------------------------------------------------------------------*/
103 
104 #ifndef heatTransferCoeffModels_ReynoldsAnalogy_H
105 #define heatTransferCoeffModels_ReynoldsAnalogy_H
106 
107 #include "heatTransferCoeffModel.H"
108 
109 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110 
111 namespace Foam
112 {
113 namespace heatTransferCoeffModels
114 {
115 
116 /*---------------------------------------------------------------------------*\
117  Class ReynoldsAnalogy Declaration
118 \*---------------------------------------------------------------------------*/
119 
120 class ReynoldsAnalogy
121 :
122  public heatTransferCoeffModel
123 {
124 protected:
125 
126  // Protected Data
127 
128  //- Name of velocity field
129  word UName_;
130 
131  //- Reference velocity
132  vector URef_;
133 
134  //- Name of fluid density field
135  word rhoName_;
136 
137  //- Reference fluid density
138  scalar rhoRef_;
139 
140  //- Name of specific heat capacity field
141  word CpName_;
142 
143  //- Reference specific heat capacity
144  scalar CpRef_;
145 
146 
147  // Protected Member Functions
148 
149  //- Return fluid density field [kg/m^3]
150  virtual tmp<scalarField> rho(const label patchi) const;
151 
152  //- Return heat capacity at constant pressure [J/kg/K]
153  virtual tmp<scalarField> Cp(const label patchi) const;
154 
155  //- Return the effective stress tensor including the laminar stress
156  virtual tmp<volSymmTensorField> devReff() const;
157 
158  //- Return skin friction coefficient field [-]
159  tmp<FieldField<Field, scalar>> Cf() const;
160 
161  //- Set the heat transfer coefficient
162  virtual void htc
163  (
165  const FieldField<Field, scalar>& q
166  );
167 
168 
169 public:
170 
171  //- Runtime type information
172  TypeName("ReynoldsAnalogy");
173 
174 
175  // Constructors
176 
177  //- Construct from components
179  (
180  const dictionary& dict,
181  const fvMesh& mesh,
182  const word& TName
183  );
184 
185  //- No copy construct
186  ReynoldsAnalogy(const ReynoldsAnalogy&) = delete;
188  //- No copy assignment
189  void operator=(const ReynoldsAnalogy&) = delete;
190 
191 
192  //- Destructor
193  virtual ~ReynoldsAnalogy() = default;
194 
195 
196  // Member Functions
197 
198  //- Read from dictionary
199  virtual bool read(const dictionary& dict);
200 };
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace heatTransferCoeffModels
206 } // End namespace Foam
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 #endif
211 
212 // ************************************************************************* //
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor including the laminar stress.
ReynoldsAnalogy(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
dictionary dict
virtual ~ReynoldsAnalogy()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
scalar rhoRef_
Reference fluid density.
const word & TName() const noexcept
Return const reference to name of temperature field.
virtual tmp< scalarField > Cp(const label patchi) const
Return heat capacity at constant pressure [J/kg/K].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual tmp< scalarField > rho(const label patchi) const
Return fluid density field [kg/m^3].
word rhoName_
Name of fluid density field.
TypeName("ReynoldsAnalogy")
Runtime type information.
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.
Heat transfer coefficient calculation based on Reynolds Analogy, which is used to relate turbulent mo...
Vector< scalar > vector
Definition: vector.H:57
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
word CpName_
Name of specific heat capacity field.
void operator=(const ReynoldsAnalogy &)=delete
No copy assignment.
scalar CpRef_
Reference specific heat capacity.
tmp< FieldField< Field, scalar > > q() const
Return boundary fields of heat-flux field.
virtual bool read(const dictionary &dict)
Read from dictionary.
Namespace for OpenFOAM.