heatTransferCoeffModel.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "heatTransferCoeffModel.H"
29 #include "fvMesh.H"
30 #include "fluidThermo.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(heatTransferCoeffModel, 0);
39  defineRunTimeSelectionTable(heatTransferCoeffModel, dictionary);
40 }
41 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
47 {
48  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
49  const volScalarField::Boundary& Tbf = T.boundaryField();
50 
51  auto tq = tmp<FieldField<Field, scalar>>::New(Tbf.size());
52  auto& q = tq.ref();
53 
54  forAll(q, patchi)
55  {
56  q.set(patchi, new Field<scalar>(Tbf[patchi].size(), Zero));
57  }
58 
59  typedef compressible::turbulenceModel cmpTurbModel;
60 
61  if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
62  {
63  const auto& turb =
64  mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
65 
66  // Note: calling he(p,T) instead of he()
67  const volScalarField he(turb.transport().he(turb.transport().p(), T));
68  const volScalarField::Boundary& hebf = he.boundaryField();
69 
70  const volScalarField alphaEff(turb.alphaEff());
71  const volScalarField::Boundary& alphaEffbf = alphaEff.boundaryField();
72 
73  for (const label patchi : patchSet_)
74  {
75  q[patchi] = alphaEffbf[patchi]*hebf[patchi].snGrad();
76  }
77  }
79  {
80  const auto& thermo =
82 
83  // Note: calling he(p,T) instead of he()
84  const volScalarField he(thermo.he(thermo.p(), T));
85  const volScalarField::Boundary& hebf = he.boundaryField();
86 
87  const volScalarField& alpha(thermo.alpha());
88  const volScalarField::Boundary& alphabf = alpha.boundaryField();
89 
90  for (const label patchi : patchSet_)
91  {
92  q[patchi] = alphabf[patchi]*hebf[patchi].snGrad();
93  }
94  }
95  else
96  {
98  << "Unable to find a valid thermo model to evaluate q. " << nl
99  << "Database contents are: " << mesh_.objectRegistry::sortedToc()
100  << exit(FatalError);
101  }
102 
103  // Add radiative heat flux contribution if present
104 
105  const auto* qrPtr = mesh_.cfindObject<volScalarField>(qrName_);
106 
107  if (qrPtr)
108  {
109  const volScalarField::Boundary& qrbf = qrPtr->boundaryField();
110 
111  for (const label patchi : patchSet_)
112  {
113  q[patchi] += qrbf[patchi];
114  }
115  }
116 
117  return tq;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
124 (
125  const dictionary& dict,
126  const fvMesh& mesh,
127  const word& TName
128 )
129 :
130  mesh_(mesh),
131  patchSet_(),
132  TName_(TName),
133  qrName_("qr")
134 {}
135 
136 
137 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
138 
140 {
141  patchSet_ = mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"));
142 
143  dict.readIfPresent("qr", qrName_);
144 
145  return true;
146 }
147 
148 
150 (
151  volScalarField& result,
153 )
154 {
155  htc(result, q);
156 
157  return true;
158 }
159 
160 
161 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::heatTransferCoeffModel::mesh_
const fvMesh & mesh_
Mesh reference.
Definition: heatTransferCoeffModel.H:111
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
fluidThermo.H
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
turbulentTransportModel.H
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
Foam::heatTransferCoeffModel::heatTransferCoeffModel
heatTransferCoeffModel(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
Definition: heatTransferCoeffModel.C:124
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::heatTransferCoeffModel::patchSet_
labelHashSet patchSet_
Optional list of (wall) patches to process.
Definition: heatTransferCoeffModel.H:114
Foam::heatTransferCoeffModel::calc
virtual bool calc(volScalarField &result, const FieldField< Field, scalar > &q)
Definition: heatTransferCoeffModel.C:150
Foam::Field< scalar >
Foam::heatTransferCoeffModel::New
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
Definition: heatTransferCoeffModelNew.C:35
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
heatTransferCoeffModel.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::objectRegistry::cfindObject
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:390
he
volScalarField & he
Definition: YEEqn.H:52
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::heatTransferCoeffModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: heatTransferCoeffModel.C:139
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::heatTransferCoeffModel::TName_
const word TName_
Temperature name.
Definition: heatTransferCoeffModel.H:117
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
turbulentFluidThermoModel.H
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::heatTransferCoeffModel::q
tmp< FieldField< Field, scalar > > q() const
Return q boundary fields.
Definition: heatTransferCoeffModel.C:46
Foam::heatTransferCoeffModel::qrName_
word qrName_
Name of radiative heat flux.
Definition: heatTransferCoeffModel.H:120
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60