convectiveHeatTransferFvPatchScalarField.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) 2011-2024 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 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44  const dictionary& dict
45 )
46 :
47  fixedValueFvPatchScalarField(p, iF, dict),
48  L_(dict.lookup<scalar>("L", dimLength))
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fieldMapper& mapper
59 )
60 :
61  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
62  L_(ptf.L_)
63 {}
64 
65 
68 (
71 )
72 :
73  fixedValueFvPatchScalarField(htcpsf, iF),
74  L_(htcpsf.L_)
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if (updated())
83  {
84  return;
85  }
86 
87  const label patchi = patch().index();
88 
90  db().lookupType<fluidThermophysicalTransportModel>
91  (
92  internalField().group()
93  );
94 
95  const compressibleMomentumTransportModel& turbModel =
96  ttm.momentumTransport();
97 
98  const tmp<scalarField> tnuw = turbModel.nu(patchi);
99  const scalarField& nuw = tnuw();
100 
101  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
102  const vectorField& Uc = turbModel.U();
103  const vectorField& Uw = turbModel.U().boundaryField()[patchi];
104  const scalarField& Cpw(ttm.thermo().Cp().boundaryField()[patchi]);
105 
106  const scalarField kappaEffw(ttm.kappaEff(patchi));
107  const scalarField Pr(rhow*nuw*Cpw/kappaEffw);
108 
109  scalarField& htc = *this;
110  forAll(htc, facei)
111  {
112  const label celli = patch().faceCells()[facei];
113 
114  const scalar Re = mag(Uc[celli] - Uw[facei])*L_/nuw[facei];
115 
116  if (Re < 5.0E+05)
117  {
118  htc[facei] = 0.664*sqrt(Re)*cbrt(Pr[facei])*kappaEffw[facei]/L_;
119  }
120  else
121  {
122  htc[facei] = 0.037*pow(Re, 0.8)*cbrt(Pr[facei])*kappaEffw[facei]/L_;
123  }
124  }
125 
126  fixedValueFvPatchScalarField::updateCoeffs();
127 }
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 {
135  writeEntry(os, "L", L_);
136  writeEntry(os, "value", *this);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
143 (
146 );
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace compressible
151 } // End namespace Foam
152 
153 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Base class for single-phase compressible turbulence models.
This boundary condition provides a convective heat transfer coefficient condition.
convectiveHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
const compressibleMomentumTransportModel & momentumTransport() const
Access function to momentum transport model.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for managing temporary objects.
Definition: tmp.H:55
label patchi
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
Namespace for OpenFOAM.
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:59
const dimensionSet dimLength
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar cbrt(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict
volScalarField & p