ReactingParcelIO.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 "ReactingParcel.H"
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 template<class ParcelType>
37 
38 
39 template<class ParcelType>
41 (
42  sizeof(scalar)
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  mass0_(0.0),
59  Y_(0)
60 {
61  if (readFields)
62  {
64 
65  if (is.format() == IOstreamOption::ASCII)
66  {
67  is >> mass0_;
68  }
69  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
70  {
71  // Non-native label or scalar size
72 
73  is.beginRawRead();
74 
75  readRawScalar(is, &mass0_);
76 
77  is.endRawRead();
78  }
79  else
80  {
81  is.read(reinterpret_cast<char*>(&mass0_), sizeofFields);
82  }
83 
84  is >> Ymix;
85 
86  Y_.transfer(Ymix);
87  }
88 
89  is.check(FUNCTION_NAME);
90 }
91 
92 
93 template<class ParcelType>
94 template<class CloudType>
96 {
98 }
99 
100 
101 template<class ParcelType>
102 template<class CloudType, class CompositionType>
104 (
105  CloudType& c,
106  const CompositionType& compModel
107 )
108 {
109  const bool readOnProc = c.size();
110 
112 
113  IOField<scalar> mass0
114  (
115  c.fieldIOobject("mass0", IOobject::MUST_READ),
116  readOnProc
117  );
118  c.checkFieldIOobject(c, mass0);
119 
120  label i = 0;
122  {
123  p.mass0_ = mass0[i];
124 
125  ++i;
126  }
127 
128  // Get names and sizes for each Y...
129  const wordList& phaseTypes = compModel.phaseTypes();
130  const label nPhases = phaseTypes.size();
131  wordList stateLabels(nPhases, "");
132  if (compModel.nPhase() == 1)
133  {
134  stateLabels = compModel.stateLabels()[0];
135  }
136 
137 
138  // Set storage for each Y... for each parcel
139  for (ReactingParcel<ParcelType>& p : c)
140  {
141  p.Y_.setSize(nPhases, 0.0);
142  }
143 
144  // Populate Y for each parcel
145  forAll(phaseTypes, j)
146  {
147  IOField<scalar> Y
148  (
149  c.fieldIOobject
150  (
151  "Y" + phaseTypes[j] + stateLabels[j],
153  ),
154  readOnProc
155  );
156 
157  label i = 0;
158  for (ReactingParcel<ParcelType>& p : c)
159  {
160  p.Y_[j] = Y[i];
161 
162  ++i;
163  }
164  }
165 }
166 
167 
168 template<class ParcelType>
169 template<class CloudType>
171 {
173 }
174 
175 
176 template<class ParcelType>
177 template<class CloudType, class CompositionType>
179 (
180  const CloudType& c,
181  const CompositionType& compModel
182 )
183 {
185 
186  const label np = c.size();
187  const bool writeOnProc = c.size();
188 
189  {
190  IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
191 
192  label i = 0;
193  for (const ReactingParcel<ParcelType>& p : c)
194  {
195  mass0[i] = p.mass0_;
196 
197  ++i;
198  }
199  mass0.write(writeOnProc);
200 
201  // Write the composition fractions
202  const wordList& phaseTypes = compModel.phaseTypes();
203  wordList stateLabels(phaseTypes.size(), "");
204  if (compModel.nPhase() == 1)
205  {
206  stateLabels = compModel.stateLabels()[0];
207  }
208 
209  forAll(phaseTypes, j)
210  {
211  IOField<scalar> Y
212  (
213  c.fieldIOobject
214  (
215  "Y" + phaseTypes[j] + stateLabels[j],
217  ),
218  np
219  );
220 
221  label i = 0;
222  for (const ReactingParcel<ParcelType>& p : c)
223  {
224  Y[i] = p.Y()[j];
225 
226  ++i;
227  }
228 
229  Y.write(writeOnProc);
230  }
231  }
232 }
233 
234 
235 template<class ParcelType>
237 (
238  Ostream& os,
239  const wordRes& filters,
240  const word& delim,
241  const bool namesOnly
242 ) const
243 {
244  ParcelType::writeProperties(os, filters, delim, namesOnly);
245 
246  #undef writeProp
247  #define writeProp(Name, Value) \
248  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
249 
250  writeProp("mass0", mass0_);
251  writeProp("Y", Y_);
252 
253  #undef writeProp
254 }
255 
256 
257 template<class ParcelType>
258 template<class CloudType>
260 (
261  CloudType& c,
262  const objectRegistry& obr
263 )
264 {
265  ParcelType::readObjects(c, obr);
266 }
267 
268 
269 template<class ParcelType>
270 template<class CloudType>
272 (
273  const CloudType& c,
274  objectRegistry& obr
275 )
276 {
277  ParcelType::writeObjects(c, obr);
278 }
279 
280 
281 template<class ParcelType>
282 template<class CloudType, class CompositionType>
284 (
285  CloudType& c,
286  const CompositionType& compModel,
287  const objectRegistry& obr
288 )
289 {
290  ParcelType::readObjects(c, obr);
291 
292  if (!c.size()) return;
293 
294 
295  auto& mass0 = cloud::lookupIOField<scalar>("mass0", obr);
296 
297  label i = 0;
299  {
300  p.mass0_ = mass0[i];
301 
302  ++i;
303  }
304 
305  // The composition fractions
306  const wordList& phaseTypes = compModel.phaseTypes();
307  wordList stateLabels(phaseTypes.size(), "");
308  if (compModel.nPhase() == 1)
309  {
310  stateLabels = compModel.stateLabels()[0];
311  }
312 
313  forAll(phaseTypes, j)
314  {
315  const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
316  auto& Y = cloud::lookupIOField<scalar>(fieldName, obr);
317 
318  label i = 0;
319  for (ReactingParcel<ParcelType>& p : c)
320  {
321  p.Y()[j] = Y[i];
322 
323  ++i;
324  }
325  }
326 }
327 
328 
329 template<class ParcelType>
330 template<class CloudType, class CompositionType>
332 (
333  const CloudType& c,
334  const CompositionType& compModel,
335  objectRegistry& obr
336 )
337 {
338  ParcelType::writeObjects(c, obr);
339 
340  const label np = c.size();
341  const bool writeOnProc = c.size();
342 
343  if (writeOnProc)
344  {
345  auto& mass0 = cloud::createIOField<scalar>("mass0", np, obr);
346 
347  label i = 0;
348  for (const ReactingParcel<ParcelType>& p : c)
349  {
350  mass0[i] = p.mass0_;
351 
352  ++i;
353  }
354 
355  // Write the composition fractions
356  const wordList& phaseTypes = compModel.phaseTypes();
357  wordList stateLabels(phaseTypes.size(), "");
358  if (compModel.nPhase() == 1)
359  {
360  stateLabels = compModel.stateLabels()[0];
361  }
362 
363  forAll(phaseTypes, j)
364  {
365  const word fieldName = "Y" + phaseTypes[j] + stateLabels[j];
366  auto& Y = cloud::createIOField<scalar>(fieldName, np, obr);
367 
368  label i = 0;
369  for (const ReactingParcel<ParcelType>& p : c)
370  {
371  Y[i] = p.Y()[j];
372 
373  ++i;
374  }
375  }
376  }
377 }
378 
379 
380 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
381 
382 template<class ParcelType>
383 Foam::Ostream& Foam::operator<<
384 (
385  Ostream& os,
386  const ReactingParcel<ParcelType>& p
387 )
388 {
390  {
391  os << static_cast<const ParcelType&>(p)
392  << token::SPACE << p.mass0()
393  << token::SPACE << p.Y();
394  }
395  else
396  {
397  os << static_cast<const ParcelType&>(p);
398  os.write
399  (
400  reinterpret_cast<const char*>(&p.mass0_),
402  );
403  os << p.Y();
404  }
405 
407  return os;
408 }
409 
410 
411 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
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
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
"ascii" (normal default)
scalarField Y_
Mass fractions of mixture [].
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
virtual bool endRawRead()=0
End of low-level raw binary read.
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual Istream & read(token &)=0
Return next token from stream.
dynamicFvMesh & mesh
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.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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
static const std::size_t sizeofFields
Size in bytes of the fields.
List< word > wordList
List of word.
Definition: fileName.H:59
virtual bool beginRawRead()=0
Start of low-level raw binary read.
PtrList< volScalarField > & Y
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
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
#define writeProp(Name, Value)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Reacting parcel class with one/two-way coupling with the continuous phase.
volScalarField & p
Registry of regIOobjects.
scalar mass0_
Initial mass [kg].
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.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.