SprayCloud.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "SprayCloud.H"
27 #include "AtomizationModel.H"
28 #include "BreakupModel.H"
29 
30 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
31 
32 template<class CloudType>
34 {
35  atomizationModel_.reset
36  (
38  (
39  this->subModelProperties(),
40  *this
41  ).ptr()
42  );
43 
44  breakupModel_.reset
45  (
47  (
48  this->subModelProperties(),
49  *this
50  ).ptr()
51  );
52 }
53 
54 
55 template<class CloudType>
57 (
59 )
60 {
61  CloudType::cloudReset(c);
62 
63  atomizationModel_.reset(c.atomizationModel_.ptr());
64  breakupModel_.reset(c.breakupModel_.ptr());
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
70 template<class CloudType>
72 (
73  const word& cloudName,
74  const volScalarField& rho,
75  const volVectorField& U,
76  const dimensionedVector& g,
77  const SLGThermo& thermo,
78  bool readFields
79 )
80 :
81  CloudType(cloudName, rho, U, g, thermo, false),
82  sprayCloud(),
83  cloudCopyPtr_(NULL),
84  averageParcelMass_(0.0),
85  atomizationModel_(NULL),
86  breakupModel_(NULL)
87 {
88  if (this->solution().active())
89  {
90  setModels();
91 
92  averageParcelMass_ = this->injectors().averageParcelMass();
93 
94  if (readFields)
95  {
96  parcelType::readFields(*this, this->composition());
97  }
98 
99  Info << "Average parcel mass: " << averageParcelMass_ << endl;
100  }
101 
102  if (this->solution().resetSourcesOnStartup())
103  {
104  CloudType::resetSourceTerms();
105  }
106 }
107 
108 
109 template<class CloudType>
111 (
113  const word& name
114 )
115 :
116  CloudType(c, name),
117  sprayCloud(),
118  cloudCopyPtr_(NULL),
119  averageParcelMass_(c.averageParcelMass_),
120  atomizationModel_(c.atomizationModel_->clone()),
121  breakupModel_(c.breakupModel_->clone())
122 {}
123 
124 
125 template<class CloudType>
127 (
128  const fvMesh& mesh,
129  const word& name,
130  const SprayCloud<CloudType>& c
131 )
132 :
133  CloudType(mesh, name, c),
134  sprayCloud(),
135  cloudCopyPtr_(NULL),
136  averageParcelMass_(0.0),
137  atomizationModel_(NULL),
138  breakupModel_(NULL)
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
143 
144 template<class CloudType>
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 template<class CloudType>
153 (
154  parcelType& parcel,
155  const scalar lagrangianDt
156 )
157 {
158  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
159 
160  const liquidMixtureProperties& liqMix = this->composition().liquids();
161 
162  const scalarField& Y(parcel.Y());
163  scalarField X(liqMix.X(Y));
164 
165  // override rho and Cp from constantProperties
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);
170 }
171 
172 
173 template<class CloudType>
175 (
176  parcelType& parcel,
177  const scalar lagrangianDt,
178  const bool fullyDescribed
179 )
180 {
181  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
182 
183  // store the injection position and initial drop size
184  parcel.position0() = parcel.position();
185  parcel.d0() = parcel.d();
186 
187  parcel.y() = breakup().y0();
188  parcel.yDot() = breakup().yDot0();
189 
190  parcel.liquidCore() = atomization().initLiquidCore();
191 }
192 
193 
194 template<class CloudType>
196 {
197  cloudCopyPtr_.reset
198  (
199  static_cast<SprayCloud<CloudType>*>
200  (
201  clone(this->name() + "Copy").ptr()
202  )
203  );
204 }
205 
206 
207 template<class CloudType>
209 {
210  cloudReset(cloudCopyPtr_());
211  cloudCopyPtr_.clear();
212 }
213 
214 
215 template<class CloudType>
217 {
218  if (this->solution().canEvolve())
219  {
220  typename parcelType::template
221  TrackingData<SprayCloud<CloudType> > td(*this);
222 
223  this->solve(td);
224  }
225 }
226 
227 
228 template<class CloudType>
230 {
231  CloudType::info();
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);
236 
237  Info << " D10, D32, Dmax (mu) = " << d10 << ", " << d32
238  << ", " << dMax << nl
239  << " Liquid penetration 95% mass (m) = " << pen << endl;
240 }
241 
242 
243 // ************************************************************************* //
Foam::SprayCloud
Templated base class for spray cloud.
Definition: SprayCloud.H:57
Foam::solution
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:62
Foam::SprayCloud::storeState
void storeState()
Store the current cloud state.
Definition: SprayCloud.C:195
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::AtomizationModel
Templated atomization model class.
Definition: SprayCloud.H:47
Foam::liquidMixtureProperties::mu
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
Definition: liquidMixtureProperties.C:420
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::SprayCloud::restoreState
void restoreState()
Reset the current cloud to the previously stored state.
Definition: SprayCloud.C:208
Foam::liquidMixtureProperties::rho
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
Definition: liquidMixtureProperties.C:275
Foam::sprayCloud
Virtual abstract base class for templated SprayCloud.
Definition: sprayCloud.H:48
U
U
Definition: pEqn.H:46
Foam::SprayCloud::info
void info()
Print cloud information.
Definition: SprayCloud.C:229
Foam::liquidMixtureProperties::X
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
Definition: liquidMixtureProperties.C:257
solve
rhoEqn solve()
Foam::SprayCloud::cloudReset
void cloudReset(SprayCloud< CloudType > &c)
Reset state of cloud.
Definition: SprayCloud.C:57
Foam::liquidMixtureProperties::sigma
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
Definition: liquidMixtureProperties.C:383
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
AtomizationModel.H
Foam::SprayCloud::checkParcelProperties
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: SprayCloud.C:175
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:36
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::BreakupModel
Templated break-up model class.
Definition: SprayCloud.H:50
Foam::SprayCloud::SprayCloud
SprayCloud(const SprayCloud &)
Disallow default bitwise copy construct.
rho
rho
Definition: pEqn.H:3
Foam::SprayCloud::~SprayCloud
virtual ~SprayCloud()
Destructor.
Definition: SprayCloud.C:145
BreakupModel.H
Foam::liquidMixtureProperties::Cp
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
Definition: liquidMixtureProperties.C:357
Foam::SprayCloud::setParcelThermoProperties
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: SprayCloud.C:153
Foam::liquidMixtureProperties
A mixture of liquids.
Definition: liquidMixtureProperties.H:74
Foam::SprayCloud::evolve
void evolve()
Evolve the spray (inject, move)
Definition: SprayCloud.C:216
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
composition
basicMultiComponentMixture & composition
Definition: createFields.H:35
cloudName
const word cloudName(propsDict.lookup("cloudName"))
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
SprayCloud.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::SprayCloud::setModels
void setModels()
Set cloud sub-models.
Definition: SprayCloud.C:33
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36