DSMCParcelIO.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-2016 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 "DSMCParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 #include "Cloud.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<class ParcelType>
35 (
36  sizeof(DSMCParcel<ParcelType>) - sizeof(ParcelType)
37 );
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class ParcelType>
44 (
45  const polyMesh& mesh,
46  Istream& is,
47  bool readFields
48 )
49 :
50  ParcelType(mesh, is, readFields),
51  U_(Zero),
52  Ei_(0.0),
53  typeId_(-1)
54 {
55  if (readFields)
56  {
57  if (is.format() == IOstream::ASCII)
58  {
59  is >> U_;
60  Ei_ = readScalar(is);
61  typeId_ = readLabel(is);
62  }
63  else
64  {
65  is.read(reinterpret_cast<char*>(&U_), sizeofFields_);
66  }
67  }
68 
69  // Check state of Istream
70  is.check
71  (
72  "DSMCParcel<ParcelType>::DSMCParcel"
73  "(const Cloud<ParcelType>&, Istream&, bool)"
74  );
75 }
76 
77 
78 template<class ParcelType>
80 {
81  if (!c.size())
82  {
83  return;
84  }
85 
87 
88  IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ));
89  c.checkFieldIOobject(c, U);
90 
91  IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::MUST_READ));
92  c.checkFieldIOobject(c, Ei);
93 
94  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::MUST_READ));
95  c.checkFieldIOobject(c, typeId);
96 
97  label i = 0;
98  forAllIter(typename Cloud<DSMCParcel<ParcelType>>, c, iter)
99  {
100  DSMCParcel<ParcelType>& p = iter();
101 
102  p.U_ = U[i];
103  p.Ei_ = Ei[i];
104  p.typeId_ = typeId[i];
105  i++;
106  }
107 }
108 
109 
110 template<class ParcelType>
112 (
114 )
115 {
116  ParcelType::writeFields(c);
117 
118  label np = c.size();
119 
120  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
121  IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::NO_READ), np);
122  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
123 
124  label i = 0;
126  {
127  const DSMCParcel<ParcelType>& p = iter();
128 
129  U[i] = p.U();
130  Ei[i] = p.Ei();
131  typeId[i] = p.typeId();
132  i++;
133  }
134 
135  U.write();
136  Ei.write();
137  typeId.write();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
142 
143 template<class ParcelType>
144 Foam::Ostream& Foam::operator<<
145 (
146  Ostream& os,
148 )
149 {
150  if (os.format() == IOstream::ASCII)
151  {
152  os << static_cast<const ParcelType& >(p)
153  << token::SPACE << p.U()
154  << token::SPACE << p.Ei()
155  << token::SPACE << p.typeId();
156  }
157  else
158  {
159  os << static_cast<const ParcelType& >(p);
160  os.write
161  (
162  reinterpret_cast<const char*>(&p.U_),
164  );
165  }
166 
167  // Check state of Ostream
168  os.check
169  (
170  "Ostream& operator<<(Ostream&, const DSMCParcel<ParcelType>&)"
171  );
172 
173  return os;
174 }
175 
176 
177 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
U
Definition: pEqn.H:83
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
#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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
label typeId_
Parcel type id.
Definition: DSMCParcel.H:163
DSMC parcel class.
Definition: DSMCParcel.H:51
DSMCParcel(const polyMesh &mesh, const vector &position, const vector &U, const scalar Ei, const label celli, const label tetFacei, const label tetPtI, const label typeId)
Construct from components.
Definition: DSMCParcelI.H:56
label typeId() const
Return type id.
Definition: DSMCParcelI.H:119
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
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
static void readFields(Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:79
static void writeFields(const Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:112
virtual Istream & read(token &)=0
Return next token from stream.
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
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
const vector & U() const
Return const access to velocity.
Definition: DSMCParcelI.H:126
scalar Ei() const
Return const access to internal energy.
Definition: DSMCParcelI.H:133
const dimensionedScalar c
Speed of light in a vacuum.
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition: DSMCParcel.H:160
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
vector U_
Velocity of Parcel [m/s].
Definition: DSMCParcel.H:156
volScalarField & p