Go to the documentation of this file.
37 namespace kineticTheoryModels
39 namespace frictionalStressModels
45 frictionalStressModel,
46 JohnsonJacksonSchaeffer,
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_)
137 "JohnsonJacksonSchaeffer:nu",
150 0.5*pf[celli]*
sin(phi_.value())
171 mag(
U.boundaryField()[patchi].snGrad())
188 coeffDict_ <<= dict_.optionalSubDict(typeName +
"Coeffs");
190 Fr_.read(coeffDict_);
191 eta_.read(coeffDict_);
194 phi_.read(coeffDict_);
197 alphaDeltaMin_.read(coeffDict_);
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
const dimensionedScalar alpha
dimensionedScalar sin(const dimensionedScalar &ds)
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=PatchField< Type >::calculatedType())
Unit conversion functions.
const Type & value() const
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
virtual ~JohnsonJacksonSchaeffer()
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label max(const labelHashSet &set, label maxValue=labelMin)
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
defineTypeNameAndDebug(JohnsonJackson, 0)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
addToRunTimeSelectionTable(frictionalStressModel, JohnsonJackson, dictionary)
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Macros for easy insertion into run-time selection tables.
Generic dimensioned Type class.
Cmpt invariantII(const SymmTensor< Cmpt > &st)
constexpr scalar degToRad(const scalar deg) noexcept
void correctBoundaryConditions()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Internal & ref(const bool updateAccessTime=true)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar & D
const polyBoundaryMesh & patches
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Generic GeometricField class.
JohnsonJacksonSchaeffer(const dictionary &dict)
const Boundary & boundaryField() const
const dimensionSet dimless