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 -------------------------------------------------------------------------------
11 License
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 
34 template<class ParcelType>
37 
38 
39 template<class ParcelType>
41 (
42  sizeof(SprayParcel<ParcelType>) - sizeof(ParcelType)
43 );
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<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),
59  position0_(Zero),
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_
80  >> liquidCore_
81  >> KHindex_
82  >> y_
83  >> yDot_
84  >> tc_
85  >> ms_
86  >> injector_
87  >> tMom_
88  >> user_;
89  }
90  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
91  {
92  // Non-native label or scalar size
93 
94  is.beginRawRead();
95 
96  readRawScalar(is, &d0_);
97  readRawScalar(is, position0_.data(), vector::nComponents);
98  readRawScalar(is, &sigma_);
99  readRawScalar(is, &mu_);
100  readRawScalar(is, &liquidCore_);
101  readRawScalar(is, &KHindex_);
102  readRawScalar(is, &y_);
103  readRawScalar(is, &yDot_);
104  readRawScalar(is, &tc_);
105  readRawScalar(is, &ms_);
106  readRawScalar(is, &injector_);
107  readRawScalar(is, &tMom_);
108  readRawScalar(is, &user_);
109 
110  is.endRawRead();
111  }
112  else
113  {
114  is.read(reinterpret_cast<char*>(&d0_), sizeofFields);
115  }
116  }
118  is.check(FUNCTION_NAME);
119 }
120 
121 
122 template<class ParcelType>
123 template<class CloudType>
125 {
127 }
128 
129 
130 template<class ParcelType>
131 template<class CloudType, class CompositionType>
133 (
134  CloudType& c,
135  const CompositionType& compModel
136 )
137 {
138  const bool readOnProc = c.size();
139 
140  ParcelType::readFields(c, compModel);
141 
142  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::MUST_READ), readOnProc);
143  c.checkFieldIOobject(c, d0);
144 
145  IOField<vector> position0
146  (
147  c.fieldIOobject("position0", IOobject::MUST_READ),
148  readOnProc
149  );
150  c.checkFieldIOobject(c, position0);
151 
153  (
154  c.fieldIOobject("sigma", IOobject::MUST_READ),
155  readOnProc
156  );
157  c.checkFieldIOobject(c, sigma);
158 
159  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::MUST_READ), readOnProc);
160  c.checkFieldIOobject(c, mu);
161 
162  IOField<scalar> liquidCore
163  (
164  c.fieldIOobject("liquidCore", IOobject::MUST_READ),
165  readOnProc
166  );
167  c.checkFieldIOobject(c, liquidCore);
168 
169  IOField<scalar> KHindex
170  (
171  c.fieldIOobject("KHindex", IOobject::MUST_READ),
172  readOnProc
173  );
174  c.checkFieldIOobject(c, KHindex);
175 
177  (
178  c.fieldIOobject("y", IOobject::MUST_READ),
179  readOnProc
180  );
181  c.checkFieldIOobject(c, y);
182 
183  IOField<scalar> yDot
184  (
185  c.fieldIOobject("yDot", IOobject::MUST_READ),
186  readOnProc
187  );
188  c.checkFieldIOobject(c, yDot);
189 
190  IOField<scalar> tc
191  (
192  c.fieldIOobject("tc", IOobject::MUST_READ),
193  readOnProc
194  );
195  c.checkFieldIOobject(c, tc);
196 
197  IOField<scalar> ms
198  (
199  c.fieldIOobject("ms", IOobject::MUST_READ),
200  readOnProc
201  );
202  c.checkFieldIOobject(c, ms);
203 
204  IOField<scalar> injector
205  (
206  c.fieldIOobject("injector", IOobject::MUST_READ),
207  readOnProc
208  );
209  c.checkFieldIOobject(c, injector);
210 
211  IOField<scalar> tMom
212  (
213  c.fieldIOobject("tMom", IOobject::MUST_READ),
214  readOnProc
215  );
216  c.checkFieldIOobject(c, tMom);
217 
218  IOField<scalar> user
219  (
220  c.fieldIOobject("user", IOobject::MUST_READ),
221  readOnProc
222  );
223  c.checkFieldIOobject(c, user);
224 
225  label i = 0;
226  for (SprayParcel<ParcelType>& p : c)
227  {
228  p.d0_ = d0[i];
229  p.position0_ = position0[i];
230  p.sigma_ = sigma[i];
231  p.mu_ = mu[i];
232  p.liquidCore_ = liquidCore[i];
233  p.KHindex_ = KHindex[i];
234  p.y_ = y[i];
235  p.yDot_ = yDot[i];
236  p.tc_ = tc[i];
237  p.ms_ = ms[i];
238  p.injector_ = injector[i];
239  p.tMom_ = tMom[i];
240  p.user_ = user[i];
241 
242  ++i;
243  }
244 }
245 
246 
247 template<class ParcelType>
248 template<class CloudType>
250 {
252 }
253 
254 
255 template<class ParcelType>
256 template<class CloudType, class CompositionType>
258 (
259  const CloudType& c,
260  const CompositionType& compModel
261 )
262 {
263  ParcelType::writeFields(c, compModel);
264 
265  const label np = c.size();
266  const bool writeOnProc = c.size();
267 
268  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::NO_READ), np);
269  IOField<vector> position0
270  (
271  c.fieldIOobject("position0", IOobject::NO_READ),
272  np
273  );
274  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::NO_READ), np);
275  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::NO_READ), np);
276  IOField<scalar> liquidCore
277  (
278  c.fieldIOobject("liquidCore", IOobject::NO_READ),
279  np
280  );
281  IOField<scalar> KHindex(c.fieldIOobject("KHindex", IOobject::NO_READ), np);
282  IOField<scalar> y(c.fieldIOobject("y", IOobject::NO_READ), np);
283  IOField<scalar> yDot(c.fieldIOobject("yDot", IOobject::NO_READ), np);
284  IOField<scalar> tc(c.fieldIOobject("tc", IOobject::NO_READ), np);
285  IOField<scalar> ms(c.fieldIOobject("ms", IOobject::NO_READ), np);
286  IOField<scalar> injector
287  (
288  c.fieldIOobject("injector", IOobject::NO_READ),
289  np
290  );
291  IOField<scalar> tMom(c.fieldIOobject("tMom", IOobject::NO_READ), np);
292  IOField<scalar> user(c.fieldIOobject("user", IOobject::NO_READ), np);
293 
294  label i = 0;
295  for (const SprayParcel<ParcelType>& p : c)
296  {
297  d0[i] = p.d0_;
298  position0[i] = p.position0_;
299  sigma[i] = p.sigma_;
300  mu[i] = p.mu_;
301  liquidCore[i] = p.liquidCore_;
302  KHindex[i] = p.KHindex_;
303  y[i] = p.y_;
304  yDot[i] = p.yDot_;
305  tc[i] = p.tc_;
306  ms[i] = p.ms_;
307  injector[i] = p.injector_;
308  tMom[i] = p.tMom_;
309  user[i] = p.user_;
310  ++i;
311  }
312 
313  d0.write(writeOnProc);
314  position0.write(writeOnProc);
315  sigma.write(writeOnProc);
316  mu.write(writeOnProc);
317  liquidCore.write(writeOnProc);
318  KHindex.write(writeOnProc);
319  y.write(writeOnProc);
320  yDot.write(writeOnProc);
321  tc.write(writeOnProc);
322  ms.write(writeOnProc);
323  injector.write(writeOnProc);
324  tMom.write(writeOnProc);
325  user.write(writeOnProc);
326 }
327 
328 
329 template<class ParcelType>
331 (
332  Ostream& os,
333  const wordRes& filters,
334  const word& delim,
335  const bool namesOnly
336 ) const
337 {
338  ParcelType::writeProperties(os, filters, delim, namesOnly);
339 
340  #undef writeProp
341  #define writeProp(Name, Value) \
342  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
343 
344  writeProp("d0", d0_);
345  writeProp("position0", position0_);
346  writeProp("sigma", sigma_);
347  writeProp("mu", mu_);
348  writeProp("liquidCore", liquidCore_);
349  writeProp("KHindex", KHindex_);
350  writeProp("y", y_);
351  writeProp("yDot", yDot_);
352  writeProp("tc", tc_);
353  writeProp("ms", ms_);
354  writeProp("injector", injector_);
355  writeProp("tMom", tMom_);
356  writeProp("user", user_);
357 
358  #undef writeProp
359 }
360 
361 
362 template<class ParcelType>
363 template<class CloudType>
365 (
366  CloudType& c,
367  const objectRegistry& obr
368 )
369 {
370  ParcelType::readObjects(c, obr);
371 }
372 
373 
374 template<class ParcelType>
375 template<class CloudType>
377 (
378  const CloudType& c,
379  objectRegistry& obr
380 )
381 {
382  ParcelType::writeObjects(c, obr);
383 }
384 
385 
386 template<class ParcelType>
387 template<class CloudType, class CompositionType>
389 (
390  CloudType& c,
391  const CompositionType& compModel,
392  const objectRegistry& obr
393 )
394 {
395  ParcelType::readObjects(c, compModel, obr);
396 
397  if (!c.size()) return;
398 
399  const auto& d0 = cloud::lookupIOField<scalar>("d0", obr);
400  const auto& position0 = cloud::lookupIOField<vector>("position0", obr);
401  const auto& sigma = cloud::lookupIOField<scalar>("sigma", obr);
402  const auto& mu = cloud::lookupIOField<scalar>("mu", obr);
403  const auto& liquidCore = cloud::lookupIOField<scalar>("liquidCore", obr);
404  const auto& KHindex = cloud::lookupIOField<scalar>("KHindex", obr);
405  const auto& y = cloud::lookupIOField<scalar>("y", obr);
406  const auto& yDot = cloud::lookupIOField<scalar>("yDot", obr);
407  const auto& tc = cloud::lookupIOField<scalar>("tc", obr);
408  const auto& ms = cloud::lookupIOField<scalar>("ms", obr);
409  const auto& injector = cloud::lookupIOField<scalar>("injector", obr);
410  const auto& tMom = cloud::lookupIOField<scalar>("tMom", obr);
411  const auto& user = cloud::lookupIOField<scalar>("user", obr);
412 
413  label i = 0;
414  for (SprayParcel<ParcelType>& p : c)
415  {
416  p.d0_ = d0[i];
417  p.position0_ = position0[i];
418  p.sigma_ = sigma[i];
419  p.mu_ = mu[i];
420  p.liquidCore_ = liquidCore[i];
421  p.KHindex_ = KHindex[i];
422  p.y_ = y[i];
423  p.yDot_ = yDot[i];
424  p.tc_ = tc[i];
425  p.ms_ = ms[i];
426  p.injector_ = injector[i];
427  p.tMom_ = tMom[i];
428  p.user_ = user[i];
429 
430  ++i;
431  }
432 }
433 
434 
435 template<class ParcelType>
436 template<class CloudType, class CompositionType>
438 (
439  const CloudType& c,
440  const CompositionType& compModel,
441  objectRegistry& obr
442 )
443 {
444  ParcelType::writeObjects(c, compModel, obr);
445 
446  const label np = c.size();
447 
448  auto& d0 = cloud::createIOField<scalar>("d0", np, obr);
449  auto& position0 = cloud::createIOField<vector>("position0", np, obr);
450  auto& sigma = cloud::createIOField<scalar>("sigma", np, obr);
451  auto& mu = cloud::createIOField<scalar>("mu", np, obr);
452  auto& liquidCore = cloud::createIOField<scalar>("liquidCore", np, obr);
453  auto& KHindex = cloud::createIOField<scalar>("KHindex", np, obr);
454  auto& y = cloud::createIOField<scalar>("y", np, obr);
455  auto& yDot= cloud::createIOField<scalar>("yDot", np, obr);
456  auto& tc = cloud::createIOField<scalar>("tc", np, obr);
457  auto& ms = cloud::createIOField<scalar>("ms", np, obr);
458  auto& injector = cloud::createIOField<scalar>("injector", np, obr);
459  auto& tMom = cloud::createIOField<scalar>("tMom", np, obr);
460  auto& user = cloud::createIOField<scalar>("user", np, obr);
461 
462  label i = 0;
463  for (const SprayParcel<ParcelType>& p : c)
464  {
465  d0[i] = p.d0_;
466  position0[i] = p.position0_;
467  sigma[i] = p.sigma_;
468  mu[i] = p.mu_;
469  liquidCore[i] = p.liquidCore_;
470  KHindex[i] = p.KHindex_;
471  y[i] = p.y_;
472  yDot[i] = p.yDot_;
473  tc[i] = p.tc_;
474  ms[i] = p.ms_;
475  injector[i] = p.injector_;
476  tMom[i] = p.tMom_;
477  user[i] = p.user_;
478 
479  ++i;
480  }
481 }
482 
483 
484 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
485 
486 template<class ParcelType>
487 Foam::Ostream& Foam::operator<<
488 (
489  Ostream& os,
490  const SprayParcel<ParcelType>& p
491 )
492 {
494  {
495  os << static_cast<const ParcelType&>(p)
496  << token::SPACE << p.d0()
497  << token::SPACE << p.position0()
498  << token::SPACE << p.sigma()
499  << token::SPACE << p.mu()
500  << token::SPACE << p.liquidCore()
501  << token::SPACE << p.KHindex()
502  << token::SPACE << p.y()
503  << token::SPACE << p.yDot()
504  << token::SPACE << p.tc()
505  << token::SPACE << p.ms()
506  << token::SPACE << p.injector()
507  << token::SPACE << p.tMom()
508  << token::SPACE << p.user();
509  }
510  else
511  {
512  os << static_cast<const ParcelType&>(p);
513  os.write
514  (
515  reinterpret_cast<const char*>(&p.d0_),
517  );
518  }
519 
521  return os;
522 }
523 
524 
525 // ************************************************************************* //
DSMCCloud< dsmcParcel > CloudType
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Check if the scalar byte-size associated with the stream is the same as the given type...
Definition: IOstream.H:379
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
scalar tMom_
Momentum relaxation time (needed for calculating parcel acc.)
Definition: SprayParcel.H:210
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
"ascii" (normal default)
virtual bool endRawRead()=0
End of low-level raw binary read.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
scalar liquidCore_
Part of liquid core ( >0.5=liquid, <0.5=droplet )
Definition: SprayParcel.H:173
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:104
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: SprayParcel.H:225
scalar y
virtual Istream & read(token &)=0
Return next token from stream.
dynamicFvMesh & mesh
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
scalar injector_
Injected from injector (needed e.g. for calculating distance.
Definition: SprayParcel.H:205
scalar user_
Passive scalar (extra variable to be defined by user)
Definition: SprayParcel.H:215
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
scalar sigma_
Liquid surface tension [N/m].
Definition: SprayParcel.H:163
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
const dimensionedScalar mu
Atomic mass unit.
vector position0_
Injection position.
Definition: SprayParcel.H:158
#define writeProp(Name, Value)
scalar yDot_
Rate of change of spherical deviation.
Definition: SprayParcel.H:188
virtual bool beginRawRead()=0
Start of low-level raw binary read.
const dimensionedScalar c
Speed of light in a vacuum.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Nothing to be read.
scalar ms_
Stripped parcel mass due to breakup.
Definition: SprayParcel.H:198
scalar d0_
Initial droplet diameter.
Definition: SprayParcel.H:153
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Check if the label byte-size associated with the stream is the same as the given type.
Definition: IOstream.H:368
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar tc_
Characteristic time (used in atomization and/or breakup model)
Definition: SprayParcel.H:193
volScalarField & p
Registry of regIOobjects.
scalar mu_
Liquid dynamic viscosity [Pa.s].
Definition: SprayParcel.H:168
scalar y_
Spherical deviation.
Definition: SprayParcel.H:183
scalar KHindex_
Index for KH Breakup.
Definition: SprayParcel.H:178
A class for handling character strings derived from std::string.
Definition: string.H:72
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
A primitive field of type <T> with automated input and output.
streamFormat format() const noexcept
Get the current stream format.
Reacting spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:40
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127