Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "SprayCloud.H"
29#include "AtomizationModel.H"
30#include "BreakupModel.H"
31
32// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33
34template<class CloudType>
36{
38 (
39 AtomizationModel<SprayCloud<CloudType>>::New
40 (
41 this->subModelProperties(),
42 *this
43 ).ptr()
44 );
45
46 breakupModel_.reset
47 (
48 BreakupModel<SprayCloud<CloudType>>::New
49 (
50 this->subModelProperties(),
51 *this
52 ).ptr()
53 );
54}
55
56
57template<class CloudType>
59(
61)
62{
64
65 atomizationModel_.reset(c.atomizationModel_.ptr());
66 breakupModel_.reset(c.breakupModel_.ptr());
67}
68
69
70// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71
72template<class CloudType>
73Foam::SprayCloud<CloudType>::SprayCloud
74(
75 const word& cloudName,
76 const volScalarField& rho,
77 const volVectorField& U,
78 const dimensionedVector& g,
79 const SLGThermo& thermo,
80 bool readFields
81)
82:
83 CloudType(cloudName, rho, U, g, thermo, false),
84 sprayCloud(),
85 cloudCopyPtr_(nullptr),
86 averageParcelMass_(0.0),
87 atomizationModel_(nullptr),
88 breakupModel_(nullptr)
89{
90 if (this->solution().active())
91 {
92 setModels();
93
94 averageParcelMass_ = this->injectors().averageParcelMass();
95
96 if (readFields)
97 {
98 parcelType::readFields(*this, this->composition());
99 this->deleteLostParticles();
100 }
101
102 Info << "Average parcel mass: " << averageParcelMass_ << endl;
103 }
104
105 if (this->solution().resetSourcesOnStartup())
107 CloudType::resetSourceTerms();
108 }
109}
110
111
112template<class CloudType>
113Foam::SprayCloud<CloudType>::SprayCloud
114(
115 SprayCloud<CloudType>& c,
116 const word& name
117)
118:
119 CloudType(c, name),
120 sprayCloud(),
121 cloudCopyPtr_(nullptr),
122 averageParcelMass_(c.averageParcelMass_),
125{}
126
127
128template<class CloudType>
129Foam::SprayCloud<CloudType>::SprayCloud
130(
131 const fvMesh& mesh,
132 const word& name,
133 const SprayCloud<CloudType>& c
134)
135:
136 CloudType(mesh, name, c),
137 sprayCloud(),
138 cloudCopyPtr_(nullptr),
139 averageParcelMass_(0.0),
140 atomizationModel_(nullptr),
142{}
143
144
145// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
146
147template<class CloudType>
149{}
150
151
152// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153
154template<class CloudType>
156(
157 parcelType& parcel,
158 const scalar lagrangianDt
159)
160{
161 CloudType::setParcelThermoProperties(parcel, lagrangianDt);
162
163 const liquidMixtureProperties& liqMix = this->composition().liquids();
164
165 const scalarField& Y(parcel.Y());
166 scalarField X(liqMix.X(Y));
167 const scalar pc = this->p()[parcel.cell()];
168
169 // override rho and Cp from constantProperties
170 parcel.Cp() = liqMix.Cp(pc, parcel.T(), X);
171 parcel.rho() = liqMix.rho(pc, parcel.T(), X);
172 parcel.sigma() = liqMix.sigma(pc, parcel.T(), X);
173 parcel.mu() = liqMix.mu(pc, parcel.T(), X);
174}
175
176
177template<class CloudType>
179(
180 parcelType& parcel,
181 const scalar lagrangianDt,
182 const bool fullyDescribed
183)
184{
185 CloudType::checkParcelProperties(parcel, lagrangianDt, fullyDescribed);
186
187 // store the injection position and initial drop size
188 parcel.position0() = parcel.position();
189 parcel.d0() = parcel.d();
190
191 parcel.y() = breakup().y0();
192 parcel.yDot() = breakup().yDot0();
193
194 parcel.liquidCore() = atomization().initLiquidCore();
195}
196
197
198template<class CloudType>
200{
201 cloudCopyPtr_.reset
202 (
203 static_cast<SprayCloud<CloudType>*>
204 (
205 clone(this->name() + "Copy").ptr()
206 )
207 );
208}
209
210
211template<class CloudType>
214 cloudReset(cloudCopyPtr_());
215 cloudCopyPtr_.clear();
216}
217
218
219template<class CloudType>
221{
222 if (this->solution().canEvolve())
223 {
224 typename parcelType::trackingData td(*this);
226 this->solve(*this, td);
227 }
228}
229
230
231template<class CloudType>
233{
235 scalar d32 = 1.0e+6*this->Dij(3, 2);
236 scalar d10 = 1.0e+6*this->Dij(1, 0);
237 scalar dMax = 1.0e+6*this->Dmax();
238 scalar pen = this->penetration(0.95);
239
240 Info << " D10, D32, Dmax (mu) = " << d10 << ", " << d32
241 << ", " << dMax << nl
242 << " Liquid penetration 95% mass (m) = " << pen << endl;
243}
244
245
246// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
Templated atomization model class.
Templated break-up model class.
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition Cloud.C:108
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition Cloud.C:125
const word & cloudName() const
const fvMesh & mesh() const
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition SLGThermo.H:63
Templated base class for spray cloud.
Definition SprayCloud.H:57
autoPtr< BreakupModel< SprayCloud< CloudType > > > breakupModel_
Break-up model.
Definition SprayCloud.H:121
virtual ~SprayCloud()
Destructor.
Definition SprayCloud.C:141
void setModels()
Set cloud sub-models.
Definition SprayCloud.C:28
const BreakupModel< SprayCloud< CloudType > > & breakup() const
Return const-access to the breakup model.
Definition SprayCloudI.H:50
void storeState()
Store the current cloud state.
Definition SprayCloud.C:192
autoPtr< AtomizationModel< SprayCloud< CloudType > > > atomizationModel_
Atomization model.
Definition SprayCloud.H:116
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
Definition SprayCloud.C:149
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
Definition SprayCloud.C:172
void cloudReset(SprayCloud< CloudType > &c)
Reset state of cloud.
Definition SprayCloud.C:52
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
Definition SprayCloud.H:70
const AtomizationModel< SprayCloud< CloudType > > & atomization() const
Return const-access to the atomization model.
Definition SprayCloudI.H:34
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
Definition SprayCloud.H:178
void evolve()
Evolve the spray (inject, move).
Definition SprayCloud.C:213
void info()
Print cloud information.
Definition SprayCloud.C:225
void restoreState()
Reset the current cloud to the previously stored state.
Definition SprayCloud.C:205
scalar penetration(const scalar fraction) const
Penetration for fraction [0-1] of the current total mass.
Definition SprayCloudI.H:73
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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].
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].
Selector class for relaxation factors, solver type and solution.
Definition solution.H:95
Virtual abstract base class for templated SprayCloud.
Definition sprayCloud.H:47
sprayCloud()=default
Null constructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
basicSpecieMixture & composition
volScalarField & p
PtrList< volScalarField > & Y
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
CEqn solve()