/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2019 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 .
\*---------------------------------------------------------------------------*/
#include "JohnsonJacksonFrictionalStress.H"
#include "addToRunTimeSelectionTable.H"
#include "unitConversion.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace kineticTheoryModels
{
namespace frictionalStressModels
{
defineTypeNameAndDebug(JohnsonJackson, 0);
addToRunTimeSelectionTable
(
frictionalStressModel,
JohnsonJackson,
dictionary
);
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::
JohnsonJackson
(
const dictionary& dict
)
:
frictionalStressModel(dict),
coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict_),
eta_("eta", dimless, coeffDict_),
p_("p", dimless, coeffDict_),
phi_("phi", dimless, coeffDict_),
alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict_)
{
phi_ *= degToRad();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::
~JohnsonJackson()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::
frictionalPressure
(
const phaseModel& phase,
const dimensionedScalar& alphaMinFriction,
const dimensionedScalar& alphaMax
) const
{
const volScalarField& alpha = phase;
return
Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
/pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
}
Foam::tmp
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::
frictionalPressurePrime
(
const phaseModel& phase,
const dimensionedScalar& alphaMinFriction,
const dimensionedScalar& alphaMax
) const
{
const volScalarField& alpha = phase;
return Fr_*
(
eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1.0)
*(alphaMax-alpha)
+ p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
)/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1.0);
}
Foam::tmp
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::nu
(
const phaseModel& phase,
const dimensionedScalar& alphaMinFriction,
const dimensionedScalar& alphaMax,
const volScalarField& pf,
const volSymmTensorField& D
) const
{
return dimensionedScalar("0.5", dimTime, 0.5)*pf*sin(phi_);
}
bool Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::read()
{
coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
Fr_.read(coeffDict_);
eta_.read(coeffDict_);
p_.read(coeffDict_);
phi_.read(coeffDict_);
phi_ *= degToRad();
alphaDeltaMin_.read(coeffDict_);
return true;
}
// ************************************************************************* //