Go to the documentation of this file.
29 #include "phaseSystem.H"
37 namespace functionObjects
60 const volScalarField::Boundary& Tbf =
T.boundaryField();
76 <<
"Unable to find a valid phaseSystem to evaluate q" <<
nl
82 for (
const label patchi : htcModelPtr_->patchSet())
88 const volScalarField::Boundary& hebf =
he.boundaryField();
102 const volScalarField::Boundary& qrbf = qrPtr->
boundaryField();
104 for (
const label patchi : htcModelPtr_->patchSet())
106 q[patchi] += qrbf[patchi];
118 htcModelPtr_->mesh().lookupObjectRef<
volScalarField>(resultName_);
120 htcModelPtr_->calc(htc, q());
136 htcModelPtr_(
nullptr)
140 setResultName(typeName,
name +
":htc:" + htcModelPtr_->type());
148 mesh_.time().timeName(),
157 mesh_.objectRegistry::store(htcPtr);
169 htcModelPtr_->read(
dict);
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< FieldField< Field, scalar > > q() const
Calculate heat flux.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A heat transfer coefficient for reactingEuler solvers.
A class for handling words, derived from Foam::string.
static autoPtr< functionObject > New(const word &name, const Time &runTime, const dictionary &dict)
Select from dictionary, based on its "type" entry.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *buf, int32_t &val)
Same as readInt32.
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field T\n"<< endl;volScalarField T(IOobject("T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Calculating field g.h\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), p_rgh);Info<< "Creating multiphaseSystem\n"<< endl;autoPtr< multiphaseSystem > fluidPtr
#define forAll(list, i)
Loop across all elements in list.
const dimensionSet dimArea(sqr(dimLength))
word name(const complex &c)
Return string representation of complex.
static autoPtr< heatTransferCoeffModel > New(const dictionary &dict, const fvMesh &mesh, const word &TName)
Return a reference to the selected heat transfer coefficient model.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimPower
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual bool calc()
Calculate the heat transfer coefficient field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
defineTypeNameAndDebug(ObukhovLength, 0)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
reactingEulerHtcModel(const reactingEulerHtcModel &)=delete
No copy construct.
Class to represent a system of phases and model interfacial transfers between them.
virtual bool read(const dictionary &dict)
Read the heatTransferCoeff data.
const Boundary & boundaryField() const
Return const-reference to the boundary field.