32template<
class CloudType>
40 SMDCalcMethod_(this->coeffDict().
lookup(
"SMDCalculationMethod"))
43 const scalar xx0 = 12.0;
45 1.0/(1.0 -
exp(-xx0)*(1.0 + xx0 +
sqr(xx0)/2.0 +
pow3(xx0)/6.0));
49 scalar xx = 0.12*(
n + 1);
51 (1.0 -
exp(-xx)*(1.0 + xx +
sqr(xx)/2.0 +
pow3(xx)/6.0))*rrd100;
54 if (SMDCalcMethod_ ==
"method1")
58 else if (SMDCalcMethod_ ==
"method2")
66 <<
"Unknown SMDCalculationMethod. Valid options are "
67 <<
"(method1 | method2). Using method2" <<
endl;
72template<
class CloudType>
76 SMDCalcMethod_(bum.SMDCalcMethod_)
82template<
class CloudType>
89template<
class CloudType>
121 scalar semiMass = nParticle*
pow3(d);
124 scalar rtd = 0.5*this->TABCmu_*
mu/(
rho*r2);
127 scalar omega2 = this->TABComega_*sigma/(
rho*r3) - rtd*rtd;
131 scalar omega =
sqrt(omega2);
132 scalar We = rhoc*
sqr(Urmag)*r/sigma;
133 scalar Wetmp = We/this->TABtwoWeCrit_;
135 scalar
y1 =
y - Wetmp;
136 scalar y2 = yDot/omega;
161 if ((Wetmp - a < -1) && (yDot < 0))
166 scalar theta =
acos((coste-Wetmp)/a);
176 tb = (theta-
phi)/omega;
182 yDot = -a*omega*
sin(omega*tb +
phi);
191 r/(1.0 + (4.0/3.0)*
sqr(
y) +
rho*r3/(8*sigma)*
sqr(yDot));
226 nParticle = semiMass/
pow3(d);
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
const uniformDimensionedVectorField & g
Templated break-up model class.
Templated base class for dsmc cloud.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
virtual ~TAB()
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool update()
Update the mesh for both mesh motion and topology change.
Lookup type of boundary radiation properties.
const volScalarField & mu
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar twoPi(2 *M_PI)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.