28 #include "surfaceInterpolate.H"
37 namespace functionObjects
47 void Foam::functionObjects::wallHeatFlux::writeFileHeader(
const label i)
62 Foam::functionObjects::wallHeatFlux::calcWallHeatFlux
67 tmp<volScalarField> twallHeatFlux
78 twallHeatFlux.ref().boundaryFieldRef();
89 if (foundObject<volScalarField>(
"qr"))
103 return twallHeatFlux;
119 phaseName_(
word::null),
144 patchSet_ = patchSet(
dict,
true);
148 if (patchSet_.empty())
152 if (isA<wallPolyPatch>(pbm[
patchi]))
158 Info<<
" processing all wall patches" <<
nl <<
endl;
162 Info<<
" processing wall patches: " <<
nl;
167 if (isA<wallPolyPatch>(pbm[
patchi]))
175 <<
"Requested wall heat-flux on non-wall boundary "
176 <<
"type patch: " << pbm[
patchi].name() <<
endl;
182 patchSet_ = filteredPatchSet;
186 resetLocalObjectName(typeName);
196 const word thermophysicalTransportModelName
203 foundObject<thermophysicalTransportModel>
205 thermophysicalTransportModelName
210 lookupObject<thermophysicalTransportModel>
212 thermophysicalTransportModelName
215 return store(fieldName, calcWallHeatFlux(ttm.
q()));
220 <<
"Unable to find thermophysicalTransportModel "
221 << thermophysicalTransportModelName
222 <<
" in the database"
254 const scalar minqp =
gMin(qp);
255 const scalar maxqp =
gMax(qp);
262 << time_.userTimeValue()
271 Log <<
" min, max, Q [W], q [W/m^2] for patch " << pp.
name() <<
" = "
272 << minqp <<
", " << maxqp <<
", " <<
Q <<
", " << q <<
endl;
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool insert(const Key &key)
Insert a new entry.
static word groupName(Name name, const word &group)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static bool master(const label communicator=0)
Am I the master process.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
Calculates and outputs the second invariant of the velocity gradient tensor [1/s^2].
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
OFstream & file()
Return access to the file (if only 1)
virtual bool write()
Write function.
Calculates the natural logarithm of the specified scalar field.
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
virtual bool read(const dictionary &)
Read optional controls.
Calculates and write the heat-flux at wall patches as the volScalarField field 'wallHeatFlux'.
wallHeatFlux(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
virtual bool execute()
Calculate the wall heat-flux.
virtual bool write()
Write the wall heat-flux.
virtual ~wallHeatFlux()
Destructor.
virtual bool read(const dictionary &)
Read the wallHeatFlux data.
void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
FunctionObject base class for managing a list of objects on behalf of the inheriting function object,...
void resetLocalObjectName(const word &name)
Reset the list of local object names from a single word.
virtual bool read(const dictionary &)
Read the list of objects to be written.
virtual bool write()
Write function.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual const word & name() const
Return name.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the gradient of the given field.
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gSum(const FieldField< Field, Type > &f)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimTime
VolField< scalar > volScalarField
const dimensionSet dimMass
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.