Go to the documentation of this file.
32 template<
class CloudType>
35 atomizationModel_.reset
39 this->subModelProperties(),
48 this->subModelProperties(),
55 template<
class CloudType>
61 CloudType::cloudReset(
c);
63 atomizationModel_.reset(
c.atomizationModel_.ptr());
64 breakupModel_.reset(
c.breakupModel_.ptr());
70 template<
class CloudType>
84 averageParcelMass_(0.0),
85 atomizationModel_(NULL),
92 averageParcelMass_ = this->injectors().averageParcelMass();
99 Info <<
"Average parcel mass: " << averageParcelMass_ <<
endl;
102 if (this->
solution().resetSourcesOnStartup())
104 CloudType::resetSourceTerms();
109 template<
class CloudType>
119 averageParcelMass_(
c.averageParcelMass_),
120 atomizationModel_(
c.atomizationModel_->clone()),
121 breakupModel_(
c.breakupModel_->clone())
125 template<
class CloudType>
136 averageParcelMass_(0.0),
137 atomizationModel_(NULL),
144 template<
class CloudType>
151 template<
class CloudType>
155 const scalar lagrangianDt
158 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
166 parcel.Cp() = liqMix.
Cp(parcel.pc(), parcel.T(), X);
167 parcel.rho() = liqMix.
rho(parcel.pc(), parcel.T(), X);
168 parcel.sigma() = liqMix.
sigma(parcel.pc(), parcel.T(), X);
169 parcel.mu() = liqMix.
mu(parcel.pc(), parcel.T(), X);
173 template<
class CloudType>
177 const scalar lagrangianDt,
178 const bool fullyDescribed
181 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
184 parcel.position0() = parcel.position();
185 parcel.d0() = parcel.d();
187 parcel.y() = breakup().y0();
188 parcel.yDot() = breakup().yDot0();
190 parcel.liquidCore() = atomization().initLiquidCore();
194 template<
class CloudType>
201 clone(this->
name() +
"Copy").ptr()
207 template<
class CloudType>
210 cloudReset(cloudCopyPtr_());
211 cloudCopyPtr_.clear();
215 template<
class CloudType>
220 typename parcelType::template
221 TrackingData<SprayCloud<CloudType> > td(*
this);
228 template<
class CloudType>
232 scalar d32 = 1.0e+6*this->Dij(3, 2);
233 scalar d10 = 1.0e+6*this->Dij(1, 0);
234 scalar dMax = 1.0e+6*this->Dmax();
235 scalar pen = this->penetration(0.95);
237 Info <<
" D10, D32, Dmax (mu) = " << d10 <<
", " << d32
238 <<
", " << dMax <<
nl
239 <<
" Liquid penetration 95% mass (m) = " << pen <<
endl;
Templated base class for spray cloud.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
void storeState()
Store the current cloud state.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Templated atomization model class.
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
This function object reads fields from the time directories and adds them to the mesh database for fu...
void restoreState()
Reset the current cloud to the previously stored state.
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
Virtual abstract base class for templated SprayCloud.
void info()
Print cloud information.
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
void cloudReset(SprayCloud< CloudType > &c)
Reset state of cloud.
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
Pre-declare SubField and related Field type.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
DSMCCloud< dsmcParcel > CloudType
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Templated break-up model class.
SprayCloud(const SprayCloud &)
Disallow default bitwise copy construct.
virtual ~SprayCloud()
Destructor.
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
void evolve()
Evolve the spray (inject, move)
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
basicMultiComponentMixture & composition
const word cloudName(propsDict.lookup("cloudName"))
const dimensionedScalar c
Speed of light in a vacuum.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Generic GeometricField class.
void setModels()
Set cloud sub-models.
word name(const complex &)
Return a string representation of a complex.
PtrList< volScalarField > & Y