ReactingParcelIO.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 | Copyright (C) 2011-2017 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 
26 #include "ReactingParcel.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
32 Foam::string Foam::ReactingParcel<ParcelType>::propertyList_ =
34 
35 template<class ParcelType>
37 (
38  sizeof(scalar)
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  mass0_(0.0),
54  Y_(0),
55  pc_(0.0)
56 {
57  if (readFields)
58  {
60 
61  if (is.format() == IOstream::ASCII)
62  {
63  is >> mass0_ >> Ymix;
64  }
65  else
66  {
67  is.read(reinterpret_cast<char*>(&mass0_), sizeofFields_);
68  is >> Ymix;
69  }
70 
71  Y_.transfer(Ymix);
72  }
73 
74  // Check state of Istream
75  is.check
76  (
77  "ReactingParcel<ParcelType>::ReactingParcel"
78  "("
79  "const polyMesh&, "
80  "Istream&, "
81  "bool"
82  ")"
83  );
84 }
85 
86 
87 template<class ParcelType>
88 template<class CloudType>
90 {
92 }
93 
94 
95 template<class ParcelType>
96 template<class CloudType, class CompositionType>
98 (
99  CloudType& c,
100  const CompositionType& compModel
101 )
102 {
103  bool valid = c.size();
104 
106 
107  IOField<scalar> mass0
108  (
109  c.fieldIOobject("mass0", IOobject::MUST_READ),
110  valid
111  );
112  c.checkFieldIOobject(c, mass0);
113 
114  label i = 0;
115  forAllIter(typename Cloud<ReactingParcel<ParcelType>>, c, iter)
116  {
117  ReactingParcel<ParcelType>& p = iter();
118  p.mass0_ = mass0[i++];
119  }
120 
121  // Get names and sizes for each Y...
122  const wordList& phaseTypes = compModel.phaseTypes();
123  const label nPhases = phaseTypes.size();
124  wordList stateLabels(nPhases, "");
125  if (compModel.nPhase() == 1)
126  {
127  stateLabels = compModel.stateLabels()[0];
128  }
129 
130 
131  // Set storage for each Y... for each parcel
132  forAllIter(typename Cloud<ReactingParcel<ParcelType>>, c, iter)
133  {
134  ReactingParcel<ParcelType>& p = iter();
135  p.Y_.setSize(nPhases, 0.0);
136  }
137 
138  // Populate Y for each parcel
139  forAll(phaseTypes, j)
140  {
142  (
143  c.fieldIOobject
144  (
145  "Y" + phaseTypes[j] + stateLabels[j],
146  IOobject::MUST_READ
147  ),
148  valid
149  );
150 
151  label i = 0;
152  forAllIter(typename Cloud<ReactingParcel<ParcelType>>, c, iter)
153  {
154  ReactingParcel<ParcelType>& p = iter();
155  p.Y_[j] = Y[i++];
156  }
157  }
158 }
159 
160 
161 template<class ParcelType>
162 template<class CloudType>
164 {
165  ParcelType::writeFields(c);
166 }
167 
168 
169 template<class ParcelType>
170 template<class CloudType, class CompositionType>
172 (
173  const CloudType& c,
174  const CompositionType& compModel
175 )
176 {
177  ParcelType::writeFields(c);
178 
179  const label np = c.size();
180 
181  {
182  IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
183 
184  label i = 0;
186  {
187  const ReactingParcel<ParcelType>& p = iter();
188  mass0[i++] = p.mass0_;
189  }
190  mass0.write(np > 0);
191 
192  // Write the composition fractions
193  const wordList& phaseTypes = compModel.phaseTypes();
194  wordList stateLabels(phaseTypes.size(), "");
195  if (compModel.nPhase() == 1)
196  {
197  stateLabels = compModel.stateLabels()[0];
198  }
199 
200  forAll(phaseTypes, j)
201  {
203  (
204  c.fieldIOobject
205  (
206  "Y" + phaseTypes[j] + stateLabels[j],
207  IOobject::NO_READ
208  ),
209  np
210  );
211  label i = 0;
213  (
215  c,
216  iter
217  )
218  {
219  const ReactingParcel<ParcelType>& p = iter();
220  Y[i++] = p.Y()[j];
221  }
222 
223  Y.write(np > 0);
224  }
225  }
226 }
227 
228 
229 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
230 
231 template<class ParcelType>
232 Foam::Ostream& Foam::operator<<
233 (
234  Ostream& os,
236 )
237 {
238  if (os.format() == IOstream::ASCII)
239  {
240  os << static_cast<const ParcelType&>(p)
241  << token::SPACE << p.mass0()
242  << token::SPACE << p.Y();
243  }
244  else
245  {
246  os << static_cast<const ParcelType&>(p);
247  os.write
248  (
249  reinterpret_cast<const char*>(&p.mass0_),
251  );
252  os << p.Y();
253  }
254 
255  // Check state of Ostream
256  os.check
257  (
258  "Ostream& operator<<(Ostream&, const ReactingParcel<ParcelType>&)"
259  );
260 
261  return os;
262 }
263 
264 
265 // ************************************************************************* //
#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
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
scalarField Y_
Mass fractions of mixture [].
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
label size() const
Definition: Cloud.H:162
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
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.
virtual Istream & read(token &)=0
Return next token from stream.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
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
PtrList< volScalarField > & Y
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:186
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Reacting parcel class with one/two-way coupling with the continuous phase.
volScalarField & p
scalar mass0_
Initial mass [kg].
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:166
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const scalarField & Y() const
Return const access to mass fractions of mixture [].