azplugins.compute¶
Overview
Compute velocity field in Cartesian coordinates. |
|
Compute velocity field in cylindrical coordinates. |
|
Center-of-mass velocity compute. |
|
Compute velocity field. |
Details
Computes.
- class hoomd.azplugins.compute.CartesianVelocityFieldCompute(num_bins, lower_bounds, upper_bounds, filter=None, include_mpcd_particles=False)¶
Bases:
VelocityFieldComputeCompute velocity field in Cartesian coordinates.
- Parameters:
num_bins (tuple[int]) – Number of bins along each of the 3 cylindrical coordinates. A value of zero indicates the coordinate is ignored.
lower_bounds (tuple[float]) – Lower bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
upper_bounds (tuple[float]) – Upper bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
filter (hoomd.filter.ParticleFilter) – HOOMD particles to include in calculation. The default value of
Nonemeans no HOOMD particles are included.include_mpcd_particles (bool) – If
True, include MPCD particles in the calculation. This argument only takes effect if HOOMD was compiled with the MPCD component.
CartesianVelocityFieldComputecalculates the mass-averaged velocity in a bin using the Cartesian coordinate system \((x, y, z)\). Particles that lie outside the lower and upper bounds are ignored.Example:
velocity_field = hoomd.azplugins.compute.CartesianVelocityFieldCompute( num_bins=(10, 8, 0), lower_bounds=(-5, 0, 0), upper_bounds=(5, 5, 0), filter=hoomd.filter.All(), )
- class hoomd.azplugins.compute.CylindricalVelocityFieldCompute(num_bins, lower_bounds, upper_bounds, filter=None, include_mpcd_particles=False)¶
Bases:
VelocityFieldComputeCompute velocity field in cylindrical coordinates.
- Parameters:
num_bins (tuple[int]) – Number of bins along each of the 3 cylindrical coordinates. A value of zero indicates the coordinate is ignored.
lower_bounds (tuple[float]) – Lower bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
upper_bounds (tuple[float]) – Upper bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
filter (hoomd.filter.ParticleFilter) – HOOMD particles to include in calculation. The default value of
Nonemeans no HOOMD particles are included.include_mpcd_particles (bool) – If
True, include MPCD particles in the calculation. This argument only takes effect if HOOMD was compiled with the MPCD component.
CylindricalVelocityFieldComputecalculates the mass-averaged velocity in a bin using the cylindrical coordinate system \((r, \theta, z)\), where \(0 \le \theta < 2\pi\). The cylindrical position coordinates are related to the Cartesian coordinates \((x, y, z)\) by\[\begin{split}r = \sqrt{x^2 + y^2} \\ \theta = \arctan\left(\frac{y}{x}\right) \\ z = z\end{split}\]Before averaging, Cartesian velocity vectors are converted to the cylindrical coordinate system using the change-of-basis matrix:
\[\begin{split}\left(\begin{array}{ccc} \cos \theta & \sin \theta & 0 \\ -\sin \theta & \cos \theta & 0 \\ 0 & 0 & 1 \end{array}\right)\end{split}\]Particles that lie outside the lower and upper bounds are ignored.
Example:
velocity_field = hoomd.azplugins.compute.CylindricalVelocityFieldCompute( num_bins=(10, 8, 0), lower_bounds=(0, 0, 0), upper_bounds=(10, 2*numpy.pi, 0), filter=hoomd.filter.All(), )
- class hoomd.azplugins.compute.VelocityCompute(filter=None, include_mpcd_particles=False)¶
Bases:
ComputeCenter-of-mass velocity compute.
- Parameters:
filter (hoomd.filter.ParticleFilter) – HOOMD particles to include. The default value of
Nonemeans no HOOMD particles are included.include_mpcd_particles (bool) – If
True, include MPCD particles in the calculation. This argument only takes effect if HOOMD was compiled with the MPCD component.
VelocityComputecalculates the center-of-mass velocity \(\mathbf{v}_{\rm cm}\) of a group of particles:\[\mathbf{v}_{\rm cm} = \dfrac{\sum_i m_i \mathbf{v}_i}{\sum_i m_i}\]where \(\mathbf{v}_i\) is the velocity and and \(m_i\) is the mass of particle i. These particles include regular HOOMD particles (specified with
filter) and MPCD particles (specified withinclude_mpcd_particles).Example:
v_com = hoomd.azplugins.compute.VelocityCompute( filter=hoomd.filter.All() ) v_com.velocity
- filter¶
HOOMD particles to include (read only).
- Type:
hoomd.filter.filter_like
- class hoomd.azplugins.compute.VelocityFieldCompute(num_bins, lower_bounds, upper_bounds, filter=None, include_mpcd_particles=False)¶
Bases:
ComputeCompute velocity field.
This class should not be instantiated directly. Use a derived type.
- Parameters:
num_bins (tuple[int]) – Number of bins along each of the 3 cylindrical coordinates. A value of zero indicates the coordinate is ignored.
lower_bounds (tuple[float]) – Lower bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
upper_bounds (tuple[float]) – Upper bounds for each coordinate. The value of this bound is ignored if the number of bins is zero.
filter (hoomd.filter.ParticleFilter) – HOOMD particles to include in calculation. The default value of
Nonemeans no HOOMD particles are included.include_mpcd_particles (bool) – If
True, include MPCD particles in the calculation. This argument only takes effect if HOOMD was compiled with the MPCD component.
- num_bins¶
Number of bins along each of the 3 cylindrical coordinates.
A value of zero indicates the coordinate is ignored.
- lower_bounds¶
Lower bounds for each coordinate.
The value of this bound is ignored if the number of bins is zero.
- upper_bounds¶
Upper bounds for each coordinate.
The value of this bound is ignored if the number of bins is zero.
- filter¶
HOOMD particles to include in calculation (read only).
Nonemeans no HOOMD particles are included.
- include_mpcd_particles¶
If
True, include MPCD particles in the calculation (read only).This property only has meaning effect if HOOMD was compiled with the MPCD component.
- Type:
- property coordinates¶
Coordinates of bin centers.
If binning is performed in more than 1 dimension, a multidimensional array is returned.
If binning is performed in 1 dimension, a 1 dimensional array is returned.
- Type:
- property velocities¶
Mass-averaged velocity vector of bin.
The velocities are returned as an array of 3-dimensional vectors matching the binning shape. This quantity is only available after the simulation has been run.
- Type: