30 template<
class CloudType>
38 weCorrCoeff_(this->coeffDict().template lookup<scalar>(
"weCorrCoeff")),
39 weBuCrit_(this->coeffDict().template lookup<scalar>(
"weBuCrit")),
40 weBuBag_(this->coeffDict().template lookup<scalar>(
"weBuBag")),
41 weBuMM_(this->coeffDict().template lookup<scalar>(
"weBuMM")),
42 ohnCoeffCrit_(this->coeffDict().template lookup<scalar>(
"ohnCoeffCrit")),
43 ohnCoeffBag_(this->coeffDict().template lookup<scalar>(
"ohnCoeffBag")),
44 ohnCoeffMM_(this->coeffDict().template lookup<scalar>(
"ohnCoeffMM")),
45 ohnExpCrit_(this->coeffDict().template lookup<scalar>(
"ohnExpCrit")),
46 ohnExpBag_(this->coeffDict().template lookup<scalar>(
"ohnExpBag")),
47 ohnExpMM_(this->coeffDict().template lookup<scalar>(
"ohnExpMM")),
48 cInit_(this->coeffDict().template lookup<scalar>(
"Cinit")),
49 c1_(this->coeffDict().template lookup<scalar>(
"C1")),
50 c2_(this->coeffDict().template lookup<scalar>(
"C2")),
51 c3_(this->coeffDict().template lookup<scalar>(
"C3")),
52 cExp1_(this->coeffDict().template lookup<scalar>(
"Cexp1")),
53 cExp2_(this->coeffDict().template lookup<scalar>(
"Cexp2")),
54 cExp3_(this->coeffDict().template lookup<scalar>(
"Cexp3")),
55 weConst_(this->coeffDict().template lookup<scalar>(
"Weconst")),
56 weCrit1_(this->coeffDict().template lookup<scalar>(
"Wecrit1")),
57 weCrit2_(this->coeffDict().template lookup<scalar>(
"Wecrit2")),
58 coeffD_(this->coeffDict().template lookup<scalar>(
"CoeffD")),
59 onExpD_(this->coeffDict().template lookup<scalar>(
"OnExpD")),
60 weExpD_(this->coeffDict().template lookup<scalar>(
"WeExpD")),
61 mu_(this->coeffDict().template lookup<scalar>(
"mu")),
62 sigma_(this->coeffDict().template lookup<scalar>(
"sigma")),
63 d32Coeff_(this->coeffDict().template lookup<scalar>(
"d32Coeff")),
64 cDmaxBM_(this->coeffDict().template lookup<scalar>(
"cDmaxBM")),
65 cDmaxS_(this->coeffDict().template lookup<scalar>(
"cDmaxS")),
66 corePerc_(this->coeffDict().template lookup<scalar>(
"corePerc"))
70 template<
class CloudType>
74 weCorrCoeff_(bum.weCorrCoeff_),
75 weBuCrit_(bum.weBuCrit_),
76 weBuBag_(bum.weBuBag_),
78 ohnCoeffCrit_(bum.ohnCoeffCrit_),
79 ohnCoeffBag_(bum.ohnCoeffBag_),
80 ohnCoeffMM_(bum.ohnCoeffMM_),
81 ohnExpCrit_(bum.ohnExpCrit_),
82 ohnExpBag_(bum.ohnExpBag_),
83 ohnExpMM_(bum.ohnExpMM_),
91 weConst_(bum.weConst_),
92 weCrit1_(bum.weCrit1_),
93 weCrit2_(bum.weCrit2_),
99 d32Coeff_(bum.d32Coeff_),
100 cDmaxBM_(bum.cDmaxBM_),
101 cDmaxS_(bum.cDmaxS_),
102 corePerc_(bum.corePerc_)
108 template<
class CloudType>
115 template<
class CloudType>
137 const label injectori,
144 bool addChild =
false;
146 scalar d03 =
pow3(d);
148 scalar mass0 = nParticle*rhopi6*d03;
151 scalar weGas = 0.5*rhoc*
sqr(Urmag)*d/
sigma;
155 scalar reLiquid = 0.5*Urmag*d/
mu;
156 scalar ohnesorge =
sqrt(weLiquid)/(reLiquid + vSmall);
158 scalar weGasCorr = weGas/(1.0 + weCorrCoeff_*ohnesorge);
164 scalar rChar = Urmag/d*
sqrt(rhoc/
rho);
168 if (tc*rChar < small)
174 scalar tChar = 1/rChar;
176 scalar tFirst = cInit_*tChar;
179 scalar tCharSecond = 0;
182 bool multimode =
false;
187 if (weGas > weConst_)
189 if (weGas < weCrit1_)
191 tCharSecond = c1_*
pow((weGas - weConst_), cExp1_);
193 else if (weGas >= weCrit1_ && weGas <= weCrit2_)
195 tCharSecond = c2_*
pow((weGas - weConst_), cExp2_);
199 tCharSecond = c3_*
pow((weGas - weConst_), cExp3_);
203 scalar weC = weBuCrit_*(1.0 + ohnCoeffCrit_*
pow(ohnesorge, ohnExpCrit_));
204 scalar weB = weBuBag_*(1.0 + ohnCoeffBag_*
pow(ohnesorge, ohnExpBag_));
205 scalar weMM = weBuMM_*(1.0 + ohnCoeffMM_*
pow(ohnesorge, ohnExpMM_));
207 if (weGas > weC && weGas < weB)
212 if (weGas >= weB && weGas <= weMM)
222 tSecond = tCharSecond*tChar;
224 scalar tBreakUP = tFirst + tSecond;
227 scalar d32 = coeffD_*d*
pow(ohnesorge, onExpD_)*
pow(weGasCorr, weExpD_);
229 if (bag || multimode)
231 scalar d05 = d32Coeff_ * d32;
246 *
exp(-0.5*
sqr((
x - mu_)/sigma_));
260 scalar dC = weConst_*
sigma/(rhoc*
sqr(Urmag));
261 scalar d32Red = 4.0*(d32*dC)/(5.0*dC - d32);
263 scalar d05 = d32Coeff_ * d32Red;
278 *
exp(-0.5*
sqr((
x - mu_)/sigma_));
288 massChild = corePerc_*mass0;
297 nParticle = mass/(rhopi6*
pow3(d));
Templated break-up model class.
Templated base class for dsmc cloud.
Secondary Breakup Model to take account of the different breakup regimes, bag, solutionmode,...
SHF(const dictionary &, CloudType &)
Construct from dictionary.
virtual bool update(const scalar dt, const vector &g, scalar &d, scalar &tc, scalar &ms, scalar &nParticle, scalar &KHindex, scalar &y, scalar &yDot, const scalar d0, const scalar rho, const scalar mu, const scalar sigma, const vector &U, const scalar rhoc, const scalar muc, const vector &Urel, const scalar Urmag, const scalar tMom, const label injectori, scalar &dChild, scalar &massChild)
Update the parcel properties.
virtual ~SHF()
Destructor.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Type sample01()
Return a type with components uniformly distributed between.
const scalar twoPi(2 *pi)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
randomGenerator rndGen(653213)