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-2022 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 "fvcGrad.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 
48 void Foam::functionObjects::wallHeatFlux::writeFileHeader(const label i)
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 [W/m^2]");
55  writeTabbed(file(), "max [W/m^2]");
56  writeTabbed(file(), "Q [W]");
57  writeTabbed(file(), "q [W/m^2]");
58  file() << endl;
59 }
60 
61 
63 Foam::functionObjects::wallHeatFlux::calcWallHeatFlux
64 (
65  const surfaceScalarField& q
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 surfaceScalarField::Boundary& qBf = q.boundaryField();
82 
83  forAllConstIter(labelHashSet, patchSet_, iter)
84  {
85  const label patchi = iter.key();
86 
87  wallHeatFluxBf[patchi] = -qBf[patchi];
88  }
89 
90  if (foundObject<volScalarField>("qr"))
91  {
92  const volScalarField& qr = lookupObject<volScalarField>("qr");
93 
94  const volScalarField::Boundary& radHeatFluxBf = qr.boundaryField();
95 
96  forAllConstIter(labelHashSet, patchSet_, iter)
97  {
98  const label patchi = iter.key();
99 
100  wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
101  }
102  }
103 
104  return twallHeatFlux;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const word& name,
113  const Time& runTime,
114  const dictionary& dict
115 )
116 :
117  fvMeshFunctionObject(name, runTime, dict),
118  logFiles(obr_, name),
119  writeLocalObjects(obr_, log),
120  phaseName_(word::null),
121  patchSet_()
122 {
123  read(dict);
124  resetLocalObjectName(IOobject::groupName(type(), phaseName_));
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 {
140 
141  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
142 
143  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
144 
145  patchSet_ =
146  mesh_.boundaryMesh().patchSet
147  (
148  wordReList(dict.lookupOrDefault("patches", wordReList()))
149  );
150 
151  Info<< type() << " " << name() << ":" << nl;
152 
153  if (patchSet_.empty())
154  {
155  forAll(pbm, patchi)
156  {
157  if (isA<wallPolyPatch>(pbm[patchi]))
158  {
159  patchSet_.insert(patchi);
160  }
161  }
162 
163  Info<< " processing all wall patches" << nl << endl;
164  }
165  else
166  {
167  Info<< " processing wall patches: " << nl;
168  labelHashSet filteredPatchSet;
169  forAllConstIter(labelHashSet, patchSet_, iter)
170  {
171  label patchi = iter.key();
172  if (isA<wallPolyPatch>(pbm[patchi]))
173  {
174  filteredPatchSet.insert(patchi);
175  Info<< " " << pbm[patchi].name() << endl;
176  }
177  else
178  {
180  << "Requested wall heat-flux on non-wall boundary "
181  << "type patch: " << pbm[patchi].name() << endl;
182  }
183  }
184 
185  Info<< endl;
186 
187  patchSet_ = filteredPatchSet;
188  }
189 
190  resetName(typeName);
191  resetLocalObjectName(typeName);
192 
193  return true;
194 }
195 
196 
198 {
199  const word fieldName(IOobject::groupName(type(), phaseName_));
200 
201  const word thermophysicalTransportModelName
202  (
203  IOobject::groupName(thermophysicalTransportModel::typeName, phaseName_)
204  );
205 
206  if
207  (
208  foundObject<thermophysicalTransportModel>
209  (
210  thermophysicalTransportModelName
211  )
212  )
213  {
214  const thermophysicalTransportModel& ttm =
215  lookupObject<thermophysicalTransportModel>
216  (
217  thermophysicalTransportModelName
218  );
219 
220  return store(fieldName, calcWallHeatFlux(ttm.q()));
221  }
222  else if (foundObject<solidThermo>(physicalProperties::typeName))
223  {
224  const solidThermo& thermo =
225  lookupObject<solidThermo>(physicalProperties::typeName);
226 
227  return store(fieldName, calcWallHeatFlux(thermo.q()));
228  }
229  else
230  {
232  << "Unable to find compressible turbulence model in the "
233  << "database" << exit(FatalError);
234  }
235 
236  return true;
237 }
238 
239 
241 {
242  Log << type() << " " << name() << " write:" << nl;
243 
245 
246  logFiles::write();
247 
249  (
250  IOobject::groupName(type(), phaseName_)
251  );
252 
253  const fvPatchList& patches = mesh_.boundary();
254 
255  const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
256 
257  forAllConstIter(labelHashSet, patchSet_, iter)
258  {
259  label patchi = iter.key();
260  const fvPatch& pp = patches[patchi];
261 
262  const scalarField& qp = wallHeatFlux.boundaryField()[patchi];
263 
264  const scalar minqp = gMin(qp);
265  const scalar maxqp = gMax(qp);
266  const scalar Q = gSum(magSf[patchi]*qp);
267  const scalar q = Q/gSum(magSf[patchi]);
268 
269  if (Pstream::master())
270  {
271  file()
272  << mesh_.time().userTimeValue()
273  << tab << pp.name()
274  << tab << minqp
275  << tab << maxqp
276  << tab << Q
277  << tab << q
278  << endl;
279  }
280 
281  Log << " min, max, Q [W], q [W/m^2] for patch " << pp.name() << " = "
282  << minqp << ", " << maxqp << ", " << Q << ", " << q << endl;
283  }
284 
285  Log << endl;
286 
287  return true;
288 }
289 
290 
291 // ************************************************************************* //
const fvPatchList & patches
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:167
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: wallHeatFlux.C:111
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:197
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:240
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:63
virtual ~wallHeatFlux()
Destructor.
Definition: wallHeatFlux.C:130
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
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:69
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
const dimensionSet dimTime
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.
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
static const word null
An empty word.
Definition: word.H:77
Calculates and outputs the second invariant of the velocity gradient tensor [1/s^2].
Definition: Q.H:59
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:53
const dimensionSet dimMass
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.
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
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W].
virtual bool read(const dictionary &)
Read the wallHeatFlux data.
Definition: wallHeatFlux.C:136
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
Specialisation 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:111
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
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57