injectedParticleIO.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) 2016-2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "injectedParticle.H"
29#include "IOstreams.H"
30#include "IOField.H"
31#include "Cloud.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
37
38
40(
41 // Does not include position_
42 sizeof(injectedParticle) - offsetof(injectedParticle, tag_)
43);
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const polyMesh& mesh,
51 Istream& is,
52 bool readFields,
53 bool newFormat
54)
55:
56 particle(mesh, is, readFields, false), // force to read old positions file
57 position_(Zero),
58 tag_(-1),
59 soi_(0.0),
60 d_(0.0),
61 U_(Zero)
62{
63 if (readFields)
64 {
65 // After the base particle class has read the fields from file and
66 // constructed the necessary barycentric coordinates we can update the
67 // particle position on this mesh
69
70 if (is.format() == IOstream::ASCII)
71 {
72 is >> tag_ >> soi_ >> d_ >> U_;
73 }
74 else if (!is.checkLabelSize<>() || !is.checkScalarSize<>())
75 {
76 // Non-native label or scalar size
77 is.beginRawRead();
78
79 readRawLabel(is, &tag_);
80 readRawScalar(is, &soi_);
81 readRawScalar(is, &d_);
82 readRawScalar(is, U_.data(), vector::nComponents);
83
84 is.endRawRead();
85 }
86 else
87 {
88 is.read(reinterpret_cast<char*>(&tag_), sizeofFields);
89 }
90 }
91
93}
94
95
97{
98 if (!c.size())
99 {
100 return;
101 }
102
103 // Note: not reading local position_ - defer to base particle class
104
106
107 IOField<label> tag(c.fieldIOobject("tag", IOobject::MUST_READ));
108 c.checkFieldIOobject(c, tag);
109
110 IOField<scalar> soi(c.fieldIOobject("soi", IOobject::MUST_READ));
111 c.checkFieldIOobject(c, soi);
112
113 IOField<scalar> d(c.fieldIOobject("d", IOobject::MUST_READ));
114 c.checkFieldIOobject(c, d);
115
116 IOField<vector> U(c.fieldIOobject("U", IOobject::MUST_READ));
117 c.checkFieldIOobject(c, U);
118
119 label i = 0;
120 for (injectedParticle& p : c)
121 {
122 p.tag_ = tag[i];
123 p.soi_ = soi[i];
124 p.d_ = d[i];
125 p.U_ = U[i];
126
127 ++i;
128 }
129}
130
131
133{
134 // Force writing positions instead of coordinates
135 const bool oldWriteCoordinates = particle::writeLagrangianCoordinates;
136 const bool oldWritePositions = particle::writeLagrangianPositions;
137
140
142
143 // Note: not writing local position_ - defer to base particle class
144
145 label np = c.size();
146
147 IOField<label> tag(c.fieldIOobject("tag", IOobject::NO_READ), np);
148 IOField<scalar> soi(c.fieldIOobject("soi", IOobject::NO_READ), np);
149 IOField<scalar> d(c.fieldIOobject("d", IOobject::NO_READ), np);
150 IOField<vector> U(c.fieldIOobject("U", IOobject::NO_READ), np);
151
152 label i = 0;
153
154 for (const injectedParticle& p : c)
155 {
156 tag[i] = p.tag();
157 soi[i] = p.soi();
158 d[i] = p.d();
159 U[i] = p.U();
160
161 ++i;
162 }
163
164 tag.write();
165 soi.write();
166 d.write();
167 U.write();
168
169 // Restore
170 particle::writeLagrangianCoordinates = oldWriteCoordinates;
171 particle::writeLagrangianPositions = oldWritePositions;
172}
173
174
176(
177 Ostream& os,
178 const wordRes& filters,
179 const word& delim,
180 const bool namesOnly
181) const
182{
183 particle::writeProperties(os, filters, delim, namesOnly);
184
185 #undef writeProp
186 #define writeProp(Name, Value) \
187 particle::writeProperty(os, Name, Value, namesOnly, delim, filters)
188
189 writeProp("tag", tag_);
190 writeProp("soi", soi_);
191 writeProp("d", d_);
192 writeProp("U", U_);
193
194 #undef writeProp
195}
196
197
199(
201 const objectRegistry& obr
202)
203{
204 particle::readObjects(c, obr);
205
206 if (!c.size()) return;
207
208 const auto& tag = cloud::lookupIOField<label>("tag", obr);
209 const auto& soi = cloud::lookupIOField<scalar>("soi", obr);
210 const auto& d = cloud::lookupIOField<scalar>("d", obr);
211 const auto& U = cloud::lookupIOField<vector>("U", obr);
212
213 label i = 0;
214
215 for (injectedParticle& p : c)
216 {
217 p.tag() = tag[i];
218 p.soi() = soi[i];
219 p.d() = d[i];
220 p.U() = U[i];
221
222 ++i;
223 }
224}
225
226
228(
230 objectRegistry& obr
231)
232{
233 // Always writes "position", not "coordinates"
235
236 const label np = c.size();
237
238 auto& tag = cloud::createIOField<label>("tag", np, obr);
239 auto& soi = cloud::createIOField<scalar>("soi", np, obr);
240 auto& d = cloud::createIOField<scalar>("d", np, obr);
241 auto& U = cloud::createIOField<vector>("U", np, obr);
242
243 label i = 0;
244
245 for (const injectedParticle& p : c)
246 {
247 tag[i] = p.tag();
248 soi[i] = p.soi();
249 d[i] = p.d();
250 U[i] = p.U();
251
252 ++i;
253 }
254}
255
256
258{
259 if (os.format() == IOstream::ASCII)
260 {
261 os << position_ << token::SPACE << cell();
262 }
263 else
264 {
266
267 const size_t s =
268 (
269 offsetof(positionsCompat1706, facei)
270 - offsetof(positionsCompat1706, position));
271
272 p.position = position_;
273 p.celli = cell();
274
275 os.write(reinterpret_cast<const char*>(&p.position), s);
276 }
277
278 // Check state of Ostream
280}
281
282
283// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
284
285Foam::Ostream& Foam::operator<<
286(
287 Ostream& os,
288 const injectedParticle& p
289)
290{
291 // Note: not writing local position_ - defer to base particle class
292
293 if (os.format() == IOstream::ASCII)
294 {
295 os << static_cast<const particle&>(p)
296 << token::SPACE << p.tag()
297 << token::SPACE << p.soi()
298 << token::SPACE << p.d()
299 << token::SPACE << p.U();
300 }
301 else
302 {
303 os << static_cast<const particle&>(p);
304 os.write
305 (
306 reinterpret_cast<const char*>(&p.tag_),
308 );
309 }
310
312 return os;
313}
314
315
316// ************************************************************************* //
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
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
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
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
static void readFields(Cloud< injectedParticle > &c)
Read fields.
static const std::size_t sizeofFields
Size in bytes of the fields.
scalar d_
Diameter [m].
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.
scalar soi_
Start of injection [s].
virtual void writePosition(Ostream &) const
Write the particle position and cell.
vector U_
Velocity [m/s].
point position_
Position.
Registry of regIOobjects.
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Base particle class.
Definition: particle.H:79
vector position() const
Return current particle position.
Definition: particleI.H:314
static string propertyList()
Definition: particle.H:372
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
static bool writeLagrangianPositions
Definition: particle.H:383
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition: particleIO.C:220
static string propertyList_
String representation of properties.
Definition: particle.H:372
static bool writeLagrangianCoordinates
Definition: particle.H:379
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
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#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.
label readRawLabel(Istream &is)
Read raw label from binary stream.
Definition: label.C:46
Old particle positions content for OpenFOAM-1706 and earlier.
Definition: particle.H:117