ReactingMultiphaseParcelIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
34 
35 template<class ParcelType>
37 (
38  sizeof(scalar)
39 );
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class ParcelType>
46 (
47  Istream& is,
48  bool readFields
49 )
50 :
51  ParcelType(is, readFields),
52  mass0_(0.0),
53  YGas_(0),
54  YLiquid_(0),
55  YSolid_(0),
56  canCombust_(0)
57 {
58  if (readFields)
59  {
63 
64  if (is.format() == IOstream::ASCII)
65  {
66  is >> mass0_ >> Yg >> Yl >> Ys;
67  }
68  else
69  {
70  is.read(reinterpret_cast<char*>(&mass0_), sizeofFields_);
71  is >> Yg >> Yl >> Ys;
72  }
73 
74  YGas_.transfer(Yg);
75  YLiquid_.transfer(Yl);
76  YSolid_.transfer(Ys);
77  }
78 
79  // Check state of Istream
80  is.check
81  (
82  "ReactingMultiphaseParcel<ParcelType>::ReactingMultiphaseParcel"
83  "("
84  "const polyMesh&, "
85  "Istream&, "
86  "bool"
87  ")"
88  );
89 }
90 
91 
92 template<class ParcelType>
93 template<class CloudType>
95 {
97 }
98 
99 
100 template<class ParcelType>
101 template<class CloudType, class CompositionType>
103 (
104  CloudType& c,
105  const CompositionType& compModel
106 )
107 {
108  bool valid = c.size();
109 
110  ParcelType::readFields(c, compModel);
111 
112  IOField<scalar> mass0
113  (
114  c.fieldIOobject("mass0", IOobject::MUST_READ),
115  valid
116  );
117  c.checkFieldIOobject(c, mass0);
118 
119  label i = 0;
120  forAllIter(typename CloudType, c, iter)
121  {
123  p.mass0_ = mass0[i++];
124  }
125 
126  // Get names and sizes for each Y...
127  const label idGas = compModel.idGas();
128  const wordList& gasNames = compModel.componentNames(idGas);
129  const label idLiquid = compModel.idLiquid();
130  const wordList& liquidNames = compModel.componentNames(idLiquid);
131  const label idSolid = compModel.idSolid();
132  const wordList& solidNames = compModel.componentNames(idSolid);
133  const wordList& stateLabels = compModel.stateLabels();
134 
135  // Set storage for each Y... for each parcel
136  forAllIter(typename CloudType, c, iter)
137  {
139  p.YGas_.setSize(gasNames.size(), 0.0);
140  p.YLiquid_.setSize(liquidNames.size(), 0.0);
141  p.YSolid_.setSize(solidNames.size(), 0.0);
142  }
143 
144  // Populate YGas for each parcel
145  forAll(gasNames, j)
146  {
147  IOField<scalar> YGas
148  (
149  c.fieldIOobject
150  (
151  "Y" + gasNames[j] + stateLabels[idGas],
153  ),
154  valid
155  );
156 
157  label i = 0;
158  forAllIter(typename CloudType, c, iter)
159  {
161  p.YGas_[j] = YGas[i++]/(p.Y()[idGas] + rootVSmall);
162  }
163  }
164  // Populate YLiquid for each parcel
165  forAll(liquidNames, j)
166  {
167  IOField<scalar> YLiquid
168  (
169  c.fieldIOobject
170  (
171  "Y" + liquidNames[j] + stateLabels[idLiquid],
173  ),
174  valid
175  );
176 
177  label i = 0;
178  forAllIter(typename CloudType, c, iter)
179  {
181  p.YLiquid_[j] = YLiquid[i++]/(p.Y()[idLiquid] + rootVSmall);
182  }
183  }
184  // Populate YSolid for each parcel
185  forAll(solidNames, j)
186  {
187  IOField<scalar> YSolid
188  (
189  c.fieldIOobject
190  (
191  "Y" + solidNames[j] + stateLabels[idSolid],
193  ),
194  valid
195  );
196 
197  label i = 0;
198  forAllIter(typename CloudType, c, iter)
199  {
201  p.YSolid_[j] = YSolid[i++]/(p.Y()[idSolid] + rootVSmall);
202  }
203  }
204 }
205 
206 
207 template<class ParcelType>
208 template<class CloudType>
210 {
211  ParcelType::writeFields(c);
212 }
213 
214 
215 template<class ParcelType>
216 template<class CloudType, class CompositionType>
218 (
219  const CloudType& c,
220  const CompositionType& compModel
221 )
222 {
223  ParcelType::writeFields(c, compModel);
224 
225  label np = c.size();
226 
227  {
228  IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
229 
230  label i = 0;
231  forAllConstIter(typename CloudType, c, iter)
232  {
233  const ReactingMultiphaseParcel<ParcelType>& p = iter();
234  mass0[i++] = p.mass0_;
235  }
236  mass0.write(np > 0);
237 
238  // Write the composition fractions
239  const wordList& stateLabels = compModel.stateLabels();
240 
241  const label idGas = compModel.idGas();
242  const wordList& gasNames = compModel.componentNames(idGas);
243  forAll(gasNames, j)
244  {
245  IOField<scalar> YGas
246  (
247  c.fieldIOobject
248  (
249  "Y" + gasNames[j] + stateLabels[idGas],
251  ),
252  np
253  );
254 
255  label i = 0;
256  forAllConstIter(typename CloudType, c, iter)
257  {
258  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
259  YGas[i++] = p0.YGas()[j]*p0.Y()[idGas];
260  }
261 
262  YGas.write(np > 0);
263  }
264 
265  const label idLiquid = compModel.idLiquid();
266  const wordList& liquidNames = compModel.componentNames(idLiquid);
267  forAll(liquidNames, j)
268  {
269  IOField<scalar> YLiquid
270  (
271  c.fieldIOobject
272  (
273  "Y" + liquidNames[j] + stateLabels[idLiquid],
275  ),
276  np
277  );
278 
279  label i = 0;
280  forAllConstIter(typename CloudType, c, iter)
281  {
282  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
283  YLiquid[i++] = p0.YLiquid()[j]*p0.Y()[idLiquid];
284  }
285 
286  YLiquid.write(np > 0);
287  }
288 
289  const label idSolid = compModel.idSolid();
290  const wordList& solidNames = compModel.componentNames(idSolid);
291  forAll(solidNames, j)
292  {
293  IOField<scalar> YSolid
294  (
295  c.fieldIOobject
296  (
297  "Y" + solidNames[j] + stateLabels[idSolid],
299  ),
300  np
301  );
302 
303  label i = 0;
304  forAllConstIter(typename CloudType, c, iter)
305  {
306  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
307  YSolid[i++] = p0.YSolid()[j]*p0.Y()[idSolid];
308  }
309 
310  YSolid.write(np > 0);
311  }
312  }
313 }
314 
315 
316 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
317 
318 template<class ParcelType>
319 Foam::Ostream& Foam::operator<<
320 (
321  Ostream& os,
323 )
324 {
325  if (os.format() == IOstream::ASCII)
326  {
327  os << static_cast<const ParcelType&>(p)
328  << token::SPACE << p.mass0()
329  << token::SPACE << p.YGas()
330  << token::SPACE << p.YLiquid()
331  << token::SPACE << p.YSolid();
332  }
333  else
334  {
335  os << static_cast<const ParcelType&>(p);
336  os.write
337  (
338  reinterpret_cast<const char*>(&p.mass0_),
339  ReactingMultiphaseParcel<ParcelType>::sizeofFields_
340  );
341  os << p.YGas() << p.YLiquid() << p.YSolid();
342  }
343 
344  // Check state of Ostream
345  os.check
346  (
347  "Ostream& operator<<"
348  "("
349  "Ostream&, "
350  "const ReactingMultiphaseParcel<ParcelType>&"
351  ")"
352  );
353 
354  return os;
355 }
356 
357 
358 // ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
virtual Istream & read(token &)=0
Return next token from stream.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
scalarField YLiquid_
Mass fractions of liquids [].
scalarField YSolid_
Mass fractions of solids [].
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
scalarField YGas_
Mass fractions of gases [].
const scalarField & YGas() const
Return const access to mass fractions of gases.
ReactingMultiphaseParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from mesh, coordinates and topology.
const scalarField & YSolid() const
Return const access to mass fractions of solids.
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
Definition: string.H:79
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
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:59
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
volScalarField & p