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-2018 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>
32 Foam::string Foam::ReactingMultiphaseParcel<ParcelType>::propertyList_ =
34 
35 template<class ParcelType>
37 (
38  0
39 );
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 template<class ParcelType>
46 (
47  const polyMesh& mesh,
48  Istream& is,
49  bool readFields
50 )
51 :
52  ParcelType(mesh, is, readFields),
53  YGas_(0),
54  YLiquid_(0),
55  YSolid_(0),
56  canCombust_(0)
57 {
58  if (readFields)
59  {
63 
64  is >> Yg >> Yl >> Ys;
65 
66  YGas_.transfer(Yg);
67  YLiquid_.transfer(Yl);
68  YSolid_.transfer(Ys);
69 
70  // scale the mass fractions
71  const scalarField& YMix = this->Y_;
72  YGas_ /= YMix[GAS] + rootVSmall;
73  YLiquid_ /= YMix[LIQ] + rootVSmall;
74  YSolid_ /= YMix[SLD] + rootVSmall;
75  }
76 
77  // Check state of Istream
78  is.check
79  (
80  "ReactingMultiphaseParcel<ParcelType>::ReactingMultiphaseParcel"
81  "("
82  "const polyMesh&, "
83  "Istream&, "
84  "bool"
85  ")"
86  );
87 }
88 
89 
90 template<class ParcelType>
91 template<class CloudType>
93 {
95 }
96 
97 
98 template<class ParcelType>
99 template<class CloudType, class CompositionType>
101 (
102  CloudType& c,
103  const CompositionType& compModel
104 )
105 {
106  bool valid = c.size();
107 
108  ParcelType::readFields(c, compModel);
109 
110  // Get names and sizes for each Y...
111  const label idGas = compModel.idGas();
112  const wordList& gasNames = compModel.componentNames(idGas);
113  const label idLiquid = compModel.idLiquid();
114  const wordList& liquidNames = compModel.componentNames(idLiquid);
115  const label idSolid = compModel.idSolid();
116  const wordList& solidNames = compModel.componentNames(idSolid);
117  const wordList& stateLabels = compModel.stateLabels();
118 
119  // Set storage for each Y... for each parcel
121  {
123  p.YGas_.setSize(gasNames.size(), 0.0);
124  p.YLiquid_.setSize(liquidNames.size(), 0.0);
125  p.YSolid_.setSize(solidNames.size(), 0.0);
126  }
127 
128  // Populate YGas for each parcel
129  forAll(gasNames, j)
130  {
131  IOField<scalar> YGas
132  (
133  c.fieldIOobject
134  (
135  "Y" + gasNames[j] + stateLabels[idGas],
136  IOobject::MUST_READ
137  ),
138  valid
139  );
140 
141  label i = 0;
142  forAllIter
143  (
145  c,
146  iter
147  )
148  {
150  p.YGas_[j] = YGas[i++]/(p.Y()[GAS] + rootVSmall);
151  }
152  }
153  // Populate YLiquid for each parcel
154  forAll(liquidNames, j)
155  {
156  IOField<scalar> YLiquid
157  (
158  c.fieldIOobject
159  (
160  "Y" + liquidNames[j] + stateLabels[idLiquid],
161  IOobject::MUST_READ
162  ),
163  valid
164  );
165 
166  label i = 0;
167  forAllIter
168  (
170  c,
171  iter
172  )
173  {
175  p.YLiquid_[j] = YLiquid[i++]/(p.Y()[LIQ] + rootVSmall);
176  }
177  }
178  // Populate YSolid for each parcel
179  forAll(solidNames, j)
180  {
181  IOField<scalar> YSolid
182  (
183  c.fieldIOobject
184  (
185  "Y" + solidNames[j] + stateLabels[idSolid],
186  IOobject::MUST_READ
187  ),
188  valid
189  );
190 
191  label i = 0;
192  forAllIter
193  (
195  c,
196  iter
197  )
198  {
200  p.YSolid_[j] = YSolid[i++]/(p.Y()[SLD] + rootVSmall);
201  }
202  }
203 }
204 
205 
206 template<class ParcelType>
207 template<class CloudType>
209 {
210  ParcelType::writeFields(c);
211 }
212 
213 
214 template<class ParcelType>
215 template<class CloudType, class CompositionType>
217 (
218  const CloudType& c,
219  const CompositionType& compModel
220 )
221 {
222  ParcelType::writeFields(c, compModel);
223 
224  label np = c.size();
225 
226  // Write the composition fractions
227  {
228  const wordList& stateLabels = compModel.stateLabels();
229 
230  const label idGas = compModel.idGas();
231  const wordList& gasNames = compModel.componentNames(idGas);
232  forAll(gasNames, j)
233  {
234  IOField<scalar> YGas
235  (
236  c.fieldIOobject
237  (
238  "Y" + gasNames[j] + stateLabels[idGas],
239  IOobject::NO_READ
240  ),
241  np
242  );
243 
244  label i = 0;
246  (
248  c,
249  iter
250  )
251  {
252  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
253  YGas[i++] = p0.YGas()[j]*p0.Y()[GAS];
254  }
255 
256  YGas.write(np > 0);
257  }
258 
259  const label idLiquid = compModel.idLiquid();
260  const wordList& liquidNames = compModel.componentNames(idLiquid);
261  forAll(liquidNames, j)
262  {
263  IOField<scalar> YLiquid
264  (
265  c.fieldIOobject
266  (
267  "Y" + liquidNames[j] + stateLabels[idLiquid],
268  IOobject::NO_READ
269  ),
270  np
271  );
272 
273  label i = 0;
275  (
277  c,
278  iter
279  )
280  {
281  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
282  YLiquid[i++] = p0.YLiquid()[j]*p0.Y()[LIQ];
283  }
284 
285  YLiquid.write(np > 0);
286  }
287 
288  const label idSolid = compModel.idSolid();
289  const wordList& solidNames = compModel.componentNames(idSolid);
290  forAll(solidNames, j)
291  {
292  IOField<scalar> YSolid
293  (
294  c.fieldIOobject
295  (
296  "Y" + solidNames[j] + stateLabels[idSolid],
297  IOobject::NO_READ
298  ),
299  np
300  );
301 
302  label i = 0;
304  (
306  c,
307  iter
308  )
309  {
310  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
311  YSolid[i++] = p0.YSolid()[j]*p0.Y()[SLD];
312  }
313 
314  YSolid.write(np > 0);
315  }
316  }
317 }
318 
319 
320 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
321 
322 template<class ParcelType>
323 Foam::Ostream& Foam::operator<<
324 (
325  Ostream& os,
327 )
328 {
329  scalarField YGasLoc(p.YGas()*p.Y()[0]);
330  scalarField YLiquidLoc(p.YLiquid()*p.Y()[1]);
331  scalarField YSolidLoc(p.YSolid()*p.Y()[2]);
332  if (os.format() == IOstream::ASCII)
333  {
334  os << static_cast<const ParcelType&>(p)
335  << token::SPACE << YGasLoc
336  << token::SPACE << YLiquidLoc
337  << token::SPACE << YSolidLoc;
338  }
339  else
340  {
341  os << static_cast<const ParcelType&>(p);
342  os << YGasLoc << YLiquidLoc << YSolidLoc;
343  }
344 
345  // Check state of Ostream
346  os.check
347  (
348  "Ostream& operator<<"
349  "("
350  "Ostream&, "
351  "const ReactingMultiphaseParcel<ParcelType>&"
352  ")"
353  );
354 
355  return os;
356 }
357 
358 
359 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:151
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.
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 [].
const scalarField & YGas() const
Return const access to mass fractions of gases.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Base cloud calls templated on particle type.
Definition: Cloud.H:52
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
scalarField YSolid_
Mass fractions of solids [].
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
scalarField YLiquid_
Mass fractions of liquids [].
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
void transfer(List< T > &)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:259
A class for handling character strings derived from std::string.
Definition: string.H:74
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:164
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const scalarField & YSolid() const
Return const access to mass fractions of solids.