Go to the documentation of this file.
43 #ifndef CollidingParcel_H
44 #define CollidingParcel_H
61 template<
class ParcelType>
66 template<
class ParcelType>
77 template<
class ParcelType>
91 public ParcelType::constantProperties
161 +
" (angularMomentumx angularMomentumy angularMomentumz)"
162 +
" (torquex torquey torquez)"
163 +
" collisionRecordsPairAccessed"
164 +
" collisionRecordsPairOrigProcOfOther"
165 +
" collisionRecordsPairOrigIdOfOther"
166 +
" (collisionRecordsPairData)"
167 +
" collisionRecordsWallAccessed"
168 +
" collisionRecordsWallPRel"
169 +
" (collisionRecordsWallData)"
182 const label tetFacei,
201 const label tetFacei,
204 const scalar nParticle0,
206 const scalar dTarget0,
209 const vector& angularMomentum0,
211 const typename ParcelType::constantProperties& constProps
220 bool newFormat =
true
269 inline const vector&
f()
const;
299 template<
class TrackCloudType>
302 TrackCloudType&
cloud,
304 const scalar trackTime
319 template<
class CloudType>
323 template<
class CloudType>
336 template<
class CloudType>
340 template<
class CloudType>
346 friend Ostream& operator<< <ParcelType>
static void writeFields(const CloudType &c)
Write.
static const std::size_t sizeofFields
Size in bytes of the fields.
A class for handling words, derived from Foam::string.
constantProperties()
Null constructor.
collisionRecordList collisionRecords_
Particle collision records.
autoPtr< CollidingParcel< ParcelType > > operator()(Istream &is) const
AddToPropertyList(ParcelType, " (fx fy fz)"+" (angularMomentumx angularMomentumy angularMomentumz)"+" (torquex torquey torquez)"+" collisionRecordsPairAccessed"+" collisionRecordsPairOrigProcOfOther"+" collisionRecordsPairOrigIdOfOther"+" (collisionRecordsPairData)"+" collisionRecordsWallAccessed"+" collisionRecordsWallPRel"+" (collisionRecordsWallData)")
String representation of properties.
static void readFields(CloudType &c)
Read.
vectorFieldCompactIOField pairDataFieldCompactIOField
A Field of objects of type <T> with automated input and output using a compact storage....
Mesh consisting of general polyhedral cells.
Class to hold thermo particle constant properties.
const vector & torque() const
Return const access to torque.
Registry of regIOobjects.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
scalar youngsModulus() const
Return const access to Young's Modulus.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
Wrapper around kinematic parcel types to add collision modelling.
vector omega() const
Particle angular velocity.
PtrList< coordinateSystem > coordinates(solidRegions.size())
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const vector & angularMomentum() const
Return const access to angular momentum.
Templated base class for dsmc cloud.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Factory class to read-construct particles used for.
CollisionRecordList< vector, vector > collisionRecordList
scalar poissonsRatio() const
Return const access to Poisson's ratio.
const collisionRecordList & collisionRecords() const
Return const access to the collision records.
TypeName("CollidingParcel")
Runtime type information.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A cloud is a registry collection of lagrangian particles.
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
iNew(const polyMesh &mesh)
ParcelType::trackingData trackingData
Use base tracking data.
A List of wordRe with additional matching capabilities.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Move the parcel.
const dimensionedScalar c
Speed of light in a vacuum.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
vector torque_
Torque on particle due to collisions in global.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
vector f_
Force on particle due to collisions [N].
virtual autoPtr< particle > clone(const polyMesh &mesh) const
Construct and return a (basic particle) clone.
vector angularMomentum_
Angular momentum of Parcel in global reference frame [kg m2/s].
const vector & f() const
Return const access to force.
CollidingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
vectorFieldCompactIOField wallDataFieldCompactIOField