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-------------------------------------------------------------------------------
11License
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"
38
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
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:
108 writeFile(obr_, name, typeName, dict),
109 patchSet_(),
110 qrName_("qr")
111{
112 volScalarField* wallHeatFluxPtr
113 (
115 (
117 (
118 scopedName(typeName),
119 mesh_.time().timeName(),
120 mesh_,
123 ),
124 mesh_,
126 )
127 );
128
129 mesh_.objectRegistry::store(wallHeatFluxPtr);
130
131 read(dict);
132
134}
135
136
137// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138
140{
143
144 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
145
146 patchSet_ =
147 mesh_.boundaryMesh().patchSet
148 (
149 dict.getOrDefault<wordRes>("patches", wordRes())
150 );
151
152 dict.readIfPresent("qr", qrName_);
153
154 Info<< type() << " " << name() << ":" << nl;
155
156 if (patchSet_.empty())
157 {
158 forAll(pbm, patchi)
159 {
160 if (isA<wallPolyPatch>(pbm[patchi]))
161 {
162 patchSet_.insert(patchi);
163 }
164 }
165
166 Info<< " processing all wall patches" << nl << endl;
167 }
168 else
169 {
170 Info<< " processing wall patches: " << nl;
171 labelHashSet filteredPatchSet;
172 for (const label patchi : patchSet_)
173 {
174 if (isA<wallPolyPatch>(pbm[patchi]))
175 {
176 filteredPatchSet.insert(patchi);
177 Info<< " " << pbm[patchi].name() << endl;
178 }
179 else
180 {
182 << "Requested wall heat-flux on non-wall boundary "
183 << "type patch: " << pbm[patchi].name() << endl;
184 }
185 }
186
187 Info<< endl;
188
189 patchSet_ = filteredPatchSet;
190 }
191
192 return true;
193}
194
195
197{
198 auto& wallHeatFlux = lookupObjectRef<volScalarField>(scopedName(typeName));
199
200 if
201 (
202 foundObject<compressible::turbulenceModel>
203 (
205 )
206 )
207 {
208 const compressible::turbulenceModel& turbModel =
209 lookupObject<compressible::turbulenceModel>
210 (
212 );
213
214 calcHeatFlux
215 (
216 turbModel.alphaEff()(),
217 turbModel.transport().he(),
219 );
220 }
221 else if (foundObject<fluidThermo>(fluidThermo::dictName))
222 {
223 const fluidThermo& thermo =
224 lookupObject<fluidThermo>(fluidThermo::dictName);
225
226 calcHeatFlux
227 (
228 thermo.alpha(),
229 thermo.he(),
231 );
232 }
233 else if (foundObject<solidThermo>(solidThermo::dictName))
234 {
235 const solidThermo& thermo =
236 lookupObject<solidThermo>(solidThermo::dictName);
237
238 calcHeatFlux(thermo.alpha(), thermo.he(), wallHeatFlux);
239 }
240 else if
241 (
242 foundObject<multiphaseInterSystem>
244 )
245 {
246 const auto& thermo = lookupObject<multiphaseInterSystem>
247 (
249 );
250
251 calcHeatFlux(thermo.kappaEff()(), thermo.T(), wallHeatFlux);
252 }
253 else
254 {
256 << "Unable to find compressible turbulence model in the "
257 << "database" << exit(FatalError);
258 }
259
260 const fvPatchList& patches = mesh_.boundary();
261
262 const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
263
264 for (const label patchi : patchSet_)
265 {
266 const fvPatch& pp = patches[patchi];
267
268 const scalarField& hfp = wallHeatFlux.boundaryField()[patchi];
269
270 const scalar minHfp = gMin(hfp);
271 const scalar maxHfp = gMax(hfp);
272 const scalar integralHfp = gSum(magSf[patchi]*hfp);
273
274 if (Pstream::master())
275 {
276 writeCurrentTime(file());
277
278 file()
279 << token::TAB << pp.name()
280 << token::TAB << minHfp
281 << token::TAB << maxHfp
282 << token::TAB << integralHfp
283 << endl;
284 }
285
286 Log << " min/max/integ(" << pp.name() << ") = "
287 << minHfp << ", " << maxHfp << ", " << integralHfp << endl;
288
289 this->setResult("min(" + pp.name() + ")", minHfp);
290 this->setResult("max(" + pp.name() + ")", maxHfp);
291 this->setResult("int(" + pp.name() + ")", integralHfp);
292 }
293
294
295 return true;
296}
297
298
300{
301 const auto& wallHeatFlux =
302 lookupObject<volScalarField>(scopedName(typeName));
303
304 Log << type() << " " << name() << " write:" << nl
305 << " writing field " << wallHeatFlux.name() << endl;
306
308
309 return true;
310}
311
312
313// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
volScalarField & he
Definition: YEEqn.H:52
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes the wall-heat flux at selected wall patches.
Definition: wallHeatFlux.H:164
virtual bool read(const dictionary &dict)
Read the wallHeatFlux data.
Definition: wallHeatFlux.C:139
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: wallHeatFlux.C:53
virtual bool execute()
Calculate the wall heat-flux.
Definition: wallHeatFlux.C:196
virtual bool write()
Write the wall heat-flux.
Definition: wallHeatFlux.C:299
void calcHeatFlux(const volScalarField &alpha, const volScalarField &he, volScalarField &wallHeatFlux)
Calculate the heat-flux.
Definition: wallHeatFlux.C:66
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:285
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:233
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
static const word phasePropertiesName
Default name of the phase properties dictionary.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Fundamental solid thermodynamic properties.
Definition: solidThermo.H:55
splitCell * master() const
Definition: splitCell.H:113
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static const word propertiesName
Default name of the turbulence properties dictionary.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const polyBoundaryMesh & patches
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Calculate the snGrad of the given volField.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333