Loading...
Searching...
No Matches
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-2022 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
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() == IOstreamOption::ASCII)
71 {
72 is >> tag_ >> soi_ >> d_ >> U_;
73 }
74 else
75 {
76 // No non-native streaming
78
79 is.read(reinterpret_cast<char*>(&tag_), sizeofFields);
80 }
81 }
82
84}
85
86
88{
89 if (!c.size())
90 {
91 return;
92 }
93
94 // Note: not reading local position_ - defer to base particle class
95
97
98 IOField<label> tag(c.newIOobject("tag", IOobject::MUST_READ));
99 c.checkFieldIOobject(c, tag);
100
101 IOField<scalar> soi(c.newIOobject("soi", IOobject::MUST_READ));
102 c.checkFieldIOobject(c, soi);
103
104 IOField<scalar> d(c.newIOobject("d", IOobject::MUST_READ));
105 c.checkFieldIOobject(c, d);
106
107 IOField<vector> U(c.newIOobject("U", IOobject::MUST_READ));
108 c.checkFieldIOobject(c, U);
109
110 label i = 0;
111 for (injectedParticle& p : c)
112 {
113 p.tag_ = tag[i];
114 p.soi_ = soi[i];
115 p.d_ = d[i];
116 p.U_ = U[i];
117
118 ++i;
119 }
120}
121
122
123void Foam::injectedParticle::writeFields(const Cloud<injectedParticle>& c)
124{
125 // Force writing positions instead of coordinates
126 const bool oldWriteCoordinates = particle::writeLagrangianCoordinates;
127 const bool oldWritePositions = particle::writeLagrangianPositions;
128
131
133
134 // Note: not writing local position_ - defer to base particle class
135
136 const label np = c.size();
137
138 IOField<label> tag(c.newIOobject("tag", IOobject::NO_READ), np);
139 IOField<scalar> soi(c.newIOobject("soi", IOobject::NO_READ), np);
140 IOField<scalar> d(c.newIOobject("d", IOobject::NO_READ), np);
141 IOField<vector> U(c.newIOobject("U", IOobject::NO_READ), np);
142
143 label i = 0;
144
145 for (const injectedParticle& p : c)
146 {
147 tag[i] = p.tag();
148 soi[i] = p.soi();
149 d[i] = p.d();
150 U[i] = p.U();
151
152 ++i;
153 }
154
155 tag.write();
156 soi.write();
157 d.write();
158 U.write();
160 // Restore
161 particle::writeLagrangianCoordinates = oldWriteCoordinates;
162 particle::writeLagrangianPositions = oldWritePositions;
163}
164
165
167(
168 Ostream& os,
169 const wordRes& filters,
170 const word& delim,
171 const bool namesOnly
172) const
173{
174 particle::writeProperties(os, filters, delim, namesOnly);
175
176 #undef writeProp
177 #define writeProp(Name, Value) \
178 particle::writeProperty(os, Name, Value, namesOnly, delim, filters)
179
180 writeProp("tag", tag_);
181 writeProp("soi", soi_);
183 writeProp("U", U_);
184
185 #undef writeProp
186}
187
188
190(
192 const objectRegistry& obr
193)
194{
195 particle::readObjects(c, obr);
196
197 if (!c.size()) return;
198
199 const auto& tag = cloud::lookupIOField<label>("tag", obr);
200 const auto& soi = cloud::lookupIOField<scalar>("soi", obr);
201 const auto& d = cloud::lookupIOField<scalar>("d", obr);
202 const auto& U = cloud::lookupIOField<vector>("U", obr);
203
204 label i = 0;
205
206 for (injectedParticle& p : c)
207 {
208 p.tag() = tag[i];
209 p.soi() = soi[i];
210 p.d() = d[i];
211 p.U() = U[i];
212
213 ++i;
214 }
215}
216
217
219(
220 const Cloud<injectedParticle>& c,
221 objectRegistry& obr
222)
223{
224 // Always writes "position", not "coordinates"
226
227 const label np = c.size();
228
229 auto& tag = cloud::createIOField<label>("tag", np, obr);
230 auto& soi = cloud::createIOField<scalar>("soi", np, obr);
231 auto& d = cloud::createIOField<scalar>("d", np, obr);
232 auto& U = cloud::createIOField<vector>("U", np, obr);
233
234 label i = 0;
235
236 for (const injectedParticle& p : c)
237 {
238 tag[i] = p.tag();
239 soi[i] = p.soi();
240 d[i] = p.d();
241 U[i] = p.U();
242
243 ++i;
244 }
245}
246
247
248void Foam::injectedParticle::writePosition(Ostream& os) const
249{
250 if (os.format() == IOstreamOption::ASCII)
251 {
252 os << position_ << token::SPACE << cell();
253 }
254 else
255 {
256 positionsCompat1706 p;
257
258 const size_t s =
259 (
260 offsetof(positionsCompat1706, facei)
261 - offsetof(positionsCompat1706, position));
262
263 p.position = position_;
264 p.celli = cell();
265
266 os.write(reinterpret_cast<const char*>(&p.position), s);
267 }
268
269 // Check state of Ostream
270 os.check(FUNCTION_NAME);
271}
272
273
274// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
275
276Foam::Ostream& Foam::operator<<
277(
278 Ostream& os,
279 const injectedParticle& p
280)
281{
282 // Note: not writing local position_ - defer to base particle class
283
284 if (os.format() == IOstreamOption::ASCII)
285 {
287 << token::SPACE << p.tag()
288 << token::SPACE << p.soi()
289 << token::SPACE << p.d()
290 << token::SPACE << p.U();
291 }
292 else
293 {
295 os.write
296 (
297 reinterpret_cast<const char*>(&p.tag_),
299 );
300 }
301
303 return os;
304}
305
306
307// ************************************************************************* //
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
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
Definition cloud.H:202
static IOField< Type > & createIOField(const word &fieldName, const label nParticle, objectRegistry &obr)
Helper to construct IOField on a supplied object registry.
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.
const vector & U() const noexcept
Return const access to velocity.
static void writeObjects(const Cloud< injectedParticle > &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
scalar d_
Diameter [m].
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
scalar d() const noexcept
Return const access to diameter.
static void writeFields(const Cloud< injectedParticle > &c)
Write fields.
label tag() const noexcept
Return const access to the tag.
scalar soi() const noexcept
Return const access to the start of injection.
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.
injectedParticle(const polyMesh &mesh, const vector &position, const label celli=-1)
Construct from a position and a cell.
vector U_
Velocity [m/s].
Registry of regIOobjects.
vector position() const
Return current particle position.
Definition particleI.H:283
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
static string propertyList()
Definition particle.H:455
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
static bool writeLagrangianPositions
Write particle positions file (v1706 format and earlier) Default is true (disable in etc/controlDict)...
Definition particle.H:472
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition particleIO.C:177
label cell() const noexcept
Return current cell particle is in.
Definition particleI.H:122
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
static string propertyList_
String representation of properties.
Definition particle.H:455
static bool writeLagrangianCoordinates
Write particle coordinates file (v1712 and later) Default is true.
Definition particle.H:466
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for handling character strings derived from std::string.
Definition string.H:76
@ SPACE
Space [isspace].
Definition token.H:144
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
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
const dimensionedScalar c
Speed of light in a vacuum.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Old particle positions content for OpenFOAM-1706 and earlier.
Definition particle.H:123