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-2020 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 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61  L_(ptf.L_)
62 {}
63 
64 
67 (
68  const fvPatch& p,
70  const dictionary& dict
71 )
72 :
73  fixedValueFvPatchScalarField(p, iF, dict),
74  L_(dict.lookup<scalar>("L"))
75 {}
76 
77 
80 (
82 )
83 :
84  fixedValueFvPatchScalarField(htcpsf),
85  L_(htcpsf.L_)
86 {}
87 
88 
91 (
94 )
95 :
96  fixedValueFvPatchScalarField(htcpsf, iF),
97  L_(htcpsf.L_)
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  if (updated())
106  {
107  return;
108  }
109 
110  const label patchi = patch().index();
111 
112  const thermophysicalTransportModel& ttm =
113  db().lookupObject<thermophysicalTransportModel>
114  (
116  (
117  thermophysicalTransportModel::typeName,
118  internalField().group()
119  )
120  );
121 
122  const compressibleMomentumTransportModel& turbModel =
123  ttm.momentumTransport();
124 
125  const scalarField alphaEffw(ttm.alphaEff(patchi));
126 
127  const tmp<scalarField> tmuw = turbModel.mu(patchi);
128  const scalarField& muw = tmuw();
129 
130  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
131  const vectorField& Uc = turbModel.U();
132  const vectorField& Uw = turbModel.U().boundaryField()[patchi];
133  const scalarField& Tw = ttm.thermo().T().boundaryField()[patchi];
134  const scalarField Cpw(ttm.thermo().Cp(Tw, patchi));
135 
136  const scalarField kappaw(Cpw*alphaEffw);
137  const scalarField Pr(muw*Cpw/kappaw);
138 
139  scalarField& htc = *this;
140  forAll(htc, facei)
141  {
142  label celli = patch().faceCells()[facei];
143 
144  scalar Re = rhow[facei]*mag(Uc[celli] - Uw[facei])*L_/muw[facei];
145 
146  if (Re < 5.0E+05)
147  {
148  htc[facei] = 0.664*sqrt(Re)*cbrt(Pr[facei])*kappaw[facei]/L_;
149  }
150  else
151  {
152  htc[facei] = 0.037*pow(Re, 0.8)*cbrt(Pr[facei])*kappaw[facei]/L_;
153  }
154  }
155 
156  fixedValueFvPatchScalarField::updateCoeffs();
157 }
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
163 {
165  writeEntry(os, "L", L_);
166  writeEntry(os, "value", *this);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
173 (
176 );
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace compressible
181 } // End namespace Foam
182 
183 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
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:158
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:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:482
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
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
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
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 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
Abstract base class for turbulence models (RAS, LES and laminar).
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:812