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-2023 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 "surfaceInterpolate.H"
29 #include "fvcGrad.H"
30 #include "wallPolyPatch.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
46 
47 void Foam::functionObjects::wallHeatFlux::writeFileHeader(const label i)
48 {
49  // Add headers to output data
50  writeHeader(file(), "Wall heat-flux");
51  writeCommented(file(), "Time");
52  writeTabbed(file(), "patch");
53  writeTabbed(file(), "min [W/m^2]");
54  writeTabbed(file(), "max [W/m^2]");
55  writeTabbed(file(), "Q [W]");
56  writeTabbed(file(), "q [W/m^2]");
57  file() << endl;
58 }
59 
60 
62 Foam::functionObjects::wallHeatFlux::calcWallHeatFlux
63 (
64  const surfaceScalarField& q
65 )
66 {
67  tmp<volScalarField> twallHeatFlux
68  (
70  (
71  type(),
72  mesh_,
74  )
75  );
76 
77  volScalarField::Boundary& wallHeatFluxBf =
78  twallHeatFlux.ref().boundaryFieldRef();
79 
80  const surfaceScalarField::Boundary& qBf = q.boundaryField();
81 
82  forAllConstIter(labelHashSet, patchSet_, iter)
83  {
84  const label patchi = iter.key();
85 
86  wallHeatFluxBf[patchi] = -qBf[patchi];
87  }
88 
89  if (foundObject<volScalarField>("qr"))
90  {
91  const volScalarField& qr = lookupObject<volScalarField>("qr");
92 
93  const volScalarField::Boundary& radHeatFluxBf = qr.boundaryField();
94 
95  forAllConstIter(labelHashSet, patchSet_, iter)
96  {
97  const label patchi = iter.key();
98 
99  wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
100  }
101  }
102 
103  return twallHeatFlux;
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
110 (
111  const word& name,
112  const Time& runTime,
113  const dictionary& dict
114 )
115 :
116  fvMeshFunctionObject(name, runTime, dict),
117  logFiles(obr_, name),
118  writeLocalObjects(obr_, log),
119  phaseName_(word::null),
120  patchSet_()
121 {
122  read(dict);
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
128 
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 {
139 
140  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
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  resetName(typeName);
190  resetLocalObjectName(typeName);
191 
192  return true;
193 }
194 
195 
197 {
198  const word fieldName(IOobject::groupName(type(), phaseName_));
199 
200  const word thermophysicalTransportModelName
201  (
202  IOobject::groupName(thermophysicalTransportModel::typeName, phaseName_)
203  );
204 
205  if
206  (
207  foundObject<thermophysicalTransportModel>
208  (
209  thermophysicalTransportModelName
210  )
211  )
212  {
213  const thermophysicalTransportModel& ttm =
214  lookupObject<thermophysicalTransportModel>
215  (
216  thermophysicalTransportModelName
217  );
218 
219  return store(fieldName, calcWallHeatFlux(ttm.q()));
220  }
221  else
222  {
224  << "Unable to find thermophysicalTransportModel in the "
225  << "database" << exit(FatalError);
226  }
227 
228  return true;
229 }
230 
231 
233 {
234  Log << type() << " " << name() << " write:" << nl;
235 
237 
238  logFiles::write();
239 
241  (
242  IOobject::groupName(type(), phaseName_)
243  );
244 
245  const fvPatchList& patches = mesh_.boundary();
246 
247  const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
248 
249  forAllConstIter(labelHashSet, patchSet_, iter)
250  {
251  label patchi = iter.key();
252  const fvPatch& pp = patches[patchi];
253 
254  const scalarField& qp = wallHeatFlux.boundaryField()[patchi];
255 
256  const scalar minqp = gMin(qp);
257  const scalar maxqp = gMax(qp);
258  const scalar Q = gSum(magSf[patchi]*qp);
259  const scalar q = Q/gSum(magSf[patchi]);
260 
261  if (Pstream::master())
262  {
263  file()
264  << time_.userTimeValue()
265  << tab << pp.name()
266  << tab << minqp
267  << tab << maxqp
268  << tab << Q
269  << tab << q
270  << endl;
271  }
272 
273  Log << " min, max, Q [W], q [W/m^2] for patch " << pp.name() << " = "
274  << minqp << ", " << maxqp << ", " << Q << ", " << q << endl;
275  }
276 
277  Log << endl;
278 
279  return true;
280 }
281 
282 
283 // ************************************************************************* //
#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.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
static word groupName(Name name, const word &group)
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
Calculates and outputs the second invariant of the velocity gradient tensor [1/s^2].
Definition: Q.H:62
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:107
virtual bool write()
Write function.
Definition: logFiles.C:167
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:106
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
virtual bool read(const dictionary &)
Read optional controls.
Calculates and write the heat-flux at wall patches as the volScalarField field 'wallHeatFlux'.
Definition: wallHeatFlux.H:116
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: wallHeatFlux.C:110
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:196
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:232
virtual ~wallHeatFlux()
Destructor.
Definition: wallHeatFlux.C:129
virtual bool read(const dictionary &)
Read the wallHeatFlux data.
Definition: wallHeatFlux.C:135
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,...
void resetLocalObjectName(const word &name)
Reset the list of local object names from a single word.
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:146
Foam::polyBoundaryMesh.
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.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Calculate the gradient of the given field.
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char tab
Definition: Ostream.H:259
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
const dimensionSet dimMass
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
error FatalError
Type gMin(const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:260
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