Go to the documentation of this file.
32 #include "surfaceInterpolate.H"
42 namespace functionObjects
71 volScalarField::Boundary& wallHeatFluxBf =
wallHeatFlux.boundaryFieldRef();
73 const volScalarField::Boundary& heBf =
he.boundaryField();
75 const volScalarField::Boundary& alphaBf =
alpha.boundaryField();
77 for (
const label patchi : patchSet_)
79 wallHeatFluxBf[patchi] = alphaBf[patchi]*heBf[patchi].snGrad();
83 const auto* qrPtr = cfindObject<volScalarField>(qrName_);
87 const volScalarField::Boundary& radHeatFluxBf = qrPtr->boundaryField();
89 for (
const label patchi : patchSet_)
91 wallHeatFluxBf[patchi] -= radHeatFluxBf[patchi];
118 mesh_.time().timeName(),
128 mesh_.objectRegistry::store(wallHeatFluxPtr);
132 writeFileHeader(file());
151 dict.readIfPresent(
"qr", qrName_);
155 if (patchSet_.empty())
159 if (isA<wallPolyPatch>(pbm[patchi]))
161 patchSet_.insert(patchi);
165 Info<<
" processing all wall patches" <<
nl <<
endl;
169 Info<<
" processing wall patches: " <<
nl;
171 for (
const label patchi : patchSet_)
173 if (isA<wallPolyPatch>(pbm[patchi]))
175 filteredPatchSet.
insert(patchi);
181 <<
"Requested wall heat-flux on non-wall boundary "
182 <<
"type patch: " << pbm[patchi].
name() <<
endl;
188 patchSet_ = filteredPatchSet;
201 foundObject<compressible::turbulenceModel>
208 lookupObject<compressible::turbulenceModel>
216 turbModel.transport().he(),
242 <<
"Unable to find compressible turbulence model in the "
248 const surfaceScalarField::Boundary& magSf = mesh_.magSf().boundaryField();
250 for (
const label patchi : patchSet_)
256 const scalar minHfp =
gMin(hfp);
257 const scalar maxHfp =
gMax(hfp);
258 const scalar integralHfp =
gSum(magSf[patchi]*hfp);
262 writeCurrentTime(file());
272 Log <<
" min/max/integ(" << pp.
name() <<
") = "
273 << minHfp <<
", " << maxHfp <<
", " << integralHfp <<
endl;
275 this->setResult(
"min(" + pp.
name() +
")", minHfp);
276 this->setResult(
"max(" + pp.
name() +
")", maxHfp);
277 this->setResult(
"int(" + pp.
name() +
")", integralHfp);
Computes the wall-heat flux at selected wall patches.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
const word & name() const
Return name.
A class for handling words, derived from Foam::string.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
static constexpr const zero Zero
Global zero (0)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void calcHeatFlux(const volScalarField &alpha, const volScalarField &he, volScalarField &wallHeatFlux)
Calculate the heat-flux.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Calculate the snGrad of the given volField.
static const word propertiesName
Default name of the turbulence properties dictionary.
Fundamental fluid thermodynamic properties.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
virtual void writeFileHeader(Ostream &os) const
File header information.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Type gSum(const FieldField< Field, Type > &f)
virtual const word & name() const
Return name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
#define forAll(list, i)
Loop across all elements in list.
virtual bool execute()
Calculate the wall heat-flux.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Fundamental solid thermodynamic properties.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual bool read(const dictionary &dict)
Read.
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (uses stdout - output is on the master only)
word name(const complex &c)
Return string representation of complex.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
word dictName() const
The local dictionary name (final part of scoped name)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool read(const dictionary &dict)
Read optional controls.
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &dict)
Read the wallHeatFlux data.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
static bool master(const label communicator=0)
Am I the master process.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
const word & name() const
Return the name of this functionObject.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
A List of wordRe with additional matching capabilities.
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.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
defineTypeNameAndDebug(ObukhovLength, 0)
const polyBoundaryMesh & patches
virtual bool write()
Write the wall heat-flux.
Base class for writing single files from the function objects.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type gMin(const FieldField< Field, Type > &f)
#define WarningInFunction
Report a warning using Foam::Warning.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Type gMax(const FieldField< Field, Type > &f)