Go to the documentation of this file.
26 #include "SchaefferFrictionalStress.H"
33 namespace kineticTheoryModels
35 namespace frictionalStressModels
41 frictionalStressModel,
54 const dictionary&
dict
57 frictionalStressModel(
dict),
58 coeffDict_(
dict.subDict(typeName +
"Coeffs")),
59 phi_(
"phi",
dimless, coeffDict_)
113 const scalar I2Dsmall = 1.0e-15;
117 tmp<volScalarField> tnu
124 alpha1.mesh().time().timeName(),
139 if (
alpha1[celli] > alphaMinFriction.value())
142 0.5*pf[celli]*
sin(phi_.value())
144 sqrt(1.0/6.0*(
sqr(D[celli].xx() - D[celli].yy())
145 +
sqr(D[celli].yy() - D[celli].zz())
146 +
sqr(D[celli].zz() - D[celli].xx()))
147 +
sqr(D[celli].xy()) +
sqr(D[celli].xz())
148 +
sqr(D[celli].yz())) + I2Dsmall
154 nuf.correctBoundaryConditions();
162 coeffDict_ <<= dict_.subDict(typeName +
"Coeffs");
164 phi_.read(coeffDict_);
dimensionedScalar phi_
Angle of internal friction.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
virtual tmp< volScalarField > nu(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
A class for managing temporary objects.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dimensionedScalar sin(const dimensionedScalar &ds)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Macros for easy insertion into run-time selection tables.
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
virtual tmp< volScalarField > frictionalPressurePrime(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Schaeffer(const dictionary &dict)
Construct from components.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
defineTypeNameAndDebug(combustionModel, 0)
virtual ~Schaeffer()
Destructor.
virtual tmp< volScalarField > frictionalPressure(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const