54 #ifndef incompressibleDriftFlux_H
55 #define incompressibleDriftFlux_H
113 virtual void correctCoNum();
182 TypeName(
"incompressibleDriftFlux");
const dictionary & alphaControls
Generic GeometricField class.
const word & name() const
Return name.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Mesh data needed to do the Finite Volume discretisation.
virtual bool addsSupToField(const word &fieldName) const
Return true if an fvModel adds a source term to the given.
Class to represent a mixture of two constant density phases.
Provides controls for the pressure reference in closed-volume simulations.
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
const fvMesh & mesh
Region mesh.
const surfaceScalarField & phi
Reference to the mass-flux field.
const volVectorField & U
Reference to the velocity field.
Solver module for 2 incompressible fluids using the mixture approach with the drift-flux approximatio...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
virtual tmp< volScalarField > psiByRho() const
Return the mixture compressibility/density.
virtual tmp< surfaceScalarField > surfaceTensionForce() const
Return the interface surface tension force for the momentum equation.
virtual const Foam::pressureReference & pressureReference() const
Return the pressure reference.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual bool incompressible() const
The flow is incompressible.
virtual bool divergent() const
Is the flow divergent?
incompressibleDriftFlux(fvMesh &mesh)
Construct from region mesh.
virtual void alphaSuSp(tmp< volScalarField::Internal > &Su, tmp< volScalarField::Internal > &Sp)
Calculate the alpha equation sources.
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
autoPtr< relativeVelocityModel > relativeVelocity
Pointer to the dispersed phase relative velocity model.
TypeName("incompressibleDriftFlux")
Runtime type information.
virtual ~incompressibleDriftFlux()
Destructor.
incompressibleDriftFluxMixture & mixture
The compressible two-phase mixture.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
virtual void setInterfaceRDeltaT(volScalarField &rDeltaT)
Adjust the rDeltaT in the vicinity of the interface.
void operator=(const incompressibleDriftFlux &)=delete
Disallow default bitwise assignment.
virtual tmp< surfaceScalarField > alphaPhi(const surfaceScalarField &phi, const volScalarField &alpha, const dictionary &alphaControls)
virtual void correctInterface()
Correct the interface properties following mesh-change.
Foam::pressureReference pressureReference_
Pressure reference.
volScalarField p
Static pressure field.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U)
Return the momentum equation stress term.
autoPtr< compressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
Solver module base-class for for 2 immiscible fluids, with optional mesh motion and mesh topology cha...
volScalarField & alpha1
Reference to the phase1-fraction.
A class for managing temporary objects.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)