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-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 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  fixedValueFvPatchScalarField(p, iF),
47  L_(1.0)
48 {}
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  fixedValueFvPatchScalarField(p, iF, dict),
60  L_(dict.lookup<scalar>("L"))
61 {}
62 
63 
66 (
68  const fvPatch& p,
70  const fvPatchFieldMapper& mapper
71 )
72 :
73  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
74  L_(ptf.L_)
75 {}
76 
77 
80 (
83 )
84 :
85  fixedValueFvPatchScalarField(htcpsf, iF),
86  L_(htcpsf.L_)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  if (updated())
95  {
96  return;
97  }
98 
99  const label patchi = patch().index();
100 
101  const thermophysicalTransportModel& ttm =
102  db().lookupObject<thermophysicalTransportModel>
103  (
105  (
106  thermophysicalTransportModel::typeName,
107  internalField().group()
108  )
109  );
110 
111  const compressibleMomentumTransportModel& turbModel =
112  ttm.momentumTransport();
113 
114  const scalarField alphaEffw(ttm.alphaEff(patchi));
115 
116  const tmp<scalarField> tnuw = turbModel.nu(patchi);
117  const scalarField& nuw = tnuw();
118 
119  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
120  const vectorField& Uc = turbModel.U();
121  const vectorField& Uw = turbModel.U().boundaryField()[patchi];
122  const scalarField& Tw = ttm.thermo().T().boundaryField()[patchi];
123  const scalarField Cpw(ttm.thermo().Cp(Tw, patchi));
124 
125  const scalarField kappaw(Cpw*alphaEffw);
126  const scalarField Pr(rhow*nuw*Cpw/kappaw);
127 
128  scalarField& htc = *this;
129  forAll(htc, facei)
130  {
131  const label celli = patch().faceCells()[facei];
132 
133  const scalar Re = mag(Uc[celli] - Uw[facei])*L_/nuw[facei];
134 
135  if (Re < 5.0E+05)
136  {
137  htc[facei] = 0.664*sqrt(Re)*cbrt(Pr[facei])*kappaw[facei]/L_;
138  }
139  else
140  {
141  htc[facei] = 0.037*pow(Re, 0.8)*cbrt(Pr[facei])*kappaw[facei]/L_;
142  }
143  }
144 
145  fixedValueFvPatchScalarField::updateCoeffs();
146 }
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
154  writeEntry(os, "L", L_);
155  writeEntry(os, "value", *this);
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
162 (
165 );
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace compressible
170 } // End namespace Foam
171 
172 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
This boundary condition provides a convective heat transfer coefficient condition.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
const compressibleMomentumTransportModel & momentumTransport() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for thermophysical transport models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
convectiveHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
virtual const volScalarField & T() const =0
Temperature [K].
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Base class for single-phase compressible turbulence models.
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864