amep.cluster.radius_of_gyration#

amep.cluster.radius_of_gyration(coords: ndarray, box_boundary: ndarray, mass: ndarray, pbc: bool = True, clusters: list | None = None) float#

Calculates the radius of gyration of the cluster containing the given coordinates. Uses the minimum image convention if pbc is True.

Notes:#

For a cluster composed of $n$ particles of masses $m_i, i=1,2, ldots, n$, located at fixed distances $s_i$ from the centre of mass, the radius of gyration is the square-root of the mass average of $s_i^2$ over all mass elements, i.e., $R_g=left(sum_{i=1}^n m_i s_i^2 / sum_{i=1}^n m_iright)^{1 / 2}$

Parameters:#

coordsnp.ndarray of shape (N,3)

Coordinates of all particles inside the cluster.

box_boundarynp.ndarray of shape (3,2)

Boundary of the simulation box in the form of np.array([[xmin, xmax], [ymin, ymax], [zmin, zmax]]).

massnp.ndarray of shape (N,)

Mass of all particles inside the cluster.

pbcbool, optional

If True, periodic boundary conditions are considered. The default is True.

clusterslist or None, optional

List of lists, where each list contains the indices of the particles that belong to the same cluster. If None, the geometric center of the given coords is calculated. If not None, the geometric center of each cluster in clusters is calculated, where coords must be the coordinates of the particles in clusters. The default is None.

returns:

Radius of gyration of the given points or each cluster in clusters.

rtype:

float or np.ndarray of shape (N,)

Examples

>>> import amep
>>> import numpy as np
>>> traj = amep.load.traj("../examples/data/lammps.h5amep")
>>> frame = traj[25]
>>> clusters, idx = amep.cluster.identify(
...     frame.coords(), frame.box, pbc=True
... )
>>> coords = frame.coords()[idx==1]
>>> mass = frame.mass()[idx==1]
>>> Rg = amep.cluster.radius_of_gyration(
...     coords, frame.box, mass, pbc=True
... )
>>> print(Rg)
8.16761958224252
>>> rgs = amep.cluster.radius_of_gyration(
...     frame.coords(), frame.box, frame.mass(),
...     pbc=True, clusters=clusters
... )
>>> centers = amep.cluster.geometric_center(
...     frame.coords(), frame.box, clusters=clusters
... )
>>> i = 1
>>> fig, axs = amep.plot.new(figsize=(3,3))
>>> amep.plot.particles(
...     axs, frame.coords()[idx==i], frame.box, color="orange",
...     radius=frame.radius()[idx==i]
... )
>>> amep.plot.particles(
...     axs, frame.coords()[idx!=i], frame.box, color="gray",
...     radius=frame.radius()[idx!=i]
... )
>>> axs.scatter(
...     centers[i,0], centers[i,1], color="k", s=50, marker='x',
...     label="geometric center"
... )
>>> angle = np.linspace(0, 2*np.pi, 150)
>>> x = rgs[i]*np.cos(angle)+centers[i,0]
>>> y = rgs[i]*np.sin(angle)+centers[i,1]
>>> axs.plot(
...     x, y, c="k", ls="--", lw=2, marker="", label="radius of gyration"
... )
>>> axs.set_xlabel(r"$x$")
>>> axs.set_ylabel(r"$y$")
>>> axs.legend(frameon=True)
>>> fig.savefig("./figures/cluster/cluster-radius_of_gyration.png")
../_images/cluster-radius_of_gyration.png