wallHeatFlux.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "wallHeatFlux.H"
31 #include "solidThermo.H"
32 #include "surfaceInterpolate.H"
33 #include "fvcSnGrad.H"
34 #include "wallPolyPatch.H"
37 #include "multiphaseInterSystem.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace functionObjects
44 {
45  defineTypeNameAndDebug(wallHeatFlux, 0);
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
52 
54 {
55  writeHeader(os, "Wall heat-flux");
56  writeCommented(os, "Time");
57  writeTabbed(os, "patch");
58  writeTabbed(os, "min");
59  writeTabbed(os, "max");
60  writeTabbed(os, "integral");
61  os << endl;
62 }
63 
64 
66 (
67  const volScalarField& alpha,
68  const volScalarField& he,
70 )
71 {
72  volScalarField::Boundary& wallHeatFluxBf = wallHeatFlux.boundaryFieldRef();
73 
74  const volScalarField::Boundary& heBf = he.boundaryField();
75 
76  const volScalarField::Boundary& alphaBf = alpha.boundaryField();
77 
78  for (const label patchi : patchSet_)
79  {
80  wallHeatFluxBf[patchi] = alphaBf[patchi]*heBf[patchi].snGrad();
81  }
82 
83 
84  const auto* qrPtr = cfindObject<volScalarField>(qrName_);
85 
86  if (qrPtr)
87  {
88  const volScalarField::Boundary& radHeatFluxBf = qrPtr->boundaryField();
89 
90  for (const label patchi : patchSet_)
91  {
92  wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
93  }
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
101 (
102  const word& name,
103  const Time& runTime,
104  const dictionary& dict
105 )
106 :
107  fvMeshFunctionObject(name, runTime, dict),
108  writeFile(obr_, name, typeName, dict),
109  patchSet_(),
110  qrName_("qr")
111 {
112  volScalarField* wallHeatFluxPtr
113  (
114  new volScalarField
115  (
116  IOobject
117  (
118  scopedName(typeName),
119  mesh_.time().timeName(),
120  mesh_,
124  ),
125  mesh_,
127  )
128  );
129 
130  mesh_.objectRegistry::store(wallHeatFluxPtr);
131 
132  read(dict);
135 }
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
144 
145  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
146 
147  patchSet_ =
148  mesh_.boundaryMesh().patchSet
149  (
150  dict.getOrDefault<wordRes>("patches", wordRes())
151  );
152 
153  dict.readIfPresent("qr", qrName_);
154 
155  Info<< type() << " " << name() << ":" << nl;
156 
157  if (patchSet_.empty())
158  {
159  forAll(pbm, patchi)
160  {
161  if (isA<wallPolyPatch>(pbm[patchi]))
162  {
163  patchSet_.insert(patchi);
164  }
165  }
166 
167  Info<< " processing all wall patches" << nl << endl;
168  }
169  else
170  {
171  Info<< " processing wall patches: " << nl;
172  labelHashSet filteredPatchSet;
173  for (const label patchi : patchSet_)
174  {
175  if (isA<wallPolyPatch>(pbm[patchi]))
176  {
177  filteredPatchSet.insert(patchi);
178  Info<< " " << pbm[patchi].name() << endl;
179  }
180  else
181  {
183  << "Requested wall heat-flux on non-wall boundary "
184  << "type patch: " << pbm[patchi].name() << endl;
185  }
186  }
187 
188  Info<< endl;
189 
190  patchSet_ = filteredPatchSet;
191  }
192 
193  return true;
194 }
195 
196 
198 {
199  auto& wallHeatFlux = lookupObjectRef<volScalarField>(scopedName(typeName));
200 
201  if
202  (
203  foundObject<compressible::turbulenceModel>
204  (
206  )
207  )
208  {
209  const compressible::turbulenceModel& turbModel =
210  lookupObject<compressible::turbulenceModel>
211  (
213  );
214 
215  calcHeatFlux
216  (
217  turbModel.alphaEff()(),
218  turbModel.transport().he(),
219  wallHeatFlux
220  );
221  }
222  else if (foundObject<fluidThermo>(fluidThermo::dictName))
223  {
224  const fluidThermo& thermo =
225  lookupObject<fluidThermo>(fluidThermo::dictName);
226 
227  calcHeatFlux
228  (
229  thermo.alpha(),
230  thermo.he(),
231  wallHeatFlux
232  );
233  }
234  else if (foundObject<solidThermo>(solidThermo::dictName))
235  {
236  const solidThermo& thermo =
237  lookupObject<solidThermo>(solidThermo::dictName);
238 
239  calcHeatFlux(thermo.alpha(), thermo.he(), wallHeatFlux);
240  }
241  else if
242  (
243  foundObject<multiphaseInterSystem>
245  )
246  {
247  const auto& thermo = lookupObject<multiphaseInterSystem>
248  (
250  );
251 
252  calcHeatFlux(thermo.kappaEff()(), thermo.T(), wallHeatFlux);
253  }
254  else
255  {
257  << "Unable to find compressible turbulence model in the "
258  << "database" << exit(FatalError);
259  }
260 
261  const fvPatchList& patches = mesh_.boundary();
262 
263  const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
264 
265  for (const label patchi : patchSet_)
266  {
267  const fvPatch& pp = patches[patchi];
268 
269  const scalarField& hfp = wallHeatFlux.boundaryField()[patchi];
270 
271  const scalar minHfp = gMin(hfp);
272  const scalar maxHfp = gMax(hfp);
273  const scalar integralHfp = gSum(magSf[patchi]*hfp);
274 
275  if (Pstream::master())
276  {
277  writeCurrentTime(file());
278 
279  file()
280  << token::TAB << pp.name()
281  << token::TAB << minHfp
282  << token::TAB << maxHfp
283  << token::TAB << integralHfp
284  << endl;
285  }
286 
287  Log << " min/max/integ(" << pp.name() << ") = "
288  << minHfp << ", " << maxHfp << ", " << integralHfp << endl;
289 
290  this->setResult("min(" + pp.name() + ")", minHfp);
291  this->setResult("max(" + pp.name() + ")", maxHfp);
292  this->setResult("int(" + pp.name() + ")", integralHfp);
293  }
294 
295 
296  return true;
297 }
298 
299 
301 {
302  const auto& wallHeatFlux =
303  lookupObject<volScalarField>(scopedName(typeName));
304 
305  Log << type() << " " << name() << " write:" << nl
306  << " writing field " << wallHeatFlux.name() << endl;
307 
308  wallHeatFlux.write();
309 
310  return true;
311 }
312 
313 
314 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const polyBoundaryMesh & pbm
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: wallHeatFlux.C:94
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:190
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static const word phasePropertiesName
Default name of the phase properties dictionary.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
Type gMin(const FieldField< Field, Type > &f)
Tab [isspace].
Definition: token.H:126
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:162
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: wallHeatFlux.C:46
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:293
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Calculate the snGrad of the given volField.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
Abstract base-class for Time/database function objects.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:227
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
psiReactionThermo & thermo
Definition: createFields.H:28
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
void calcHeatFlux(const volScalarField &alpha, const volScalarField &he, volScalarField &wallHeatFlux)
Calculate the heat-flux.
Definition: wallHeatFlux.C:59
dimensionedScalar pow3(const dimensionedScalar &ds)
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1037
const polyBoundaryMesh & patches
Nothing to be read.
#define Log
Definition: PDRblock.C:28
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Computes the wall-heat flux at selected wall patches.
Definition: wallHeatFlux.H:156
virtual bool read(const dictionary &dict)
Read the wallHeatFlux data.
Definition: wallHeatFlux.C:133
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
Request registration (bool: true)
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329