Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / src / lagrangian / dsmc / parcels / Templates / DsmcParcel / DsmcParcelIO.C
blob7ec34aad080ab6de1c8c769ee05b61c82cfb2666
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 "DsmcParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 #include "Cloud.H"
31 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
33 template <class ParcelType>
34 Foam::DsmcParcel<ParcelType>::DsmcParcel
36     const polyMesh& mesh,
37     Istream& is,
38     bool readFields
41     ParcelType(mesh, is, readFields),
42     U_(vector::zero),
43     Ei_(0.0),
44     typeId_(-1)
46     if (readFields)
47     {
48         if (is.format() == IOstream::ASCII)
49         {
50             is >> U_;
51             Ei_ = readScalar(is);
52             typeId_ = readLabel(is);
53         }
54         else
55         {
56             is.read
57             (
58                 reinterpret_cast<char*>(&U_),
59                 sizeof(U_)
60               + sizeof(Ei_)
61               + sizeof(typeId_)
62             );
63         }
64     }
66     // Check state of Istream
67     is.check
68     (
69         "DsmcParcel<ParcelType>::DsmcParcel"
70         "(const Cloud<ParcelType>&, Istream&, bool)"
71     );
75 template<class ParcelType>
76 void Foam::DsmcParcel<ParcelType>::readFields(Cloud<DsmcParcel<ParcelType> >& c)
78     if (!c.size())
79     {
80         return;
81     }
83     ParcelType::readFields(c);
85     IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ));
86     c.checkFieldIOobject(c, U);
88     IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::MUST_READ));
89     c.checkFieldIOobject(c, Ei);
91     IOField<label> typeId(c.fieldIOobject("typeId", IOobject::MUST_READ));
92     c.checkFieldIOobject(c, typeId);
94     label i = 0;
95     forAllIter(typename Cloud<DsmcParcel<ParcelType> >, c, iter)
96     {
97         DsmcParcel<ParcelType>& p = iter();
99         p.U_ = U[i];
100         p.Ei_ = Ei[i];
101         p.typeId_ = typeId[i];
102         i++;
103     }
107 template<class ParcelType>
108 void Foam::DsmcParcel<ParcelType>::writeFields
110     const Cloud<DsmcParcel<ParcelType> >& c
113     ParcelType::writeFields(c);
115     label np =  c.size();
117     IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
118     IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::NO_READ), np);
119     IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
121     label i = 0;
122     forAllConstIter(typename Cloud<DsmcParcel<ParcelType> >, c, iter)
123     {
124         const DsmcParcel<ParcelType>& p = iter();
126         U[i] = p.U();
127         Ei[i] = p.Ei();
128         typeId[i] = p.typeId();
129         i++;
130     }
132     U.write();
133     Ei.write();
134     typeId.write();
138 // * * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * //
140 template<class ParcelType>
141 Foam::Ostream& Foam::operator<<
143     Ostream& os,
144     const DsmcParcel<ParcelType>& p
147     if (os.format() == IOstream::ASCII)
148     {
149         os  << static_cast<const ParcelType& >(p)
150             << token::SPACE << p.U()
151             << token::SPACE << p.Ei()
152             << token::SPACE << p.typeId();
153     }
154     else
155     {
156         os  << static_cast<const ParcelType& >(p);
157         os.write
158         (
159             reinterpret_cast<const char*>(&p.U_),
160             sizeof(p.U())
161           + sizeof(p.Ei())
162           + sizeof(p.typeId())
163         );
164     }
166     // Check state of Ostream
167     os.check
168     (
169         "Ostream& operator<<(Ostream&, const DsmcParcel<ParcelType>&)"
170     );
172     return os;
176 // ************************************************************************* //