|
|
|
Particle< ParticleType > Class Template Reference
Inheritance diagram for Particle< ParticleType >:
![]()
Collaboration diagram for Particle< ParticleType >:
![]()
Detailed Descriptiontemplate<class ParticleType>
| |||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
| virtual ~Particle | ( | ) | [inline, virtual] |
| Foam::scalar lambda | ( | const vector & | from, | |
| const vector & | to, | |||
| const label | facei, | |||
| const scalar | stepFraction | |||
| ) | const [inline, protected] |
Return the 'lambda' value for the position, p, on the face,.
where, p = from + lamda*(to - from) for non-static meshes
Definition at line 24 of file ParticleI.H.
References b, beta(), primitiveMesh::cellCentres(), cp, primitiveMesh::faceAreas(), primitiveMesh::faceCentres(), polyMesh::faces(), Foam::mag(), mesh, polyMesh::moving(), polyMesh::oldPoints(), and Foam::sqrt().

| Foam::scalar lambda | ( | const vector & | from, | |
| const vector & | to, | |||
| const label | facei | |||
| ) | const [inline, protected] |
Return the 'lambda' value for the position, p, on the face,.
where, p = from + lamda*(to - from) for static meshes
Definition at line 172 of file ParticleI.H.
| void findFaces | ( | const vector & | position, | |
| DynamicList< label > & | faceList | |||
| ) | const [inline, protected] |
Find the faces between position and cell centre.
Definition at line 31 of file Particle.C.
References primitiveMesh::cellCentres(), primitiveMesh::cells(), DynamicList< T, SizeInc, SizeMult, SizeDiv >::clear(), forAll, and mesh.

| void findFaces | ( | const vector & | position, | |
| const label | celli, | |||
| const scalar | stepFraction, | |||
| DynamicList< label > & | faceList | |||
| ) | const [inline, protected] |
Find the faces between position and cell centre.
Definition at line 56 of file Particle.C.
References DynamicList< T, SizeInc, SizeMult, SizeDiv >::append().

