SprayParcelI.H
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-2020 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
28 template<class ParcelType>
30 :
31  ParcelType::constantProperties(),
32  sigma0_(this->dict_, 0.0),
33  mu0_(this->dict_, 0.0)
34 {}
35 
36 
37 template<class ParcelType>
39 (
40  const constantProperties& cp
41 )
42 :
43  ParcelType::constantProperties(cp),
44  sigma0_(cp.sigma0_),
45  mu0_(cp.mu0_)
46 {}
47 
48 
49 template<class ParcelType>
51 (
52  const dictionary& parentDict
53 )
54 :
55  ParcelType::constantProperties(parentDict),
56  sigma0_(this->dict_, "sigma0"),
57  mu0_(this->dict_, "mu0")
58 {}
59 
60 
61 template<class ParcelType>
63 (
64  const label parcelTypeId,
65  const scalar rhoMin,
66  const scalar rho0,
67  const scalar minParcelMass,
68  const scalar youngsModulus,
69  const scalar poissonsRatio,
70  const scalar T0,
71  const scalar TMin,
72  const scalar TMax,
73  const scalar Cp0,
74  const scalar epsilon0,
75  const scalar f0,
76  const scalar Pr,
77  const scalar pMin,
78  const Switch& constantVolume,
79  const scalar sigma0,
80  const scalar mu0
81 )
82 :
83  ParcelType::constantProperties
84  (
85  parcelTypeId,
86  rhoMin,
87  rho0,
88  minParcelMass,
89  youngsModulus,
90  poissonsRatio,
91  T0,
92  TMin,
93  TMax,
94  Cp0,
95  epsilon0,
96  f0,
97  Pr,
98  pMin,
99  constantVolume
100  ),
101  sigma0_(this->dict_, sigma0),
102  mu0_(this->dict_, mu0)
103 {}
104 
105 
106 template<class ParcelType>
108 (
109  const polyMesh& mesh,
110  const barycentric& coordinates,
111  const label celli,
112  const label tetFacei,
113  const label tetPti
114 )
115 :
116  ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
117  d0_(this->d()),
118  position0_(this->position()),
119  sigma_(0.0),
120  mu_(0.0),
121  liquidCore_(0.0),
122  KHindex_(0.0),
123  y_(0.0),
124  yDot_(0.0),
125  tc_(0.0),
126  ms_(0.0),
127  injector_(1.0),
128  tMom_(great),
129  user_(0.0)
130 {}
131 
132 
133 template<class ParcelType>
135 (
136  const polyMesh& mesh,
137  const vector& position,
138  const label celli
139 )
140 :
141  ParcelType(mesh, position, celli),
142  d0_(this->d()),
143  position0_(this->position()),
144  sigma_(0.0),
145  mu_(0.0),
146  liquidCore_(0.0),
147  KHindex_(0.0),
148  y_(0.0),
149  yDot_(0.0),
150  tc_(0.0),
151  ms_(0.0),
152  injector_(1.0),
153  tMom_(great),
154  user_(0.0)
155 {}
156 
157 
158 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
159 
160 template<class ParcelType>
161 inline Foam::scalar
163 {
164  return sigma0_.value();
165 }
166 
167 
168 template<class ParcelType>
169 inline Foam::scalar
171 {
172  return mu0_.value();
173 }
174 
175 
176 // * * * * * * * * * * SprayParcel Member Functions * * * * * * * * * * * * //
177 
178 template<class ParcelType>
179 inline Foam::scalar Foam::SprayParcel<ParcelType>::d0() const
180 {
181  return d0_;
182 }
183 
184 
185 template<class ParcelType>
187 {
188  return position0_;
189 }
190 
191 
192 template<class ParcelType>
193 inline Foam::scalar Foam::SprayParcel<ParcelType>::sigma() const
194 {
195  return sigma_;
196 }
197 
198 
199 template<class ParcelType>
200 inline Foam::scalar Foam::SprayParcel<ParcelType>::mu() const
201 {
202  return mu_;
203 }
204 
205 
206 template<class ParcelType>
207 inline Foam::scalar Foam::SprayParcel<ParcelType>::liquidCore() const
208 {
209  return liquidCore_;
210 }
211 
212 
213 template<class ParcelType>
214 inline Foam::scalar Foam::SprayParcel<ParcelType>::KHindex() const
215 {
216  return KHindex_;
217 }
218 
219 
220 template<class ParcelType>
221 inline Foam::scalar Foam::SprayParcel<ParcelType>::y() const
222 {
223  return y_;
224 }
225 
226 
227 template<class ParcelType>
228 inline Foam::scalar Foam::SprayParcel<ParcelType>::yDot() const
229 {
230  return yDot_;
231 }
232 
233 
234 template<class ParcelType>
235 inline Foam::scalar Foam::SprayParcel<ParcelType>::tc() const
236 {
237  return tc_;
238 }
239 
240 
241 template<class ParcelType>
242 inline Foam::scalar Foam::SprayParcel<ParcelType>::ms() const
243 {
244  return ms_;
245 }
246 
247 
248 template<class ParcelType>
249 inline Foam::scalar Foam::SprayParcel<ParcelType>::injector() const
250 {
251  return injector_;
252 }
253 
254 
255 template<class ParcelType>
256 inline Foam::scalar Foam::SprayParcel<ParcelType>::tMom() const
257 {
258  return tMom_;
259 }
260 
261 
262 template<class ParcelType>
263 inline Foam::scalar Foam::SprayParcel<ParcelType>::user() const
264 {
265  return user_;
266 }
267 
268 
269 template<class ParcelType>
270 inline Foam::scalar& Foam::SprayParcel<ParcelType>::d0()
271 {
272  return d0_;
273 }
274 
275 
276 template<class ParcelType>
278 {
279  return position0_;
280 }
281 
282 
283 template<class ParcelType>
285 {
286  return sigma_;
287 }
288 
289 
290 template<class ParcelType>
291 inline Foam::scalar& Foam::SprayParcel<ParcelType>::mu()
292 {
293  return mu_;
294 }
295 
296 
297 template<class ParcelType>
299 {
300  return liquidCore_;
301 }
302 
303 
304 template<class ParcelType>
306 {
307  return KHindex_;
308 }
309 
310 
311 template<class ParcelType>
312 inline Foam::scalar& Foam::SprayParcel<ParcelType>::y()
313 {
314  return y_;
315 }
316 
317 
318 template<class ParcelType>
320 {
321  return yDot_;
322 }
323 
324 
325 template<class ParcelType>
326 inline Foam::scalar& Foam::SprayParcel<ParcelType>::tc()
327 {
328  return tc_;
329 }
330 
331 
332 template<class ParcelType>
333 inline Foam::scalar& Foam::SprayParcel<ParcelType>::ms()
334 {
335  return ms_;
336 }
337 
338 
339 template<class ParcelType>
341 {
342  return injector_;
343 }
344 
345 
346 template<class ParcelType>
348 {
349  return tMom_;
350 }
351 
352 
353 template<class ParcelType>
355 {
356  return user_;
357 }
358 
359 
360 // ************************************************************************* //
scalar mu0() const
Return const access to the initial dynamic viscosity.
Definition: SprayParcelI.H:170
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:183
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
Definition: SprayParcelI.H:214
scalar tMom() const
Return const access to momentum relaxation time.
Definition: SprayParcelI.H:256
scalar sigma0() const
Return const access to the initial surface tension.
Definition: SprayParcelI.H:162
scalar liquidCore() const
Return const access to liquid core.
Definition: SprayParcelI.H:207
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
scalar ms() const
Return const access to stripped parcel mass.
Definition: SprayParcelI.H:242
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:161
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: SprayParcelI.H:108
scalar y() const
Return const access to spherical deviation.
Definition: SprayParcelI.H:221
const Type & value() const
Return the value.
scalar d0() const
Return const access to initial droplet diameter.
Definition: SprayParcelI.H:179
scalar mu() const
Return const access to the liquid dynamic viscosity.
Definition: SprayParcelI.H:200
scalar tc() const
Return const access to atomisation characteristic time.
Definition: SprayParcelI.H:235
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:180
Class to hold reacting particle constant properties.
Definition: SprayParcel.H:79
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:186
scalar yDot() const
Return const access to rate of change of spherical deviation.
Definition: SprayParcelI.H:228
scalar user() const
Return const access to passive user scalar.
Definition: SprayParcelI.H:263
constantProperties()
Null constructor.
Definition: SprayParcelI.H:29
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:155
vector position0_
Injection position.
Definition: SprayParcel.H:152
scalar injector() const
Return const access to injector id.
Definition: SprayParcelI.H:249
const vector & position0() const
Return const access to initial droplet position.
Definition: SprayParcelI.H:186
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:170
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:176
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:149
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
scalar tc_
Characteristic time (used in atomisation and/or breakup model)
Definition: SprayParcel.H:173
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:158
scalar y_
Spherical deviation.
Definition: SprayParcel.H:167
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:164
scalar sigma() const
Return const access to the liquid surface tension.
Definition: SprayParcelI.H:193