MPPICParcelIO.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) 2013-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 "MPPICParcel.H"
30#include "IOstreams.H"
31#include "IOField.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35template<class ParcelType>
38
39
40template<class ParcelType>
42(
43 sizeof(MPPICParcel<ParcelType>) - sizeof(ParcelType)
44);
45
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
49template<class ParcelType>
51(
52 const polyMesh& mesh,
53 Istream& is,
54 bool readFields,
55 bool newFormat
56)
57:
58 ParcelType(mesh, is, readFields, newFormat),
59 UCorrect_(Zero)
60{
61 if (readFields)
62 {
63 if (is.format() == IOstream::ASCII)
64 {
65 is >> UCorrect_;
66 }
67 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
68 {
69 // Non-native label or scalar size
70
71 is.beginRawRead();
72
73 readRawScalar(is, UCorrect_.data(), vector::nComponents);
74
75 is.endRawRead();
76 }
77 else
78 {
79 is.read(reinterpret_cast<char*>(&UCorrect_), sizeofFields);
80 }
81 }
82
84}
85
86
87template<class ParcelType>
88template<class CloudType>
90{
91 bool valid = c.size();
92
94
95 IOField<vector> UCorrect
96 (
97 c.fieldIOobject("UCorrect", IOobject::MUST_READ),
98 valid
99 );
100 c.checkFieldIOobject(c, UCorrect);
101
102 label i = 0;
103 for (MPPICParcel<ParcelType>& p : c)
104 {
105 p.UCorrect_ = UCorrect[i];
106
107 ++i;
108 }
109}
110
111
112template<class ParcelType>
113template<class CloudType>
115{
117
118 const label np = c.size();
119
121 UCorrect(c.fieldIOobject("UCorrect", IOobject::NO_READ), np);
122
123 label i = 0;
124
125 for (const MPPICParcel<ParcelType>& p : c)
126 {
127 UCorrect[i] = p.UCorrect();
128
129 ++i;
130 }
131
132 UCorrect.write(np > 0);
133}
134
135
136template<class ParcelType>
138(
139 Ostream& os,
140 const wordRes& filters,
141 const word& delim,
142 const bool namesOnly
143) const
144{
145 ParcelType::writeProperties(os, filters, delim, namesOnly);
146
147 #undef writeProp
148 #define writeProp(Name, Value) \
149 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
150
151 writeProp("UCorrect", UCorrect_);
152
153 #undef writeProp
154}
155
156
157template<class ParcelType>
158template<class CloudType>
160(
161 CloudType& c,
162 const objectRegistry& obr
163)
164{
166
167 if (!c.size()) return;
168
169 const auto& UCorrect = cloud::lookupIOField<vector>("UCorrect", obr);
170
171 label i = 0;
172 for (MPPICParcel<ParcelType>& p : c)
173 {
174 p.UCorrect() = UCorrect[i];
175
176 ++i;
177 }
178}
179
180
181template<class ParcelType>
182template<class CloudType>
184(
185 const CloudType& c,
186 objectRegistry& obr
187)
188{
190
191 const label np = c.size();
192
193 auto& UCorrect = cloud::createIOField<vector>("UCorrect", np, obr);
194
195 label i = 0;
196 for (const MPPICParcel<ParcelType>& p : c)
197 {
198 UCorrect[i] = p.UCorrect();
199
200 ++i;
201 }
202}
203
204
205// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
206
207template<class ParcelType>
208Foam::Ostream& Foam::operator<<
209(
210 Ostream& os,
212)
213{
214 if (os.format() == IOstream::ASCII)
215 {
216 os << static_cast<const ParcelType&>(p)
217 << token::SPACE << p.UCorrect();
218 }
219 else
220 {
221 os << static_cast<const ParcelType&>(p);
222 os.write
223 (
224 reinterpret_cast<const char*>(&p.UCorrect_),
226 );
227 }
228
230 return os;
231}
232
233
234// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
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.
Wrapper around kinematic parcel types to add MPPIC modelling.
Definition: MPPICParcel.H:80
static const std::size_t sizeofFields
Size in bytes of the fields.
Definition: MPPICParcel.H:84
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
vector UCorrect_
Velocity correction due to collisions [m/s].
Definition: MPPICParcel.H:168
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
static void readFields(CloudType &c)
Read.
Definition: MPPICParcelIO.C:89
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
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Registry of regIOobjects.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
static string propertyList()
Definition: particle.H:372
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling character strings derived from std::string.
Definition: string.H:79
@ SPACE
Space [isspace].
Definition: token.H:125
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#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.