Loading...
Searching...
No Matches
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-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 "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() == IOstreamOption::ASCII)
62 {
63 is >> U_ >> Ei_ >> typeId_;
64 }
65 else
66 {
67 // No non-native streaming
69
70 is.read(reinterpret_cast<char*>(&U_), sizeofFields);
71 }
72 }
73
75}
76
77
78template<class ParcelType>
79void Foam::DSMCParcel<ParcelType>::readFields(Cloud<DSMCParcel<ParcelType>>& c)
80{
81 const bool readOnProc = c.size();
82
83 ParcelType::readFields(c);
84
85 IOField<vector> U(c.newIOobject("U", IOobject::MUST_READ), readOnProc);
86 c.checkFieldIOobject(c, U);
87
88 IOField<scalar> Ei(c.newIOobject("Ei", IOobject::MUST_READ), readOnProc);
89 c.checkFieldIOobject(c, Ei);
90
91 IOField<label> typeId
92 (
93 c.newIOobject("typeId", IOobject::MUST_READ),
94 readOnProc
95 );
96 c.checkFieldIOobject(c, typeId);
97
98 label i = 0;
99 for (DSMCParcel<ParcelType>& p : c)
100 {
101 p.U_ = U[i];
102 p.Ei_ = Ei[i];
103 p.typeId_ = typeId[i];
104 ++i;
105 }
106}
107
108
109template<class ParcelType>
111(
112 const Cloud<DSMCParcel<ParcelType>>& c
113)
114{
115 ParcelType::writeFields(c);
116
117 const label np = c.size();
118 const bool writeOnProc = c.size();
119
120 IOField<vector> U(c.newIOobject("U", IOobject::NO_READ), np);
121 IOField<scalar> Ei(c.newIOobject("Ei", IOobject::NO_READ), np);
122 IOField<label> typeId(c.newIOobject("typeId", IOobject::NO_READ), np);
123
124 label i = 0;
125 for (const DSMCParcel<ParcelType>& p : c)
126 {
127 U[i] = p.U();
128 Ei[i] = p.Ei();
129 typeId[i] = p.typeId();
130 ++i;
131 }
132
133 U.write(writeOnProc);
134 Ei.write(writeOnProc);
135 typeId.write(writeOnProc);
136}
137
138
139// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
140
141template<class ParcelType>
142Foam::Ostream& Foam::operator<<
143(
144 Ostream& os,
146)
147{
149 {
151 << token::SPACE << p.U()
152 << token::SPACE << p.Ei()
153 << token::SPACE << p.typeId();
154 }
155 else
156 {
158 os.write
159 (
160 reinterpret_cast<const char*>(&p.U_),
162 );
163 }
164
166 return os;
167}
168
169
170// ************************************************************************* //
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
DSMC parcel class.
Definition DSMCParcel.H:68
label typeId_
Parcel type id.
Definition DSMCParcel.H:179
label typeId() const
Return type id.
const vector & U() const
Return const access to velocity.
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition DSMCParcel.H:74
scalar Ei() const
Return const access to internal energy.
static void readFields(Cloud< DSMCParcel< ParcelType > > &c)
static void writeFields(const Cloud< DSMCParcel< ParcelType > > &c)
DSMCParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const vector &U, const scalar Ei, const label typeId)
Construct from components.
Definition DSMCParcelI.H:48
vector U_
Velocity of Parcel [m/s].
Definition DSMCParcel.H:167
scalar Ei_
Internal energy of the Parcel, covering all non-translational.
Definition DSMCParcel.H:174
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
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
@ SPACE
Space [isspace].
Definition token.H:144
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127