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-2023 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 : patchIDs_)
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 : patchIDs_)
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  qrName_("qr")
110 {
111  volScalarField* wallHeatFluxPtr
112  (
113  new volScalarField
114  (
115  IOobject
116  (
117  scopedName(typeName),
118  mesh_.time().timeName(),
119  mesh_,
123  ),
124  mesh_,
126  )
127  );
128 
129  mesh_.objectRegistry::store(wallHeatFluxPtr);
130 
131  read(dict);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
140 {
141  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
142 
145 
146  dict.readIfPresent("qr", qrName_);
147 
149  labelHashSet patchSet(0);
150  if (dict.readIfPresent("patches", patchNames) && !patchNames.empty())
151  {
152  patchSet = pbm.patchSet(patchNames);
153  }
154 
155  labelHashSet allWalls(pbm.findPatchIDs<wallPolyPatch>());
156 
157  Info<< type() << ' ' << name() << ':' << nl;
158 
159  if (patchSet.empty())
160  {
161  patchIDs_ = allWalls.sortedToc();
162 
163  Info<< " processing all (" << patchIDs_.size()
164  << ") wall patches" << nl << endl;
165  }
166  else
167  {
168  allWalls &= patchSet;
169  patchSet -= allWalls;
170  patchIDs_ = allWalls.sortedToc();
171 
172  if (!patchSet.empty())
173  {
175  << "Requested wall heat-flux on ("
176  << patchSet.size() << ") non-wall patches:" << nl;
177 
178  for (const label patchi : patchSet.sortedToc())
179  {
180  Info<< " " << pbm[patchi].name() << nl;
181  }
182  Info<< nl;
183  }
184 
185  Info<< " processing (" << patchIDs_.size()
186  << ") wall patches:" << nl;
187 
188  for (const label patchi : patchIDs_)
189  {
190  Info<< " " << pbm[patchi].name() << nl;
191  }
193  }
194 
195  return true;
196 }
197 
198 
200 {
201  auto& wallHeatFlux = lookupObjectRef<volScalarField>(scopedName(typeName));
202 
203  if
204  (
205  foundObject<compressible::turbulenceModel>
206  (
208  )
209  )
210  {
211  const compressible::turbulenceModel& turbModel =
212  lookupObject<compressible::turbulenceModel>
213  (
215  );
216 
217  calcHeatFlux
218  (
219  turbModel.alphaEff()(),
220  turbModel.transport().he(),
221  wallHeatFlux
222  );
223  }
224  else if (foundObject<fluidThermo>(fluidThermo::dictName))
225  {
226  const fluidThermo& thermo =
227  lookupObject<fluidThermo>(fluidThermo::dictName);
228 
229  calcHeatFlux
230  (
231  thermo.alpha(),
232  thermo.he(),
233  wallHeatFlux
234  );
235  }
236  else if (foundObject<solidThermo>(solidThermo::dictName))
237  {
238  const solidThermo& thermo =
239  lookupObject<solidThermo>(solidThermo::dictName);
240 
241  calcHeatFlux(thermo.alpha(), thermo.he(), wallHeatFlux);
242  }
243  else if
244  (
245  foundObject<multiphaseInterSystem>
247  )
248  {
249  const auto& thermo = lookupObject<multiphaseInterSystem>
250  (
252  );
253 
254  calcHeatFlux(thermo.kappaEff()(), thermo.T(), wallHeatFlux);
255  }
256  else
257  {
259  << "Unable to find compressible turbulence model in the "
260  << "database" << exit(FatalError);
261  }
262 
263 
264  const fvPatchList& patches = mesh_.boundary();
265 
266  const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
267 
268  for (const label patchi : patchIDs_)
269  {
270  const fvPatch& pp = patches[patchi];
271 
272  const scalarField& hfp = wallHeatFlux.boundaryField()[patchi];
273 
274  const MinMax<scalar> limits = gMinMax(hfp);
275  const scalar integralHfp = gSum(magSf[patchi]*hfp);
276 
277  if (Pstream::master())
278  {
279  writeCurrentTime(file());
280 
281  file()
282  << token::TAB << pp.name()
283  << token::TAB << limits.min()
284  << token::TAB << limits.max()
285  << token::TAB << integralHfp
286  << endl;
287  }
288 
289  Log << " min/max/integ(" << pp.name() << ") = "
290  << limits.min() << ", " << limits.max()
291  << ", " << integralHfp << endl;
292 
293  this->setResult("min(" + pp.name() + ")", limits.min());
294  this->setResult("max(" + pp.name() + ")", limits.max());
295  this->setResult("int(" + pp.name() + ")", integralHfp);
296  }
297 
298  return true;
299 }
300 
301 
303 {
304  const auto& wallHeatFlux =
305  lookupObject<volScalarField>(scopedName(typeName));
306 
307  Log << type() << ' ' << name() << " write:" << nl
308  << " writing field " << wallHeatFlux.name() << endl;
309 
310  wallHeatFlux.write();
311 
312  return true;
313 }
314 
315 
316 // ************************************************************************* //
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:192
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Tab [isspace].
Definition: token.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: wallHeatFlux.C:46
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:295
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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:360
Abstract base-class for Time/database function objects.
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.
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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.
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:46
wordList patchNames(nPatches)
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
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
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:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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:1082
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:132
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:172
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:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329