solidParticleIO.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 "solidParticle.H"
27 #include "solidParticleCloud.H"
28 #include "IOstreams.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 const std::size_t Foam::solidParticle::sizeofFields_
33 (
34  sizeof(solidParticle) - sizeof(particle)
35 );
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 :
42  particle(is, readFields)
43 {
44  if (readFields)
45  {
46  if (is.format() == IOstream::ASCII)
47  {
48  d_ = readScalar(is);
49  is >> U_;
50  }
51  else
52  {
53  is.read(reinterpret_cast<char*>(&d_), sizeofFields_);
54  }
55  }
56 
57  // Check state of Istream
58  is.check("solidParticle::solidParticle(Istream&)");
59 }
60 
61 
63 {
64  bool valid = c.size();
65 
67 
68  IOField<scalar> d(c.fieldIOobject("d", IOobject::MUST_READ), valid);
69  c.checkFieldIOobject(c, d);
70 
71  IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ), valid);
72  c.checkFieldIOobject(c, U);
73 
74  label i = 0;
76  {
77  solidParticle& p = iter();
78 
79  p.d_ = d[i];
80  p.U_ = U[i];
81  i++;
82  }
83 }
84 
85 
87 {
89 
90  label np = c.size();
91 
92  IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
93  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
94 
95  label i = 0;
97  {
98  const solidParticle& p = iter();
99 
100  d[i] = p.d_;
101  U[i] = p.U_;
102  i++;
103  }
104 
105  d.write(np > 0);
106  U.write(np > 0);
107 }
108 
109 
110 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
111 
113 {
114  if (os.format() == IOstream::ASCII)
115  {
116  os << static_cast<const particle&>(p)
117  << token::SPACE << p.d_
118  << token::SPACE << p.U_;
119  }
120  else
121  {
122  os << static_cast<const particle&>(p);
123  os.write
124  (
125  reinterpret_cast<const char*>(&p.d_),
126  solidParticle::sizeofFields_
127  );
128  }
129 
130  // Check state of Ostream
131  os.check("Ostream& operator<<(Ostream&, const solidParticle&)");
132 
133  return os;
134 }
135 
136 
137 // ************************************************************************* //
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
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 Ostream & write(const char)=0
Write character.
Base particle class.
Definition: particle.H:83
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
virtual bool write(const bool write=true) const
Write using setting from DB.
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:69
solidParticle(Istream &is, bool readFields=true)
Construct from Istream.
static void readFields(Cloud< solidParticle > &c)
static void writeFields(const Cloud< solidParticle > &c)
U
Definition: pEqn.H:72
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
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
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
volScalarField & p