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-2019 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() == IOstream::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  }
117 
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 valid = c.size();
139 
140  ParcelType::readFields(c, compModel);
141 
142  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::MUST_READ), valid);
143  c.checkFieldIOobject(c, d0);
144 
145  IOField<vector> position0
146  (
147  c.fieldIOobject("position0", IOobject::MUST_READ),
148  valid
149  );
150  c.checkFieldIOobject(c, position0);
151 
152  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::MUST_READ), valid);
153  c.checkFieldIOobject(c, sigma);
154 
155  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::MUST_READ), valid);
156  c.checkFieldIOobject(c, mu);
157 
158  IOField<scalar> liquidCore
159  (
160  c.fieldIOobject("liquidCore", IOobject::MUST_READ),
161  valid
162  );
163  c.checkFieldIOobject(c, liquidCore);
164 
165  IOField<scalar> KHindex
166  (
167  c.fieldIOobject("KHindex", IOobject::MUST_READ),
168  valid
169  );
170  c.checkFieldIOobject(c, KHindex);
171 
173  (
174  c.fieldIOobject("y", IOobject::MUST_READ),
175  valid
176  );
177  c.checkFieldIOobject(c, y);
178 
179  IOField<scalar> yDot
180  (
181  c.fieldIOobject("yDot", IOobject::MUST_READ),
182  valid
183  );
184  c.checkFieldIOobject(c, yDot);
185 
186  IOField<scalar> tc
187  (
188  c.fieldIOobject("tc", IOobject::MUST_READ),
189  valid
190  );
191  c.checkFieldIOobject(c, tc);
192 
193  IOField<scalar> ms
194  (
195  c.fieldIOobject("ms", IOobject::MUST_READ),
196  valid
197  );
198  c.checkFieldIOobject(c, ms);
199 
200  IOField<scalar> injector
201  (
202  c.fieldIOobject("injector", IOobject::MUST_READ),
203  valid
204  );
205  c.checkFieldIOobject(c, injector);
206 
207  IOField<scalar> tMom
208  (
209  c.fieldIOobject("tMom", IOobject::MUST_READ),
210  valid
211  );
212  c.checkFieldIOobject(c, tMom);
213 
214  IOField<scalar> user
215  (
216  c.fieldIOobject("user", IOobject::MUST_READ),
217  valid
218  );
219  c.checkFieldIOobject(c, user);
220 
221  label i = 0;
222  for (SprayParcel<ParcelType>& p : c)
223  {
224  p.d0_ = d0[i];
225  p.position0_ = position0[i];
226  p.sigma_ = sigma[i];
227  p.mu_ = mu[i];
228  p.liquidCore_ = liquidCore[i];
229  p.KHindex_ = KHindex[i];
230  p.y_ = y[i];
231  p.yDot_ = yDot[i];
232  p.tc_ = tc[i];
233  p.ms_ = ms[i];
234  p.injector_ = injector[i];
235  p.tMom_ = tMom[i];
236  p.user_ = user[i];
237 
238  ++i;
239  }
240 }
241 
242 
243 template<class ParcelType>
244 template<class CloudType>
246 {
248 }
249 
250 
251 template<class ParcelType>
252 template<class CloudType, class CompositionType>
254 (
255  const CloudType& c,
256  const CompositionType& compModel
257 )
258 {
259  ParcelType::writeFields(c, compModel);
260 
261  const label np = c.size();
262  const bool valid = np;
263 
264 
265  IOField<scalar> d0(c.fieldIOobject("d0", IOobject::NO_READ), np);
266  IOField<vector> position0
267  (
268  c.fieldIOobject("position0", IOobject::NO_READ),
269  np
270  );
271  IOField<scalar> sigma(c.fieldIOobject("sigma", IOobject::NO_READ), np);
272  IOField<scalar> mu(c.fieldIOobject("mu", IOobject::NO_READ), np);
273  IOField<scalar> liquidCore
274  (
275  c.fieldIOobject("liquidCore", IOobject::NO_READ),
276  np
277  );
278  IOField<scalar> KHindex(c.fieldIOobject("KHindex", IOobject::NO_READ), np);
279  IOField<scalar> y(c.fieldIOobject("y", IOobject::NO_READ), np);
280  IOField<scalar> yDot(c.fieldIOobject("yDot", IOobject::NO_READ), np);
281  IOField<scalar> tc(c.fieldIOobject("tc", IOobject::NO_READ), np);
282  IOField<scalar> ms(c.fieldIOobject("ms", IOobject::NO_READ), np);
283  IOField<scalar> injector
284  (
285  c.fieldIOobject("injector", IOobject::NO_READ),
286  np
287  );
288  IOField<scalar> tMom(c.fieldIOobject("tMom", IOobject::NO_READ), np);
289  IOField<scalar> user(c.fieldIOobject("user", IOobject::NO_READ), np);
290 
291  label i = 0;
292  for (const SprayParcel<ParcelType>& p : c)
293  {
294  d0[i] = p.d0_;
295  position0[i] = p.position0_;
296  sigma[i] = p.sigma_;
297  mu[i] = p.mu_;
298  liquidCore[i] = p.liquidCore_;
299  KHindex[i] = p.KHindex_;
300  y[i] = p.y_;
301  yDot[i] = p.yDot_;
302  tc[i] = p.tc_;
303  ms[i] = p.ms_;
304  injector[i] = p.injector_;
305  tMom[i] = p.tMom_;
306  user[i] = p.user_;
307  ++i;
308  }
309 
310  d0.write(valid);
311  position0.write(valid);
312  sigma.write(valid);
313  mu.write(valid);
314  liquidCore.write(valid);
315  KHindex.write(valid);
316  y.write(valid);
317  yDot.write(valid);
318  tc.write(valid);
319  ms.write(valid);
320  injector.write(valid);
321  tMom.write(valid);
322  user.write(valid);
323 }
324 
325 
326 template<class ParcelType>
328 (
329  Ostream& os,
330  const wordRes& filters,
331  const word& delim,
332  const bool namesOnly
333 ) const
334 {
335  ParcelType::writeProperties(os, filters, delim, namesOnly);
336 
337  #undef writeProp
338  #define writeProp(Name, Value) \
339  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
340 
341  writeProp("d0", d0_);
342  writeProp("position0", position0_);
343  writeProp("sigma", sigma_);
344  writeProp("mu", mu_);
345  writeProp("liquidCore", liquidCore_);
346  writeProp("KHindex", KHindex_);
347  writeProp("y", y_);
348  writeProp("yDot", yDot_);
349  writeProp("tc", tc_);
350  writeProp("ms", ms_);
351  writeProp("injector", injector_);
352  writeProp("tMom", tMom_);
353  writeProp("user", user_);
354 
355  #undef writeProp
356 }
357 
358 
359 template<class ParcelType>
360 template<class CloudType>
362 (
363  CloudType& c,
364  const objectRegistry& obr
365 )
366 {
367  ParcelType::readObjects(c, obr);
368 }
369 
370 
371 template<class ParcelType>
372 template<class CloudType>
374 (
375  const CloudType& c,
376  objectRegistry& obr
377 )
378 {
379  ParcelType::writeObjects(c, obr);
380 }
381 
382 
383 template<class ParcelType>
384 template<class CloudType, class CompositionType>
386 (
387  CloudType& c,
388  const CompositionType& compModel,
389  const objectRegistry& obr
390 )
391 {
392  ParcelType::readObjects(c, compModel, obr);
393 
394  if (!c.size()) return;
395 
396  const auto& d0 = cloud::lookupIOField<scalar>("d0", obr);
397  const auto& position0 = cloud::lookupIOField<vector>("position0", obr);
398  const auto& sigma = cloud::lookupIOField<scalar>("sigma", obr);
399  const auto& mu = cloud::lookupIOField<scalar>("mu", obr);
400  const auto& liquidCore = cloud::lookupIOField<scalar>("liquidCore", obr);
401  const auto& KHindex = cloud::lookupIOField<scalar>("KHindex", obr);
402  const auto& y = cloud::lookupIOField<scalar>("y", obr);
403  const auto& yDot = cloud::lookupIOField<scalar>("yDot", obr);
404  const auto& tc = cloud::lookupIOField<scalar>("tc", obr);
405  const auto& ms = cloud::lookupIOField<scalar>("ms", obr);
406  const auto& injector = cloud::lookupIOField<scalar>("injector", obr);
407  const auto& tMom = cloud::lookupIOField<scalar>("tMom", obr);
408  const auto& user = cloud::lookupIOField<scalar>("user", obr);
409 
410  label i = 0;
411  for (SprayParcel<ParcelType>& p : c)
412  {
413  p.d0_ = d0[i];
414  p.position0_ = position0[i];
415  p.sigma_ = sigma[i];
416  p.mu_ = mu[i];
417  p.liquidCore_ = liquidCore[i];
418  p.KHindex_ = KHindex[i];
419  p.y_ = y[i];
420  p.yDot_ = yDot[i];
421  p.tc_ = tc[i];
422  p.ms_ = ms[i];
423  p.injector_ = injector[i];
424  p.tMom_ = tMom[i];
425  p.user_ = user[i];
426 
427  ++i;
428  }
429 }
430 
431 
432 template<class ParcelType>
433 template<class CloudType, class CompositionType>
435 (
436  const CloudType& c,
437  const CompositionType& compModel,
438  objectRegistry& obr
439 )
440 {
441  ParcelType::writeObjects(c, compModel, obr);
442 
443  const label np = c.size();
444 
445  auto& d0 = cloud::createIOField<scalar>("d0", np, obr);
446  auto& position0 = cloud::createIOField<vector>("position0", np, obr);
447  auto& sigma = cloud::createIOField<scalar>("sigma", np, obr);
448  auto& mu = cloud::createIOField<scalar>("mu", np, obr);
449  auto& liquidCore = cloud::createIOField<scalar>("liquidCore", np, obr);
450  auto& KHindex = cloud::createIOField<scalar>("KHindex", np, obr);
451  auto& y = cloud::createIOField<scalar>("y", np, obr);
452  auto& yDot= cloud::createIOField<scalar>("yDot", np, obr);
453  auto& tc = cloud::createIOField<scalar>("tc", np, obr);
454  auto& ms = cloud::createIOField<scalar>("ms", np, obr);
455  auto& injector = cloud::createIOField<scalar>("injector", np, obr);
456  auto& tMom = cloud::createIOField<scalar>("tMom", np, obr);
457  auto& user = cloud::createIOField<scalar>("user", np, obr);
458 
459  label i = 0;
460  for (const SprayParcel<ParcelType>& p : c)
461  {
462  d0[i] = p.d0_;
463  position0[i] = p.position0_;
464  sigma[i] = p.sigma_;
465  mu[i] = p.mu_;
466  liquidCore[i] = p.liquidCore_;
467  KHindex[i] = p.KHindex_;
468  y[i] = p.y_;
469  yDot[i] = p.yDot_;
470  tc[i] = p.tc_;
471  ms[i] = p.ms_;
472  injector[i] = p.injector_;
473  tMom[i] = p.tMom_;
474  user[i] = p.user_;
475 
476  ++i;
477  }
478 }
479 
480 
481 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
482 
483 template<class ParcelType>
484 Foam::Ostream& Foam::operator<<
485 (
486  Ostream& os,
488 )
489 {
490  if (os.format() == IOstream::ASCII)
491  {
492  os << static_cast<const ParcelType&>(p)
493  << token::SPACE << p.d0()
494  << token::SPACE << p.position0()
495  << token::SPACE << p.sigma()
496  << token::SPACE << p.mu()
497  << token::SPACE << p.liquidCore()
498  << token::SPACE << p.KHindex()
499  << token::SPACE << p.y()
500  << token::SPACE << p.yDot()
501  << token::SPACE << p.tc()
502  << token::SPACE << p.ms()
503  << token::SPACE << p.injector()
504  << token::SPACE << p.tMom()
505  << token::SPACE << p.user();
506  }
507  else
508  {
509  os << static_cast<const ParcelType&>(p);
510  os.write
511  (
512  reinterpret_cast<const char*>(&p.d0_),
513  SprayParcel<ParcelType>::sizeofFields
514  );
515  }
516 
517  os.check(FUNCTION_NAME);
518  return os;
519 }
520 
521 
522 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::SprayParcel::writeProperties
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Definition: SprayParcelIO.C:328
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::SprayParcel::writeObjects
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
Definition: SprayParcelIO.C:374
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::SprayParcel::sizeofFields
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: SprayParcel.H:182
Foam::IOstreamOption::format
streamFormat format() const noexcept
Get the current stream format.
Definition: IOstreamOption.H:286
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:76
Foam::SprayParcel::readFields
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
Definition: SprayParcelIO.C:133
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::IOstream::checkLabelSize
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Definition: IOstream.H:300
Foam::Istream::endRawRead
virtual bool endRawRead()=0
End of low-level raw binary read.
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::SprayParcel
Reacting spray parcel, with added functionality for atomization and breakup.
Definition: SprayParcel.H:46
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::IOstream::checkScalarSize
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Definition: IOstream.H:309
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::SprayParcel::writeFields
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
Definition: SprayParcelIO.C:254
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
writeProp
#define writeProp(Name, Value)
Foam::Istream::beginRawRead
virtual bool beginRawRead()=0
Start of low-level raw binary read.
SprayParcel.H
Foam::SprayParcel::SprayParcel
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:111
Foam::SprayParcel::readObjects
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
Definition: SprayParcelIO.C:362
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::writeFields
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields, const bool writeFaceFields)
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::Istream::read
virtual Istream & read(token &)=0
Return next token from stream.