Loading...
Searching...
No Matches
KinematicParcelIO.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 "KinematicParcel.H"
30#include "IOstreams.H"
31#include "IOField.H"
32#include "Cloud.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36template<class ParcelType>
39
40
41template<class ParcelType>
43(
44 sizeof(KinematicParcel<ParcelType>)
46);
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51template<class ParcelType>
53(
54 const polyMesh& mesh,
55 Istream& is,
56 bool readFields,
57 bool newFormat
58)
59:
60 ParcelType(mesh, is, readFields, newFormat),
61 active_(false),
62 typeId_(0),
63 nParticle_(0.0),
64 d_(0.0),
65 dTarget_(0.0),
66 U_(Zero),
67 rho_(0.0),
68 age_(0.0),
69 tTurb_(0.0),
70 UTurb_(Zero),
72{
73 if (readFields)
74 {
75 if (is.format() == IOstreamOption::ASCII)
76 {
77 is >> active_
78 >> typeId_
79 >> nParticle_
80 >> d_
81 >> dTarget_
82 >> U_
83 >> rho_
84 >> age_
85 >> tTurb_
86 >> UTurb_
87 >> UCorrect_;
88 }
89 else
90 {
91 // No non-native streaming
93
94 is.read(reinterpret_cast<char*>(&active_), sizeofFields);
95 }
96 }
99}
100
101
102template<class ParcelType>
103template<class CloudType>
105{
106 const bool readOnProc = c.size();
107
108 ParcelType::readFields(c);
109
110 IOField<label> active
111 (
112 c.newIOobject("active", IOobject::MUST_READ),
113 readOnProc
114 );
115 c.checkFieldIOobject(c, active);
116
117 IOField<label> typeId
118 (
119 c.newIOobject("typeId", IOobject::MUST_READ),
120 readOnProc
121 );
122 c.checkFieldIOobject(c, typeId);
123
124 IOField<scalar> nParticle
125 (
126 c.newIOobject("nParticle", IOobject::MUST_READ),
127 readOnProc
128 );
129 c.checkFieldIOobject(c, nParticle);
130
132 (
133 c.newIOobject("d", IOobject::MUST_READ),
134 readOnProc
135 );
136 c.checkFieldIOobject(c, d);
137
138 IOField<scalar> dTarget
139 (
140 c.newIOobject("dTarget", IOobject::MUST_READ),
141 readOnProc
142 );
143 c.checkFieldIOobject(c, dTarget);
144
146 (
147 c.newIOobject("U", IOobject::MUST_READ),
148 readOnProc
149 );
150 c.checkFieldIOobject(c, U);
151
153 (
154 c.newIOobject("rho", IOobject::MUST_READ),
155 readOnProc
156 );
157 c.checkFieldIOobject(c, rho);
158
160 (
161 c.newIOobject("age", IOobject::MUST_READ),
162 readOnProc
163 );
164 c.checkFieldIOobject(c, age);
165
166 IOField<scalar> tTurb
167 (
168 c.newIOobject("tTurb", IOobject::MUST_READ),
169 readOnProc
170 );
171 c.checkFieldIOobject(c, tTurb);
172
173 IOField<vector> UTurb
174 (
175 c.newIOobject("UTurb", IOobject::MUST_READ),
176 readOnProc
177 );
178 c.checkFieldIOobject(c, UTurb);
179
180 IOField<vector> UCorrect
181 (
182 c.newIOobject("UCorrect", IOobject::MUST_READ),
183 readOnProc
184 );
185 c.checkFieldIOobject(c, UCorrect);
186
187 label i = 0;
188
190 {
191 p.active_ = active[i];
192 p.typeId_ = typeId[i];
193 p.nParticle_ = nParticle[i];
194 p.d_ = d[i];
195 p.dTarget_ = dTarget[i];
196 p.U_ = U[i];
197 p.rho_ = rho[i];
198 p.age_ = age[i];
199 p.tTurb_ = tTurb[i];
200 p.UTurb_ = UTurb[i];
201 p.UCorrect_ = UCorrect[i];
202
203 ++i;
204 }
205}
206
207
208template<class ParcelType>
209template<class CloudType>
211{
212 ParcelType::writeFields(c);
213
214 const label np = c.size();
215 const bool writeOnProc = c.size();
216
217 IOField<label> active(c.newIOobject("active", IOobject::NO_READ), np);
218 IOField<label> typeId(c.newIOobject("typeId", IOobject::NO_READ), np);
219 IOField<scalar> nParticle
220 (
221 c.newIOobject("nParticle", IOobject::NO_READ),
222 np
223 );
224 IOField<scalar> d(c.newIOobject("d", IOobject::NO_READ), np);
225 IOField<scalar> dTarget(c.newIOobject("dTarget", IOobject::NO_READ), np);
226 IOField<vector> U(c.newIOobject("U", IOobject::NO_READ), np);
227 IOField<scalar> rho(c.newIOobject("rho", IOobject::NO_READ), np);
228 IOField<scalar> age(c.newIOobject("age", IOobject::NO_READ), np);
229 IOField<scalar> tTurb(c.newIOobject("tTurb", IOobject::NO_READ), np);
230 IOField<vector> UTurb(c.newIOobject("UTurb", IOobject::NO_READ), np);
231 IOField<vector> UCorrect(c.newIOobject("UCorrect", IOobject::NO_READ), np);
232
233 label i = 0;
234
235 for (const KinematicParcel<ParcelType>& p : c)
236 {
237 active[i] = p.active();
238 typeId[i] = p.typeId();
239 nParticle[i] = p.nParticle();
240 d[i] = p.d();
241 dTarget[i] = p.dTarget();
242 U[i] = p.U();
243 rho[i] = p.rho();
244 age[i] = p.age();
245 tTurb[i] = p.tTurb();
246 UTurb[i] = p.UTurb();
247 UCorrect[i] = p.UCorrect();
248
249 ++i;
250 }
251
252 active.write(writeOnProc);
253 typeId.write(writeOnProc);
254 nParticle.write(writeOnProc);
255 d.write(writeOnProc);
256 dTarget.write(writeOnProc);
257 U.write(writeOnProc);
258 rho.write(writeOnProc);
259 age.write(writeOnProc);
260 tTurb.write(writeOnProc);
261 UTurb.write(writeOnProc);
262 UCorrect.write(writeOnProc);
263}
264
265
266template<class ParcelType>
268(
269 Ostream& os,
270 const wordRes& filters,
271 const word& delim,
272 const bool namesOnly
273) const
274{
275 ParcelType::writeProperties(os, filters, delim, namesOnly);
276
277 #undef writeProp
278 #define writeProp(Name, Value) \
279 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
280
281 writeProp("active", active_);
282 writeProp("typeId", typeId_);
283 writeProp("nParticle", nParticle_);
284 writeProp("d", d_);
285 writeProp("dTarget", dTarget_);
286 writeProp("U", U_);
287 writeProp("rho", rho_);
288 writeProp("age", age_);
289 writeProp("tTurb", tTurb_);
290 writeProp("UTurb", UTurb_);
291 writeProp("UCorrect", UCorrect_);
293 #undef writeProp
294}
295
296
297template<class ParcelType>
298template<class CloudType>
300(
301 CloudType& c,
302 const objectRegistry& obr
303)
304{
305 ParcelType::readObjects(c, obr);
306
307 if (!c.size()) return;
308
309 const auto& active = cloud::lookupIOField<label>("active", obr);
310 const auto& typeId = cloud::lookupIOField<label>("typeId", obr);
311 const auto& nParticle = cloud::lookupIOField<scalar>("nParticle", obr);
312 const auto& d = cloud::lookupIOField<scalar>("d", obr);
313 const auto& dTarget = cloud::lookupIOField<scalar>("dTarget", obr);
314 const auto& U = cloud::lookupIOField<vector>("U", obr);
315 const auto& rho = cloud::lookupIOField<scalar>("rho", obr);
316 const auto& age = cloud::lookupIOField<scalar>("age", obr);
317 const auto& tTurb = cloud::lookupIOField<scalar>("tTurb", obr);
318 const auto& UTurb = cloud::lookupIOField<vector>("UTurb", obr);
319 const auto& UCorrect = cloud::lookupIOField<vector>("UCorrect", obr);
320
321 label i = 0;
322
324 {
325 p.active_ = active[i];
326 p.typeId_ = typeId[i];
327 p.nParticle_ = nParticle[i];
328 p.d_ = d[i];
329 p.dTarget_ = dTarget[i];
330 p.U_ = U[i];
331 p.rho_ = rho[i];
332 p.age_ = age[i];
333 p.tTurb_ = tTurb[i];
334 p.UTurb_ = UTurb[i];
335 p.UCorrect_ = UCorrect[i];
336
337 ++i;
338 }
339}
340
341
342template<class ParcelType>
343template<class CloudType>
345(
346 const CloudType& c,
347 objectRegistry& obr
348)
349{
350 ParcelType::writeObjects(c, obr);
351
352 const label np = c.size();
353
354 auto& active = cloud::createIOField<label>("active", np, obr);
355 auto& typeId = cloud::createIOField<label>("typeId", np, obr);
356 auto& nParticle = cloud::createIOField<scalar>("nParticle", np, obr);
357 auto& d = cloud::createIOField<scalar>("d", np, obr);
358 auto& dTarget = cloud::createIOField<scalar>("dTarget", np, obr);
359 auto& U = cloud::createIOField<vector>("U", np, obr);
360 auto& rho = cloud::createIOField<scalar>("rho", np, obr);
361 auto& age = cloud::createIOField<scalar>("age", np, obr);
362 auto& tTurb = cloud::createIOField<scalar>("tTurb", np, obr);
363 auto&& UTurb = cloud::createIOField<vector>("UTurb", np, obr);
364 auto&& UCorrect = cloud::createIOField<vector>("UCorrect", np, obr);
365
366 label i = 0;
367
368 for (const KinematicParcel<ParcelType>& p : c)
369 {
370 active[i] = p.active();
371 typeId[i] = p.typeId();
372 nParticle[i] = p.nParticle();
373 d[i] = p.d();
374 dTarget[i] = p.dTarget();
375 U[i] = p.U();
376 rho[i] = p.rho();
377 age[i] = p.age();
378 tTurb[i] = p.tTurb();
379 UTurb[i] = p.UTurb();
380 UCorrect[i] = p.UCorrect();
381
382 ++i;
383 }
384}
385
386
387// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
388
389template<class ParcelType>
390Foam::Ostream& Foam::operator<<
391(
392 Ostream& os,
394)
395{
397 {
399 << token::SPACE << bool(p.active())
400 << token::SPACE << p.typeId()
401 << token::SPACE << p.nParticle()
402 << token::SPACE << p.d()
403 << token::SPACE << p.dTarget()
404 << token::SPACE << p.U()
405 << token::SPACE << p.rho()
406 << token::SPACE << p.age()
407 << token::SPACE << p.tTurb()
408 << token::SPACE << p.UTurb()
409 << token::SPACE << p.UCorrect();
410 }
411 else
412 {
414 os.write
415 (
416 reinterpret_cast<const char*>(&p.active_),
418 );
419 }
420
422 return os;
423}
424
425
426// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
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.
Kinematic parcel class with rotational motion (as spherical particles only) and one/two-way coupling ...
const vector & UCorrect() const
label typeId_
Parcel type id.
vector UTurb_
Turbulent velocity fluctuation [m/s].
scalar rho_
Density [kg/m3].
const vector & U() const
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
label active_
Active flag - tracking inactive when active = false.
static void readFields(TrackCloudType &c)
Read.
scalar tTurb_
Time spent in turbulent eddy [s].
scalar dTarget_
Target diameter [m].
static const std::size_t sizeofFields
Size in bytes of the fields.
scalar d_
Diameter [m].
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
static void writeFields(const TrackCloudType &c)
Write.
KinematicParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
vector UCorrect_
Velocity correction due to collisions MPPIC [m/s].
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
Write individual parcel properties to stream.
vector U_
Velocity of Parcel [m/s].
const vector & UTurb() const
scalar nParticle_
Number of particles in Parcel.
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
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.
Registry of regIOobjects.
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
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
#define FUNCTION_NAME
DSMCCloud< dsmcParcel > CloudType
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127