wallHeatTransferCoeff.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) 2020-2022 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 
26 #include "wallHeatTransferCoeff.H"
29 #include "fvsPatchField.H"
30 #include "basicThermo.H"
31 #include "wallPolyPatch.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(wallHeatTransferCoeff, 0);
42  (
43  functionObject,
44  wallHeatTransferCoeff,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 (
55  const label i
56 )
57 {
58  // Add headers to output data
59  writeHeader(file(), "Wall heat transfer coefficient");
60  writeCommented(file(), "Time");
61  writeTabbed(file(), "patch");
62  writeTabbed(file(), "min");
63  writeTabbed(file(), "max");
64  writeTabbed(file(), "average");
65  file() << endl;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const word& name,
74  const Time& runTime,
75  const dictionary& dict
76 )
77 :
78  fvMeshFunctionObject(name, runTime, dict),
79  logFiles(obr_, name),
80  writeLocalObjects(obr_, log),
81  coeffModel_(wallHeatTransferCoeffModel::New(dict.name(), mesh_, dict)),
82  rho_("rho", dimDensity, Zero),
83  Cp_("Cp", dimArea/sqr(dimTime)/dimTemperature, Zero),
84  runTime_(runTime),
85  patchSet_()
86 {
87  read(dict);
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
103 
104  if (!foundObject<basicThermo>(physicalProperties::typeName))
105  {
106  rho_.read(dict);
107  Cp_.read(dict);
108  }
109 
110  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
111 
112  patchSet_ =
113  pbm.patchSet
114  (
115  wordReList(dict.lookupOrDefault("patches", wordReList()))
116  );
117 
118  Info<< type() << ":" << nl;
119 
120  if (patchSet_.empty())
121  {
122  forAll(pbm, patchi)
123  {
124  if (isA<wallPolyPatch>(pbm[patchi]))
125  {
126  patchSet_.insert(patchi);
127  }
128  }
129 
130  Info<< " processing all wall patches" << nl << endl;
131  }
132  else
133  {
134  Info<< " processing wall patches: " << nl;
135  labelHashSet filteredPatchSet;
136  forAllConstIter(labelHashSet, patchSet_, iter)
137  {
138  label patchi = iter.key();
139  if (isA<wallPolyPatch>(pbm[patchi]))
140  {
141  filteredPatchSet.insert(patchi);
142  Info<< " " << pbm[patchi].name() << endl;
143  }
144  else
145  {
147  << "Requested wall heat-transferCoeff on non-wall boundary"
148  << " type patch: " << pbm[patchi].name() << nl << endl;
149  }
150  }
151 
152  Info<< endl;
153 
154  patchSet_ = filteredPatchSet;
155  }
156 
157  coeffModel_->read(dict);
158 
159  resetName(typeName);
160  resetLocalObjectName(typeName);
161 
162  return true;
163 }
164 
165 
167 {
168  const momentumTransportModel& mmtm =
169  lookupObject<momentumTransportModel>
170  (
171  momentumTransportModel::typeName
172  );
173 
174  tmp<volScalarField> thtc;
175  thtc = coeffModel_->htcByRhoCp(mmtm, patchSet_);
176 
177  if (!foundObject<basicThermo>(physicalProperties::typeName))
178  {
179  thtc.ref() *= rho_*Cp_;
180  }
181  else
182  {
183  const basicThermo& thermo =
184  lookupObject<basicThermo>(physicalProperties::typeName);
185 
186  thtc.ref() *= thermo.rho()*thermo.Cp();
187  }
188 
189  store("wallHeatTransferCoeff", thtc);
190 
191  return true;
192 }
193 
194 
196 {
197  Log << name() << " write:" << nl;
198 
200  logFiles::write();
201 
202  const volScalarField& htc = obr_.lookupObject<volScalarField>(type());
203 
204  const fvPatchList& patches = mesh_.boundary();
205 
206  const surfaceScalarField::Boundary& magSf =
207  mesh_.magSf().boundaryField();
208 
209  forAllConstIter(labelHashSet, patchSet_, iter)
210  {
211  label patchi = iter.key();
212  const fvPatch& pp = patches[patchi];
213 
214  const scalarField& hfp = htc.boundaryField()[patchi];
215 
216  const scalar minHtcp = gMin(hfp);
217  const scalar maxHtcp = gMax(hfp);
218  const scalar averageHtcp =
219  gSum(magSf[patchi]*hfp)/gSum(magSf[patchi]);
220 
221  if (Pstream::master())
222  {
223  file()
224  << mesh_.time().userTimeValue()
225  << tab << pp.name()
226  << tab << minHtcp
227  << tab << maxHtcp
228  << tab << averageHtcp
229  << endl;
230  }
231 
232  Log << " min/max/average(" << pp.name() << ") = "
233  << minHtcp << ", " << maxHtcp << ", " << averageHtcp << endl;
234  }
235 
236  Log << endl;
237 
238  return true;
239 }
240 
241 
242 // ************************************************************************* //
const fvPatchList & patches
virtual bool write()
Write function.
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:103
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write function.
Definition: logFiles.C:167
const dimensionSet dimArea
const word & name() const
Return name.
Definition: IOobject.H:315
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
static const char tab
Definition: Ostream.H:259
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
virtual bool write()
Write the wall heat transfer coefficient.
Type gMin(const FieldField< Field, Type > &f)
virtual bool read(const dictionary &)
Read the wallHeatTransferCoeffs data.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual void writeFileHeader(const label i)
File header information.
static autoPtr< wallHeatTransferCoeffModel > New(const word &name, const fvMesh &mesh, const dictionary &)
Return a reference to the selected subset.
Macros for easy insertion into run-time selection tables.
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
const dimensionSet dimTime
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:109
Type gSum(const FieldField< Field, Type > &f)
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimDensity
static const zero Zero
Definition: zero.H:97
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
wallHeatTransferCoeff(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, mesh and dict.
Abstract base class for turbulence models (RAS, LES and laminar).
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
FunctionObject base class for managing a list of objects on behalf of the inheriting function object...
virtual bool read(const dictionary &)
Read the list of objects to be written.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual bool execute()
Calculate the wall heat transfer coefficient.
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionSet dimTemperature
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57