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 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 {
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 cmpTurbModel& turb =
64  mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
65 
66  const volScalarField& he = turb.transport().he();
67  const volScalarField::Boundary& hebf = he.boundaryField();
68 
69  const volScalarField alphaEff(turb.alphaEff());
70  const volScalarField::Boundary& alphaEffbf = alphaEff.boundaryField();
71 
72  for (const label patchi : patchSet_)
73  {
74  q[patchi] = alphaEffbf[patchi]*hebf[patchi].snGrad();
75  }
76  }
78  {
79  const fluidThermo& thermo =
81 
82  const volScalarField& he = thermo.he();
83  const volScalarField::Boundary& hebf = he.boundaryField();
84 
85  const volScalarField& alpha(thermo.alpha());
86  const volScalarField::Boundary& alphabf = alpha.boundaryField();
87 
88  for (label patchi : patchSet_)
89  {
90  q[patchi] = alphabf[patchi]*hebf[patchi].snGrad();
91  }
92  }
93  else
94  {
96  << "Unable to find a valid thermo model to evaluate q"
97  << exit(FatalError);
98  }
99 
100  // Add radiative heat flux contribution if present
102  {
104  const volScalarField::Boundary& qrbf = qr.boundaryField();
105 
106  for (const label patchi : patchSet_)
107  {
108  q[patchi] += qrbf[patchi];
109  }
110  }
111 
112  return tq;
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
118 Foam::heatTransferCoeffModel::heatTransferCoeffModel
119 (
120  const dictionary& dict,
121  const fvMesh& mesh,
122  const word& TName
123 )
124 :
125  mesh_(mesh),
126  TName_(TName),
127  patchSet_(),
128  qrName_("qr")
129 {}
130 
131 
132 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
133 
135 {
136  patchSet_ = mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"));
137 
138  dict.readIfPresent("qr", qrName_);
139 
140  return true;
141 }
142 
143 
145 {
146  htc(result);
147 
148  return true;
149 }
150 
151 
152 // ************************************************************************* //
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::heatTransferCoeffModel::mesh_
const fvMesh & mesh_
Definition: heatTransferCoeffModel.H:79
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:81
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::heatTransferCoeffModel::patchSet_
labelHashSet patchSet_
Definition: heatTransferCoeffModel.H:83
Foam::label
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:62
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::heatTransferCoeffModel::calc
virtual bool calc(volScalarField &result)
Definition: heatTransferCoeffModel.C:144
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
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
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:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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:355
he
volScalarField & he
Definition: YEEqn.H:52
Foam::heatTransferCoeffModel::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: heatTransferCoeffModel.C:134
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_
Definition: heatTransferCoeffModel.H:81
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:417
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
Definition: heatTransferCoeffModel.C:46
Foam::heatTransferCoeffModel::qrName_
word qrName_
Definition: heatTransferCoeffModel.H:85