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-2023 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);
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  const fvMesh& mesh,
49  const word& TName
50 )
51 :
52  mesh_(mesh),
53  TName_(TName),
54  qrName_("qr")
55 {}
56 
57 
58 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
59 
62 {
63  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
64  const volScalarField::Boundary& Tbf = T.boundaryField();
65 
66  auto tq = tmp<FieldField<Field, scalar>>::New(Tbf.size());
67  auto& q = tq.ref();
68 
69  forAll(q, patchi)
70  {
71  q.set(patchi, new Field<scalar>(Tbf[patchi].size(), Zero));
72  }
73 
74  typedef compressible::turbulenceModel cmpTurbModel;
75 
76  if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
77  {
78  const auto& turb =
79  mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
80 
81  // Note: calling he(p,T) instead of he()
82  const volScalarField he(turb.transport().he(turb.transport().p(), T));
83  const volScalarField::Boundary& hebf = he.boundaryField();
84 
85  const volScalarField alphaEff(turb.alphaEff());
86  const volScalarField::Boundary& alphaEffbf = alphaEff.boundaryField();
87 
88  for (const label patchi : patchIDs_)
89  {
90  q[patchi] = alphaEffbf[patchi]*hebf[patchi].snGrad();
91  }
92  }
93  else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
94  {
95  const auto& thermo =
97 
98  // Note: calling he(p,T) instead of he()
99  const volScalarField he(thermo.he(thermo.p(), T));
100  const volScalarField::Boundary& hebf = he.boundaryField();
101 
102  const volScalarField& alpha(thermo.alpha());
103  const volScalarField::Boundary& alphabf = alpha.boundaryField();
104 
105  for (const label patchi : patchIDs_)
106  {
107  q[patchi] = alphabf[patchi]*hebf[patchi].snGrad();
108  }
109  }
110  else
111  {
113  << "Unable to find a valid thermo model to evaluate q. " << nl
114  << "Database contents are: " << mesh_.objectRegistry::sortedToc()
115  << exit(FatalError);
116  }
117 
118  // Add radiative heat flux contribution if present
119 
120  const auto* qrPtr = mesh_.cfindObject<volScalarField>(qrName_);
121 
122  if (qrPtr)
123  {
124  const volScalarField::Boundary& qrbf = qrPtr->boundaryField();
125 
126  for (const label patchi : patchIDs_)
127  {
128  q[patchi] += qrbf[patchi];
129  }
130  }
131 
132  return tq;
133 }
134 
135 
137 (
138  volScalarField& result,
139  const FieldField<Field, scalar>& q
140 )
141 {
142  htc(result, q);
143 
144  return true;
145 }
146 
147 
149 {
150  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
151 
152  dict.readIfPresent("qr", qrName_);
153 
154  patchIDs_ = pbm.patchSet(dict.get<wordRes>("patches")).sortedToc();
155 
156  return true;
157 }
158 
159 
160 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const polyBoundaryMesh & pbm
dictionary dict
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
word qrName_
Name of radiative heat flux field.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
labelHashSet patchSet(const UList< wordRe > &select, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
virtual bool read(const dictionary &dict)
Read from dictionary.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A base class for heat transfer coefficient models.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
dynamicFvMesh & mesh
virtual bool calc(volScalarField &result, const FieldField< Field, scalar > &q)
Calculate the heat transfer coefficient field and return true if successful.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
const fvMesh & mesh_
Const reference to the mesh.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const word TName_
Name of temperature field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
heatTransferCoeffModel(const dictionary &dict, const fvMesh &mesh, const word &TName)
Construct from components.
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
tmp< FieldField< Field, scalar > > q() const
Return boundary fields of heat-flux field.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
labelList patchIDs_
List of (wall) patches to process (selected by name)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127