MomentumParcelIO.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-2020 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 "MomentumParcel.H"
27 #include "IOstreams.H"
28 #include "IOField.H"
29 #include "Cloud.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<class ParcelType>
34 Foam::string Foam::MomentumParcel<ParcelType>::propertyList_ =
36 
37 template<class ParcelType>
39 (
40  sizeof(MomentumParcel<ParcelType>)
41  - offsetof(MomentumParcel<ParcelType>, moving_)
42 );
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class ParcelType>
49 (
50  const polyMesh& mesh,
51  Istream& is,
52  bool readFields
53 )
54 :
55  ParcelType(mesh, is, readFields),
56  moving_(false),
57  typeId_(0),
58  nParticle_(0.0),
59  d_(0.0),
60  dTarget_(0.0),
61  U_(Zero),
62  rho_(0.0),
63  age_(0.0),
64  tTurb_(0.0),
65  UTurb_(Zero)
66 {
67  if (readFields)
68  {
69  if (is.format() == IOstream::ASCII)
70  {
71  moving_ = readBool(is);
72  typeId_ = readLabel(is);
73  nParticle_ = readScalar(is);
74  d_ = readScalar(is);
75  dTarget_ = readScalar(is);
76  is >> U_;
77  rho_ = readScalar(is);
78  age_ = readScalar(is);
79  tTurb_ = readScalar(is);
80  is >> UTurb_;
81  }
82  else
83  {
84  is.read(reinterpret_cast<char*>(&moving_), sizeofFields_);
85  }
86  }
87 
88  // Check state of Istream
89  is.check
90  (
91  "MomentumParcel<ParcelType>::MomentumParcel"
92  "(const polyMesh&, Istream&, bool)"
93  );
94 }
95 
96 
97 template<class ParcelType>
98 template<class CloudType>
100 {
101  bool write = c.size();
102 
104 
105  IOField<label> moving
106  (
107  c.fieldIOobject("active", IOobject::MUST_READ),
108  write
109  );
110  c.checkFieldIOobject(c, moving);
111 
112  IOField<label> typeId
113  (
114  c.fieldIOobject("typeId", IOobject::MUST_READ),
115  write
116  );
117  c.checkFieldIOobject(c, typeId);
118 
119  IOField<scalar> nParticle
120  (
121  c.fieldIOobject("nParticle", IOobject::MUST_READ),
122  write
123  );
124  c.checkFieldIOobject(c, nParticle);
125 
127  (
128  c.fieldIOobject("d", IOobject::MUST_READ),
129  write
130  );
131  c.checkFieldIOobject(c, d);
132 
133  IOField<scalar> dTarget
134  (
135  c.fieldIOobject("dTarget", IOobject::MUST_READ),
136  write
137  );
138  c.checkFieldIOobject(c, dTarget);
139 
141  (
142  c.fieldIOobject("U", IOobject::MUST_READ),
143  write
144  );
145  c.checkFieldIOobject(c, U);
146 
148  (
149  c.fieldIOobject("rho", IOobject::MUST_READ),
150  write
151  );
152  c.checkFieldIOobject(c, rho);
153 
154  IOField<scalar> age
155  (
156  c.fieldIOobject("age", IOobject::MUST_READ),
157  write
158  );
159  c.checkFieldIOobject(c, age);
160 
161  IOField<scalar> tTurb
162  (
163  c.fieldIOobject("tTurb", IOobject::MUST_READ),
164  write
165  );
166  c.checkFieldIOobject(c, tTurb);
167 
168  IOField<vector> UTurb
169  (
170  c.fieldIOobject("UTurb", IOobject::MUST_READ),
171  write
172  );
173  c.checkFieldIOobject(c, UTurb);
174 
175  label i = 0;
176 
177  forAllIter(typename CloudType, c, iter)
178  {
179  MomentumParcel<ParcelType>& p = iter();
180 
181  p.moving_ = moving[i];
182  p.typeId_ = typeId[i];
183  p.nParticle_ = nParticle[i];
184  p.d_ = d[i];
185  p.dTarget_ = dTarget[i];
186  p.U_ = U[i];
187  p.rho_ = rho[i];
188  p.age_ = age[i];
189  p.tTurb_ = tTurb[i];
190  p.UTurb_ = UTurb[i];
191 
192  i++;
193  }
194 }
195 
196 
197 template<class ParcelType>
198 template<class CloudType>
200 {
201  ParcelType::writeFields(c);
202 
203  label np = c.size();
204 
205  IOField<label> moving(c.fieldIOobject("active", IOobject::NO_READ), np);
206  IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
207  IOField<scalar> nParticle
208  (
209  c.fieldIOobject("nParticle", IOobject::NO_READ),
210  np
211  );
212  IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
213  IOField<scalar> dTarget(c.fieldIOobject("dTarget", IOobject::NO_READ), np);
214  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
215  IOField<scalar> rho(c.fieldIOobject("rho", IOobject::NO_READ), np);
216  IOField<scalar> age(c.fieldIOobject("age", IOobject::NO_READ), np);
217  IOField<scalar> tTurb(c.fieldIOobject("tTurb", IOobject::NO_READ), np);
218  IOField<vector> UTurb(c.fieldIOobject("UTurb", IOobject::NO_READ), np);
219 
220  label i = 0;
221 
222  forAllConstIter(typename CloudType, c, iter)
223  {
224  const MomentumParcel<ParcelType>& p = iter();
225 
226  moving[i] = p.moving();
227  typeId[i] = p.typeId();
228  nParticle[i] = p.nParticle();
229  d[i] = p.d();
230  dTarget[i] = p.dTarget();
231  U[i] = p.U();
232  rho[i] = p.rho();
233  age[i] = p.age();
234  tTurb[i] = p.tTurb();
235  UTurb[i] = p.UTurb();
236 
237  i++;
238  }
239 
240  const bool write = np > 0;
241 
242  moving.write(write);
243  typeId.write(write);
244  nParticle.write(write);
245  d.write(write);
246  dTarget.write(write);
247  U.write(write);
248  rho.write(write);
249  age.write(write);
250  tTurb.write(write);
251  UTurb.write(write);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
256 
257 template<class ParcelType>
258 Foam::Ostream& Foam::operator<<
259 (
260  Ostream& os,
262 )
263 {
264  if (os.format() == IOstream::ASCII)
265  {
266  os << static_cast<const ParcelType&>(p)
267  << token::SPACE << p.moving()
268  << token::SPACE << p.typeId()
269  << token::SPACE << p.nParticle()
270  << token::SPACE << p.d()
271  << token::SPACE << p.dTarget()
272  << token::SPACE << p.U()
273  << token::SPACE << p.rho()
274  << token::SPACE << p.age()
275  << token::SPACE << p.tTurb()
276  << token::SPACE << p.UTurb();
277  }
278  else
279  {
280  os << static_cast<const ParcelType&>(p);
281  os.write
282  (
283  reinterpret_cast<const char*>(&p.moving_),
285  );
286  }
287 
288  // Check state of Ostream
289  os.check
290  (
291  "Ostream& operator<<(Ostream&, const MomentumParcel<ParcelType>&)"
292  );
293 
294  return os;
295 }
296 
297 
298 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
const vector & UTurb() const
Return const access to turbulent velocity fluctuation.
bool moving_
Moving flag - tracking stopped when moving = false.
scalar d_
Diameter [m].
scalar dTarget_
Target diameter [m].
U
Definition: pEqn.H:72
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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
scalar age() const
Return const access to the age.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label typeId_
Parcel type id.
label size() const
Return the number of particles in the cloud.
Definition: Cloud.H:190
bool readBool(Istream &)
Definition: boolIO.C:60
scalar dTarget() const
Return const access to target diameter.
bool moving() const
Return const access to moving flag.
const dimensionedScalar c
Speed of light in a vacuum.
scalar d() const
Return const access to diameter.
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
scalar rho_
Density [kg/m^3].
scalar nParticle_
Number of particles in Parcel.
static void readFields(TrackCloudType &c)
Read.
scalar nParticle() const
Return const access to number of particles.
virtual Istream & read(token &)=0
Return next token from stream.
scalar age_
Age [s].
scalar tTurb_
Time spent in turbulent eddy [s].
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
label readLabel(Istream &is)
Definition: label.H:64
vector U_
Velocity of Parcel [m/s].
Momentum parcel class with rotational motion (as spherical particles only) and one/two-way coupling w...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
scalar rho() const
Return const access to density.
MomentumParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
static void writeFields(const TrackCloudType &c)
Write.
const vector & U() const
Return const access to velocity.
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:190
label typeId() const
Return const access to type id.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
volScalarField & p
scalar tTurb() const
Return const access to time spent in turbulent eddy.
vector UTurb_
Turbulent velocity fluctuation [m/s].
A class for handling character strings derived from std::string.
Definition: string.H:76
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:170
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50