Go to the documentation of this file.
37 namespace kineticTheoryModels
39 namespace frictionalStressModels
63 coeffDict_(
dict.optionalSubDict(typeName +
"Coeffs")),
65 eta_(
"eta",
dimless, coeffDict_),
67 phi_(
"phi",
dimless, coeffDict_),
68 alphaDeltaMin_(
"alphaDeltaMin",
dimless, coeffDict_)
95 Fr_*
pow(
max(
alpha - alphaMinFriction, scalar(0)), eta_)
113 eta_*
pow(
max(
alpha - alphaMinFriction, scalar(0)), eta_ - 1)
115 + p_*
pow(
max(
alpha - alphaMinFriction, scalar(0)), eta_)
136 coeffDict_ <<= dict_.optionalSubDict(typeName +
"Coeffs");
138 Fr_.read(coeffDict_);
139 eta_.read(coeffDict_);
142 phi_.read(coeffDict_);
145 alphaDeltaMin_.read(coeffDict_);
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensionedScalar sin(const dimensionedScalar &ds)
Unit conversion functions.
Dimension set for the base types.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
virtual ~JohnsonJackson()
Destructor.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
JohnsonJackson(const dictionary &dict)
Construct from components.
const dimensionedScalar & D
addToRunTimeSelectionTable(frictionalStressModel, JohnsonJackson, dictionary)
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
defineTypeNameAndDebug(JohnsonJackson, 0)