wallHeatFlux.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) 2016-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 "wallHeatFlux.H"
28 #include "solidThermo.H"
29 #include "surfaceInterpolate.H"
30 #include "fvcSnGrad.H"
31 #include "wallPolyPatch.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(wallHeatFlux, 0);
41  addToRunTimeSelectionTable(functionObject, wallHeatFlux, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
47 
49 {
50  // Add headers to output data
51  writeHeader(file(), "Wall heat-flux");
52  writeCommented(file(), "Time");
53  writeTabbed(file(), "patch");
54  writeTabbed(file(), "min");
55  writeTabbed(file(), "max");
56  writeTabbed(file(), "integral");
57  file() << endl;
58 }
59 
60 
63 (
64  const volScalarField& alpha,
65  const volScalarField& he
66 )
67 {
68  tmp<volScalarField> twallHeatFlux
69  (
70  new volScalarField
71  (
72  IOobject
73  (
74  type(),
75  mesh_.time().timeName(),
76  mesh_
77  ),
78  mesh_,
80  )
81  );
82 
83  volScalarField::Boundary& wallHeatFluxBf =
84  twallHeatFlux.ref().boundaryFieldRef();
85 
86  const volScalarField::Boundary& heBf =
87  he.boundaryField();
88 
89  const volScalarField::Boundary& alphaBf =
90  alpha.boundaryField();
91 
92  forAll(wallHeatFluxBf, patchi)
93  {
94  if (!wallHeatFluxBf[patchi].coupled())
95  {
96  wallHeatFluxBf[patchi] = alphaBf[patchi]*heBf[patchi].snGrad();
97  }
98  }
99 
100  if (foundObject<volScalarField>("qr"))
101  {
102  const volScalarField& qr = lookupObject<volScalarField>("qr");
103 
104  const volScalarField::Boundary& radHeatFluxBf =
105  qr.boundaryField();
106 
107  forAll(wallHeatFluxBf, patchi)
108  {
109  if (!wallHeatFluxBf[patchi].coupled())
110  {
111  wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
112  }
113  }
114  }
115 
116  return twallHeatFlux;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
122 Foam::functionObjects::wallHeatFlux::wallHeatFlux
123 (
124  const word& name,
125  const Time& runTime,
126  const dictionary& dict
127 )
128 :
129  fvMeshFunctionObject(name, runTime, dict),
130  logFiles(obr_, name),
131  writeLocalObjects(obr_, log),
132  patchSet_()
133 {
134  read(dict);
135  resetName(typeName);
136  resetLocalObjectName(typeName);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
141 
143 {}
144 
145 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
149 {
152 
153  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
154 
155  patchSet_ =
156  mesh_.boundaryMesh().patchSet
157  (
158  wordReList(dict.lookupOrDefault("patches", wordReList()))
159  );
160 
161  Info<< type() << " " << name() << ":" << nl;
162 
163  if (patchSet_.empty())
164  {
165  forAll(pbm, patchi)
166  {
167  if (isA<wallPolyPatch>(pbm[patchi]))
168  {
169  patchSet_.insert(patchi);
170  }
171  }
172 
173  Info<< " processing all wall patches" << nl << endl;
174  }
175  else
176  {
177  Info<< " processing wall patches: " << nl;
178  labelHashSet filteredPatchSet;
179  forAllConstIter(labelHashSet, patchSet_, iter)
180  {
181  label patchi = iter.key();
182  if (isA<wallPolyPatch>(pbm[patchi]))
183  {
184  filteredPatchSet.insert(patchi);
185  Info<< " " << pbm[patchi].name() << endl;
186  }
187  else
188  {
190  << "Requested wall heat-flux on non-wall boundary "
191  << "type patch: " << pbm[patchi].name() << endl;
192  }
193  }
194 
195  Info<< endl;
196 
197  patchSet_ = filteredPatchSet;
198  }
199 
200  return true;
201 }
202 
203 
205 {
206  word name(type());
207 
208  if
209  (
210  foundObject<compressible::turbulenceModel>
211  (
213  )
214  )
215  {
216  const compressible::turbulenceModel& turbModel =
217  lookupObject<compressible::turbulenceModel>
218  (
220  );
221 
222  return store
223  (
224  name,
225  calcWallHeatFlux(turbModel.alphaEff(), turbModel.transport().he())
226  );
227  }
228  else if (foundObject<solidThermo>(solidThermo::dictName))
229  {
230  const solidThermo& thermo =
231  lookupObject<solidThermo>(solidThermo::dictName);
232 
233  return store(name, calcWallHeatFlux(thermo.alpha(), thermo.he()));
234  }
235  else
236  {
238  << "Unable to find compressible turbulence model in the "
239  << "database" << exit(FatalError);
240  }
241 
242  return true;
243 }
244 
245 
247 {
248  Log << type() << " " << name() << " write:" << nl;
249 
251 
252  logFiles::write();
253 
256 
257  const fvPatchList& patches = mesh_.boundary();
258 
259  const surfaceScalarField::Boundary& magSf =
260  mesh_.magSf().boundaryField();
261 
262  forAllConstIter(labelHashSet, patchSet_, iter)
263  {
264  label patchi = iter.key();
265  const fvPatch& pp = patches[patchi];
266 
267  const scalarField& hfp = wallHeatFlux.boundaryField()[patchi];
268 
269  const scalar minHfp = gMin(hfp);
270  const scalar maxHfp = gMax(hfp);
271  const scalar integralHfp = gSum(magSf[patchi]*hfp);
272 
273  if (Pstream::master())
274  {
275  file()
276  << mesh_.time().value()
277  << tab << pp.name()
278  << tab << minHfp
279  << tab << maxHfp
280  << tab << integralHfp
281  << endl;
282  }
283 
284  Log << " min/max/integ(" << pp.name() << ") = "
285  << minHfp << ", " << maxHfp << ", " << integralHfp << endl;
286  }
287 
288  Log << endl;
289 
290  return true;
291 }
292 
293 
294 // ************************************************************************* //
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
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:204
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
Type gMin(const FieldField< Field, Type > &f)
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
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:246
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
virtual ~wallHeatFlux()
Destructor.
Definition: wallHeatFlux.C:142
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:421
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:511
Calculate the snGrad of the given volField.
patches[0]
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
rhoReactionThermo & thermo
Definition: createFields.H:28
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
const word & name() const
Return name.
Definition: fvPatch.H:149
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
virtual tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
tmp< volScalarField > calcWallHeatFlux(const volScalarField &alpha, const volScalarField &he)
Calculate the heat-flux.
Definition: wallHeatFlux.C:63
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
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)
virtual void writeFileHeader(const label i)
File header information.
Definition: wallHeatFlux.C:48
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:48
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.
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
virtual bool read(const dictionary &)
Read the wallHeatFlux data.
Definition: wallHeatFlux.C:148
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Calculates and write the heat-flux at wall patches as the volScalarField field &#39;wallHeatFlux&#39;.
Definition: wallHeatFlux.H:109
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
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.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57