solidParticleIO.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-2015 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  if (!c.size())
69  {
70  return;
71  }
72 
74 
76  c.checkFieldIOobject(c, d);
77 
79  c.checkFieldIOobject(c, U);
80 
81  label i = 0;
83  {
84  solidParticle& p = iter();
85 
86  p.d_ = d[i];
87  p.U_ = U[i];
88  i++;
89  }
90 }
91 
92 
94 {
96 
97  label np = c.size();
98 
101 
102  label i = 0;
104  {
105  const solidParticle& p = iter();
106 
107  d[i] = p.d_;
108  U[i] = p.U_;
109  i++;
110  }
111 
112  d.write();
113  U.write();
114 }
115 
116 
117 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
118 
120 {
121  if (os.format() == IOstream::ASCII)
122  {
123  os << static_cast<const particle&>(p)
124  << token::SPACE << p.d_
125  << token::SPACE << p.U_;
126  }
127  else
128  {
129  os << static_cast<const particle&>(p);
130  os.write
131  (
132  reinterpret_cast<const char*>(&p.d_),
133  solidParticle::sizeofFields_
134  );
135  }
136 
137  // Check state of Ostream
138  os.check("Ostream& operator<<(Ostream&, const solidParticle&)");
139 
140  return os;
141 }
142 
143 
144 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
streamFormat format() const
Return current stream format.
Definition: IOstream.H:377
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 checkFieldIOobject(const Cloud< ParticleType > &c, const IOField< DataType > &data) const
Check lagrangian data field.
Definition: CloudIO.C:215
#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
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
Base particle class.
Definition: particle.H:78
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
virtual Istream & read(token &)=0
Return next token from stream.
static void writeFields(const CloudType &c)
Write the fields associated with the owner cloud.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
solidParticle(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI, const scalar d, const vector &U)
Construct from components.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static void writeFields(const Cloud< solidParticle > &c)
label size() const
Definition: Cloud.H:175
scalar d() const
Return diameter.
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const ensightPart &)
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject fieldIOobject(const word &fieldName, const IOobject::readOption r) const
Helper to construct IOobject for field and current time.
Definition: CloudIO.C:195
volScalarField & p
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const vector & U() const
Return velocity.
static void readFields(Cloud< solidParticle > &c)
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:65