DSMCParcelIO.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) 2016-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 "DSMCParcel.H"
30#include "IOstreams.H"
31#include "IOField.H"
32#include "Cloud.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36template<class ParcelType>
38(
39 sizeof(DSMCParcel<ParcelType>) - sizeof(ParcelType)
40);
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class ParcelType>
47(
48 const polyMesh& mesh,
49 Istream& is,
50 bool readFields,
51 bool newFormat
52)
53:
54 ParcelType(mesh, is, readFields, newFormat),
55 U_(Zero),
56 Ei_(0.0),
57 typeId_(-1)
58{
59 if (readFields)
60 {
61 if (is.format() == IOstream::ASCII)
62 {
63 is >> U_ >> Ei_ >> typeId_;
64 }
65 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
66 {
67 // Non-native label or scalar size
68
69 is.beginRawRead();
70
71 readRawScalar(is, U_.data(), vector::nComponents);
72 readRawScalar(is, &Ei_);
74
75 is.endRawRead();
76 }
77 else
78 {
79 is.read(reinterpret_cast<char*>(&U_), sizeofFields);
80 }
81 }
82
84}
85
86
87template<class ParcelType>
89{
90 bool valid = c.size();
91
93
94 IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ), valid);
95 c.checkFieldIOobject(c, U);
96
97 IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::MUST_READ), valid);
98 c.checkFieldIOobject(c, Ei);
99
100 IOField<label> typeId
101 (
102 c.fieldIOobject("typeId", IOobject::MUST_READ),
103 valid
104 );
105 c.checkFieldIOobject(c, typeId);
106
107 label i = 0;
108 for (DSMCParcel<ParcelType>& p : c)
109 {
110 p.U_ = U[i];
111 p.Ei_ = Ei[i];
112 p.typeId_ = typeId[i];
113 ++i;
114 }
115}
116
117
118template<class ParcelType>
120(
122)
123{
125
126 label np = c.size();
127
128 IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
129 IOField<scalar> Ei(c.fieldIOobject("Ei", IOobject::NO_READ), np);
130 IOField<label> typeId(c.fieldIOobject("typeId", IOobject::NO_READ), np);
131
132 label i = 0;
133 for (const DSMCParcel<ParcelType>& p : c)
134 {
135 U[i] = p.U();
136 Ei[i] = p.Ei();
137 typeId[i] = p.typeId();
138 ++i;
139 }
140
141 U.write(np > 0);
142 Ei.write(np > 0);
143 typeId.write(np > 0);
144}
145
146
147// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
148
149template<class ParcelType>
150Foam::Ostream& Foam::operator<<
151(
152 Ostream& os,
154)
155{
156 if (os.format() == IOstream::ASCII)
157 {
158 os << static_cast<const ParcelType& >(p)
159 << token::SPACE << p.U()
160 << token::SPACE << p.Ei()
161 << token::SPACE << p.typeId();
162 }
163 else
164 {
165 os << static_cast<const ParcelType& >(p);
166 os.write
167 (
168 reinterpret_cast<const char*>(&p.U_),
170 );
171 }
172
174 return os;
175}
176
177
178// ************************************************************************* //
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
DSMC parcel class.
Definition: DSMCParcel.H:73
label typeId_
Parcel type id.
Definition: DSMCParcel.H:148
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: DSMCParcel.H:77
static void readFields(Cloud< DSMCParcel< ParcelType > > &c)
Definition: DSMCParcelIO.C:88
vector U_
Velocity of Parcel [m/s].
Definition: DSMCParcel.H:141
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition: DSMCParcel.H:145
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
Cmpt * data() noexcept
Return pointer to the first data element.
Definition: VectorSpaceI.H:185
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
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
@ 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.
label readRawLabel(Istream &is)
Read raw label from binary stream.
Definition: label.C:46