Go to the documentation of this file.
29 #include "phasePair.H"
31 #include "surfaceInterpolate.H"
51 const dictionary&
dict,
52 const phasePair& pair,
53 const bool registerObject
56 dragModel(
dict, pair, registerObject),
74 <<
"Drag coefficient not defined for the segregated model."
77 return pair_.phase1();
83 const fvMesh&
mesh(pair_.phase1().mesh());
91 tmp<volScalarField> tnu1(pair_.phase1().nu());
92 tmp<volScalarField> tnu2(pair_.phase2().nu());
102 mesh.time().timeName(),
107 zeroGradientFvPatchField<scalar>::typeName
110 L.correctBoundaryConditions();
118 pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha()
126 (pair_.phase1().residualAlpha() + pair_.phase2().residualAlpha())/
L
const vector L(dict.get< vector >("L"))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
defineTypeNameAndDebug(blended, 0)
A class for managing temporary objects.
static constexpr const zero Zero
const volScalarField & alpha2
const volScalarField & alpha1
CGAL::Exact_predicates_exact_constructions_kernel K
virtual tmp< volScalarField > K() const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< surfaceScalarField > Kf() const
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
addToRunTimeSelectionTable(dragModel, blended, dictionary)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual tmp< volScalarField > CdRe() const
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Calculate the gradient of the given field.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar cbrt(const dimensionedScalar &ds)
segregated(const dictionary &dict, const phasePair &pair, const bool registerObject)
static const Identity< scalar > I
const dimensionSet dimless