Loading...
Searching...
No Matches
SprayParcelIO.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 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "SprayParcel.H"
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34template<class ParcelType>
37
38
39template<class ParcelType>
41(
42 sizeof(SprayParcel<ParcelType>) - sizeof(ParcelType)
43);
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
48template<class ParcelType>
50(
51 const polyMesh& mesh,
52 Istream& is,
53 bool readFields,
54 bool newFormat
55)
56:
57 ParcelType(mesh, is, readFields, newFormat),
58 d0_(0.0),
60 sigma_(0.0),
61 mu_(0.0),
62 liquidCore_(0.0),
63 KHindex_(0.0),
64 y_(0.0),
65 yDot_(0.0),
66 tc_(0.0),
67 ms_(0.0),
68 injector_(1.0),
69 tMom_(GREAT),
70 user_(0.0)
71{
72 if (readFields)
73 {
74 if (is.format() == IOstreamOption::ASCII)
75 {
76 is >> d0_
77 >> position0_
78 >> sigma_
79 >> mu_
81 >> KHindex_
82 >> y_
83 >> yDot_
84 >> tc_
85 >> ms_
86 >> injector_
87 >> tMom_
88 >> user_;
89 }
90 else
91 {
92 // No non-native streaming
94
95 is.read(reinterpret_cast<char*>(&d0_), sizeofFields);
96 }
97 }
100}
101
102
103template<class ParcelType>
104template<class CloudType>
107 ParcelType::readFields(c);
108}
109
110
111template<class ParcelType>
112template<class CloudType, class CompositionType>
114(
115 CloudType& c,
116 const CompositionType& compModel
117)
118{
119 const bool readOnProc = c.size();
120
121 ParcelType::readFields(c, compModel);
122
123 IOField<scalar> d0(c.newIOobject("d0", IOobject::MUST_READ), readOnProc);
124 c.checkFieldIOobject(c, d0);
125
126 IOField<vector> position0
127 (
128 c.newIOobject("position0", IOobject::MUST_READ),
129 readOnProc
130 );
131 c.checkFieldIOobject(c, position0);
132
134 (
135 c.newIOobject("sigma", IOobject::MUST_READ),
136 readOnProc
137 );
138 c.checkFieldIOobject(c, sigma);
139
140 IOField<scalar> mu(c.newIOobject("mu", IOobject::MUST_READ), readOnProc);
141 c.checkFieldIOobject(c, mu);
142
143 IOField<scalar> liquidCore
144 (
145 c.newIOobject("liquidCore", IOobject::MUST_READ),
146 readOnProc
147 );
148 c.checkFieldIOobject(c, liquidCore);
149
150 IOField<scalar> KHindex
151 (
152 c.newIOobject("KHindex", IOobject::MUST_READ),
153 readOnProc
154 );
155 c.checkFieldIOobject(c, KHindex);
156
158 (
159 c.newIOobject("y", IOobject::MUST_READ),
160 readOnProc
161 );
162 c.checkFieldIOobject(c, y);
163
164 IOField<scalar> yDot
165 (
166 c.newIOobject("yDot", IOobject::MUST_READ),
167 readOnProc
168 );
169 c.checkFieldIOobject(c, yDot);
170
172 (
173 c.newIOobject("tc", IOobject::MUST_READ),
174 readOnProc
175 );
176 c.checkFieldIOobject(c, tc);
177
179 (
180 c.newIOobject("ms", IOobject::MUST_READ),
181 readOnProc
182 );
183 c.checkFieldIOobject(c, ms);
184
185 IOField<scalar> injector
186 (
187 c.newIOobject("injector", IOobject::MUST_READ),
188 readOnProc
189 );
190 c.checkFieldIOobject(c, injector);
191
192 IOField<scalar> tMom
193 (
194 c.newIOobject("tMom", IOobject::MUST_READ),
195 readOnProc
196 );
197 c.checkFieldIOobject(c, tMom);
198
199 IOField<scalar> user
200 (
201 c.newIOobject("user", IOobject::MUST_READ),
202 readOnProc
203 );
204 c.checkFieldIOobject(c, user);
205
206 label i = 0;
207 for (SprayParcel<ParcelType>& p : c)
208 {
209 p.d0_ = d0[i];
210 p.position0_ = position0[i];
211 p.sigma_ = sigma[i];
212 p.mu_ = mu[i];
213 p.liquidCore_ = liquidCore[i];
214 p.KHindex_ = KHindex[i];
215 p.y_ = y[i];
216 p.yDot_ = yDot[i];
217 p.tc_ = tc[i];
218 p.ms_ = ms[i];
219 p.injector_ = injector[i];
220 p.tMom_ = tMom[i];
221 p.user_ = user[i];
222
223 ++i;
224 }
225}
226
227
228template<class ParcelType>
229template<class CloudType>
232 ParcelType::writeFields(c);
233}
234
235
236template<class ParcelType>
237template<class CloudType, class CompositionType>
239(
240 const CloudType& c,
241 const CompositionType& compModel
242)
243{
244 ParcelType::writeFields(c, compModel);
245
246 const label np = c.size();
247 const bool writeOnProc = c.size();
248
249 IOField<scalar> d0(c.newIOobject("d0", IOobject::NO_READ), np);
250 IOField<vector> position0
251 (
252 c.newIOobject("position0", IOobject::NO_READ),
253 np
254 );
255 IOField<scalar> sigma(c.newIOobject("sigma", IOobject::NO_READ), np);
256 IOField<scalar> mu(c.newIOobject("mu", IOobject::NO_READ), np);
257 IOField<scalar> liquidCore
258 (
259 c.newIOobject("liquidCore", IOobject::NO_READ),
260 np
261 );
262 IOField<scalar> KHindex(c.newIOobject("KHindex", IOobject::NO_READ), np);
263 IOField<scalar> y(c.newIOobject("y", IOobject::NO_READ), np);
264 IOField<scalar> yDot(c.newIOobject("yDot", IOobject::NO_READ), np);
265 IOField<scalar> tc(c.newIOobject("tc", IOobject::NO_READ), np);
266 IOField<scalar> ms(c.newIOobject("ms", IOobject::NO_READ), np);
267 IOField<scalar> injector
268 (
269 c.newIOobject("injector", IOobject::NO_READ),
270 np
271 );
272 IOField<scalar> tMom(c.newIOobject("tMom", IOobject::NO_READ), np);
273 IOField<scalar> user(c.newIOobject("user", IOobject::NO_READ), np);
274
275 label i = 0;
276 for (const SprayParcel<ParcelType>& p : c)
277 {
278 d0[i] = p.d0_;
279 position0[i] = p.position0_;
280 sigma[i] = p.sigma_;
281 mu[i] = p.mu_;
282 liquidCore[i] = p.liquidCore_;
283 KHindex[i] = p.KHindex_;
284 y[i] = p.y_;
285 yDot[i] = p.yDot_;
286 tc[i] = p.tc_;
287 ms[i] = p.ms_;
288 injector[i] = p.injector_;
289 tMom[i] = p.tMom_;
290 user[i] = p.user_;
291 ++i;
292 }
293
294 d0.write(writeOnProc);
295 position0.write(writeOnProc);
296 sigma.write(writeOnProc);
297 mu.write(writeOnProc);
298 liquidCore.write(writeOnProc);
299 KHindex.write(writeOnProc);
300 y.write(writeOnProc);
301 yDot.write(writeOnProc);
302 tc.write(writeOnProc);
303 ms.write(writeOnProc);
304 injector.write(writeOnProc);
305 tMom.write(writeOnProc);
306 user.write(writeOnProc);
307}
308
309
310template<class ParcelType>
312(
313 Ostream& os,
314 const wordRes& filters,
315 const word& delim,
316 const bool namesOnly
317) const
318{
319 ParcelType::writeProperties(os, filters, delim, namesOnly);
320
321 #undef writeProp
322 #define writeProp(Name, Value) \
323 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
324
325 writeProp("d0", d0_);
326 writeProp("position0", position0_);
327 writeProp("sigma", sigma_);
328 writeProp("mu", mu_);
329 writeProp("liquidCore", liquidCore_);
330 writeProp("KHindex", KHindex_);
331 writeProp("y", y_);
332 writeProp("yDot", yDot_);
333 writeProp("tc", tc_);
334 writeProp("ms", ms_);
335 writeProp("injector", injector_);
336 writeProp("tMom", tMom_);
337 writeProp("user", user_);
339 #undef writeProp
340}
341
342
343template<class ParcelType>
344template<class CloudType>
346(
347 CloudType& c,
348 const objectRegistry& obr
349)
351 ParcelType::readObjects(c, obr);
352}
353
354
355template<class ParcelType>
356template<class CloudType>
358(
359 const CloudType& c,
360 objectRegistry& obr
361)
363 ParcelType::writeObjects(c, obr);
364}
365
366
367template<class ParcelType>
368template<class CloudType, class CompositionType>
370(
371 CloudType& c,
372 const CompositionType& compModel,
373 const objectRegistry& obr
374)
375{
376 ParcelType::readObjects(c, compModel, obr);
377
378 if (!c.size()) return;
379
380 const auto& d0 = cloud::lookupIOField<scalar>("d0", obr);
381 const auto& position0 = cloud::lookupIOField<vector>("position0", obr);
382 const auto& sigma = cloud::lookupIOField<scalar>("sigma", obr);
383 const auto& mu = cloud::lookupIOField<scalar>("mu", obr);
384 const auto& liquidCore = cloud::lookupIOField<scalar>("liquidCore", obr);
385 const auto& KHindex = cloud::lookupIOField<scalar>("KHindex", obr);
386 const auto& y = cloud::lookupIOField<scalar>("y", obr);
387 const auto& yDot = cloud::lookupIOField<scalar>("yDot", obr);
388 const auto& tc = cloud::lookupIOField<scalar>("tc", obr);
389 const auto& ms = cloud::lookupIOField<scalar>("ms", obr);
390 const auto& injector = cloud::lookupIOField<scalar>("injector", obr);
391 const auto& tMom = cloud::lookupIOField<scalar>("tMom", obr);
392 const auto& user = cloud::lookupIOField<scalar>("user", obr);
393
394 label i = 0;
395 for (SprayParcel<ParcelType>& p : c)
396 {
397 p.d0_ = d0[i];
398 p.position0_ = position0[i];
399 p.sigma_ = sigma[i];
400 p.mu_ = mu[i];
401 p.liquidCore_ = liquidCore[i];
402 p.KHindex_ = KHindex[i];
403 p.y_ = y[i];
404 p.yDot_ = yDot[i];
405 p.tc_ = tc[i];
406 p.ms_ = ms[i];
407 p.injector_ = injector[i];
408 p.tMom_ = tMom[i];
409 p.user_ = user[i];
410
411 ++i;
412 }
413}
414
415
416template<class ParcelType>
417template<class CloudType, class CompositionType>
419(
420 const CloudType& c,
421 const CompositionType& compModel,
422 objectRegistry& obr
423)
424{
425 ParcelType::writeObjects(c, compModel, obr);
426
427 const label np = c.size();
428
429 auto& d0 = cloud::createIOField<scalar>("d0", np, obr);
430 auto& position0 = cloud::createIOField<vector>("position0", np, obr);
431 auto& sigma = cloud::createIOField<scalar>("sigma", np, obr);
432 auto& mu = cloud::createIOField<scalar>("mu", np, obr);
433 auto& liquidCore = cloud::createIOField<scalar>("liquidCore", np, obr);
434 auto& KHindex = cloud::createIOField<scalar>("KHindex", np, obr);
435 auto& y = cloud::createIOField<scalar>("y", np, obr);
436 auto& yDot= cloud::createIOField<scalar>("yDot", np, obr);
437 auto& tc = cloud::createIOField<scalar>("tc", np, obr);
438 auto& ms = cloud::createIOField<scalar>("ms", np, obr);
439 auto& injector = cloud::createIOField<scalar>("injector", np, obr);
440 auto& tMom = cloud::createIOField<scalar>("tMom", np, obr);
441 auto& user = cloud::createIOField<scalar>("user", np, obr);
442
443 label i = 0;
444 for (const SprayParcel<ParcelType>& p : c)
445 {
446 d0[i] = p.d0_;
447 position0[i] = p.position0_;
448 sigma[i] = p.sigma_;
449 mu[i] = p.mu_;
450 liquidCore[i] = p.liquidCore_;
451 KHindex[i] = p.KHindex_;
452 y[i] = p.y_;
453 yDot[i] = p.yDot_;
454 tc[i] = p.tc_;
455 ms[i] = p.ms_;
456 injector[i] = p.injector_;
457 tMom[i] = p.tMom_;
458 user[i] = p.user_;
459
460 ++i;
461 }
462}
463
464
465// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
466
467template<class ParcelType>
468Foam::Ostream& Foam::operator<<
469(
470 Ostream& os,
472)
473{
475 {
477 << token::SPACE << p.d0()
478 << token::SPACE << p.position0()
479 << token::SPACE << p.sigma()
480 << token::SPACE << p.mu()
481 << token::SPACE << p.liquidCore()
482 << token::SPACE << p.KHindex()
483 << token::SPACE << p.y()
484 << token::SPACE << p.yDot()
485 << token::SPACE << p.tc()
486 << token::SPACE << p.ms()
487 << token::SPACE << p.injector()
488 << token::SPACE << p.tMom()
489 << token::SPACE << p.user();
490 }
491 else
492 {
494 os.write
495 (
496 reinterpret_cast<const char*>(&p.d0_),
498 );
499 }
500
502 return os;
503}
504
505
506// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
scalar y
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
bool fatalCheckNativeSizes(const char *operation) const
Assert that the label/scalar byte-size associated with the stream are the native label/scalar sizes.
Definition IOstream.C:67
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
virtual Istream & read(token &)=0
Return next token from stream.
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Reacting spray parcel, with added functionality for atomization and breakup.
Definition SprayParcel.H:57
scalar injector() const
Return const access to injector id.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
scalar liquidCore() const
Return const access to liquid core.
scalar d0_
Initial droplet diameter.
scalar tMom() const
Return const access to momentum relaxation time.
vector position0_
Injection position.
scalar mu() const
Return const access to the liquid dynamic viscosity.
scalar mu_
Liquid dynamic viscosity [Pa.s].
static const std::size_t sizeofFields
Size in bytes of the fields.
scalar sigma_
Liquid surface tension [N/m].
scalar user_
Passive scalar (extra variable to be defined by user).
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
scalar yDot() const
Return const access to rate of change of spherical deviation.
scalar d0() const
Return const access to initial droplet diameter.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.).
scalar KHindex() const
Return const access to Kelvin-Helmholtz breakup index.
scalar tc() const
Return const access to atomization characteristic time.
scalar y_
Spherical deviation.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
scalar tc_
Characteristic time (used in atomization and/or breakup model).
SprayParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
scalar injector_
Injected from injector (needed e.g. for calculating distance.
scalar KHindex_
Index for KH Breakup.
scalar ms() const
Return const access to stripped parcel mass.
scalar y() const
Return const access to spherical deviation.
scalar sigma() const
Return const access to the liquid surface tension.
scalar ms_
Stripped parcel mass due to breakup.
scalar yDot_
Rate of change of spherical deviation.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
scalar user() const
Return const access to passive user scalar.
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet ).
const vector & position0() const
Return const access to initial droplet position.
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
Definition cloud.H:202
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for handling character strings derived from std::string.
Definition string.H:76
@ SPACE
Space [isspace].
Definition token.H:144
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
DSMCCloud< dsmcParcel > CloudType
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)