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-2018 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 "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 const std::size_t Foam::solidParticle::sizeofFields_
32 (
33  sizeof(solidParticle) - sizeof(particle)
34 );
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const polyMesh& mesh,
42  Istream& is,
43  bool readFields
44 )
45 :
46  particle(mesh, is, readFields)
47 {
48  if (readFields)
49  {
50  if (is.format() == IOstream::ASCII)
51  {
52  d_ = readScalar(is);
53  is >> U_;
54  }
55  else
56  {
57  is.read(reinterpret_cast<char*>(&d_), sizeofFields_);
58  }
59  }
60 
61  // Check state of Istream
62  is.check("solidParticle::solidParticle(Istream&)");
63 }
64 
65 
67 {
68  bool valid = c.size();
69 
71 
73  c.checkFieldIOobject(c, d);
74 
76  c.checkFieldIOobject(c, U);
77 
78  label i = 0;
80  {
81  solidParticle& p = iter();
82 
83  p.d_ = d[i];
84  p.U_ = U[i];
85  i++;
86  }
87 }
88 
89 
91 {
93 
94  label np = c.size();
95 
98 
99  label i = 0;
101  {
102  const solidParticle& p = iter();
103 
104  d[i] = p.d_;
105  U[i] = p.U_;
106  i++;
107  }
108 
109  d.write(np > 0);
110  U.write(np > 0);
111 }
112 
113 
114 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
115 
117 {
118  if (os.format() == IOstream::ASCII)
119  {
120  os << static_cast<const particle&>(p)
121  << token::SPACE << p.d_
122  << token::SPACE << p.U_;
123  }
124  else
125  {
126  os << static_cast<const particle&>(p);
127  os.write
128  (
129  reinterpret_cast<const char*>(&p.d_),
130  solidParticle::sizeofFields_
131  );
132  }
133 
134  // Check state of Ostream
135  os.check("Ostream& operator<<(Ostream&, const solidParticle&)");
136 
137  return os;
138 }
139 
140 
141 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
virtual Ostream & write(const char)=0
Write character.
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
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
#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 size() const
Return the number of particles in the cloud.
Definition: Cloud.H:148
scalar d() const
Return diameter.
const dimensionedScalar c
Speed of light in a vacuum.
Base particle class.
Definition: particle.H:83
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
const vector & U() const
Return velocity.
virtual Istream & read(token &)=0
Return next token from stream.
streamFormat format() const
Return current stream format.
Definition: IOstream.H:374
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static void writeFields(const Cloud< solidParticle > &c)
solidParticle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const scalar d, const vector &U)
Construct from components.
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.
Ostream & operator<<(Ostream &, const ensightPart &)
void checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:187
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:167
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
static void readFields(Cloud< solidParticle > &c)
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:66