OpenMPCD
ImplementationDetails/Simulation.hpp
Go to the documentation of this file.
1 /**
2  * @file
3  * Declares various functions that are needed in MPCD.
4  */
5 
6 #ifndef OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
7 #define OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
8 
9 #include <OpenMPCD/Types.hpp>
10 #include <OpenMPCD/Vector3D.hpp>
11 
12 namespace OpenMPCD
13 {
14 namespace ImplementationDetails
15 {
16 
17 /**
18  * Returns the collision cell index for the given position.
19  *
20  * The collision cell index is determined as follows:
21  * The three cartesian coordinates of the given particle position are rounded
22  * down towards the nearest integer. This triple of numbers then determine the
23  * coordinates of the collision cell the particle lies in. If they are called
24  * `cellX`, `cellY`, and `cellZ`, respectively, what is returned is
25  * `cellX + cellY * boxSizeX` + cellZ * boxSizeX * boxSizeY`.
26  *
27  * @param[in] position The position to consider, which must lie in the primary
28  * simulation volume, i.e. the `x`, `y`, and `z` coordinates
29  * must be non-negative and smaller than the simulation box
30  * size along the respective direction.
31  * @param[in] boxSizeX The number of collision cells in the primary simulation
32  * volume along the `x` axis.
33  * @param[in] boxSizeY The number of collision cells in the primary simulation
34  * volume along the `y` axis.
35  * @param[in] boxSizeZ The number of collision cells in the primary simulation
36  * volume along the `z` axis.
37  */
38 unsigned int getCollisionCellIndex(
39  const Vector3D<MPCParticlePositionType>& position,
40  const unsigned int boxSizeX,
41  const unsigned int boxSizeY,
42  const unsigned int boxSizeZ);
43 
44 /**
45  * Sorts the given particle into the collision cells, temporarily applying
46  * Lees-Edwards boundary conditions.
47  *
48  * @see LeesEdwardsBoundaryConditions
49  *
50  * @param[in] particleID
51  * The particle ID.
52  * @param[in] gridShift_
53  * The grid shift vector to temporarily add to the position.
54  * @param[in] shearRate
55  * The Lees-Edwards shear rate \f$ \dot{\gamma} \f$.
56  * @param[in] mpcTime
57  * The MPC time that has passed since the start of the
58  * simulation.
59  * @param[in] positions
60  * The particle positions.
61  * @param[in,out] velocities
62  * The particle velocities, which may be changed due to
63  * Lees-Edwards boundary conditions.
64  * @param[out] velocityCorrections
65  * An array, long enough to store one element per particle; at
66  * position `particleID`, the velocity corrections along the `x`
67  * axis applied due to Lees-Edwards boundary conditions will be
68  * saved.
69  * @param[out] collisionCellIndices
70  * The collision cell indices for the particles.
71  * @param[in] boxSizeX
72  * The number of collision cells in the primary simulation volume
73  * along the `x` axis.
74  * @param[in] boxSizeY
75  * The number of collision cells in the primary simulation volume
76  * along the `y` axis.
77  * @param[in] boxSizeZ
78  * The number of collision cells in the primary simulation volume
79  * along the `z` axis.
80  */
82  const unsigned int particleID,
83  const MPCParticlePositionType* const gridShift_,
84  const FP shearRate,
85  const FP mpcTime,
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);
93 
94 /**
95  * Computes the collision cell mass and momentum contributions by the given
96  * particles.
97  *
98  * @param[in] particleCount The number of particles there are.
99  * @param[in] velocities The particle velocities.
100  * @param[in] collisionCellIndices The collision cell indices for the
101  * particles.
102  * @param[in,out] collisionCellMomenta The collision cell momenta. This array is
103  * expected to be long enough to store `3`
104  * coordinates for each collision cell, and
105  * is assumed to have each element set to
106  * `0` prior to calling this function.
107  * @param[in,out] collisionCellMasses The masses of the collision cells. This
108  * array is expected to be long enough to
109  * store `1` element for each collision
110  * cell, and is assumed to have each element
111  * set to `0` prior to calling this
112  * function.
113  * @param[in] particleMass The mass of any one particle.
114  */
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);
122 
123 /**
124  * Returns the center-of-mass momentum of a collision cell.
125  *
126  * @param[in] collisionCellIndex The index of the collision cell.
127  * @param[in] collisionCellMomenta The momenta of the collision cells.
128  */
129 Vector3D<MPCParticleVelocityType> getCollisionCellCenterOfMassMomentum(
130  const unsigned int collisionCellIndex,
131  const MPCParticleVelocityType* const collisionCellMomenta);
132 
133 /**
134  * Returns the center-of-mass velocity of a collision cell.
135  *
136  * @param[in] collisionCellIndex The index of the collision cell.
137  * @param[in] collisionCellMomenta The momenta of the collision cells.
138  * @param[in] collisionCellMasses The masses of the collision cells.
139  */
140 Vector3D<MPCParticleVelocityType> getCollisionCellCenterOfMassVelocity(
141  const unsigned int collisionCellIndex,
142  const MPCParticleVelocityType* const collisionCellMomenta,
143  const FP* const collisionCellMasses);
144 
145 /**
146  * Applies the first step of the SRD rotation to the given particles.
147  *
148  * @param[in] particleCount
149  * The number of particles there are.
150  * @param[in,out] velocities
151  * The particle velocities as input; after the function call, the
152  * values in this array will contain the velocities in the
153  * collision cell's center-of-mass frame.
154  * @param[in] particleMass
155  * The mass of any one particle.
156  * @param[in] collisionCellIndices
157  * An array which holds, for each particle, the index of the
158  * collision cell it is currently in.
159  * @param[in] collisionCellMomenta
160  * An array holding, for each collision cell, its three cartesian
161  * momentum coordinates.
162  * @param[in] collisionCellMasses
163  * An array holding, for each collision cell, its total mass
164  * content.
165  * @param[in] collisionCellRotationAxes
166  * An array holding, for each collision cell, the random rotation
167  * axis it is to use during the collision. Each axis consists of
168  * the `x`, `y`, and `z` coordinate of a normalized vector.
169  * @param[in] collisionAngle
170  * The rotation angle to use for the SRD collision.
171  * @param[out] collisionCellFrameInternalKineticEnergies
172  * An array, long enough to hold an element for each collision
173  * cell. The respective elements are set to the sum of the
174  * kinetic energies of the particles in that collision
175  * cell, as measured in that collision cell's center-of-mass
176  * frame. The array is not set to `0` by this function prior to
177  * execution of the algorithm.
178  * @param[out] collisionCellParticleCounts
179  * An array, long enough to hold an element for each collision
180  * cell. The respective elements are set to the number of
181  * particles in that collision cell. The array is assumed to
182  * contain only values `0` prior to calling this function.
183  */
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);
195 
196 /**
197  * Applies the first step of the SRD rotation to the given particles.
198  *
199  * @param[in] particleCount
200  * The number of particles there are.
201  * @param[in,out] velocities
202  * The particle new velocities as output; as input, the
203  * values in this array are to contain the rotated velocities in
204  * the collision cell's center-of-mass frame.
205  * @param[in] collisionCellIndices
206  * An array which holds, for each particle, the index of the
207  * collision cell it is currently in.
208  * @param[in] collisionCellMomenta
209  * An array holding, for each collision cell, its three cartesian
210  * momentum coordinates.
211  * @param[in] collisionCellMasses
212  * An array holding, for each collision cell, its total mass
213  * content.
214  * @param[in] collisionCellVelocityScalings
215  * For each collision cell, stores the factor by which to scale
216  * the relative velocities.
217  */
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);
225 
226 /**
227  * Undoes the velocity corrections applied by
228  * `sortIntoCollisionCellsLeesEdwards`.
229  *
230  * @param[in] particleCount
231  * The number of particles there are.
232  * @param[in,out] velocities
233  * The particle velocities.
234  * @param[in] velocityCorrections
235  * The velocity corrections as returned by
236  * `sortIntoCollisionCellsLeesEdwards`.
237  */
239  const unsigned int particleCount,
240  MPCParticleVelocityType* const velocities,
241  const MPCParticleVelocityType* const velocityCorrections);
242 
243 } //namespace ImplementationDetails
244 } //namespace OpenMPCD
245 
246 #endif //OPENMPCD_IMPLEMENTATIONDETAILS_SIMULATION_HPP
OpenMPCD::ImplementationDetails::getCollisionCellCenterOfMassMomentum
Vector3D< MPCParticleVelocityType > getCollisionCellCenterOfMassMomentum(const unsigned int collisionCellIndex, const MPCParticleVelocityType *const collisionCellMomenta)
Returns the center-of-mass momentum of a collision cell.
Definition: Simulation.cpp:102
Vector3D.hpp
OpenMPCD::ImplementationDetails::getCollisionCellIndex
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.
Definition: Simulation.cpp:12
OpenMPCD::ImplementationDetails::collisionCellStochasticRotationStep1
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.
Definition: Simulation.cpp:131
Types.hpp
OpenMPCD::ImplementationDetails::undoLeesEdwardsVelocityCorrections
void undoLeesEdwardsVelocityCorrections(const unsigned int particleCount, MPCParticleVelocityType *const velocities, const MPCParticleVelocityType *const velocityCorrections)
Undoes the velocity corrections applied by sortIntoCollisionCellsLeesEdwards.
Definition: Simulation.cpp:206
OpenMPCD::ImplementationDetails::collisionCellContributions
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.
Definition: Simulation.cpp:74
OpenMPCD::ImplementationDetails::sortIntoCollisionCellsLeesEdwards
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...
Definition: Simulation.cpp:36
OpenMPCD::ImplementationDetails::getCollisionCellCenterOfMassVelocity
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.
Definition: Simulation.cpp:113
OpenMPCD::ImplementationDetails::collisionCellStochasticRotationStep2
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.
Definition: Simulation.cpp:176