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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "solidParticle.H"
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 const std::size_t Foam::solidParticle::sizeofFields
35 (
36  sizeof(solidParticle) - sizeof(particle)
37 );
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const polyMesh& mesh,
45  Istream& is,
46  bool readFields,
47  bool newFormat
48 )
49 :
50  particle(mesh, is, readFields, newFormat)
51 {
52  if (readFields)
53  {
54  if (is.format() == IOstreamOption::ASCII)
55  {
56  is >> d_ >> U_;
57  }
58  else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
59  {
60  // Non-native label or scalar size
61 
62  is.beginRawRead();
63 
64  readRawScalar(is, &d_);
65  readRawScalar(is, U_.data(), vector::nComponents);
66 
67  is.endRawRead();
68  }
69  else
70  {
71  is.read(reinterpret_cast<char*>(&d_), sizeofFields);
72  }
73  }
74 
75  is.check(FUNCTION_NAME);
76 }
77 
78 
79 void Foam::solidParticle::readFields(Cloud<solidParticle>& c)
80 {
81  const bool readOnProc = c.size();
82 
84 
85  IOField<scalar> d(c.fieldIOobject("d", IOobject::MUST_READ), readOnProc);
86  c.checkFieldIOobject(c, d);
87 
88  IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ), readOnProc);
89  c.checkFieldIOobject(c, U);
90 
91  label i = 0;
92  for (solidParticle& p : c)
93  {
94  p.d_ = d[i];
95  p.U_ = U[i];
96  ++i;
97  }
98 }
99 
100 
101 void Foam::solidParticle::writeFields(const Cloud<solidParticle>& c)
102 {
104 
105  const label np = c.size();
106  const bool writeOnProc = c.size();
107 
108  IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
109  IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
110 
111  label i = 0;
112  for (const solidParticle& p : c)
113  {
114  d[i] = p.d_;
115  U[i] = p.U_;
116  ++i;
117  }
118 
119  d.write(writeOnProc);
120  U.write(writeOnProc);
121 }
122 
123 
124 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
125 
126 Foam::Ostream& Foam::operator<<(Ostream& os, const solidParticle& p)
127 {
129  {
130  os << static_cast<const particle&>(p)
131  << token::SPACE << p.d_
132  << token::SPACE << p.U_;
133  }
134  else
135  {
136  os << static_cast<const particle&>(p);
137  os.write
138  (
139  reinterpret_cast<const char*>(&p.d_),
141  );
142  }
143 
145  return os;
146 }
147 
148 
149 // ************************************************************************* //
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Check if the scalar byte-size associated with the stream is the same as the given type...
Definition: IOstream.H:379
virtual Ostream & write(const char c) override
Write character.
Definition: OBJstream.C:69
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
"ascii" (normal default)
virtual bool endRawRead()=0
End of low-level raw binary read.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Base particle class.
Definition: particle.H:69
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.
dynamicFvMesh & mesh
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
Space [isspace].
Definition: token.H:131
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
static void writeFields(const Cloud< solidParticle > &c)
solidParticle(const polyMesh &mesh, const vector &position, const label celli=-1)
Construct from a position and a cell.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
U
Definition: pEqn.H:72
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
const dimensionedScalar c
Speed of light in a vacuum.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Nothing to be read.
static const std::size_t sizeofFields
Size in bytes of the fields.
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Check if the label byte-size associated with the stream is the same as the given type.
Definition: IOstream.H:368
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
volScalarField & p
streamFormat format() const noexcept
Get the current stream format.
static void readFields(Cloud< solidParticle > &c)
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:60