multiphaseEulerFoam::kineticTheoryModels, phasePressureModel: now use the phase.alphaMax()

rather than read and use a potentially inconsistent local value for alphaMax.
This commit is contained in:
Henry Weller
2022-05-25 13:47:18 +01:00
parent 51307a5808
commit 6d3e31f8e0
13 changed files with 34 additions and 73 deletions

View File

@ -192,21 +192,6 @@ void Foam::JohnsonJacksonParticleSlipFvPatchVectorField::updateCoeffs()
: alpha
);
// lookup the packed volume fraction
dimensionedScalar alphaMax
(
"alphaMax",
dimless,
db()
.lookupObject<IOdictionary>
(
IOobject::groupName("momentumTransport", phase.name())
)
.subDict("RAS")
.subDict("kineticTheoryCoeffs")
.lookup("alphaMax")
);
// calculate the slip value fraction
scalarField c
(
@ -215,7 +200,7 @@ void Foam::JohnsonJacksonParticleSlipFvPatchVectorField::updateCoeffs()
*gs0
*specularityCoefficient_.value()
*sqrt(3*Theta)
/max(6*(nu - nuFric)*alphaMax.value(), small)
/max(6*(nu - nuFric)*phase.alphaMax(), small)
);
this->valueFraction() = c/(c + patch().deltaCoeffs());

View File

@ -205,21 +205,6 @@ void Foam::JohnsonJacksonParticleThetaFvPatchScalarField::updateCoeffs()
const scalarField Theta(patchInternalField());
// lookup the packed volume fraction
dimensionedScalar alphaMax
(
"alphaMax",
dimless,
db()
.lookupObject<IOdictionary>
(
IOobject::groupName("momentumTransport", phase.name())
)
.subDict("RAS")
.subDict("kineticTheoryCoeffs")
.lookup("alphaMax")
);
// calculate the reference value and the value fraction
if (restitutionCoefficient_.value() != 1.0)
{
@ -238,7 +223,7 @@ void Foam::JohnsonJacksonParticleThetaFvPatchScalarField::updateCoeffs()
*gs0
*(scalar(1) - sqr(restitutionCoefficient_.value()))
*sqrt(3*Theta)
/max(4*kappa*alphaMax.value(), small)
/max(4*kappa*phase.alphaMax(), small)
);
this->valueFraction() = c/(c + patch().deltaCoeffs());
@ -258,7 +243,7 @@ void Foam::JohnsonJacksonParticleThetaFvPatchScalarField::updateCoeffs()
*gs0
*sqrt(3*Theta)
*magSqr(U)
/max(6*kappa*alphaMax.value(), small);
/max(6*kappa*phase.alphaMax(), small);
this->valueFraction() = 0;
}

View File

@ -130,7 +130,6 @@ Foam::RASModels::kineticTheoryModel::kineticTheoryModel
equilibrium_(coeffDict_.lookup("equilibrium")),
e_("e", dimless, coeffDict_),
alphaMax_("alphaMax", dimless, coeffDict_),
alphaMinFriction_
(
"alphaMinFriction",
@ -245,7 +244,6 @@ bool Foam::RASModels::kineticTheoryModel::read()
{
coeffDict().lookup("equilibrium") >> equilibrium_;
e_.readIfPresent(coeffDict());
alphaMax_.readIfPresent(coeffDict());
alphaMinFriction_.readIfPresent(coeffDict());
viscosityModel_->read();
@ -316,8 +314,18 @@ Foam::RASModels::kineticTheoryModel::pPrime() const
*granularPressureModel_->granularPressureCoeffPrime
(
alpha_,
radialModel_->g0(alpha_, alphaMinFriction_, alphaMax_),
radialModel_->g0prime(alpha_, alphaMinFriction_, alphaMax_),
radialModel_->g0
(
alpha_,
alphaMinFriction_,
phase_.alphaMax()
),
radialModel_->g0prime
(
alpha_,
alphaMinFriction_,
phase_.alphaMax()
),
rho,
e_
)
@ -325,7 +333,7 @@ Foam::RASModels::kineticTheoryModel::pPrime() const
(
phase_,
alphaMinFriction_,
alphaMax_
phase_.alphaMax()
)
)
);
@ -413,7 +421,7 @@ void Foam::RASModels::kineticTheoryModel::correct()
const volSymmTensorField D(symm(gradU));
// Calculating the radial distribution function
gs0_ = radialModel_->g0(alpha, alphaMinFriction_, alphaMax_);
gs0_ = radialModel_->g0(alpha, alphaMinFriction_, phase_.alphaMax());
if (!equilibrium_)
{
@ -587,7 +595,7 @@ void Foam::RASModels::kineticTheoryModel::correct()
(
phase_,
alphaMinFriction_,
alphaMax_
phase_.alphaMax()
)
);
@ -600,7 +608,7 @@ void Foam::RASModels::kineticTheoryModel::correct()
(
phase_,
alphaMinFriction_,
alphaMax_,
phase_.alphaMax(),
pf/rho,
D
),

