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