ReactingMultiphaseParcelIO.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 
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 template<class ParcelType>
37 
38 
39 template<class ParcelType>
41 (
42  0
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  YGas_(0),
59  YLiquid_(0),
60  YSolid_(0),
61  canCombust_(0)
62 {
63  if (readFields)
64  {
68 
69  is >> Yg >> Yl >> Ys;
70 
71  YGas_.transfer(Yg);
72  YLiquid_.transfer(Yl);
73  YSolid_.transfer(Ys);
74 
75  // scale the mass fractions
76  const scalarField& YMix = this->Y_;
77  YGas_ /= YMix[GAS] + ROOTVSMALL;
78  YLiquid_ /= YMix[LIQ] + ROOTVSMALL;
79  YSolid_ /= YMix[SLD] + ROOTVSMALL;
80  }
81 
82  is.check(FUNCTION_NAME);
83 }
84 
85 
86 template<class ParcelType>
87 template<class CloudType>
89 {
91 }
92 
93 
94 template<class ParcelType>
95 template<class CloudType, class CompositionType>
97 (
98  CloudType& c,
99  const CompositionType& compModel
100 )
101 {
102  bool valid = c.size();
103 
104  ParcelType::readFields(c, compModel);
105 
106  // Get names and sizes for each Y...
107  const label idGas = compModel.idGas();
108  const wordList& gasNames = compModel.componentNames(idGas);
109  const label idLiquid = compModel.idLiquid();
110  const wordList& liquidNames = compModel.componentNames(idLiquid);
111  const label idSolid = compModel.idSolid();
112  const wordList& solidNames = compModel.componentNames(idSolid);
113  const wordList& stateLabels = compModel.stateLabels();
114 
115  // Set storage for each Y... for each parcel
117  {
118  p.YGas_.setSize(gasNames.size(), 0.0);
119  p.YLiquid_.setSize(liquidNames.size(), 0.0);
120  p.YSolid_.setSize(solidNames.size(), 0.0);
121  }
122 
123  // Populate YGas for each parcel
124  forAll(gasNames, j)
125  {
126  IOField<scalar> YGas
127  (
128  c.fieldIOobject
129  (
130  "Y" + gasNames[j] + stateLabels[idGas],
131  IOobject::MUST_READ
132  ),
133  valid
134  );
135 
136  label i = 0;
138  {
139  p.YGas_[j] = YGas[i]/(p.Y()[GAS] + ROOTVSMALL);
140  ++i;
141  }
142  }
143  // Populate YLiquid for each parcel
144  forAll(liquidNames, j)
145  {
146  IOField<scalar> YLiquid
147  (
148  c.fieldIOobject
149  (
150  "Y" + liquidNames[j] + stateLabels[idLiquid],
151  IOobject::MUST_READ
152  ),
153  valid
154  );
155 
156  label i = 0;
158  {
159  p.YLiquid_[j] = YLiquid[i]/(p.Y()[LIQ] + ROOTVSMALL);
160  ++i;
161  }
162  }
163  // Populate YSolid for each parcel
164  forAll(solidNames, j)
165  {
166  IOField<scalar> YSolid
167  (
168  c.fieldIOobject
169  (
170  "Y" + solidNames[j] + stateLabels[idSolid],
171  IOobject::MUST_READ
172  ),
173  valid
174  );
175 
176  label i = 0;
178  {
179  p.YSolid_[j] = YSolid[i]/(p.Y()[SLD] + ROOTVSMALL);
180  ++i;
181  }
182  }
183 }
184 
185 
186 template<class ParcelType>
187 template<class CloudType>
189 {
191 }
192 
193 
194 template<class ParcelType>
195 template<class CloudType, class CompositionType>
197 (
198  const CloudType& c,
199  const CompositionType& compModel
200 )
201 {
202  ParcelType::writeFields(c, compModel);
203 
204  label np = c.size();
205 
206  // Write the composition fractions
207  {
208  const wordList& stateLabels = compModel.stateLabels();
209 
210  const label idGas = compModel.idGas();
211  const wordList& gasNames = compModel.componentNames(idGas);
212  forAll(gasNames, j)
213  {
214  IOField<scalar> YGas
215  (
216  c.fieldIOobject
217  (
218  "Y" + gasNames[j] + stateLabels[idGas],
219  IOobject::NO_READ
220  ),
221  np
222  );
223 
224  label i = 0;
226  {
227  YGas[i] = p0.YGas()[j]*p0.Y()[GAS];
228  ++i;
229  }
230 
231  YGas.write(np > 0);
232  }
233 
234  const label idLiquid = compModel.idLiquid();
235  const wordList& liquidNames = compModel.componentNames(idLiquid);
236  forAll(liquidNames, j)
237  {
238  IOField<scalar> YLiquid
239  (
240  c.fieldIOobject
241  (
242  "Y" + liquidNames[j] + stateLabels[idLiquid],
243  IOobject::NO_READ
244  ),
245  np
246  );
247 
248  label i = 0;
250  {
251  YLiquid[i] = p0.YLiquid()[j]*p0.Y()[LIQ];
252  ++i;
253  }
254 
255  YLiquid.write(np > 0);
256  }
257 
258  const label idSolid = compModel.idSolid();
259  const wordList& solidNames = compModel.componentNames(idSolid);
260  forAll(solidNames, j)
261  {
262  IOField<scalar> YSolid
263  (
264  c.fieldIOobject
265  (
266  "Y" + solidNames[j] + stateLabels[idSolid],
267  IOobject::NO_READ
268  ),
269  np
270  );
271 
272  label i = 0;
274  {
275  YSolid[i] = p0.YSolid()[j]*p0.Y()[SLD];
276  ++i;
277  }
278 
279  YSolid.write(np > 0);
280  }
281  }
282 }
283 
284 
285 template<class ParcelType>
287 (
288  Ostream& os,
289  const wordRes& filters,
290  const word& delim,
291  const bool namesOnly
292 ) const
293 {
294  ParcelType::writeProperties(os, filters, delim, namesOnly);
295 
296  #undef writeProp
297  #define writeProp(Name, Value) \
298  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
299 
300  writeProp("YGas", YGas_);
301  writeProp("YLiquid", YLiquid_);
302  writeProp("YSolid", YSolid_);
303  writeProp("canCombust", canCombust_);
304 
305  #undef writeProp
306 }
307 
308 
309 template<class ParcelType>
310 template<class CloudType>
312 (
313  CloudType& c,
314  const objectRegistry& obr
315 )
316 {
317  ParcelType::readObjects(c, obr);
318 }
319 
320 
321 template<class ParcelType>
322 template<class CloudType>
324 (
325  const CloudType& c,
326  objectRegistry& obr
327 )
328 {
329  ParcelType::writeObjects(c, obr);
330 }
331 
332 
333 template<class ParcelType>
334 template<class CloudType, class CompositionType>
336 (
337  CloudType& c,
338  const CompositionType& compModel,
339  const objectRegistry& obr
340 )
341 {
342  ParcelType::readObjects(c, obr);
343 
344  const label np = c.size();
345 
346  // The composition fractions
347  if (np > 0)
348  {
349  const wordList& stateLabels = compModel.stateLabels();
350 
351  const label idGas = compModel.idGas();
352  const wordList& gasNames = compModel.componentNames(idGas);
353  forAll(gasNames, j)
354  {
355  const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
356  const auto& YGas = cloud::lookupIOField<scalar>(fieldName, obr);
357 
358  label i = 0;
360  {
361  p0.YGas()[j]*p0.Y()[GAS] = YGas[i];
362  ++i;
363  }
364  }
365 
366  const label idLiquid = compModel.idLiquid();
367  const wordList& liquidNames = compModel.componentNames(idLiquid);
368  forAll(liquidNames, j)
369  {
370  const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
371  const auto& YLiquid = cloud::lookupIOField<scalar>(fieldName, obr);
372 
373  label i = 0;
375  {
376  p0.YLiquid()[j]*p0.Y()[LIQ] = YLiquid[i];
377  ++i;
378  }
379  }
380 
381  const label idSolid = compModel.idSolid();
382  const wordList& solidNames = compModel.componentNames(idSolid);
383  forAll(solidNames, j)
384  {
385  const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
386  const auto& YSolid = cloud::lookupIOField<scalar>(fieldName, obr);
387 
388  label i = 0;
390  {
391  p0.YSolid()[j]*p0.Y()[SLD] = YSolid[i];
392  ++i;
393  }
394  }
395  }
396 }
397 
398 
399 template<class ParcelType>
400 template<class CloudType, class CompositionType>
402 (
403  const CloudType& c,
404  const CompositionType& compModel,
405  objectRegistry& obr
406 )
407 {
408  ParcelType::writeObjects(c, obr);
409 
410  const label np = c.size();
411 
412  // Write the composition fractions
413  if (np > 0)
414  {
415  const wordList& stateLabels = compModel.stateLabels();
416 
417  const label idGas = compModel.idGas();
418  const wordList& gasNames = compModel.componentNames(idGas);
419  forAll(gasNames, j)
420  {
421  const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
422  auto& YGas = cloud::createIOField<scalar>(fieldName, np, obr);
423 
424  label i = 0;
426  {
427  YGas[i] = p0.YGas()[j]*p0.Y()[GAS];
428  ++i;
429  }
430  }
431 
432  const label idLiquid = compModel.idLiquid();
433  const wordList& liquidNames = compModel.componentNames(idLiquid);
434  forAll(liquidNames, j)
435  {
436  const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
437  auto& YLiquid = cloud::createIOField<scalar>(fieldName, np, obr);
438 
439  label i = 0;
441  {
442  YLiquid[i] = p0.YLiquid()[j]*p0.Y()[LIQ];
443  ++i;
444  }
445  }
446 
447  const label idSolid = compModel.idSolid();
448  const wordList& solidNames = compModel.componentNames(idSolid);
449  forAll(solidNames, j)
450  {
451  const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
452  auto& YSolid = cloud::createIOField<scalar>(fieldName, np, obr);
453 
454  label i = 0;
456  {
457  YSolid[i] = p0.YSolid()[j]*p0.Y()[SLD];
458  ++i;
459  }
460  }
461  }
462 }
463 
464 
465 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
466 
467 template<class ParcelType>
468 Foam::Ostream& Foam::operator<<
469 (
470  Ostream& os,
472 )
473 {
474  scalarField YGasLoc(p.YGas()*p.Y()[0]);
475  scalarField YLiquidLoc(p.YLiquid()*p.Y()[1]);
476  scalarField YSolidLoc(p.YSolid()*p.Y()[2]);
477 
478  if (os.format() == IOstream::ASCII)
479  {
480  os << static_cast<const ParcelType&>(p)
481  << token::SPACE << YGasLoc
482  << token::SPACE << YLiquidLoc
483  << token::SPACE << YSolidLoc;
484  }
485  else
486  {
487  os << static_cast<const ParcelType&>(p);
488  os << YGasLoc << YLiquidLoc << YSolidLoc;
489  }
490 
491  os.check(FUNCTION_NAME);
492  return os;
493 }
494 
495 
496 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
writeProp
#define writeProp(Name, Value)
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::ReactingMultiphaseParcel::sizeofFields
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: ReactingMultiphaseParcel.H:72
Foam::ReactingMultiphaseParcel::writeFields
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
Definition: ReactingMultiphaseParcelIO.C:197
Foam::DynamicList< scalar >
Foam::ReactingMultiphaseParcel::writeProperties
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
Definition: ReactingMultiphaseParcelIO.C:287
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:73
Foam::ReactingMultiphaseParcel::readFields
static void readFields(CloudType &c, const CompositionType &compModel)
Read - composition supplied.
Definition: ReactingMultiphaseParcelIO.C:97
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::label
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:62
Foam::Field< scalar >
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::ReactingMultiphaseParcel
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
Definition: ReactingMultiphaseParcel.H:55
ReactingMultiphaseParcel.H
Foam::writeFields
void writeFields(const fvMesh &mesh, const wordHashSet &selectedFields)
Foam::ReactingMultiphaseParcel::readObjects
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
Definition: ReactingMultiphaseParcelIO.C:312
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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
Foam::List< word >
Foam::ReactingMultiphaseParcel::writeObjects
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
Definition: ReactingMultiphaseParcelIO.C:324
Foam::DynamicList::transfer
void transfer(List< T > &lst)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:426
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::ReactingMultiphaseParcel::ReactingMultiphaseParcel
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
Definition: ReactingMultiphaseParcelI.H:71
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
p0
const volScalarField & p0
Definition: EEqn.H:36
solidNames
const wordList solidNames(rp["solid"])