Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / lagrangian / intermediate / parcels / Templates / ReactingParcel / ReactingParcelIO.C
blob3bd35f709675e9940ed21b3602bcc581b862c2ed
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2011 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "ReactingParcel.H"
27 #include "IOstreams.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 template<class ParcelType>
32 Foam::string Foam::ReactingParcel<ParcelType>::propHeader =
33     ParcelType::propHeader
34   + " mass0"
35   + " nPhases(Y1..YN)";
38 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
40 template<class ParcelType>
41 Foam::ReactingParcel<ParcelType>::ReactingParcel
43     const polyMesh& mesh,
44     Istream& is,
45     bool readFields
48     ParcelType(mesh, is, readFields),
49     mass0_(0.0),
50     Y_(0),
51     pc_(0.0)
53     if (readFields)
54     {
55         DynamicList<scalar> Ymix;
57         if (is.format() == IOstream::ASCII)
58         {
59             is >> mass0_ >> Ymix;
60         }
61         else
62         {
63             is.read
64             (
65                 reinterpret_cast<char*>(&mass0_),
66               + sizeof(mass0_)
67             );
68             is >> Ymix;
69         }
71         Y_.transfer(Ymix);
72     }
74     // Check state of Istream
75     is.check
76     (
77         "ReactingParcel<ParcelType>::ReactingParcel"
78         "("
79             "const polyMesh&, "
80             "Istream&, "
81             "bool"
82         ")"
83     );
87 template<class ParcelType>
88 template<class CloudType>
89 void Foam::ReactingParcel<ParcelType>::readFields(CloudType& c)
91     if (!c.size())
92     {
93         return;
94     }
96     ParcelType::readFields(c);
100 template<class ParcelType>
101 template<class CloudType, class CompositionType>
102 void Foam::ReactingParcel<ParcelType>::readFields
104     CloudType& c,
105     const CompositionType& compModel
108     if (!c.size())
109     {
110         return;
111     }
113     ParcelType::readFields(c);
115     IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::MUST_READ));
116     c.checkFieldIOobject(c, mass0);
118     label i = 0;
119     forAllIter(typename Cloud<ReactingParcel<ParcelType> >, c, iter)
120     {
121         ReactingParcel<ParcelType>& p = iter();
122         p.mass0_ = mass0[i++];
123     }
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();
132     }
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     }
142     // Populate Y for each parcel
143     forAll(phaseTypes, j)
144     {
145         IOField<scalar> Y
146         (
147             c.fieldIOobject
148             (
149                 "Y" + phaseTypes[j] + stateLabels[j],
150                  IOobject::MUST_READ
151             )
152         );
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     }
164 template<class ParcelType>
165 template<class CloudType>
166 void Foam::ReactingParcel<ParcelType>::writeFields(const CloudType& c)
168     ParcelType::writeFields(c);
172 template<class ParcelType>
173 template<class CloudType, class CompositionType>
174 void Foam::ReactingParcel<ParcelType>::writeFields
176     const CloudType& c,
177     const CompositionType& compModel
180     ParcelType::writeFields(c);
182     const label np = c.size();
184     if (np > 0)
185     {
186         IOField<scalar> mass0(c.fieldIOobject("mass0", IOobject::NO_READ), np);
188         label i = 0;
189         forAllConstIter(typename Cloud<ReactingParcel<ParcelType> >, c, iter)
190         {
191             const ReactingParcel<ParcelType>& p = iter();
192             mass0[i++] = p.mass0_;
193         }
194         mass0.write();
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();
202         }
204         forAll(phaseTypes, j)
205         {
206             IOField<scalar> Y
207             (
208                 c.fieldIOobject
209                 (
210                     "Y" + phaseTypes[j] + stateLabels[j],
211                     IOobject::NO_READ
212                 ),
213                 np
214             );
215             label i = 0;
216             forAllConstIter
217             (
218                 typename Cloud<ReactingParcel<ParcelType> >,
219                 c,
220                 iter
221             )
222             {
223                 const ReactingParcel<ParcelType>& p = iter();
224                 Y[i++] = p.Y()[j];
225             }
227             Y.write();
228         }
229     }
233 // * * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * //
235 template<class ParcelType>
236 Foam::Ostream& Foam::operator<<
238     Ostream& os,
239     const ReactingParcel<ParcelType>& p
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_),
254             sizeof(p.mass0())
255         );
256         os  << p.Y();
257     }
259     // Check state of Ostream
260     os.check
261     (
262         "Ostream& operator<<(Ostream&, const ReactingParcel<ParcelType>&)"
263     );
265     return os;
269 // ************************************************************************* //