Overridable function to handle the particle hitting a patch.
Executed before other patch-hitting functions
Reimplemented in DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 421 of file Particle.C.
| void hitWedgePatch | ( | const wedgePolyPatch & | wpp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a wedgePatch.
Definition at line 434 of file Particle.C.
| void hitSymmetryPatch | ( | const symmetryPolyPatch & | spp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a.
symmetryPatch
Definition at line 449 of file Particle.C.
| void hitCyclicPatch | ( | const cyclicPolyPatch & | cpp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a cyclicPatch.
Definition at line 464 of file Particle.C.
| void hitProcessorPatch | ( | const processorPolyPatch & | spp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a.
processorPatch
Reimplemented in DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 496 of file Particle.C.
| void hitWallPatch | ( | const wallPolyPatch & | spp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a wallPatch.
Reimplemented in DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 506 of file Particle.C.
| void hitPatch | ( | const polyPatch & | spp, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Overridable function to handle the particle hitting a.
general patch
Reimplemented in DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 516 of file Particle.C.
| void transformPosition | ( | const tensor & | T | ) | [inline, protected, virtual] |
Transform the position the particle.
according to the given transformation tensor
Definition at line 402 of file Particle.C.
| void transformProperties | ( | const tensor & | T | ) | [inline, protected, virtual] |
Transform the physical properties of the particle.
according to the given transformation tensor
Reimplemented in parcel, DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, molecule, solidParticle, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 409 of file Particle.C.
References Particle< ParticleType >::position_, and Foam::transform().

| void transformProperties | ( | const vector & | separation | ) | [inline, protected, virtual] |
Transform the physical properties of the particle.
according to the given separation vector
Reimplemented in parcel, DsmcParcel< ParcelType >, KinematicParcel< ParcelType >, molecule, solidParticle, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 414 of file Particle.C.
| void prepareForParallelTransfer | ( | const label | patchi, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Convert global addressing to the processor patch.
local equivalents
Definition at line 84 of file Particle.C.
| void correctAfterParallelTransfer | ( | const label | patchi, | |
| TrackData & | td | |||
| ) | [inline, protected] |
Convert processor patch addressing to the global equivalents.
and set the celli to the face-neighbour
Definition at line 97 of file Particle.C.
| TypeName | ( | "Particle< ParticleType >" | ) |
Runtime type information.
| autoPtr<ParticleType> clone | ( | ) | const [inline] |
Construct a clone.
Reimplemented in trackedParticle, indexedParticle, passiveParticle, CoalParcel< ThermoType >, dsmcParcel, DsmcParcel< ParcelType >, basicKinematicParcel, BasicReactingMultiphaseParcel< ThermoType >, BasicReactingParcel< ThermoType >, basicThermoParcel, KinematicParcel< ParcelType >, ReactingMultiphaseParcel< ParcelType >, ReactingParcel< ParcelType >, ThermoParcel< ParcelType >, molecule, solidParticle, DsmcParcel< dsmcParcel >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, KinematicParcel< basicThermoParcel >, ReactingMultiphaseParcel< CoalParcel< ThermoType > >, ReactingMultiphaseParcel< BasicReactingMultiphaseParcel< ThermoType > >, ReactingParcel< CoalParcel< ThermoType > >, ReactingParcel< BasicReactingParcel< ThermoType > >, ReactingParcel< BasicReactingMultiphaseParcel< ThermoType > >, ReactingParcel< ParcelType >, ThermoParcel< CoalParcel< ThermoType > >, ThermoParcel< BasicReactingMultiphaseParcel< ThermoType > >, ThermoParcel< BasicReactingParcel< ThermoType > >, ThermoParcel< basicThermoParcel >, and ThermoParcel< ParcelType >.
Definition at line 318 of file Particle.H.
| bool inCell | ( | ) | const [inline] |
| const Foam::vector & position | ( | ) | const [inline] |
Return current particle position.
Definition at line 278 of file ParticleI.H.
References Particle< ParticleType >::cloud_.
| Foam::vector & position | ( | ) | [inline] |
Return current particle position.
Definition at line 285 of file ParticleI.H.
References Particle< ParticleType >::position_.
| Foam::label & cell | ( | ) | [inline] |
Return current cell particle is in.
Definition at line 299 of file ParticleI.H.
References Particle< ParticleType >::celli_.
| Foam::label cell | ( | ) | const [inline] |
Return current cell particle is in.
Definition at line 292 of file ParticleI.H.
References Particle< ParticleType >::position_.
| Foam::label face | ( | ) | const [inline] |
Return current face particle is on otherwise -1.
Definition at line 306 of file ParticleI.H.
References Particle< ParticleType >::celli_.
| const Foam::Cloud< ParticleType > & cloud | ( | ) | const [inline] |
| bool softImpact | ( | ) | const [inline] |
Return the impact model to be used, soft or hard (default).
Definition at line 348 of file ParticleI.H.
References Particle< ParticleType >::origId_.
| Foam::scalar currentTime | ( | ) | const [inline] |
| bool onBoundary | ( | ) | const [inline] |
Is the particle on the boundary/(or outside the domain)?
Definition at line 313 of file ParticleI.H.
References Particle< ParticleType >::facei_.
| Foam::label patch | ( | const label | facei | ) | const [inline] |
Which patch is particle on.
Definition at line 364 of file ParticleI.H.
References Particle< ParticleType >::cloud_, and Particle< ParticleType >::stepFraction_.
| Foam::label patchFace | ( | const label | patchi, | |
| const label | facei | |||
| ) | const [inline] |
Which face of this patch is this particle on.
Definition at line 372 of file ParticleI.H.
References Particle< ParticleType >::cloud_.
| Foam::scalar wallImpactDistance | ( | const vector & | n | ) | const [inline] |
The nearest distance to a wall that.
the particle can be in the n direction
Reimplemented in KinematicParcel< ParcelType >, solidParticle, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 383 of file ParticleI.H.
References Particle< ParticleType >::cloud_.
| Foam::scalar & stepFraction | ( | ) | [inline] |
Return the fraction of time-step completed.
Definition at line 320 of file ParticleI.H.
References Particle< ParticleType >::cloud_, and Particle< ParticleType >::facei_.
Referenced by spray::inject().

| Foam::scalar stepFraction | ( | ) | const [inline] |
Return the fraction of time-step completed.
Definition at line 327 of file ParticleI.H.
References Particle< ParticleType >::stepFraction_.
| Foam::label origProc | ( | ) | const [inline] |
Return the originating processor id.
Definition at line 334 of file ParticleI.H.
References Particle< ParticleType >::stepFraction_.
| Foam::label origId | ( | ) | const [inline] |
Return the particle id on originating processor.
Definition at line 341 of file ParticleI.H.
References Particle< ParticleType >::origProc_.
| Foam::label track | ( | const vector & | endPosition, | |
| TrackData & | td | |||
| ) | [inline] |
Track particle to end of trajectory.
or until it hits the boundary. On entry 'stepFraction()' should be set to the fraction of the time-step at which the tracking starts and on exit it contains the fraction of the time-step completed. Returns the boundary face index if the track stops at the boundary, -1 otherwise.
Definition at line 194 of file Particle.C.
Referenced by Particle< ParticleType >::trackToFace().

| Foam::label track | ( | const vector & | endPosition | ) | [inline] |
Calls the templated track with dummy TrackData.
Reimplemented in ExactParticle< ParticleType >, and ExactParticle< trackedParticle >.
Definition at line 213 of file Particle.C.
| Foam::scalar trackToFace | ( | const vector & | endPosition, | |
| TrackData & | td | |||
| ) | [inline] |
Track particle to a given position and returns 1.0 if the.
trajectory is completed without hitting a face otherwise stops at the face and returns the fraction of the trajectory completed. on entry 'stepFraction()' should be set to the fraction of the time-step at which the tracking starts.
Definition at line 222 of file Particle.C.
References Particle< ParticleType >::track().

| Foam::scalar trackToFace | ( | const vector & | endPosition | ) | [inline] |
Calls the templated trackToFace with dummy TrackData.
Reimplemented in ExactParticle< ParticleType >, and ExactParticle< trackedParticle >.
Definition at line 393 of file Particle.C.
References primitiveMesh::cellCentres(), and Particle< ParticleType >::position_.

| Foam::label faceInterpolation | ( | ) | const [inline] |
Return the index of the face to be used in the interpolation.
routine
Reimplemented in KinematicParcel< ParcelType >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, and KinematicParcel< basicThermoParcel >.
Definition at line 390 of file ParticleI.H.
| void readFields | ( | Cloud< ParticleType > & | c | ) | [inline, static] |
Read the fields associated with the owner cloud.
Reimplemented in parcel, and solidParticle.
Definition at line 101 of file ParticleIO.C.
| void writeFields | ( | const Cloud< ParticleType > & | c | ) | [inline, static] |
Write the fields associated with the owner cloud.
Reimplemented in parcel, and solidParticle.
Definition at line 134 of file ParticleIO.C.
| void write | ( | Ostream & | os, | |
| bool | writeFields | |||
| ) | const [inline] |
friend class Cloud< ParticleType > [friend] |
Reimplemented in ExactParticle< ParticleType >, and ExactParticle< trackedParticle >.
Definition at line 284 of file Particle.H.
Reference to the particle cloud.
Definition at line 107 of file Particle.H.
Referenced by Particle< ParticleType >::Particle(), Particle< ParticleType >::patch(), Particle< ParticleType >::patchFace(), Particle< ParticleType >::position(), Particle< ParticleType >::stepFraction(), and Particle< ParticleType >::wallImpactDistance().
Position of particle.
Definition at line 110 of file Particle.H.
Referenced by Particle< ParticleType >::cell(), Particle< ParticleType >::Particle(), Particle< ParticleType >::position(), Particle< ParticleType >::trackToFace(), and Particle< ParticleType >::transformProperties().
Index of the cell it is in.
Definition at line 113 of file Particle.H.
Referenced by Particle< ParticleType >::cell(), Particle< ParticleType >::face(), and Particle< ParticleType >::Particle().
Face index if the particle is on a face otherwise -1.
Definition at line 116 of file Particle.H.
Referenced by Particle< ParticleType >::onBoundary(), Particle< ParticleType >::Particle(), and Particle< ParticleType >::stepFraction().
scalar stepFraction_ [protected] |
Fraction of time-step completed.
Definition at line 119 of file Particle.H.
Referenced by Particle< ParticleType >::origProc(), Particle< ParticleType >::Particle(), Particle< ParticleType >::patch(), and Particle< ParticleType >::stepFraction().
Originating processor id.
Definition at line 122 of file Particle.H.
Referenced by Particle< ParticleType >::origId(), and Particle< ParticleType >::Particle().
Local particle id on originating processor.
Definition at line 125 of file Particle.H.
Referenced by Particle< ParticleType >::Particle(), and Particle< ParticleType >::softImpact().
Foam::string propHeader [inline, static] |
Initial value:
"(Px Py Pz) cellI origProc origId"
Reimplemented in KinematicParcel< ParcelType >, ReactingMultiphaseParcel< ParcelType >, ReactingParcel< ParcelType >, ThermoParcel< ParcelType >, KinematicParcel< CoalParcel< ThermoType > >, KinematicParcel< basicKinematicParcel >, KinematicParcel< ParcelType >, KinematicParcel< BasicReactingMultiphaseParcel< ThermoType > >, KinematicParcel< BasicReactingParcel< ThermoType > >, KinematicParcel< basicThermoParcel >, ReactingMultiphaseParcel< CoalParcel< ThermoType > >, ReactingMultiphaseParcel< BasicReactingMultiphaseParcel< ThermoType > >, ReactingParcel< CoalParcel< ThermoType > >, ReactingParcel< BasicReactingParcel< ThermoType > >, ReactingParcel< BasicReactingMultiphaseParcel< ThermoType > >, ReactingParcel< ParcelType >, ThermoParcel< CoalParcel< ThermoType > >, ThermoParcel< BasicReactingMultiphaseParcel< ThermoType > >, ThermoParcel< BasicReactingParcel< ThermoType > >, ThermoParcel< basicThermoParcel >, and ThermoParcel< ParcelType >.
Definition at line 293 of file Particle.H.