Go to the documentation of this file.
29 #include "diameterModel.H"
33 #include "surfaceInterpolate.H"
40 const word& phaseName,
49 IOobject::groupName(
"alpha", phaseName),
50 mesh.time().timeName(),
58 phaseDict_(phaseDict),
87 IOobject::groupName(
"U", phaseName),
88 mesh.time().timeName(),
99 IOobject::groupName(
"DDtU", phaseName),
100 mesh.time().timeName(),
110 IOobject::groupName(
"alphaPhi", phaseName),
111 mesh.time().timeName(),
118 alphaPhi_.setOriented();
120 const word phiName = IOobject::groupName(
"phi", name_);
125 mesh.time().timeName(),
132 Info<<
"Reading face flux field " << phiName <<
endl;
141 mesh.time().timeName(),
152 Info<<
"Calculating face flux field " << phiName <<
endl;
156 U_.boundaryField().size(),
157 calculatedFvPatchScalarField::typeName
160 forAll(U_.boundaryField(), i)
164 isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
165 || isA<slipFvPatchVectorField>(U_.boundaryField()[i])
166 || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
169 phiTypes[i] = fixedValueFvPatchScalarField::typeName;
180 mesh.time().timeName(),
223 phaseDict_ = phaseDict;
228 phaseDict_.readEntry(
"kappa", kappa_.value());
229 phaseDict_.readEntry(
"Cp", Cp_.value());
230 phaseDict_.readEntry(
"rho", rho_.value());
242 const volScalarField::Boundary& alphaBf = boundaryField();
245 forAll(alphaPhiBf, patchi)
251 alphaPhip = phiBf[patchi]*alphaBf[patchi];
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimVelocity
const dimensionSet dimDensity
tmp< volScalarField > d() const
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
void correct()
Correct the phase properties.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
#define forAll(list, i)
Loop across all elements in list.
virtual ~phaseModel()
Destructor.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionSet dimSpecificHeatCapacity(dimGasConstant)
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimViscosity
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
virtual bool read()
Read phase properties dictionary.
Calculate the face-flux of the given field.
autoPtr< phaseModel > clone() const
Return clone.
virtual bool coupled() const
Return true if this patch field is coupled.
const Boundary & boundaryField() const
Return const-reference to the boundary field.