View File

@ -56,7 +56,6 @@ SourceFiles
#include "granularPressureModel.H"
#include "frictionalStressModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
@ -74,12 +73,10 @@ class kineticTheoryModel
{
// Private Data
// Input Fields
const phaseModel& phase_;
const phaseModel& phase_;
//- Name of the continuous phase
word continuousPhaseName_;
//- Name of the continuous phase
word continuousPhaseName_;
// Sub-models
@ -110,9 +107,6 @@ class kineticTheoryModel
//- Coefficient of restitution
dimensionedScalar e_;
//- Maximum packing phase-fraction
dimensionedScalar alphaMax_;
//- Min value for which the frictional stresses are zero
dimensionedScalar alphaMinFriction_;

View File

@ -49,7 +49,8 @@ Foam::RASModels::phasePressureModel::phasePressureModel
viscosity
),
alphaMax_(coeffDict_.lookup<scalar>("alphaMax")),
phase_(refCast<const phaseModel>(viscosity)),
preAlphaExp_(coeffDict_.lookup<scalar>("preAlphaExp")),
expMax_(coeffDict_.lookup<scalar>("expMax")),
g0_
@ -84,7 +85,6 @@ bool Foam::RASModels::phasePressureModel::read()
read()
)
{
coeffDict().lookup("alphaMax") >> alphaMax_;
coeffDict().lookup("preAlphaExp") >> preAlphaExp_;
coeffDict().lookup("expMax") >> expMax_;
g0_.readIfPresent(coeffDict());
@ -153,7 +153,7 @@ Foam::RASModels::phasePressureModel::pPrime() const
g0_
*min
(
exp(preAlphaExp_*(alpha_ - alphaMax_)),
exp(preAlphaExp_*(alpha_ - phase_.alphaMax())),
expMax_
)
)
@ -185,7 +185,8 @@ Foam::RASModels::phasePressureModel::pPrimef() const
g0_
*min
(
exp(preAlphaExp_*(fvc::interpolate(alpha_) - alphaMax_)),
exp(preAlphaExp_
*(fvc::interpolate(alpha_) - phase_.alphaMax())),
expMax_
)
)

View File

@ -38,7 +38,6 @@ Description
{
preAlphaExp 500;
expMax 1000;
alphaMax 0.62;
g0 1000;
}
\endverbatim
@ -73,10 +72,9 @@ class phasePressureModel
{
// Private Data
// Kinetic Theory Model coefficients
const phaseModel& phase_;
//- Maximum packing phase-fraction
scalar alphaMax_;
// Phase pressure coefficients
//- Pre-exponential factor
scalar preAlphaExp_;

View File

@ -28,7 +28,6 @@ RAS
equilibrium off;
e 0.8;
alphaMax 0.65;
alphaMinFriction 0.5;
residualAlpha 1e-6;

View File

@ -29,7 +29,8 @@ solids
d 462e-6;
}
residualAlpha 1e-5;
alphaMax 0.65;
residualAlpha 1e-5;
}
gas

View File

@ -28,7 +28,6 @@ RAS
equilibrium off;
e 0.8;
alphaMax 0.62;
alphaMinFriction 0.5;
residualAlpha 1e-4;
@ -52,7 +51,6 @@ RAS
{
preAlphaExp 500;
expMax 1000;
alphaMax 0.62;
g0 1000;
}
}

View File

@ -53,13 +53,7 @@ blending
}
surfaceTension
{
air_particles
{
type constant;
sigma 0;
}
}
{}
interfaceCompression
{}

View File

@ -26,7 +26,6 @@ RAS
{
preAlphaExp 500;
expMax 1000;
alphaMax 0.62;
g0 1000;
}
}

View File

@ -66,6 +66,7 @@ particles
);
}
alphaMax 0.62;
residualAlpha 1e-8;
}
@ -176,4 +177,4 @@ turbulentDispersion
}
}
// ************************************************************************* //
// ************************************************************************* //

View File

@ -28,7 +28,6 @@ RAS
equilibrium on;
e 0.8;
alphaMax 0.62;
alphaMinFriction 0.5;
granularViscosityModel Gidaspow;
@ -51,7 +50,6 @@ RAS
{
preAlphaExp 500;
expMax 1000;
alphaMax 0.62;
g0 1000;
}
}