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) 2017-2018 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"
28 #include "wallPolyPatch.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
37  defineTypeNameAndDebug(wallHeatTransferCoeff, 0);
39  (
40  functionObject,
41  wallHeatTransferCoeff,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
49 
51 (
52  const label i
53 )
54 {
55  // Add headers to output data
56  writeHeader(file(), "Wall heat transfer coefficient");
57  writeCommented(file(), "Time");
58  writeTabbed(file(), "patch");
59  writeTabbed(file(), "min");
60  writeTabbed(file(), "max");
61  writeTabbed(file(), "integral");
62  file() << endl;
63 }
64 
65 
68 (
69  const volScalarField& nu,
70  const volScalarField& nut
71 )
72 {
73  tmp<volScalarField> twallHeatTransferCoeff
74  (
75  new volScalarField
76  (
77  IOobject
78  (
79  type(),
80  mesh_.time().timeName(),
81  mesh_
82  ),
83  mesh_,
85  (
86  "0",
88  0
89  )
90  )
91  );
92 
93  volScalarField::Boundary& wallHeatTransferCoeffBf =
94  twallHeatTransferCoeff.ref().boundaryFieldRef();
95 
96  const volScalarField::Boundary& nuBf = nu.boundaryField();
97  const volScalarField::Boundary& nutBf = nut.boundaryField();
98 
99  forAll(wallHeatTransferCoeffBf, patchi)
100  {
101  if (!wallHeatTransferCoeffBf[patchi].coupled())
102  {
103  wallHeatTransferCoeffBf[patchi] =
104  rho_*Cp_*(nuBf[patchi]/Prl_ + nutBf[patchi]/Prt_);
105  }
106  }
107 
108  return twallHeatTransferCoeff;
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
114 Foam::functionObjects::wallHeatTransferCoeff::wallHeatTransferCoeff
115 (
116  const word& name,
117  const Time& runTime,
118  const dictionary& dict
119 )
120 :
121  fvMeshFunctionObject(name, runTime, dict),
122  logFiles(obr_, name),
123  writeLocalObjects(obr_, log),
124  patchSet_()
125 {
126  read(dict);
127  resetName(typeName);
128  resetLocalObjectName(typeName);
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
133 
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
144 
145  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
146 
147  patchSet_ =
148  mesh_.boundaryMesh().patchSet
149  (
150  wordReList(dict.lookupOrDefault("patches", wordReList()))
151  );
152 
153  Info<< type() << " " << name() << ":" << nl;
154 
155  if (patchSet_.empty())
156  {
157  forAll(pbm, patchi)
158  {
159  if (isA<wallPolyPatch>(pbm[patchi]))
160  {
161  patchSet_.insert(patchi);
162  }
163  }
164 
165  Info<< " processing all wall patches" << nl << endl;
166  }
167  else
168  {
169  Info<< " processing wall patches: " << nl;
170  labelHashSet filteredPatchSet;
171  forAllConstIter(labelHashSet, patchSet_, iter)
172  {
173  label patchi = iter.key();
174  if (isA<wallPolyPatch>(pbm[patchi]))
175  {
176  filteredPatchSet.insert(patchi);
177  Info<< " " << pbm[patchi].name() << endl;
178  }
179  else
180  {
182  << "Requested wall heat-transferCoeff on non-wall boundary "
183  << "type patch: " << pbm[patchi].name() << endl;
184  }
185  }
186 
187  Info<< endl;
188 
189  patchSet_ = filteredPatchSet;
190  }
191 
192  dict.lookup("rho") >> rho_;
193  dict.lookup("Cp") >> Cp_;
194  dict.lookup("Prl") >> Prl_;
195  dict.lookup("Prt") >> Prt_;
196 
197  return true;
198 }
199 
200 
202 {
203  word name(type());
204 
205  if
206  (
207  foundObject<incompressible::turbulenceModel>
208  (
210  )
211  )
212  {
213  const incompressible::turbulenceModel& turbModel =
214  lookupObject<incompressible::turbulenceModel>
215  (
217  );
218 
219  return store
220  (
221  name,
222  calcHeatTransferCoeff(turbModel.nu(), turbModel.nut())
223  );
224  }
225  else
226  {
228  << "Unable to find incompressible turbulence model in the "
229  << "database" << exit(FatalError);
230 
231  return false;
232  }
233 }
234 
235 
237 {
238  Log << type() << " " << name() << " write:" << nl;
239 
241 
242  logFiles::write();
243 
246 
247  const fvPatchList& patches = mesh_.boundary();
248 
249  const surfaceScalarField::Boundary& magSf =
250  mesh_.magSf().boundaryField();
251 
252  forAllConstIter(labelHashSet, patchSet_, iter)
253  {
254  label patchi = iter.key();
255  const fvPatch& pp = patches[patchi];
256 
257  const scalarField& hfp = wallHeatTransferCoeff.boundaryField()[patchi];
258 
259  const scalar minHtcp = gMin(hfp);
260  const scalar maxHtcp = gMax(hfp);
261  const scalar integralHtcp = gSum(magSf[patchi]*hfp);
262 
263  if (Pstream::master())
264  {
265  file()
266  << mesh_.time().value()
267  << tab << pp.name()
268  << tab << minHtcp
269  << tab << maxHtcp
270  << tab << integralHtcp
271  << endl;
272  }
273 
274  Log << " min/max/integ(" << pp.name() << ") = "
275  << minHtcp << ", " << maxHtcp << ", " << integralHtcp << endl;
276  }
277 
278  Log << endl;
279 
280  return true;
281 }
282 
283 
284 // ************************************************************************* //
virtual bool write()
Write function.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool write()
Write function.
Definition: logFiles.C:180
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 word & name() const
Return name.
Definition: IOobject.H:297
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:264
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool write()
Write the wall heat transfer coefficient.
Type gMin(const FieldField< Field, Type > &f)
virtual bool read(const dictionary &)
Read the wallHeatTransferCoeff 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:174
Calculates and write the estimated incompressible flow heat transfer coefficient at wall patches as t...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
tmp< volScalarField > calcHeatTransferCoeff(const volScalarField &nu, const volScalarField &nut)
Calculate the heat transfer coefficient.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:421
patches[0]
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
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.
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
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.
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Type gSum(const FieldField< Field, Type > &f)
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:265
Type gMax(const FieldField< Field, Type > &f)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
dimensionedScalar pow3(const dimensionedScalar &ds)
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
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 execute()
Calculate the wall heat transfer coefficient.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.