30#include "diameterModel.H"
34#include "surfaceInterpolate.H"
41 const word& phaseName,
59 phaseDict_(phaseDict),
119 alphaPhi_.setOriented();
124 mesh.time().timeName(),
133 Info<<
"Reading face flux field " <<
io.name() <<
endl;
139 Info<<
"Calculating face flux field " <<
io.name() <<
endl;
145 U_.boundaryField().size(),
149 forAll(U_.boundaryField(), i)
158 phiTypes[i] = fixedValueFvPatchScalarField::typeName;
205 phaseDict_ = phaseDict;
209 phaseDict_.readEntry(
"nu", nu_.value());
210 phaseDict_.readEntry(
"kappa", kappa_.value());
211 phaseDict_.readEntry(
"Cp", Cp_.value());
212 phaseDict_.readEntry(
"rho", rho_.value());
225 const tmp<surfaceScalarField> tphi(
phi());
226 const auto& phiBf = tphi().boundaryField();
228 forAll(alphaPhiBf, patchi)
232 if (!alphaPhip.coupled())
234 alphaPhip = phiBf[patchi]*alphaBf[patchi];
const Mesh & mesh() const noexcept
Return const reference to mesh.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
virtual bool coupled() const
True if the patch field is coupled.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
static autoPtr< diameterModel > New(const dictionary &dict, const phaseModel &phase)
void correct()
Correct the phase properties.
virtual ~phaseModel()
Destructor.
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
tmp< volScalarField > d() const
const surfaceScalarField & alphaPhi() const
autoPtr< phaseModel > clone() const
Return clone.
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
virtual bool read()
Read phase properties dictionary.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Calculate the face-flux of the given field.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
const dimensionSet dimViscosity
List< word > wordList
List of word.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimSpecificHeatCapacity(dimGasConstant)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))