convectiveHeatTransferFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 
28 #include "fvPatchFieldMapper.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace compressible
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedValueFvPatchScalarField(p, iF),
48  L_(1.0)
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
62  L_(ptf.L_)
63 {}
64 
65 
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
74  fixedValueFvPatchScalarField(p, iF, dict),
75  L_(readScalar(dict.lookup("L")))
76 {}
77 
78 
81 (
83 )
84 :
85  fixedValueFvPatchScalarField(htcpsf),
86  L_(htcpsf.L_)
87 {}
88 
89 
92 (
95 )
96 :
97  fixedValueFvPatchScalarField(htcpsf, iF),
98  L_(htcpsf.L_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  if (updated())
107  {
108  return;
109  }
110 
111  const label patchi = patch().index();
112 
113  const compressible::turbulenceModel& turbModel =
114  db().lookupObject<compressible::turbulenceModel>
115  (
117  (
118  compressible::turbulenceModel::propertiesName,
119  internalField().group()
120  )
121  );
122 
123  const scalarField alphaEffw(turbModel.alphaEff(patchi));
124 
125  const tmp<scalarField> tmuw = turbModel.mu(patchi);
126  const scalarField& muw = tmuw();
127 
128  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
129  const vectorField& Uc = turbModel.U();
130  const vectorField& Uw = turbModel.U().boundaryField()[patchi];
131  const scalarField& Tw = turbModel.transport().T().boundaryField()[patchi];
132  const scalarField& pw = turbModel.transport().p().boundaryField()[patchi];
133  const scalarField Cpw(turbModel.transport().Cp(pw, Tw, patchi));
134 
135  const scalarField kappaw(Cpw*alphaEffw);
136  const scalarField Pr(muw*Cpw/kappaw);
137 
138  scalarField& htc = *this;
139  forAll(htc, facei)
140  {
141  label faceCelli = patch().faceCells()[facei];
142 
143  scalar Re = rhow[facei]*mag(Uc[faceCelli] - Uw[facei])*L_/muw[facei];
144 
145  if (Re < 5.0E+05)
146  {
147  htc[facei] = 0.664*sqrt(Re)*cbrt(Pr[facei])*kappaw[facei]/L_;
148  }
149  else
150  {
151  htc[facei] = 0.037*pow(Re, 0.8)*cbrt(Pr[facei])*kappaw[facei]/L_;
152  }
153  }
154 
155  fixedValueFvPatchScalarField::updateCoeffs();
156 }
157 
158 
159 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
160 
162 {
164  os.writeKeyword("L") << L_ << token::END_STATEMENT << nl;
165  writeEntry("value", os);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
172 (
175 );
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace compressible
180 } // End namespace Foam
181 
182 // ************************************************************************* //
const char *const group
Group name for atomic constants.
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:137
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
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:65
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
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.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
convectiveHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
label patchi
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:54
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
bool compressible
Definition: pEqn.H:30