SprayCloud.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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_(nullptr),
84  averageParcelMass_(0.0),
85  atomizationModel_(nullptr),
86  breakupModel_(nullptr)
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  this->deleteLostParticles();
98  }
99 
100  Info << "Average parcel mass: " << averageParcelMass_ << endl;
101  }
102 
103  if (this->solution().resetSourcesOnStartup())
104  {
105  CloudType::resetSourceTerms();
106  }
107 }
108 
109 
110 template<class CloudType>
112 (
114  const word& name
115 )
116 :
117  CloudType(c, name),
118  sprayCloud(),
119  cloudCopyPtr_(nullptr),
120  averageParcelMass_(c.averageParcelMass_),
121  atomizationModel_(c.atomizationModel_->clone()),
122  breakupModel_(c.breakupModel_->clone())
123 {}
124 
125 
126 template<class CloudType>
128 (
129  const fvMesh& mesh,
130  const word& name,
131  const SprayCloud<CloudType>& c
132 )
133 :
134  CloudType(mesh, name, c),
135  sprayCloud(),
136  cloudCopyPtr_(nullptr),
137  averageParcelMass_(0.0),
138  atomizationModel_(nullptr),
139  breakupModel_(nullptr)
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
145 template<class CloudType>
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
152 template<class CloudType>
154 (
155  parcelType& parcel,
156  const scalar lagrangianDt
157 )
158 {
159  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
160 
161  const liquidMixtureProperties& liqMix = this->composition().liquids();
162 
163  const scalarField& Y(parcel.Y());
164  scalarField X(liqMix.X(Y));
165  const scalar pc = this->p()[parcel.cell()];
166 
167  // override rho and Cp from constantProperties
168  parcel.Cp() = liqMix.Cp(pc, parcel.T(), X);
169  parcel.rho() = liqMix.rho(pc, parcel.T(), X);
170  parcel.sigma() = liqMix.sigma(pc, parcel.T(), X);
171  parcel.mu() = liqMix.mu(pc, parcel.T(), X);
172 }
173 
174 
175 template<class CloudType>
177 (
178  parcelType& parcel,
179  const scalar lagrangianDt,
180  const bool fullyDescribed
181 )
182 {
183  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
184 
185  // store the injection position and initial drop size
186  parcel.position0() = parcel.position();
187  parcel.d0() = parcel.d();
188 
189  parcel.y() = breakup().y0();
190  parcel.yDot() = breakup().yDot0();
191 
192  parcel.liquidCore() = atomization().initLiquidCore();
193 }
194 
195 
196 template<class CloudType>
198 {
199  cloudCopyPtr_.reset
200  (
201  static_cast<SprayCloud<CloudType>*>
202  (
203  clone(this->name() + "Copy").ptr()
204  )
205  );
206 }
207 
208 
209 template<class CloudType>
211 {
212  cloudReset(cloudCopyPtr_());
213  cloudCopyPtr_.clear();
214 }
215 
216 
217 template<class CloudType>
219 {
220  if (this->solution().canEvolve())
221  {
222  typename parcelType::trackingData td(*this);
223 
224  this->solve(*this, td);
225  }
226 }
227 
228 
229 template<class CloudType>
231 {
232  CloudType::info();
233  scalar d32 = 1.0e+6*this->Dij(3, 2);
234  scalar d10 = 1.0e+6*this->Dij(1, 0);
235  scalar dMax = 1.0e+6*this->Dmax();
236  scalar pen = this->penetration(0.95);
237 
238  Info << " D10, D32, Dmax (mu) = " << d10 << ", " << d32
239  << ", " << dMax << nl
240  << " Liquid penetration 95% mass (m) = " << pen << endl;
241 }
242 
243 
244 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture heat capacity [J/(kg K)].
DSMCCloud< dsmcParcel > CloudType
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
void storeState()
Store the current cloud state.
Definition: SprayCloud.C:197
basicSpecieMixture & composition
void restoreState()
Reset the current cloud to the previously stored state.
Definition: SprayCloud.C:210
Templated atomization model class.
Definition: SprayCloud.H:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
void setModels()
Set cloud sub-models.
Definition: SprayCloud.C:33
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: SprayCloud.C:154
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual ~SprayCloud()
Destructor.
Definition: SprayCloud.C:146
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: SprayCloud.C:177
scalarField X(const scalarField &Y) const
Returns the mole fractions corresponding to the given mass fractions.
A class for handling words, derived from string.
Definition: word.H:59
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
rhoEqn solve()
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:218
Templated base class for spray cloud.
Definition: SprayCloud.H:57
autoPtr< BreakupModel< SprayCloud< CloudType > > > breakupModel_
Break-up model.
Definition: SprayCloud.H:107
static const char nl
Definition: Ostream.H:265
void evolve()
Evolve the spray (inject, move)
Definition: SprayCloud.C:218
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void info()
Print cloud information.
Definition: SprayCloud.C:230
autoPtr< AtomizationModel< SprayCloud< CloudType > > > atomizationModel_
Atomization model.
Definition: SprayCloud.H:104
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void cloudReset(SprayCloud< CloudType > &c)
Reset state of cloud.
Definition: SprayCloud.C:57
messageStream Info
Templated break-up model class.
Definition: SprayCloud.H:50
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
volScalarField & p
Virtual abstract base class for templated SprayCloud.
Definition: sprayCloud.H:48