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-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 
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 
76  is.check(FUNCTION_NAME);
77 }
78 
79 
80 template<class ParcelType>
81 template<class CloudType>
83 {
85 }
86 
87 
88 template<class ParcelType>
89 template<class CloudType, class CompositionType>
91 (
92  CloudType& c,
93  const CompositionType& compModel
94 )
95 {
96  const bool readOnProc = c.size();
97 
98  ParcelType::readFields(c, compModel);
99 
100  // Get names and sizes for each Y...
101  const label idGas = compModel.idGas();
102  const wordList& gasNames = compModel.componentNames(idGas);
103  const label idLiquid = compModel.idLiquid();
104  const wordList& liquidNames = compModel.componentNames(idLiquid);
105  const label idSolid = compModel.idSolid();
106  const wordList& solidNames = compModel.componentNames(idSolid);
107  const wordList& stateLabels = compModel.stateLabels();
108 
109  // Set storage for each Y... for each parcel
111  {
112  p.YGas_.setSize(gasNames.size(), 0.0);
113  p.YLiquid_.setSize(liquidNames.size(), 0.0);
114  p.YSolid_.setSize(solidNames.size(), 0.0);
115  }
116 
117  // Populate YGas for each parcel
118  forAll(gasNames, j)
119  {
120  IOField<scalar> YGas
121  (
122  c.fieldIOobject
123  (
124  "Y" + gasNames[j] + stateLabels[idGas],
126  ),
127  readOnProc
128  );
129 
130  label i = 0;
131  for (ReactingMultiphaseParcel<ParcelType>& p : c)
132  {
133  p.YGas_[j] = YGas[i]/(max(p.Y()[GAS], SMALL));
134  ++i;
135  }
136  }
137  // Populate YLiquid for each parcel
138  forAll(liquidNames, j)
139  {
140  IOField<scalar> YLiquid
141  (
142  c.fieldIOobject
143  (
144  "Y" + liquidNames[j] + stateLabels[idLiquid],
146  ),
147  readOnProc
148  );
149 
150  label i = 0;
151  for (ReactingMultiphaseParcel<ParcelType>& p : c)
152  {
153  p.YLiquid_[j] = YLiquid[i]/(max(p.Y()[LIQ], SMALL));
154  ++i;
155  }
156  }
157  // Populate YSolid for each parcel
158  forAll(solidNames, j)
159  {
160  IOField<scalar> YSolid
161  (
162  c.fieldIOobject
163  (
164  "Y" + solidNames[j] + stateLabels[idSolid],
166  ),
167  readOnProc
168  );
169 
170  label i = 0;
171  for (ReactingMultiphaseParcel<ParcelType>& p : c)
172  {
173  p.YSolid_[j] = YSolid[i]/(max(p.Y()[SLD], SMALL));
174  ++i;
175  }
176  }
177 }
178 
179 
180 template<class ParcelType>
181 template<class CloudType>
183 {
185 }
186 
187 
188 template<class ParcelType>
189 template<class CloudType, class CompositionType>
191 (
192  const CloudType& c,
193  const CompositionType& compModel
194 )
195 {
196  ParcelType::writeFields(c, compModel);
197 
198  const label np = c.size();
199  const bool writeOnProc = c.size();
200 
201  // Write the composition fractions
202  {
203  const wordList& stateLabels = compModel.stateLabels();
204 
205  const label idGas = compModel.idGas();
206  const wordList& gasNames = compModel.componentNames(idGas);
207  forAll(gasNames, j)
208  {
209  IOField<scalar> YGas
210  (
211  c.fieldIOobject
212  (
213  "Y" + gasNames[j] + stateLabels[idGas],
215  ),
216  np
217  );
218 
219  label i = 0;
221  {
222  YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
223  ++i;
224  }
225 
226  YGas.write(writeOnProc);
227  }
228 
229  const label idLiquid = compModel.idLiquid();
230  const wordList& liquidNames = compModel.componentNames(idLiquid);
231  forAll(liquidNames, j)
232  {
233  IOField<scalar> YLiquid
234  (
235  c.fieldIOobject
236  (
237  "Y" + liquidNames[j] + stateLabels[idLiquid],
239  ),
240  np
241  );
242 
243  label i = 0;
244  for (const ReactingMultiphaseParcel<ParcelType>& p0 : c)
245  {
246  YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
247  ++i;
248  }
249 
250  YLiquid.write(writeOnProc);
251  }
252 
253  const label idSolid = compModel.idSolid();
254  const wordList& solidNames = compModel.componentNames(idSolid);
255  forAll(solidNames, j)
256  {
257  IOField<scalar> YSolid
258  (
259  c.fieldIOobject
260  (
261  "Y" + solidNames[j] + stateLabels[idSolid],
263  ),
264  np
265  );
266 
267  label i = 0;
268  for (const ReactingMultiphaseParcel<ParcelType>& p0 : c)
269  {
270  YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
271  ++i;
272  }
273 
274  YSolid.write(writeOnProc);
275  }
276  }
277 }
278 
279 
280 template<class ParcelType>
282 (
283  Ostream& os,
284  const wordRes& filters,
285  const word& delim,
286  const bool namesOnly
287 ) const
288 {
289  ParcelType::writeProperties(os, filters, delim, namesOnly);
290 
291  #undef writeProp
292  #define writeProp(Name, Value) \
293  ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
294 
295  writeProp("YGas", YGas_);
296  writeProp("YLiquid", YLiquid_);
297  writeProp("YSolid", YSolid_);
298  writeProp("canCombust", canCombust_);
299 
300  #undef writeProp
301 }
302 
303 
304 template<class ParcelType>
305 template<class CloudType>
307 (
308  CloudType& c,
309  const objectRegistry& obr
310 )
311 {
312  ParcelType::readObjects(c, obr);
313 }
314 
315 
316 template<class ParcelType>
317 template<class CloudType>
319 (
320  const CloudType& c,
321  objectRegistry& obr
322 )
323 {
324  ParcelType::writeObjects(c, obr);
325 }
326 
327 
328 template<class ParcelType>
329 template<class CloudType, class CompositionType>
331 (
332  CloudType& c,
333  const CompositionType& compModel,
334  const objectRegistry& obr
335 )
336 {
337  ParcelType::readObjects(c, obr);
338 
339  // const label np = c.size();
340  const bool readOnProc = c.size();
341 
342  // The composition fractions
343  if (readOnProc)
344  {
345  const wordList& stateLabels = compModel.stateLabels();
346 
347  const label idGas = compModel.idGas();
348  const wordList& gasNames = compModel.componentNames(idGas);
349  forAll(gasNames, j)
350  {
351  const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
352  const auto& YGas = cloud::lookupIOField<scalar>(fieldName, obr);
353 
354  label i = 0;
356  {
357  p0.YGas()[j]*max(p0.Y()[GAS], SMALL) = YGas[i];
358  ++i;
359  }
360  }
361 
362  const label idLiquid = compModel.idLiquid();
363  const wordList& liquidNames = compModel.componentNames(idLiquid);
364  forAll(liquidNames, j)
365  {
366  const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
367  const auto& YLiquid = cloud::lookupIOField<scalar>(fieldName, obr);
368 
369  label i = 0;
370  for (ReactingMultiphaseParcel<ParcelType>& p0 : c)
371  {
372  p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL) = YLiquid[i];
373  ++i;
374  }
375  }
376 
377  const label idSolid = compModel.idSolid();
378  const wordList& solidNames = compModel.componentNames(idSolid);
379  forAll(solidNames, j)
380  {
381  const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
382  const auto& YSolid = cloud::lookupIOField<scalar>(fieldName, obr);
383 
384  label i = 0;
385  for (ReactingMultiphaseParcel<ParcelType>& p0 : c)
386  {
387  p0.YSolid()[j]*max(p0.Y()[SLD], SMALL) = YSolid[i];
388  ++i;
389  }
390  }
391  }
392 }
393 
394 
395 template<class ParcelType>
396 template<class CloudType, class CompositionType>
398 (
399  const CloudType& c,
400  const CompositionType& compModel,
401  objectRegistry& obr
402 )
403 {
404  ParcelType::writeObjects(c, obr);
405 
406  const label np = c.size();
407  const bool writeOnProc = c.size();
408 
409  // Write the composition fractions
410  if (writeOnProc)
411  {
412  const wordList& stateLabels = compModel.stateLabels();
413 
414  const label idGas = compModel.idGas();
415  const wordList& gasNames = compModel.componentNames(idGas);
416  forAll(gasNames, j)
417  {
418  const word fieldName = "Y" + gasNames[j] + stateLabels[idGas];
419  auto& YGas = cloud::createIOField<scalar>(fieldName, np, obr);
420 
421  label i = 0;
422  for (const ReactingMultiphaseParcel<ParcelType>& p0 : c)
423  {
424  YGas[i] = p0.YGas()[j]*max(p0.Y()[GAS], SMALL);
425  ++i;
426  }
427  }
428 
429  const label idLiquid = compModel.idLiquid();
430  const wordList& liquidNames = compModel.componentNames(idLiquid);
431  forAll(liquidNames, j)
432  {
433  const word fieldName = "Y" + liquidNames[j] + stateLabels[idLiquid];
434  auto& YLiquid = cloud::createIOField<scalar>(fieldName, np, obr);
435 
436  label i = 0;
437  for (const ReactingMultiphaseParcel<ParcelType>& p0 : c)
438  {
439  YLiquid[i] = p0.YLiquid()[j]*max(p0.Y()[LIQ], SMALL);
440  ++i;
441  }
442  }
443 
444  const label idSolid = compModel.idSolid();
445  const wordList& solidNames = compModel.componentNames(idSolid);
446  forAll(solidNames, j)
447  {
448  const word fieldName = "Y" + solidNames[j] + stateLabels[idSolid];
449  auto& YSolid = cloud::createIOField<scalar>(fieldName, np, obr);
450 
451  label i = 0;
452  for (const ReactingMultiphaseParcel<ParcelType>& p0 : c)
453  {
454  YSolid[i] = p0.YSolid()[j]*max(p0.Y()[SLD], SMALL);
455  ++i;
456  }
457  }
458  }
459 }
460 
461 
462 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
463 
464 template<class ParcelType>
465 Foam::Ostream& Foam::operator<<
466 (
467  Ostream& os,
468  const ReactingMultiphaseParcel<ParcelType>& p
469 )
470 {
471  scalarField YGasLoc(p.YGas());
472  scalarField YLiquidLoc(p.YLiquid());
473  scalarField YSolidLoc(p.YSolid());
474 
476  {
477  os << static_cast<const ParcelType&>(p)
478  << token::SPACE << YGasLoc
479  << token::SPACE << YLiquidLoc
480  << token::SPACE << YSolidLoc;
481  }
482  else
483  {
484  os << static_cast<const ParcelType&>(p);
485  os << YGasLoc << YLiquidLoc << YSolidLoc;
486  }
487 
489  return os;
490 }
491 
492 
493 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
DSMCCloud< dsmcParcel > CloudType
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write - composition supplied.
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
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
"ascii" (normal default)
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
const wordList solidNames(rp["solid"])
#define writeProp(Name, Value)
#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.
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, position and topology.
scalarField YGas_
Mass fractions of gases [].
void setSize(const label n)
Alias for resize()
Definition: List.H:316
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.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
scalarField YSolid_
Mass fractions of solids [].
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
static const std::size_t sizeofFields
Size in bytes of the fields.
scalarField YLiquid_
Mass fractions of liquids [].
List< word > wordList
List of word.
Definition: fileName.H:59
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase...
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
Registry of regIOobjects.
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.
const volScalarField & p0
Definition: EEqn.H:36
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.