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-2015 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>
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  {
59  DynamicList<scalar> Ymix;
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 {
91  if (!c.size())
92  {
93  return;
94  }
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  if (!c.size())
109  {
110  return;
111  }
112 
114 
115  IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::MUST_READ));
116  c.checkFieldIOobject(c, mass0);
117 
118  label i = 0;
119  forAllIter(typename Cloud<ReactingParcel<ParcelType> >, c, iter)
120  {
121  ReactingParcel<ParcelType>& p = iter();
122  p.mass0_ = mass0[i++];
123  }
124 
125  // Get names and sizes for each Y...
126  const wordList& phaseTypes = compModel.phaseTypes();
127  const label nPhases = phaseTypes.size();
128  wordList stateLabels(nPhases, "");
129  if (compModel.nPhase() == 1)
130  {
131  stateLabels = compModel.stateLabels()[0];
132  }
133 
134 
135  // Set storage for each Y... for each parcel
136  forAllIter(typename Cloud<ReactingParcel<ParcelType> >, c, iter)
137  {
138  ReactingParcel<ParcelType>& p = iter();
139  p.Y_.setSize(nPhases, 0.0);
140  }
141 
142  // Populate Y for each parcel
143  forAll(phaseTypes, j)
144  {
146  (
147  c.fieldIOobject
148  (
149  "Y" + phaseTypes[j] + stateLabels[j],
151  )
152  );
153 
154  label i = 0;
155  forAllIter(typename Cloud<ReactingParcel<ParcelType> >, c, iter)
156  {
157  ReactingParcel<ParcelType>& p = iter();
158  p.Y_[j] = Y[i++];
159  }
160  }
161 }
162 
163 
164 template<class ParcelType>
165 template<class CloudType>
167 {
168  ParcelType::writeFields(c);
169 }
170 
171 
172 template<class ParcelType>
173 template<class CloudType, class CompositionType>
175 (
176  const CloudType& c,
177  const CompositionType& compModel
178 )
179 {
180  ParcelType::writeFields(c);
181 
182  const label np = c.size();
183 
184  if (np > 0)
185  {
186  IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
187 
188  label i = 0;
190  {
191  const ReactingParcel<ParcelType>& p = iter();
192  mass0[i++] = p.mass0_;
193  }
194  mass0.write();
195 
196  // Write the composition fractions
197  const wordList& phaseTypes = compModel.phaseTypes();
198  wordList stateLabels(phaseTypes.size(), "");
199  if (compModel.nPhase() == 1)
200  {
201  stateLabels = compModel.stateLabels()[0];
202  }
203 
204  forAll(phaseTypes, j)
205  {
207  (
208  c.fieldIOobject
209  (
210  "Y" + phaseTypes[j] + stateLabels[j],
212  ),
213  np
214  );
215  label i = 0;
217  (
219  c,
220  iter
221  )
222  {
223  const ReactingParcel<ParcelType>& p = iter();
224  Y[i++] = p.Y()[j];
225  }
226 
227  Y.write();
228  }
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
234 
235 template<class ParcelType>
236 Foam::Ostream& Foam::operator<<
237 (
238  Ostream& os,
240 )
241 {
242  if (os.format() == IOstream::ASCII)
243  {
244  os << static_cast<const ParcelType&>(p)
245  << token::SPACE << p.mass0()
246  << token::SPACE << p.Y();
247  }
248  else
249  {
250  os << static_cast<const ParcelType&>(p);
251  os.write
252  (
253  reinterpret_cast<const char*>(&p.mass0_),
255  );
256  os << p.Y();
257  }
258 
259  // Check state of Ostream
260  os.check
261  (
262  "Ostream& operator<<(Ostream&, const ReactingParcel<ParcelType>&)"
263  );
264 
265  return os;
266 }
267 
268 
269 // ************************************************************************* //
ReactingParcel.H
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
p
p
Definition: pEqn.H:62
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::ReactingParcel::readFields
static void readFields(CloudType &c, const CompositionType &compModel)
Read.
Definition: ReactingParcelIO.C:103
Foam::ReactingParcel::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: ReactingParcel.H:72
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:74
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::IOstream::ASCII
@ ASCII
Definition: IOstream.H:88
Foam::label
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
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::ReactingParcel::ReactingParcel
ReactingParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
Definition: ReactingParcelI.H:64
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::ReactingParcel::mass0_
scalar mass0_
Initial mass [kg].
Definition: ReactingParcel.H:161
Foam::ReactingParcel::writeFields
static void writeFields(const CloudType &c, const CompositionType &compModel)
Write.
Definition: ReactingParcelIO.C:175
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
readFields
void readFields(const boolList &haveMesh, const fvMesh &mesh, const autoPtr< fvMeshSubset > &subsetterPtr, IOobjectList &allObjects, PtrList< GeoField > &fields)
Definition: redistributePar.C:589
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:52
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::ReactingParcel::Y_
scalarField Y_
Mass fractions of mixture [].
Definition: ReactingParcel.H:164
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::ReactingParcel
Reacting parcel class with one/two-way coupling with the continuous phase.
Definition: ReactingParcel.H:50
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36
Foam::token::SPACE
@ SPACE
Definition: token.H:95