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-2020 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
120  forAllIter(typename CloudType, c, iter)
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(typename CloudType, c, iter)
143  {
145  p.YGas_[j] = YGas[i++]/(p.Y()[GAS] + rootVSmall);
146  }
147  }
148  // Populate YLiquid for each parcel
149  forAll(liquidNames, j)
150  {
151  IOField<scalar> YLiquid
152  (
153  c.fieldIOobject
154  (
155  "Y" + liquidNames[j] + stateLabels[idLiquid],
156  IOobject::MUST_READ
157  ),
158  valid
159  );
160 
161  label i = 0;
162  forAllIter(typename CloudType, c, iter)
163  {
165  p.YLiquid_[j] = YLiquid[i++]/(p.Y()[LIQ] + rootVSmall);
166  }
167  }
168  // Populate YSolid for each parcel
169  forAll(solidNames, j)
170  {
171  IOField<scalar> YSolid
172  (
173  c.fieldIOobject
174  (
175  "Y" + solidNames[j] + stateLabels[idSolid],
176  IOobject::MUST_READ
177  ),
178  valid
179  );
180 
181  label i = 0;
182  forAllIter(typename CloudType, c, iter)
183  {
185  p.YSolid_[j] = YSolid[i++]/(p.Y()[SLD] + rootVSmall);
186  }
187  }
188 }
189 
190 
191 template<class ParcelType>
192 template<class CloudType>
194 {
195  ParcelType::writeFields(c);
196 }
197 
198 
199 template<class ParcelType>
200 template<class CloudType, class CompositionType>
202 (
203  const CloudType& c,
204  const CompositionType& compModel
205 )
206 {
207  ParcelType::writeFields(c, compModel);
208 
209  label np = c.size();
210 
211  // Write the composition fractions
212  {
213  const wordList& stateLabels = compModel.stateLabels();
214 
215  const label idGas = compModel.idGas();
216  const wordList& gasNames = compModel.componentNames(idGas);
217  forAll(gasNames, j)
218  {
219  IOField<scalar> YGas
220  (
221  c.fieldIOobject
222  (
223  "Y" + gasNames[j] + stateLabels[idGas],
224  IOobject::NO_READ
225  ),
226  np
227  );
228 
229  label i = 0;
230  forAllConstIter(typename CloudType, c, iter)
231  {
232  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
233  YGas[i++] = p0.YGas()[j]*p0.Y()[GAS];
234  }
235 
236  YGas.write(np > 0);
237  }
238 
239  const label idLiquid = compModel.idLiquid();
240  const wordList& liquidNames = compModel.componentNames(idLiquid);
241  forAll(liquidNames, j)
242  {
243  IOField<scalar> YLiquid
244  (
245  c.fieldIOobject
246  (
247  "Y" + liquidNames[j] + stateLabels[idLiquid],
248  IOobject::NO_READ
249  ),
250  np
251  );
252 
253  label i = 0;
254  forAllConstIter(typename CloudType, c, iter)
255  {
256  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
257  YLiquid[i++] = p0.YLiquid()[j]*p0.Y()[LIQ];
258  }
259 
260  YLiquid.write(np > 0);
261  }
262 
263  const label idSolid = compModel.idSolid();
264  const wordList& solidNames = compModel.componentNames(idSolid);
265  forAll(solidNames, j)
266  {
267  IOField<scalar> YSolid
268  (
269  c.fieldIOobject
270  (
271  "Y" + solidNames[j] + stateLabels[idSolid],
272  IOobject::NO_READ
273  ),
274  np
275  );
276 
277  label i = 0;
278  forAllConstIter(typename CloudType, c, iter)
279  {
280  const ReactingMultiphaseParcel<ParcelType>& p0 = iter();
281  YSolid[i++] = p0.YSolid()[j]*p0.Y()[SLD];
282  }
283 
284  YSolid.write(np > 0);
285  }
286  }
287 }
288 
289 
290 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
291 
292 template<class ParcelType>
293 Foam::Ostream& Foam::operator<<
294 (
295  Ostream& os,
297 )
298 {
299  scalarField YGasLoc(p.YGas()*p.Y()[0]);
300  scalarField YLiquidLoc(p.YLiquid()*p.Y()[1]);
301  scalarField YSolidLoc(p.YSolid()*p.Y()[2]);
302  if (os.format() == IOstream::ASCII)
303  {
304  os << static_cast<const ParcelType&>(p)
305  << token::SPACE << YGasLoc
306  << token::SPACE << YLiquidLoc
307  << token::SPACE << YSolidLoc;
308  }
309  else
310  {
311  os << static_cast<const ParcelType&>(p);
312  os << YGasLoc << YLiquidLoc << YSolidLoc;
313  }
314 
315  // Check state of Ostream
316  os.check
317  (
318  "Ostream& operator<<"
319  "("
320  "Ostream&, "
321  "const ReactingMultiphaseParcel<ParcelType>&"
322  ")"
323  );
324 
325  return os;
326 }
327 
328 
329 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
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:164
const wordList solidNames(rp.found("solid") ? rp["solid"] :wordList(0))
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:190
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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:76
volScalarField & p
void transfer(List< T > &)
Transfer contents of the argument List into this.
Definition: DynamicListI.H:285
A class for handling character strings derived from std::string.
Definition: string.H:76
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:170
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.