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-2020 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  (
71  (
72  type(),
73  mesh_,
75  )
76  );
77 
78  volScalarField::Boundary& wallHeatFluxBf =
79  twallHeatFlux.ref().boundaryFieldRef();
80 
81  const volScalarField::Boundary& heBf = he.boundaryField();
82  const volScalarField::Boundary& alphaBf = alpha.boundaryField();
83 
84  forAllConstIter(labelHashSet, patchSet_, iter)
85  {
86  const label patchi = iter.key();
87 
88  wallHeatFluxBf[patchi] = alphaBf[patchi]*heBf[patchi].snGrad();
89  }
90 
91  if (foundObject<volScalarField>("qr"))
92  {
93  const volScalarField& qr = lookupObject<volScalarField>("qr");
94 
95  const volScalarField::Boundary& radHeatFluxBf = qr.boundaryField();
96 
97  forAllConstIter(labelHashSet, patchSet_, iter)
98  {
99  const label patchi = iter.key();
100 
101  wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
102  }
103  }
104 
105  return twallHeatFlux;
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
110 
112 (
113  const word& name,
114  const Time& runTime,
115  const dictionary& dict
116 )
117 :
118  fvMeshFunctionObject(name, runTime, dict),
119  logFiles(obr_, name),
120  writeLocalObjects(obr_, log),
121  patchSet_()
122 {
123  read(dict);
124  resetName(typeName);
125  resetLocalObjectName(typeName);
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
130 
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
141 
142  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
143 
144  patchSet_ =
145  mesh_.boundaryMesh().patchSet
146  (
147  wordReList(dict.lookupOrDefault("patches", wordReList()))
148  );
149 
150  Info<< type() << " " << name() << ":" << nl;
151 
152  if (patchSet_.empty())
153  {
154  forAll(pbm, patchi)
155  {
156  if (isA<wallPolyPatch>(pbm[patchi]))
157  {
158  patchSet_.insert(patchi);
159  }
160  }
161 
162  Info<< " processing all wall patches" << nl << endl;
163  }
164  else
165  {
166  Info<< " processing wall patches: " << nl;
167  labelHashSet filteredPatchSet;
168  forAllConstIter(labelHashSet, patchSet_, iter)
169  {
170  label patchi = iter.key();
171  if (isA<wallPolyPatch>(pbm[patchi]))
172  {
173  filteredPatchSet.insert(patchi);
174  Info<< " " << pbm[patchi].name() << endl;
175  }
176  else
177  {
179  << "Requested wall heat-flux on non-wall boundary "
180  << "type patch: " << pbm[patchi].name() << endl;
181  }
182  }
183 
184  Info<< endl;
185 
186  patchSet_ = filteredPatchSet;
187  }
188 
189  return true;
190 }
191 
192 
194 {
195  word name(type());
196 
197  if
198  (
199  foundObject<thermophysicalTransportModel>
200  (
201  thermophysicalTransportModel::typeName
202  )
203  )
204  {
205  const thermophysicalTransportModel& ttm =
206  lookupObject<thermophysicalTransportModel>
207  (
208  thermophysicalTransportModel::typeName
209  );
210 
211  return store
212  (
213  name,
214  calcWallHeatFlux(ttm.alphaEff(), ttm.thermo().he())
215  );
216  }
217  else if (foundObject<solidThermo>(solidThermo::dictName))
218  {
219  const solidThermo& thermo =
220  lookupObject<solidThermo>(solidThermo::dictName);
221 
222  return store(name, calcWallHeatFlux(thermo.alpha(), thermo.he()));
223  }
224  else
225  {
227  << "Unable to find compressible turbulence model in the "
228  << "database" << exit(FatalError);
229  }
230 
231  return true;
232 }
233 
234 
236 {
237  Log << type() << " " << name() << " write:" << nl;
238 
240 
241  logFiles::write();
242 
245 
246  const fvPatchList& patches = mesh_.boundary();
247 
248  const surfaceScalarField::Boundary& magSf =
249  mesh_.magSf().boundaryField();
250 
251  forAllConstIter(labelHashSet, patchSet_, iter)
252  {
253  label patchi = iter.key();
254  const fvPatch& pp = patches[patchi];
255 
256  const scalarField& hfp = wallHeatFlux.boundaryField()[patchi];
257 
258  const scalar minHfp = gMin(hfp);
259  const scalar maxHfp = gMax(hfp);
260  const scalar integralHfp = gSum(magSf[patchi]*hfp);
261 
262  if (Pstream::master())
263  {
264  file()
265  << mesh_.time().value()
266  << tab << pp.name()
267  << tab << minHfp
268  << tab << maxHfp
269  << tab << integralHfp
270  << endl;
271  }
272 
273  Log << " min/max/integ(" << pp.name() << ") = "
274  << minHfp << ", " << maxHfp << ", " << integralHfp << endl;
275  }
276 
277  Log << endl;
278 
279  return true;
280 }
281 
282 
283 // ************************************************************************* //
virtual bool write()
Write function.
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:103
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:303
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: wallHeatFlux.C:112
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:193
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:259
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
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:235
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
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:61
virtual ~wallHeatFlux()
Destructor.
Definition: wallHeatFlux.C:131
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:494
Calculate the snGrad of the given volField.
patches[0]
rhoReactionThermo & thermo
Definition: createFields.H:28
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.
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
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:123
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)
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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:260
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
defineTypeNameAndDebug(Qdot, 0)
Abstract base class for thermophysical transport models (RAS, LES and laminar).
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.
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
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
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:137
virtual const word & name() const
Return name.
Definition: fvPatch.H:143
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
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
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57