36template<
class ParcelType>
41template<
class ParcelType>
44 sizeof(KinematicParcel<ParcelType>)
51template<
class ParcelType>
102template<
class ParcelType>
103template<
class CloudType>
106 const bool readOnProc = c.size();
108 ParcelType::readFields(c);
115 c.checkFieldIOobject(c, active);
122 c.checkFieldIOobject(c, typeId);
129 c.checkFieldIOobject(c, nParticle);
136 c.checkFieldIOobject(c, d);
143 c.checkFieldIOobject(c, dTarget);
150 c.checkFieldIOobject(c,
U);
157 c.checkFieldIOobject(c,
rho);
164 c.checkFieldIOobject(c, age);
171 c.checkFieldIOobject(c, tTurb);
178 c.checkFieldIOobject(c, UTurb);
185 c.checkFieldIOobject(c, UCorrect);
191 p.active_ = active[i];
192 p.typeId_ = typeId[i];
193 p.nParticle_ = nParticle[i];
195 p.dTarget_ = dTarget[i];
201 p.UCorrect_ = UCorrect[i];
208template<
class ParcelType>
209template<
class CloudType>
212 ParcelType::writeFields(c);
214 const label np = c.size();
215 const bool writeOnProc = c.size();
219 IOField<scalar> nParticle
235 for (
const KinematicParcel<ParcelType>&
p : c)
237 active[i] =
p.active();
238 typeId[i] =
p.typeId();
239 nParticle[i] =
p.nParticle();
241 dTarget[i] =
p.dTarget();
245 tTurb[i] =
p.tTurb();
246 UTurb[i] =
p.UTurb();
247 UCorrect[i] =
p.UCorrect();
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);
261 UTurb.write(writeOnProc);
266template<
class ParcelType>
275 ParcelType::writeProperties(
os, filters, delim, namesOnly);
278 #define writeProp(Name, Value) \
279 ParcelType::writeProperty(os, Name, Value, namesOnly, delim, filters)
297template<
class ParcelType>
298template<
class CloudType>
305 ParcelType::readObjects(c, obr);
307 if (!c.size())
return;
325 p.active_ = active[i];
326 p.typeId_ = typeId[i];
327 p.nParticle_ = nParticle[i];
329 p.dTarget_ = dTarget[i];
335 p.UCorrect_ = UCorrect[i];
342template<
class ParcelType>
343template<
class CloudType>
350 ParcelType::writeObjects(c, obr);
352 const label np = c.size();
368 for (
const KinematicParcel<ParcelType>&
p : c)
370 active[i] =
p.active();
371 typeId[i] =
p.typeId();
372 nParticle[i] =
p.nParticle();
374 dTarget[i] =
p.dTarget();
378 tTurb[i] =
p.tTurb();
379 UTurb[i] =
p.UTurb();
380 UCorrect[i] =
p.UCorrect();
389template<
class ParcelType>
390Foam::Ostream& Foam::operator<<
416 reinterpret_cast<const char*
>(&
p.active_),
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.
@ 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.
bool fatalCheckNativeSizes(const char *operation) const
Assert that the label/scalar byte-size associated with the stream are the native label/scalar sizes.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
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].
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.
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static const IOField< Type > & lookupIOField(const word &fieldName, const objectRegistry &obr)
Lookup an IOField within object registry.
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.
A class for handling character strings derived from std::string.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
OBJstream os(runTime.globalPath()/outputName)
#define writeProp(Name, Value)
DSMCCloud< dsmcParcel > CloudType
static constexpr const zero Zero
Global zero (0).