33 template<
class ParcelType>
34 template<
class TrackCloudType>
37 TrackCloudType&
cloud,
41 ParcelType::setCellValues(
cloud, td);
45 template<
class ParcelType>
46 template<
class TrackCloudType>
49 TrackCloudType&
cloud,
54 ParcelType::cellValueSourceCorrection(
cloud, td, dt);
58 template<
class ParcelType>
59 template<
class TrackCloudType>
62 TrackCloudType&
cloud,
67 typedef typename TrackCloudType::thermoCloudType thermoCloudType;
72 if (liquidCore() > 0.5)
75 cloud.forces().setCalcCoupled(
false);
83 const scalar
T0 = this->
T();
84 const scalar pc0 = td.pc();
92 cloud.constProps().setTMax(TMax);
99 const scalar mass0 = this->mass();
102 ParcelType::calc(
cloud,td, dt);
108 this->ms() -= this->ms()*(mass0 - this->mass())/mass0;
112 scalar T1 = this->
T();
115 this->
Cp() = composition.
liquids().
Cp(td.pc(), T1, X1);
119 scalar rho1 = composition.
liquids().
rho(td.pc(), T1, X1);
122 mu_ = composition.
liquids().
mu(td.pc(), T1, X1);
124 scalar d1 = this->d()*
cbrt(
rho0/rho1);
127 if (liquidCore() > 0.5)
129 calcAtomisation(
cloud, td, dt);
133 scalar d2 = this->d();
134 this->nParticle() *=
pow3(d1/d2);
138 calcBreakup(
cloud, td, dt);
143 cloud.forces().setCalcCoupled(
true);
147 template<
class ParcelType>
148 template<
class TrackCloudType>
151 TrackCloudType&
cloud,
156 if (injector_ == -1)
return;
158 typedef typename TrackCloudType::thermoCloudType thermoCloudType;
162 typedef typename TrackCloudType::sprayCloudType sprayCloudType;
167 scalar Wc = td.rhoc()*
RR*td.Tc()/td.pc();
169 scalar Tav = atomisation.
Taverage(this->
T(), td.Tc());
172 scalar rhoAv = td.pc()/(
R*Tav);
174 scalar soi =
cloud.injectors()[injector_].timeStart();
176 const vector&
pos = this->position(td.mesh);
177 const vector& injectionPos = this->position0();
181 scalar Urel =
mag(this->
U());
183 scalar t0 =
max(0.0, currentTime - this->age() - soi);
184 scalar t1 =
min(t0 + dt,
cloud.injectors()[injector_].timeEnd() - soi);
186 scalar
rho0 = mass0_/this->volume(d0_);
189 scalar volFlowRate =
cloud.injectors()[injector_].massToInject(t0, t1)/
rho0;
194 chi = this->chi(
cloud, td, composition.
liquids().
X(this->Y()));
218 template<
class ParcelType>
219 template<
class TrackCloudType>
222 TrackCloudType&
cloud,
227 const typename TrackCloudType::parcelType&
p =
228 static_cast<const typename TrackCloudType::parcelType&
>(*this);
229 typename TrackCloudType::parcelType::trackingData& ttd =
230 static_cast<typename TrackCloudType::parcelType::trackingData&
>(td);
232 const typename TrackCloudType::forceType& forces =
cloud.forces();
234 if (
cloud.breakup().solveOscillationEq())
236 solveTABEq(
cloud, td, dt);
240 scalar Wc = td.rhoc()*
RR*td.Tc()/td.pc();
242 scalar Tav =
cloud.atomisation().Taverage(this->
T(), td.Tc());
245 scalar rhoAv = td.pc()/(
R*Tav);
246 scalar muAv = td.muc();
247 vector Urel = this->
U() - td.Uc();
248 scalar Urmag =
mag(Urel);
249 scalar
Re = this->
Re(rhoAv, this->
U(), td.Uc(), this->d(), muAv);
251 const scalar mass =
p.mass();
252 const forceSuSp Fcp = forces.calcCoupled(
p, ttd, dt, mass,
Re, muAv);
253 const forceSuSp Fncp = forces.calcNonCoupled(
p, ttd, dt, mass,
Re, muAv);
254 this->tMom() = mass/(Fcp.
Sp() + Fncp.
Sp());
258 scalar parcelMassChild = 0.0;
262 cloud.breakup().update
289 scalar
Re = rhoAv*Urmag*dChild/muAv;
293 child->origId() = this->getNewParticleIndex();
295 child->
d0() = dChild;
296 const scalar massChild = child->mass();
297 child->
mass0() = massChild;
298 child->nParticle() = parcelMassChild/massChild;
301 forces.calcCoupled(*child, ttd, dt, massChild,
Re, muAv);
303 forces.calcNonCoupled(*child, ttd, dt, massChild,
Re, muAv);
308 child->
y() =
cloud.breakup().y0();
311 child->
ms() = -great;
312 child->
injector() = this->injector();
313 child->
tMom() = massChild/(Fcp.
Sp() + Fncp.
Sp());
314 child->calcDispersion(
cloud, td, dt);
316 cloud.addParticle(child);
321 template<
class ParcelType>
322 template<
class TrackCloudType>
325 TrackCloudType&
cloud,
332 typedef typename TrackCloudType::thermoCloudType thermoCloudType;
337 scalar
T0 = this->
T();
339 scalar pAmb =
cloud.pAmbient();
345 if (pv >= 0.999*pAmb)
352 scalar hl = liq.
hl(pAmb, TBoil);
353 scalar iTp = liq.
ha(pAmb,
T0) - pAmb/liq.
rho(pAmb,
T0);
354 scalar iTb = liq.
ha(pAmb, TBoil) - pAmb/liq.
rho(pAmb, TBoil);
356 chi += X[i]*(iTp - iTb)/hl;
360 chi =
min(1.0,
max(chi, 0.0));
366 template<
class ParcelType>
367 template<
class TrackCloudType>
370 TrackCloudType&
cloud,
375 const scalar& TABCmu =
cloud.breakup().TABCmu();
376 const scalar& TABtwoWeCrit =
cloud.breakup().TABtwoWeCrit();
377 const scalar& TABComega =
cloud.breakup().TABComega();
379 scalar r = 0.5*this->d();
384 scalar rtd = 0.5*TABCmu*mu_/(this->
rho()*r2);
387 scalar omega2 = TABComega*sigma_/(this->
rho()*r3) - rtd*rtd;
393 this->We(td.rhoc(), this->U(), td.Uc(), r, sigma_)/TABtwoWeCrit;
396 scalar
y0 = this->
y() - We;
397 scalar yDot0 = this->yDot() +
y0*rtd;
402 scalar
e =
exp(-rtd*dt);
405 this->yDot() = (We - this->
y())*rtd +
e*(yDot0*
c -
omega*
y0*
s);
418 template<
class ParcelType>
423 position0_(
p.position0_),
426 liquidCore_(
p.liquidCore_),
427 KHindex_(
p.KHindex_),
432 injector_(
p.injector_),
scalar Cp(const scalar p, const scalar T) const
#define R(A, B, C, D, E, F, K, M)
#define forAll(list, i)
Loop across all elements in list.
Templated atomisation model class.
virtual void update(const scalar dt, scalar &d, scalar &liquidCore, scalar &tc, const scalar rho, const scalar mu, const scalar sigma, const scalar volFlowRate, const scalar rhoAv, const scalar Urel, const vector &pos, const vector &injectionPos, const scalar pAmbient, const scalar chi, randomGenerator &rndGen) const =0
virtual bool calcChi() const =0
Flag to indicate if chi needs to be calculated.
scalar Taverage(const scalar &Tliq, const scalar &Tc) const
Average temperature calculation.
Templated reacting parcel composition model class Consists of carrier species (via thermo package),...
const liquidMixtureProperties & liquids() const
Return the global (additional) liquids.
const objectRegistry & db() const
Return the local objectRegistry.
Reaching spray parcel, with added functionality for atomisation and breakup.
scalar chi(TrackCloudType &cloud, trackingData &td, const scalarField &X) const
Calculate the chi-factor for flash-boiling for the.
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
scalar liquidCore() const
Return const access to liquid core.
label injector() const
Return const access to injector id.
scalar tMom() const
Return const access to momentum relaxation time.
void cellValueSourceCorrection(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct cell values using latest transfer information.
void calcAtomisation(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to atomisation model.
scalar yDot() const
Return const access to rate of change of spherical deviation.
scalar d0() const
Return const access to initial droplet diameter.
ParcelType::trackingData trackingData
Use base tracking data.
void solveTABEq(TrackCloudType &cloud, trackingData &td, const scalar dt)
Solve the TAB equation.
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
scalar tc() const
Return const access to atomisation characteristic time.
scalar mass0() const
Return const access to initial mass [kg].
void setCellValues(TrackCloudType &cloud, trackingData &td)
Set cell values.
void calcBreakup(TrackCloudType &cloud, trackingData &td, const scalar dt)
Correct parcel properties according to breakup model.
scalar ms() const
Return const access to stripped parcel mass.
scalar y() const
Return const access to spherical deviation.
void calc(TrackCloudType &cloud, trackingData &td, const scalar dt)
Update parcel properties over the time interval.
A cloud is a collection of lagrangian particles.
const Type & value() const
Return const reference to value.
Helper container for force Su and Sp terms.
scalar Sp() const
Return const access to the implicit coefficient [kg/s].
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
const PtrList< liquidProperties > & properties() const
Return the liquid properties.
scalar Tc(const scalarField &X) const
Calculate the critical temperature of mixture.
scalar pv(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture vapour pressure [Pa].
scalar pvInvert(const scalar p, const scalarField &X) const
Invert the vapour pressure relationship to retrieve the boiling.
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/kg/K].
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
The thermophysical properties of a liquid.
virtual scalar hl(scalar p, scalar T) const =0
Heat of vapourisation [J/kg].
virtual scalar ha(scalar p, scalar T) const =0
Liquid absolute enthalpy [J/kg].
virtual scalar pvInvert(scalar p) const
Invert the vapour pressure relationship to retrieve the.
const Time & time() const
Return time.
Convenience class to handle the input of constant rotational speed. Reads an omega entry with default...
virtual scalar rho(scalar p, scalar T) const =0
Density [kg/m^3].
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar e
Elementary charge.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar cbrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar cos(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
scalarList X0(nSpecie, 0.0)