Go to the documentation of this file.
6 #ifndef OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
7 #define OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
14 namespace ImplementationDetails
39 const Vector3D<MPCParticlePositionType>& position,
40 const unsigned int boxSizeX,
41 const unsigned int boxSizeY,
42 const unsigned int boxSizeZ);
82 const unsigned int particleID,
83 const MPCParticlePositionType*
const gridShift_,
86 const MPCParticlePositionType*
const positions,
87 MPCParticleVelocityType*
const velocities,
88 MPCParticleVelocityType*
const velocityCorrections,
89 unsigned int*
const collisionCellIndices,
90 const unsigned int boxSizeX,
91 const unsigned int boxSizeY,
92 const unsigned int boxSizeZ);
116 const unsigned int particleCount,
117 const MPCParticleVelocityType*
const velocities,
118 const unsigned int*
const collisionCellIndices,
119 MPCParticleVelocityType*
const collisionCellMomenta,
120 FP*
const collisionCellMasses,
121 const FP particleMass);
130 const unsigned int collisionCellIndex,
131 const MPCParticleVelocityType*
const collisionCellMomenta);
141 const unsigned int collisionCellIndex,
142 const MPCParticleVelocityType*
const collisionCellMomenta,
143 const FP*
const collisionCellMasses);
185 const unsigned int particleCount,
186 MPCParticleVelocityType*
const velocities,
187 const FP particleMass,
188 const unsigned int*
const collisionCellIndices,
189 const MPCParticleVelocityType*
const collisionCellMomenta,
190 const FP*
const collisionCellMasses,
191 const MPCParticlePositionType*
const collisionCellRotationAxes,
192 const FP collisionAngle,
193 FP*
const collisionCellFrameInternalKineticEnergies,
194 unsigned int*
const collisionCellParticleCounts);
219 const unsigned int particleCount,
220 MPCParticleVelocityType*
const velocities,
221 const unsigned int*
const collisionCellIndices,
222 const MPCParticleVelocityType*
const collisionCellMomenta,
223 const FP*
const collisionCellMasses,
224 const FP*
const collisionCellVelocityScalings);
239 const unsigned int particleCount,
240 MPCParticleVelocityType*
const velocities,
241 const MPCParticleVelocityType*
const velocityCorrections);
246 #endif //OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
Vector3D< MPCParticleVelocityType > getCollisionCellCenterOfMassMomentum(const unsigned int collisionCellIndex, const MPCParticleVelocityType *const collisionCellMomenta)
Returns the center-of-mass momentum of a collision cell.
unsigned int getCollisionCellIndex(const Vector3D< MPCParticlePositionType > &position, const unsigned int boxSizeX, const unsigned int boxSizeY, const unsigned int boxSizeZ)
Returns the collision cell index for the given position.
void collisionCellStochasticRotationStep1(const unsigned int particleCount, MPCParticleVelocityType *const velocities, const FP particleMass, const unsigned int *const collisionCellIndices, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses, const MPCParticlePositionType *const collisionCellRotationAxes, const FP collisionAngle, FP *const collisionCellFrameInternalKineticEnergies, unsigned int *const collisionCellParticleCounts)
Applies the first step of the SRD rotation to the given particles.
void undoLeesEdwardsVelocityCorrections(const unsigned int particleCount, MPCParticleVelocityType *const velocities, const MPCParticleVelocityType *const velocityCorrections)
Undoes the velocity corrections applied by sortIntoCollisionCellsLeesEdwards.
void collisionCellContributions(const unsigned int particleCount, const MPCParticleVelocityType *const velocities, const unsigned int *const collisionCellIndices, MPCParticleVelocityType *const collisionCellMomenta, FP *const collisionCellMasses, const FP particleMass)
Computes the collision cell mass and momentum contributions by the given particles.
void sortIntoCollisionCellsLeesEdwards(const unsigned int particleID, const MPCParticlePositionType *const gridShift_, const FP shearRate, const FP mpcTime, const MPCParticlePositionType *const positions, MPCParticleVelocityType *const velocities, MPCParticleVelocityType *const velocityCorrections, unsigned int *const collisionCellIndices, const unsigned int boxSizeX, const unsigned int boxSizeY, const unsigned int boxSizeZ)
Sorts the given particle into the collision cells, temporarily applying Lees-Edwards boundary conditi...
Vector3D< MPCParticleVelocityType > getCollisionCellCenterOfMassVelocity(const unsigned int collisionCellIndex, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses)
Returns the center-of-mass velocity of a collision cell.
void collisionCellStochasticRotationStep2(const unsigned int particleCount, MPCParticleVelocityType *const velocities, const unsigned int *const collisionCellIndices, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses, const FP *const collisionCellVelocityScalings)
Applies the first step of the SRD rotation to the given particles.