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-2021 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 "basicThermo.H"
30 #include "wallPolyPatch.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(wallHeatTransferCoeff, 0);
41  (
42  functionObject,
43  wallHeatTransferCoeff,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51 
53 (
54  const label i
55 )
56 {
57  // Add headers to output data
58  writeHeader(file(), "Wall heat transfer coefficient");
59  writeCommented(file(), "Time");
60  writeTabbed(file(), "patch");
61  writeTabbed(file(), "min");
62  writeTabbed(file(), "max");
63  writeTabbed(file(), "average");
64  file() << endl;
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const word& name,
73  const Time& runTime,
74  const dictionary& dict
75 )
76 :
77  fvMeshFunctionObject(name, runTime, dict),
78  logFiles(obr_, name),
79  writeLocalObjects(obr_, log),
80  coeffModel_(wallHeatTransferCoeffModel::New(dict.name(), mesh_, dict)),
81  rho_("rho", dimDensity, Zero),
82  Cp_("Cp", dimArea/sqr(dimTime)/dimTemperature, Zero),
83  runTime_(runTime),
84  patchSet_()
85 {
86  read(dict);
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
102 
103  if (!foundObject<basicThermo>(basicThermo::dictName))
104  {
105  rho_.read(dict);
106  Cp_.read(dict);
107  }
108 
109  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
110 
111  patchSet_ =
112  pbm.patchSet
113  (
114  wordReList(dict.lookupOrDefault("patches", wordReList()))
115  );
116 
117  Info<< type() << ":" << nl;
118 
119  if (patchSet_.empty())
120  {
121  forAll(pbm, patchi)
122  {
123  if (isA<wallPolyPatch>(pbm[patchi]))
124  {
125  patchSet_.insert(patchi);
126  }
127  }
128 
129  Info<< " processing all wall patches" << nl << endl;
130  }
131  else
132  {
133  Info<< " processing wall patches: " << nl;
134  labelHashSet filteredPatchSet;
135  forAllConstIter(labelHashSet, patchSet_, iter)
136  {
137  label patchi = iter.key();
138  if (isA<wallPolyPatch>(pbm[patchi]))
139  {
140  filteredPatchSet.insert(patchi);
141  Info<< " " << pbm[patchi].name() << endl;
142  }
143  else
144  {
146  << "Requested wall heat-transferCoeff on non-wall boundary"
147  << " type patch: " << pbm[patchi].name() << nl << endl;
148  }
149  }
150 
151  Info<< endl;
152 
153  patchSet_ = filteredPatchSet;
154  }
155 
156  coeffModel_->read(dict);
157 
158  resetName(typeName);
159  resetLocalObjectName(typeName);
160 
161  return true;
162 }
163 
164 
166 {
167  const momentumTransportModel& mmtm =
168  lookupObject<momentumTransportModel>
169  (
170  momentumTransportModel::typeName
171  );
172 
173  tmp<volScalarField> thtc;
174  thtc = coeffModel_->htcByRhoCp(mmtm, patchSet_);
175 
176  if (!foundObject<basicThermo>(basicThermo::dictName))
177  {
178  thtc.ref() *= rho_*Cp_;
179  }
180  else
181  {
182  const basicThermo& thermo =
183  lookupObject<basicThermo>(basicThermo::dictName);
184 
185  thtc.ref() *= thermo.rho()*thermo.Cp();
186  }
187 
188  store("wallHeatTransferCoeff", thtc);
189 
190  return true;
191 }
192 
193 
195 {
196  Log << name() << " write:" << nl;
197 
199  logFiles::write();
200 
201  const volScalarField& htc = obr_.lookupObject<volScalarField>(type());
202 
203  const fvPatchList& patches = mesh_.boundary();
204 
205  const surfaceScalarField::Boundary& magSf =
206  mesh_.magSf().boundaryField();
207 
208  forAllConstIter(labelHashSet, patchSet_, iter)
209  {
210  label patchi = iter.key();
211  const fvPatch& pp = patches[patchi];
212 
213  const scalarField& hfp = htc.boundaryField()[patchi];
214 
215  const scalar minHtcp = gMin(hfp);
216  const scalar maxHtcp = gMax(hfp);
217  const scalar averageHtcp =
218  gSum(magSf[patchi]*hfp)/gSum(magSf[patchi]);
219 
220  if (Pstream::master())
221  {
222  file()
223  << mesh_.time().value()
224  << tab << pp.name()
225  << tab << minHtcp
226  << tab << maxHtcp
227  << tab << averageHtcp
228  << endl;
229  }
230 
231  Log << " min/max/average(" << pp.name() << ") = "
232  << minHtcp << ", " << maxHtcp << ", " << averageHtcp << endl;
233  }
234 
235  Log << endl;
236 
237  return true;
238 }
239 
240 
241 // ************************************************************************* //
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
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
const dimensionSet dimArea
const word & name() const
Return name.
Definition: IOobject.H:303
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
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:62
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
patches[0]
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:68
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.
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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
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:144
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