ReynoldsAnalogy.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "ReynoldsAnalogy.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace wallHeatTransferCoeffModels
36 {
39  (
42  word
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const fvMesh& mesh,
54  const dictionary& dict
55 )
56 :
58  mesh_(mesh),
59  Uref_("Uref", dimVelocity, dict)
60 {
61  read(dict);
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
74 (
75  const dictionary& dict
76 )
77 {
78  Uref_.read(dict);
79  return true;
80 }
81 
82 
85 (
86  const momentumTransportModel& mmtm,
87  const labelHashSet& patches
88 ) const
89 {
90  tmp<volSymmTensorField> ttau(this->tau(mmtm, mesh_));
91  const volSymmTensorField::Boundary& tauBf = ttau.ref().boundaryField();
92 
93  // Create temporary field for heat transfer coefficient
94  tmp<volScalarField> thtcByRhoCp
95  (
97  (
98  type(),
99  mesh_,
101  )
102  );
103 
104  volScalarField::Boundary& thtcByRhoCpBf =
105  thtcByRhoCp.ref().boundaryFieldRef();
106 
107  forAll(thtcByRhoCpBf, patchi)
108  {
109  if (!thtcByRhoCpBf[patchi].coupled())
110  {
111  const vectorField nf(tauBf[patchi].patch().nf());
112 
113  // Wall shear stress in [m^2/s^2]
114  tmp<vectorField> tauwp(-nf&tauBf[patchi]);
115 
116  // Non-dimensional skin friction coefficient [-]
117  const scalarField Cf(2*mag(tauwp)/sqr(Uref_.value()));
118 
119  thtcByRhoCpBf[patchi] = 0.5*Uref_.value()*Cf;
120  }
121  }
122 
123  return thtcByRhoCp;
124 }
125 
126 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base class for turbulence models (RAS, LES and laminar).
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for run time selection of heat transfer coefficient models.
Calculates and writes the estimated flow heat transfer coefficient at wall patches as the volScalarFi...
virtual tmp< volScalarField > htcByRhoCp(const momentumTransportModel &mmtm, const labelHashSet &patches) const
Calculate the heat transfer coefficient.
ReynoldsAnalogy(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from name, mesh and dict.
virtual bool read(const dictionary &)
Read the ReynoldsAnalogy data.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const fvPatchList & patches
addToRunTimeSelectionTable(wallHeatTransferCoeffModel, kappaEff, word)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict