Go to the documentation of this file.
37 #ifndef CollidingParcel_H
38 #define CollidingParcel_H
54 template<
class ParcelType>
59 template<
class ParcelType>
70 template<
class ParcelType>
86 public ParcelType::constantProperties
152 +
" (angularMomentumx angularMomentumy angularMomentumz)"
153 +
" (torquex torquey torquez)"
154 +
" collisionRecordsPairAccessed"
155 +
" collisionRecordsPairOrigProcOfOther"
156 +
" collisionRecordsPairOrigIdOfOther"
157 +
" (collisionRecordsPairData)"
158 +
" collisionRecordsWallAccessed"
159 +
" collisionRecordsWallPRel"
160 +
" (collisionRecordsWallData)"
173 const label tetFaceI,
183 const label tetFaceI,
186 const scalar nParticle0,
188 const scalar dTarget0,
191 const vector& angularMomentum0,
193 const typename ParcelType::constantProperties& constProps
250 inline const vector&
f()
const;
280 template<
class TrackData>
281 bool move(TrackData& td,
const scalar trackTime);
295 template<
class CloudType>
299 template<
class CloudType>
305 friend Ostream& operator<< <ParcelType>
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
static void writeFields(const CloudType &c)
Write.
constantProperties()
Null constructor.
collisionRecordList collisionRecords_
Particle collision records.
demandDrivenEntry< scalar > youngsModulus_
Young's modulus [N/m2].
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....
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
This function object reads fields from the time directories and adds them to the mesh database for fu...
Mesh consisting of general polyhedral cells.
Class to hold thermo particle constant properties.
const vector & torque() const
Return const access to torque.
static const std::size_t sizeofFields_
Size in bytes of the fields.
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)....
Wrapper around kinematic parcel types to add collision modelling.
vector omega() const
Particle angular velocity.
demandDrivenEntry< scalar > poissonsRatio_
Poisson's ratio.
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 any number of values (e....
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.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
virtual autoPtr< particle > clone() const
Construct and return a (basic particle) clone.
iNew(const polyMesh &mesh)
const dimensionedScalar c
Speed of light in a vacuum.
vector torque_
Torque on particle due to collisions in global.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
CollidingParcel(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from owner, position, and cloud owner.
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.
bool move(TrackData &td, const scalar trackTime)
Move the parcel.
vectorFieldCompactIOField wallDataFieldCompactIOField