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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
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() == IOstream::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
76}
77
78
80{
81 bool valid = c.size();
82
84
85 IOField<scalar> d(c.fieldIOobject("d", IOobject::MUST_READ), valid);
86 c.checkFieldIOobject(c, d);
87
88 IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ), valid);
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
102{
104
105 const label np = c.size();
106
107 IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
108 IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
109
110 label i = 0;
111 for (const solidParticle& p : c)
112 {
113 d[i] = p.d_;
114 U[i] = p.U_;
115 ++i;
116 }
117
118 d.write(np > 0);
119 U.write(np > 0);
120}
121
122
123// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
124
126{
127 if (os.format() == IOstream::ASCII)
128 {
129 os << static_cast<const particle&>(p)
130 << token::SPACE << p.d_
131 << token::SPACE << p.U_;
132 }
133 else
134 {
135 os << static_cast<const particle&>(p);
136 os.write
137 (
138 reinterpret_cast<const char*>(&p.d_),
140 );
141 }
142
144 return os;
145}
146
147
148// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Base cloud calls templated on particle type.
Definition: Cloud.H:68
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
streamFormat format() const noexcept
Get the current stream format.
std::enable_if< std::is_integral< T >::value, bool >::type checkLabelSize() const noexcept
Definition: IOstream.H:300
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
std::enable_if< std::is_floating_point< T >::value, bool >::type checkScalarSize() const noexcept
Definition: IOstream.H:309
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
virtual bool endRawRead()=0
End of low-level raw binary read.
virtual bool beginRawRead()=0
Start of low-level raw binary read.
virtual Istream & read(token &)=0
Return next token from stream.
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Base particle class.
Definition: particle.H:79
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:68
static const std::size_t sizeofFields
Size in bytes of the fields.
static void readFields(Cloud< solidParticle > &c)
@ SPACE
Space [isspace].
Definition: token.H:125
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83