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-2021 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 "AtomisationModel.H"
28 #include "BreakupModel.H"
29 #include "parcelThermo.H"
30 
31 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
36  atomisationModel_.reset
37  (
39  (
40  this->subModelProperties(),
41  *this
42  ).ptr()
43  );
44 
45  breakupModel_.reset
46  (
48  (
49  this->subModelProperties(),
50  *this
51  ).ptr()
52  );
53 }
54 
55 
56 template<class CloudType>
58 (
60 )
61 {
62  CloudType::cloudReset(c);
63 
64  atomisationModel_.reset(c.atomisationModel_.ptr());
65  breakupModel_.reset(c.breakupModel_.ptr());
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class CloudType>
73 (
74  const word& cloudName,
75  const volScalarField& rho,
76  const volVectorField& U,
77  const dimensionedVector& g,
78  const fluidThermo& carrierThermo,
79  bool readFields
80 )
81 :
82  CloudType(cloudName, rho, U, g, carrierThermo, false),
83  cloudCopyPtr_(nullptr),
84  averageParcelMass_(0.0),
85  atomisationModel_(nullptr),
86  breakupModel_(nullptr)
87 {
88  setModels();
89 
90  averageParcelMass_ = this->injectors().averageParcelMass();
91 
92  if (readFields)
93  {
94  parcelType::readFields(*this, this->composition());
95  this->deleteLostParticles();
96  }
97 
98  Info << "Average parcel mass: " << averageParcelMass_ << endl;
99 
100  if (this->solution().resetSourcesOnStartup())
101  {
102  CloudType::resetSourceTerms();
103  }
104 }
105 
106 
107 template<class CloudType>
109 (
111  const word& name
112 )
113 :
114  CloudType(c, name),
115  cloudCopyPtr_(nullptr),
116  averageParcelMass_(c.averageParcelMass_),
117  atomisationModel_(c.atomisationModel_->clone()),
118  breakupModel_(c.breakupModel_->clone())
119 {}
120 
121 
122 template<class CloudType>
124 (
125  const fvMesh& mesh,
126  const word& name,
127  const SprayCloud<CloudType>& c
128 )
129 :
130  CloudType(mesh, name, c),
131  cloudCopyPtr_(nullptr),
132  averageParcelMass_(0.0),
133  atomisationModel_(nullptr),
134  breakupModel_(nullptr)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
139 
140 template<class CloudType>
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class CloudType>
149 (
150  parcelType& parcel,
151  const scalar lagrangianDt
152 )
153 {
154  CloudType::setParcelThermoProperties(parcel, lagrangianDt);
155 
156  const liquidMixtureProperties& liqMix = this->composition().liquids();
157 
158  const scalarField& Y(parcel.Y());
159  scalarField X(liqMix.X(Y));
160  const scalar pc = this->p()[parcel.cell()];
161 
162  // override rho and Cp from constantProperties
163  parcel.Cp() = liqMix.Cp(pc, parcel.T(), X);
164  parcel.rho() = liqMix.rho(pc, parcel.T(), X);
165  parcel.sigma() = liqMix.sigma(pc, parcel.T(), X);
166  parcel.mu() = liqMix.mu(pc, parcel.T(), X);
167 }
168 
169 
170 template<class CloudType>
172 (
173  parcelType& parcel,
174  const scalar lagrangianDt,
175  const bool fullyDescribed
176 )
177 {
178  CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
179 
180  // store the injection position and initial drop size
181  parcel.position0() = parcel.position();
182  parcel.d0() = parcel.d();
183 
184  parcel.y() = breakup().y0();
185  parcel.yDot() = breakup().yDot0();
186 
187  parcel.liquidCore() = atomisation().initLiquidCore();
188 }
189 
190 
191 template<class CloudType>
193 {
194  cloudCopyPtr_.reset
195  (
196  static_cast<SprayCloud<CloudType>*>
197  (
198  clone(this->name() + "Copy").ptr()
199  )
200  );
201 }
202 
203 
204 template<class CloudType>
206 {
207  cloudReset(cloudCopyPtr_());
208  cloudCopyPtr_.clear();
209 }
210 
211 
212 template<class CloudType>
214 {
215  if (this->solution().canEvolve())
216  {
217  typename parcelType::trackingData td(*this);
218 
219  this->solve(*this, td);
220  }
221 }
222 
223 
224 template<class CloudType>
226 {
227  CloudType::info();
228  scalar d32 = 1.0e+6*this->Dij(3, 2);
229  scalar d10 = 1.0e+6*this->Dij(1, 0);
230  scalar dMax = 1.0e+6*this->Dmax();
231  scalar pen = this->penetration(0.95);
232 
233  Info << " D10, D32, Dmax (mu) = " << d10 << ", " << d32
234  << ", " << dMax << nl
235  << " Liquid penetration 95% mass (m) = " << pen << endl;
236 }
237 
238 
239 // ************************************************************************* //
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:192
basicSpecieMixture & composition
void restoreState()
Reset the current cloud to the previously stored state.
Definition: SprayCloud.C:205
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void setModels()
Set cloud sub-models.
Definition: SprayCloud.C:34
Templated atomisation model class.
Definition: SprayCloud.H:52
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition: SprayCloud.C:149
SprayCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const dimensionedVector &g, const fluidThermo &carrierThermo, bool readFields=true)
Construct given carrier gas fields.
Definition: SprayCloud.C:73
scalar mu(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture viscosity [Pa s].
virtual ~SprayCloud()
Destructor.
Definition: SprayCloud.C:141
scalar rho(const scalar p, const scalar T, const scalarField &X) const
Calculate the mixture density [kg/m^3].
T clone(const T &t)
Definition: List.H:54
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition: SprayCloud.C:172
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
scalar sigma(const scalar p, const scalar T, const scalarField &X) const
Estimate mixture surface tension [N/m].
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:221
Templated base class for spray cloud.
Definition: SprayCloud.H:69
autoPtr< BreakupModel< SprayCloud< CloudType > > > breakupModel_
Break-up model.
Definition: SprayCloud.H:110
static const char nl
Definition: Ostream.H:260
void evolve()
Evolve the spray (inject, move)
Definition: SprayCloud.C:213
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void info()
Print cloud information.
Definition: SprayCloud.C:225
rhoEqn solve()
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:58
messageStream Info
Templated break-up model class.
Definition: SprayCloud.H:55
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:48
volScalarField & p
autoPtr< AtomisationModel< SprayCloud< CloudType > > > atomisationModel_
Atomisation model.
Definition: SprayCloud.H:107