DSMCParcelIO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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  ParcelType(is, readFields),
46  U_(Zero),
47  Ei_(0.0),
48  typeId_(-1)
49 {
50  if (readFields)
51  {
52  if (is.format() == IOstream::ASCII)
53  {
54  is >> U_;
55  Ei_ = readScalar(is);
56  typeId_ = readLabel(is);
57  }
58  else
59  {
60  is.read(reinterpret_cast<char*>(&U_), sizeofFields_);
61  }
62  }
63 
64  // Check state of Istream
65  is.check
66  (
67  "DSMCParcel<ParcelType>::DSMCParcel"
68  "(const Cloud<ParcelType>&, Istream&, bool)"
69  );
70 }
71 
72 
73 template<class ParcelType>
75 {
76  bool valid = c.size();
77 
79 
80  IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ), valid);
81  c.checkFieldIOobject(c, U);
82 
83  IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::MUST_READ), valid);
84  c.checkFieldIOobject(c, Ei);
85 
86  IOField<label> typeId
87  (
88  c.fieldIOobject("typeId", IOobject::MUST_READ),
89  valid
90  );
91  c.checkFieldIOobject(c, typeId);
92 
93  label i = 0;
94  forAllIter(typename Cloud<DSMCParcel<ParcelType>>, c, iter)
95  {
96  DSMCParcel<ParcelType>& p = iter();
97 
98  p.U_ = U[i];
99  p.Ei_ = Ei[i];
100  p.typeId_ = typeId[i];
101  i++;
102  }
103 }
104 
105 
106 template<class ParcelType>
108 (
110 )
111 {
112  ParcelType::writeFields(c);
113 
114  label np = c.size();
115 
116  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
117  IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::NO_READ), np);
118  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
119 
120  label i = 0;
122  {
123  const DSMCParcel<ParcelType>& p = iter();
124 
125  U[i] = p.U();
126  Ei[i] = p.Ei();
127  typeId[i] = p.typeId();
128  i++;
129  }
130 
131  U.write(np > 0);
132  Ei.write(np > 0);
133  typeId.write(np > 0);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
138 
139 template<class ParcelType>
140 Foam::Ostream& Foam::operator<<
141 (
142  Ostream& os,
144 )
145 {
146  if (os.format() == IOstream::ASCII)
147  {
148  os << static_cast<const ParcelType& >(p)
149  << token::SPACE << p.U()
150  << token::SPACE << p.Ei()
151  << token::SPACE << p.typeId();
152  }
153  else
154  {
155  os << static_cast<const ParcelType& >(p);
156  os.write
157  (
158  reinterpret_cast<const char*>(&p.U_),
159  DSMCParcel<ParcelType>::sizeofFields_
160  );
161  }
162 
163  // Check state of Ostream
164  os.check
165  (
166  "Ostream& operator<<(Ostream&, const DSMCParcel<ParcelType>&)"
167  );
168 
169  return os;
170 }
171 
172 
173 // ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Base cloud calls templated on particle type.
Definition: Cloud.H:74
DSMC parcel class.
Definition: DSMCParcel.H:70
static void readFields(Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:74
label typeId_
Parcel type id.
Definition: DSMCParcel.H:147
DSMCParcel(const polyMesh &mesh, const vector &position, const label celli, const vector &U, const scalar Ei, const label typeId)
Construct from a position and a cell, searching for the rest of the.
Definition: DSMCParcelI.H:57
static void writeFields(const Cloud< DSMCParcel< ParcelType >> &c)
Definition: DSMCParcelIO.C:108
vector U_
Velocity of Parcel [m/s].
Definition: DSMCParcel.H:140
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition: DSMCParcel.H:144
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
virtual Istream & read(token &)=0
Return next token from stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual bool write(const bool write=true) const
Write using setting from DB.
U
Definition: pEqn.H:72
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
static const zero Zero
Definition: zero.H:97
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
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
label readLabel(Istream &is)
Definition: label.H:64
volScalarField & p