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-2024 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 {
42  (
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_ = patchSet(dict, true);
113 
114  Info<< type() << ":" << nl;
115 
116  if (patchSet_.empty())
117  {
118  forAll(pbm, patchi)
119  {
120  if (isA<wallPolyPatch>(pbm[patchi]))
121  {
122  patchSet_.insert(patchi);
123  }
124  }
125 
126  Info<< " processing all wall patches" << nl << endl;
127  }
128  else
129  {
130  Info<< " processing wall patches: " << nl;
131  labelHashSet filteredPatchSet;
132  forAllConstIter(labelHashSet, patchSet_, iter)
133  {
134  label patchi = iter.key();
135  if (isA<wallPolyPatch>(pbm[patchi]))
136  {
137  filteredPatchSet.insert(patchi);
138  Info<< " " << pbm[patchi].name() << endl;
139  }
140  else
141  {
143  << "Requested wall heat-transferCoeff on non-wall boundary"
144  << " type patch: " << pbm[patchi].name() << nl << endl;
145  }
146  }
147 
148  Info<< endl;
149 
150  patchSet_ = filteredPatchSet;
151  }
152 
153  coeffModel_->read(dict);
154 
155  resetName(typeName);
156  resetLocalObjectName(typeName);
157 
158  return true;
159 }
160 
161 
163 {
164  const momentumTransportModel& mmtm =
165  lookupObject<momentumTransportModel>
166  (
167  momentumTransportModel::typeName
168  );
169 
170  tmp<volScalarField> thtc;
171  thtc = coeffModel_->htcByRhoCp(mmtm, patchSet_);
172 
173  if (!foundObject<basicThermo>(physicalProperties::typeName))
174  {
175  thtc.ref() *= rho_*Cp_;
176  }
177  else
178  {
179  const basicThermo& thermo =
180  lookupObject<basicThermo>(physicalProperties::typeName);
181 
182  thtc.ref() *= thermo.rho()*thermo.Cp();
183  }
184 
185  store("wallHeatTransferCoeff", thtc);
186 
187  return true;
188 }
189 
190 
192 {
193  Log << name() << " write:" << nl;
194 
196  logFiles::write();
197 
198  const volScalarField& htc = obr_.lookupObject<volScalarField>(type());
199 
200  const fvPatchList& patches = mesh_.boundary();
201 
202  const surfaceScalarField::Boundary& magSf =
203  mesh_.magSf().boundaryField();
204 
205  forAllConstIter(labelHashSet, patchSet_, iter)
206  {
207  label patchi = iter.key();
208  const fvPatch& pp = patches[patchi];
209 
210  const scalarField& hfp = htc.boundaryField()[patchi];
211 
212  const scalar minHtcp = gMin(hfp);
213  const scalar maxHtcp = gMax(hfp);
214  const scalar averageHtcp =
215  gSum(magSf[patchi]*hfp)/gSum(magSf[patchi]);
216 
217  if (Pstream::master())
218  {
219  file()
220  << time_.userTimeValue()
221  << tab << pp.name()
222  << tab << minHtcp
223  << tab << maxHtcp
224  << tab << averageHtcp
225  << endl;
226  }
227 
228  Log << " min/max/average(" << pp.name() << ") = "
229  << minHtcp << ", " << maxHtcp << ", " << averageHtcp << endl;
230  }
231 
232  Log << endl;
233 
234  return true;
235 }
236 
237 
238 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
const word & name() const
Return name.
Definition: IOobject.H:310
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:60
OFstream & file()
Return access to the file (if only 1)
Definition: logFiles.C:113
virtual bool write()
Write function.
Definition: logFiles.C:173
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:106
virtual bool read(const dictionary &)
Read optional controls.
Calculates and writes the estimated heat transfer coefficient at wall patches as the volScalarField f...
wallHeatTransferCoeff(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, mesh and dict.
virtual void writeFileHeader(const label i)
File header information.
virtual bool execute()
Calculate the wall heat transfer coefficient.
virtual bool write()
Write the wall heat transfer coefficient.
virtual bool read(const dictionary &)
Read the wallHeatTransferCoeffs data.
void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:121
void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:131
void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:110
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.
virtual bool write()
Write function.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
Abstract base class for turbulence models (RAS, LES and laminar).
Foam::polyBoundaryMesh.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for run time selection of heat transfer coefficient models.
A class for handling words, derived from string.
Definition: word.H:62
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char tab
Definition: Ostream.H:265
messageStream Info
const dimensionSet dimTemperature
const dimensionSet dimTime
const dimensionSet dimDensity
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const dimensionSet dimArea
Type gMin(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:266
Type gMax(const FieldField< Field, Type > &f)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31