driftFluxFoam::relativeVelocityModels: Added MichaelsBolger model for flocculated suspensions

Description
    Michaels & Bolger relative velocity model

    Reference:
    \verbatim
        Michaels, A. S., & Bolger, J. C. (1962).
        Settling rates and sediment volumes
        of flocculated kaolin suspensions.
        Industrial & Engineering Chemistry Fundamentals, 1(1), 24-33.
    \endverbatim

Usage
    Example usage:
    \verbatim
        relativeVelocityModel MichaelsBolger;

        MichaelsBolgerCoeffs
        {
            a0          0;    // Extended Michaels & Bolger coefficient,
            a1          4.65; // Exponent

            alphaMax    0.6;  // Maximum dispersed phase-fraction
                              // (packing fraction)
        }
    \endverbatim
This commit is contained in:
Henry Weller
2022-04-22 16:50:41 +01:00
parent 929c69340f
commit 64a2a3bf75
22 changed files with 326 additions and 155 deletions

View File

@ -24,6 +24,7 @@ License
\*---------------------------------------------------------------------------*/
#include "incompressibleTwoPhaseInteractingMixture.H"
#include "mixtureViscosityModel.H"
#include "relativeVelocityModel.H"
#include "fvcDiv.H"
#include "addToRunTimeSelectionTable.H"
@ -48,7 +49,9 @@ incompressibleTwoPhaseInteractingMixture
:
twoPhaseMixture(U.mesh()),
muModel_(mixtureViscosityModel::New(U.mesh(), phase1Name())),
U_(U),
muModel_(mixtureViscosityModel::New(*this)),
nucModel_(viscosityModel::New(U.mesh(), phase2Name())),
rhod_("rho", dimDensity, muModel_()),
@ -59,9 +62,7 @@ incompressibleTwoPhaseInteractingMixture
dimLength,
muModel_->lookupOrDefault("d", 0.0)
),
alphaMax_(muModel_->lookupOrDefault("alphaMax", 1.0)),
U_(U),
alphaMax_(lookupOrDefault("alphaMax", 1.0)),
g_(g),
@ -95,6 +96,20 @@ Foam::incompressibleTwoPhaseInteractingMixture::
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
const Foam::volScalarField&
Foam::incompressibleTwoPhaseInteractingMixture::alphad() const
{
return alpha1();
}
const Foam::volScalarField&
Foam::incompressibleTwoPhaseInteractingMixture::alphac() const
{
return alpha2();
}
const Foam::mixtureViscosityModel&
Foam::incompressibleTwoPhaseInteractingMixture::muModel() const
{

View File

@ -39,7 +39,6 @@ SourceFiles
#include "twoPhaseMixture.H"
#include "viscosityModel.H"
#include "mixtureViscosityModel.H"
#include "uniformDimensionedFields.H"
#include "IOMRFZoneList.H"
@ -48,6 +47,7 @@ SourceFiles
namespace Foam
{
class mixtureViscosityModel;
class relativeVelocityModel;
/*---------------------------------------------------------------------------*\
@ -61,6 +61,8 @@ class incompressibleTwoPhaseInteractingMixture
{
// Private data
volVectorField& U_;
autoPtr<mixtureViscosityModel> muModel_;
autoPtr<viscosityModel> nucModel_;
@ -73,8 +75,6 @@ class incompressibleTwoPhaseInteractingMixture
//- Optional maximum dispersed phase-fraction (e.g. packing limit)
scalar alphaMax_;
volVectorField& U_;
//- Acceleration due to gravity
const uniformDimensionedVectorField& g_;
@ -109,6 +109,12 @@ public:
// Member Functions
//- Return const-access to the continuous phase-fraction
const volScalarField& alphac() const;
//- Return const-access to the dispersed phase-fraction
const volScalarField& alphad() const;
//- Return const-access to the mixture viscosityModel
const mixtureViscosityModel& muModel() const;

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -49,11 +49,10 @@ namespace mixtureViscosityModels
Foam::mixtureViscosityModels::BinghamPlastic::BinghamPlastic
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
:
plastic(mesh, group),
plastic(mixture),
yieldStressCoeff_
(
"BinghamCoeff",
@ -95,7 +94,7 @@ Foam::mixtureViscosityModels::BinghamPlastic::mu
(
log10(vGreat),
yieldStressExponent_
*(max(alpha_, scalar(0)) + yieldStressOffset_)
*(max(mixture_.alphad(), scalar(0)) + yieldStressOffset_)
)
)
- pow

View File

@ -74,8 +74,11 @@ public:
// Constructors
//- Construct from components
BinghamPlastic(const fvMesh& mesh, const word& group);
//- Construct from mixture
BinghamPlastic
(
const incompressibleTwoPhaseInteractingMixture& mixture
);
//- Destructor

View File

@ -2,9 +2,12 @@ EXE_INC = \
-I../incompressibleTwoPhaseInteractingMixture \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/physicalProperties/lnInclude \
-I$(LIB_SRC)/twoPhaseModels/twoPhaseMixture/lnInclude
-I$(LIB_SRC)/twoPhaseModels/twoPhaseMixture/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
LIB_LIBS = \
-ltwoPhaseMixture \
-lphysicalProperties \
-lfiniteVolume
-lfiniteVolume \
-lmeshTools

View File

@ -49,34 +49,16 @@ namespace mixtureViscosityModels
Foam::mixtureViscosityModels::Quemada::Quemada
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
:
mixtureViscosityModel(mesh, group),
alphaMax_
(
"alphaMax",
dimless,
optionalSubDict(typeName + "Coeffs").lookup("alphaMax")
),
mixtureViscosityModel(mixture),
q_(optionalSubDict(typeName + "Coeffs").lookupOrDefault("q", scalar(2))),
muMax_
(
"muMax",
dimDynamicViscosity,
optionalSubDict(typeName + "Coeffs").lookup("muMax")
),
alpha_
(
mesh.lookupObject<volScalarField>
(
IOobject::groupName
(
lookupOrDefault<word>("alpha", "alpha"),
group
)
)
)
{}
@ -90,7 +72,11 @@ Foam::mixtureViscosityModels::Quemada::mu
const volVectorField& U
) const
{
return min(muc*pow(1.0 - alpha_/alphaMax_, -q_), muMax_);
return min
(
muc*pow(max(1 - mixture_.alphad()/mixture_.alphaMax(), 0.0), -q_),
muMax_
);
}
@ -100,8 +86,8 @@ bool Foam::mixtureViscosityModels::Quemada::read()
{
const dictionary& dict = optionalSubDict(typeName + "Coeffs");
dict.lookup("alphaMax") >> alphaMax_;
dict.lookup("q") >> q_;
dict.lookup("muMax") >> muMax_;
return true;
}

View File

@ -40,10 +40,11 @@ Usage
\verbatim
viscosityModel Quemada;
alphaMax 0.6; // Maximum dispersed phase-fraction (packing fraction)
q 2; // Exponent, defaults to 2
rho 1996; // Dispersed phase density
q 2; // Exponent, defaults to 2
muMax 1e-2; // Maximum viscosity (for numerical stability)
rho 1996;
\endverbatim
SourceFiles
@ -75,18 +76,12 @@ class Quemada
{
// Private data
//- Maximum dispersed phase-fraction (packing fraction)
dimensionedScalar alphaMax_;
//- Exponent (defaults to 2)
scalar q_;
//- Maximum viscosity
dimensionedScalar muMax_;
//- Dispersed phase fraction
const volScalarField& alpha_;
public:
@ -96,8 +91,8 @@ public:
// Constructors
//- Construct from components
Quemada(const fvMesh& mesh, const word& group);
//- Construct from mixture
Quemada(const incompressibleTwoPhaseInteractingMixture& mixture);
//- Destructor

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -39,11 +39,11 @@ namespace Foam
Foam::mixtureViscosityModel::mixtureViscosityModel
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
:
viscosityModel(mesh, group)
viscosityModel(mixture.U().mesh(), mixture.phase1Name()),
mixture_(mixture)
{}

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -42,7 +42,7 @@ SourceFiles
#ifndef mixtureViscosityModel_H
#define mixtureViscosityModel_H
#include "dictionary.H"
#include "incompressibleTwoPhaseInteractingMixture.H"
#include "viscosityModel.H"
#include "runTimeSelectionTables.H"
@ -60,6 +60,14 @@ class mixtureViscosityModel
public viscosityModel
{
protected:
// Protected data
//- Mixture properties
const incompressibleTwoPhaseInteractingMixture& mixture_;
public:
//- Runtime type information
@ -74,20 +82,22 @@ public:
mixtureViscosityModel,
dictionary,
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
),
(mesh, group)
(mixture)
);
// Constructors
//- Construct from components
mixtureViscosityModel(const fvMesh& U, const word& group);
//- Construct from mixture
mixtureViscosityModel
(
const incompressibleTwoPhaseInteractingMixture& mixture
);
//- Disallow default bitwise copy construction
mixtureViscosityModel(const mixtureViscosityModel&);
mixtureViscosityModel(const mixtureViscosityModel&) = delete;
// Selectors
@ -95,8 +105,7 @@ public:
//- Return a reference to the selected viscosity model
static autoPtr<mixtureViscosityModel> New
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
);

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -31,15 +31,18 @@ License
Foam::autoPtr<Foam::mixtureViscosityModel> Foam::mixtureViscosityModel::New
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
{
const word modelType
(
IOdictionary
(
viscosityModel::findModelDict(mesh, group)
viscosityModel::findModelDict
(
mixture.U().mesh(),
mixture.phase1Name()
)
).lookup("viscosityModel")
);
@ -58,7 +61,7 @@ Foam::autoPtr<Foam::mixtureViscosityModel> Foam::mixtureViscosityModel::New
<< exit(FatalError);
}
return autoPtr<mixtureViscosityModel>(cstrIter()(mesh, group));
return autoPtr<mixtureViscosityModel>(cstrIter()(mixture));
}

View File

@ -48,11 +48,10 @@ namespace mixtureViscosityModels
Foam::mixtureViscosityModels::plastic::plastic
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
:
mixtureViscosityModel(mesh, group),
mixtureViscosityModel(mixture),
plasticCoeffs_(optionalSubDict(typeName + "Coeffs")),
plasticViscosityCoeff_
(
@ -66,18 +65,7 @@ Foam::mixtureViscosityModels::plastic::plastic
dimless,
plasticCoeffs_.lookup("exponent")
),
muMax_("muMax", dimDynamicViscosity, plasticCoeffs_.lookup("muMax")),
alpha_
(
mesh.lookupObject<volScalarField>
(
IOobject::groupName
(
lookupOrDefault<word>("alpha", "alpha"),
group
)
)
)
muMax_("muMax", dimDynamicViscosity, plasticCoeffs_.lookup("muMax"))
{}
@ -98,7 +86,7 @@ Foam::mixtureViscosityModels::plastic::mu
pow
(
scalar(10),
plasticViscosityExponent_*alpha_
plasticViscosityExponent_*mixture_.alphad()
) - scalar(1)
),
muMax_

View File

@ -73,9 +73,6 @@ protected:
//- Maximum viscosity
dimensionedScalar muMax_;
//- Plastic phase fraction
const volScalarField& alpha_;
public:
@ -85,8 +82,8 @@ public:
// Constructors
//- Construct from components
plastic(const fvMesh& mesh, const word& group);
//- Construct from mixture
plastic(const incompressibleTwoPhaseInteractingMixture& mixture);
//- Destructor

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2014-2021 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -48,22 +48,10 @@ namespace mixtureViscosityModels
Foam::mixtureViscosityModels::slurry::slurry
(
const fvMesh& mesh,
const word& group
const incompressibleTwoPhaseInteractingMixture& mixture
)
:
mixtureViscosityModel(mesh, group),
alpha_
(
mesh.lookupObject<volScalarField>
(
IOobject::groupName
(
lookupOrDefault<word>("alpha", "alpha"),
group
)
)
)
mixtureViscosityModel(mixture)
{}
@ -76,9 +64,11 @@ Foam::mixtureViscosityModels::slurry::mu
const volVectorField& U
) const
{
const volScalarField& alphad = mixture_.alphad();
return
(
muc*(1.0 + 2.5*alpha_ + 10.05*sqr(alpha_) + 0.00273*exp(16.6*alpha_))
muc*(1.0 + 2.5*alphad + 10.05*sqr(alphad) + 0.00273*exp(16.6*alphad))
);
}

View File

@ -66,11 +66,6 @@ class slurry
:
public mixtureViscosityModel
{
// Private data
//- Slurry phase fraction
const volScalarField& alpha_;
public:
@ -80,8 +75,8 @@ public:
// Constructors
//- Construct from components
slurry(const fvMesh& mesh, const word& group);
//- Construct from mixture
slurry(const incompressibleTwoPhaseInteractingMixture& mixture);
//- Destructor

View File

@ -1,5 +1,6 @@
relativeVelocityModel/relativeVelocityModel.C
simple/simple.C
general/general.C
MichaelsBolger/MichaelsBolger.C
LIB = $(FOAM_LIBBIN)/libdriftFluxRelativeVelocityModels

View File

@ -0,0 +1,80 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "MichaelsBolger.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace relativeVelocityModels
{
defineTypeNameAndDebug(MichaelsBolger, 0);
addToRunTimeSelectionTable
(
relativeVelocityModel,
MichaelsBolger,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::relativeVelocityModels::MichaelsBolger::MichaelsBolger
(
const dictionary& dict,
const incompressibleTwoPhaseInteractingMixture& mixture,
const uniformDimensionedVectorField& g
)
:
relativeVelocityModel(dict, mixture, g),
a0_("a0", dimless, dict),
a1_("a1", dimless, dict),
alphaMax_("alphaMax", dimless, dict),
Vc_("Vc", dimTime, dict)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::relativeVelocityModels::MichaelsBolger::~MichaelsBolger()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void Foam::relativeVelocityModels::MichaelsBolger::correct()
{
Udm_ =
(mixture_.rhoc()/mixture_.rho())*Vc_
*(a0_ + pow(max(1 - mixture_.alphad()/mixture().alphaMax(), 0.0), a1_))
*acceleration();
}
// ************************************************************************* //

View File

@ -0,0 +1,130 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::relativeVelocityModels::MichaelsBolger
Description
Michaels & Bolger relative velocity model
Reference:
\verbatim
Michaels, A. S., & Bolger, J. C. (1962).
Settling rates and sediment volumes
of flocculated kaolin suspensions.
Industrial & Engineering Chemistry Fundamentals, 1(1), 24-33.
\endverbatim
Usage
Example usage:
\verbatim
relativeVelocityModel MichaelsBolger;
MichaelsBolgerCoeffs
{
a0 0; // Extended Michaels & Bolger coefficient,
a1 4.65; // Exponent
alphaMax 0.6; // Maximum dispersed phase-fraction
// (packing fraction)
}
\endverbatim
SourceFiles
MichaelsBolger.C
\*---------------------------------------------------------------------------*/
#ifndef MichaelsBolger_H
#define MichaelsBolger_H
#include "relativeVelocityModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace relativeVelocityModels
{
/*---------------------------------------------------------------------------*\
Class MichaelsBolger Declaration
\*---------------------------------------------------------------------------*/
class MichaelsBolger
:
public relativeVelocityModel
{
// Private Data
//- Extended Michaels & Bolger coefficient, defaults to 0
dimensionedScalar a0_;
//- Exponent, defaults to 4.65
dimensionedScalar a1_;
//- alphaMax coefficient
dimensionedScalar alphaMax_;
//- Drift velocity coefficient
dimensionedScalar Vc_;
public:
//- Runtime type information
TypeName("MichaelsBolger");
// Constructors
//- Construct from components
MichaelsBolger
(
const dictionary& dict,
const incompressibleTwoPhaseInteractingMixture& mixture,
const uniformDimensionedVectorField& g
);
//- Destructor
~MichaelsBolger();
// Member Functions
//- Update the diffusion velocity
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace relativeVelocityModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -65,11 +65,13 @@ Foam::relativeVelocityModels::general::~general()
void Foam::relativeVelocityModels::general::correct()
{
const volScalarField& alphad = mixture_.alphad();
Udm_ =
(rhoc_/rho())*Vc_*acceleration()
(mixture_.rhoc()/mixture_.rho())*Vc_*acceleration()
*(
exp(-a_*max(alphad_ - residualAlpha_, scalar(0)))
- exp(-a1_*max(alphad_ - residualAlpha_, scalar(0)))
exp(-a_*max(alphad - residualAlpha_, scalar(0)))
- exp(-a1_*max(alphad - residualAlpha_, scalar(0)))
);
}

View File

@ -76,22 +76,18 @@ Foam::relativeVelocityModel::relativeVelocityModel
)
:
mixture_(mixture),
alphac_(mixture.alpha2()),
alphad_(mixture.alpha1()),
rhoc_(mixture.rhoc()),
rhod_(mixture.rhod()),
g_(g),
Udm_
(
IOobject
(
"Udm",
alphac_.time().timeName(),
alphac_.mesh(),
mixture.U().time().timeName(),
mixture.U().mesh(),
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
alphac_.mesh(),
mixture_.U().mesh(),
dimensionedVector(dimVelocity, Zero),
UdmPatchFieldTypes()
)
@ -145,12 +141,6 @@ Foam::relativeVelocityModel::~relativeVelocityModel()
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
Foam::tmp<Foam::volScalarField> Foam::relativeVelocityModel::rho() const
{
return alphac_*rhoc_ + alphad_*rhod_;
}
Foam::tmp<Foam::volVectorField>
Foam::relativeVelocityModel::acceleration() const
{
@ -162,8 +152,8 @@ Foam::relativeVelocityModel::acceleration() const
Foam::tmp<Foam::volSymmTensorField> Foam::relativeVelocityModel::tauDm() const
{
const volScalarField betac(alphac_*rhoc_);
const volScalarField betad(alphad_*rhod_);
const volScalarField betac(mixture_.alphac()*mixture_.rhoc());
const volScalarField betad(mixture_.alphad()*mixture_.rhod());
// Calculate the relative velocity of the continuous phase w.r.t the mean
const volVectorField Ucm(betad*Udm_/betac);

View File

@ -62,27 +62,13 @@ protected:
//- Mixture properties
const incompressibleTwoPhaseInteractingMixture& mixture_;
//- Name of the continuous phase
const word continuousPhaseName_;
//- Continuous phase fraction
const volScalarField& alphac_;
//- Dispersed phase fraction
const volScalarField& alphad_;
//- Continuous density
const dimensionedScalar& rhoc_;
//- Dispersed density
const dimensionedScalar& rhod_;
//- Acceleration due to gravity
const uniformDimensionedVectorField& g_;
//- Dispersed diffusion velocity
mutable volVectorField Udm_;
public:
//- Runtime type information
@ -138,9 +124,6 @@ public:
return mixture_;
}
//- Return the mixture mean density
tmp<volScalarField> rho() const;
//- Return the diffusion velocity of the dispersed phase
const volVectorField& Udm() const
{

View File

@ -49,8 +49,7 @@ Foam::relativeVelocityModels::simple::simple
:
relativeVelocityModel(dict, mixture, g),
a_("a", dimless, dict),
Vc_("Vc", dimTime, dict),
residualAlpha_("residualAlpha", dimless, dict)
Vc_("Vc", dimTime, dict)
{}
@ -65,8 +64,8 @@ Foam::relativeVelocityModels::simple::~simple()
void Foam::relativeVelocityModels::simple::correct()
{
Udm_ =
(rhoc_/rho())*Vc_*acceleration()
*pow(scalar(10), -a_*max(alphad_, scalar(0)));
(mixture_.rhoc()/mixture_.rho())*Vc_*acceleration()
*pow(scalar(10), -a_*max(mixture_.alphad(), scalar(0)));
}

View File

@ -60,9 +60,6 @@ class simple
//- Drift velocity coefficient
dimensionedScalar Vc_;
//- Residual phase fraction
dimensionedScalar residualAlpha_;
public: