Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
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
59 {
60 // No non-native streaming
62
63 is.read(reinterpret_cast<char*>(&d_), sizeofFields);
64 }
65 }
66
68}
69
70
72{
73 const bool readOnProc = c.size();
74
76
77 IOField<scalar> d(c.newIOobject("d", IOobject::MUST_READ), readOnProc);
78 c.checkFieldIOobject(c, d);
79
80 IOField<vector> U(c.newIOobject("U", IOobject::MUST_READ), readOnProc);
81 c.checkFieldIOobject(c, U);
82
83 label i = 0;
84 for (solidParticle& p : c)
85 {
86 p.d_ = d[i];
87 p.U_ = U[i];
88 ++i;
89 }
90}
91
92
93void Foam::solidParticle::writeFields(const Cloud<solidParticle>& c)
94{
96
97 const label np = c.size();
98 const bool writeOnProc = c.size();
99
100 IOField<scalar> d(c.newIOobject("d", IOobject::NO_READ), np);
101 IOField<vector> U(c.newIOobject("U", IOobject::NO_READ), np);
102
103 label i = 0;
104 for (const solidParticle& p : c)
105 {
106 d[i] = p.d_;
107 U[i] = p.U_;
108 ++i;
109 }
110
111 d.write(writeOnProc);
112 U.write(writeOnProc);
113}
114
115
116// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
117
118Foam::Ostream& Foam::operator<<(Ostream& os, const solidParticle& p)
119{
120 if (os.format() == IOstreamOption::ASCII)
121 {
123 << token::SPACE << p.d_
124 << token::SPACE << p.U_;
125 }
126 else
127 {
129 os.write
130 (
131 reinterpret_cast<const char*>(&p.d_),
133 );
134 }
135
137 return os;
138}
139
140
141// ************************************************************************* //
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:64
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
streamFormat format() const noexcept
Get the current stream format.
@ ASCII
"ascii" (normal default)
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
bool fatalCheckNativeSizes(const char *operation) const
Assert that the label/scalar byte-size associated with the stream are the native label/scalar sizes.
Definition IOstream.C:67
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.
virtual Ostream & write(const char c) override
Write character.
Definition OBJstream.C:69
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Simple solid spherical particle class with one-way coupling with the continuous phase.
static const std::size_t sizeofFields
Size in bytes of the fields.
const vector & U() const noexcept
Return velocity.
scalar d() const noexcept
Return diameter.
solidParticle(const polyMesh &mesh, const vector &position, const label celli=-1)
Construct from a position and a cell.
static void readFields(Cloud< solidParticle > &c)
static void writeFields(const Cloud< solidParticle > &c)
@ SPACE
Space [isspace].
Definition token.H:144
U
Definition pEqn.H:72
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
const dimensionedScalar c
Speed of light in a vacuum.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).