vasppy.rdf module
- class RadialDistributionFunction(structures: list[Structure], indices_i: list[int], indices_j: list[int] | None = None, nbins: int = 500, r_min: float = 0.0, r_max: float = 10.0, weights: list[float] | None = None)[source]
Bases:
objectClass for computing radial distribution functions.
- nbins
Number of bins.
- Type:
int
- range
Minimum and maximum values of r.
- Type:
(float, float)
- intervals
r values of the bin edges.
- Type:
np.array(float)
- dr
bin width.
- Type:
float
- r
mid-points of each bin.
- Type:
float
- rdf
RDF values.
- Type:
np.array(float)
- coordination_number
Volume integral of the RDF.
- Type:
np.array(float)
Initialise a RadialDistributionFunction instance.
- Parameters:
structures (list(pymatgen.Structure)) – List of pymatgen Structure objects.
indices_i (list(int)) – List of indices for species i.
indices_j (
list(int), optional) – List of indices for species j. Optional, default is None.nbins (
int, optional) – Number of bins used for the RDF. Optional, default is 500.rmin (
float, optional) – Minimum r value. Optional, default is 0.0.rmax (
float, optional) – Maximum r value. Optional, default is 10.0.weights (
list(float), optional) – List of weights for each structure. Optional, default is None.
- Returns:
None
- classmethod from_species_strings(structures: list[Structure], species_i: str, species_j: str | None = None, **kwargs) RadialDistributionFunction[source]
Initialise a RadialDistributionFunction instance by specifying species strings.
- Parameters:
structures (list(pymatgen.Structure)) – List of pymatgen Structure objects.
species_i (str) – String for species i, e.g.
"Na".species_j (
str, optional) – String for species j, e.g."Cl". Optional default is None.**kwargs – Variable length keyword argument list. See
vasppy.rdf.RadialDistributionFunction()for the full list of accepted arguments.
- Returns:
(RadialDistributionFunction)
- class VanHoveAnalysis(structures: list[Structure], indices: list[int], d_steps: int, nbins: int = 500, r_min: float = 0.0, r_max: float = 10.0)[source]
Bases:
objectClass for computing Van Hove correlation functions.
- nbins
Number of bins.
- Type:
int
- range
Minimum and maximum values of r.
- Type:
(float, float)
- intervals
r values of the bin edges.
- Type:
np.array(float)
- dr
bin width.
- Type:
float
- r
mid-points of each bin.
- Type:
float
- gsrt
Self part of the Van Hove correlation function.
- Type:
np.array(float)
- gdrt
Distinct part of the Van Hove correlation function.
- Type:
np.array(float)
Initialise a VanHoveCorrelationFunction instance.
- Parameters:
structures (list(pymatgen.Structure)) – List of pymatgen Structure objects.
indices (list(int)) – List of indices for species to consider.
d_steps (int) – number of steps between structures at dt=0 and dt=t.
nbins (
int, optional) – Number of bins used for the RDF. Optional, default is 500.rmin (
float, optional) – Minimum r value. Optional, default is 0.0.rmax (
float, optional) – Maximum r value. Optional, default is 10.0.
- Returns:
None
- distinct(sigma: float | None = None) ndarray[source]
Returns the distinct part of the Van Hove correlation function.
- Parameters:
sigma (
float, optional) – Optional smearing width.- Returns:
(np.ndarray)
- self(sigma: float | None = None) ndarray[source]
Returns the self part of the Van Hove correlation function.
- Parameters:
sigma (
float, optional) – Optional smearing width.- Returns:
(np.ndarray)
- smeared_gdrt(sigma: float = 0.1) ndarray[tuple[int, ...], dtype[floating[Any]]][source]
Smear the distinct part of the Van Hove correlation function with a Gaussian kernel.
- Parameters:
sigma (
float, optional) – Standard deviation for Gaussian kernel. Optional, default is 0.1.- Returns:
Smeared data.
- Return type:
(np.array)
- smeared_gsrt(sigma: float = 0.1) ndarray[tuple[int, ...], dtype[floating[Any]]][source]
Smear the self part of the Van Hove correlation function with a Gaussian kernel.
- Parameters:
sigma (
float, optional) – Standard deviation for Gaussian kernel. Optional, default is 0.1.- Returns:
Smeared data.
- Return type:
(np.array)
- shell_volumes(intervals: ndarray[tuple[int, ...], dtype[floating[Any]]]) ndarray[tuple[int, ...], dtype[floating[Any]]][source]
Volumes of concentric spherical shells.
- Parameters:
intervals (np.array) – N radial boundaries used to define the set of N-1 shells.
- Returns:
Volumes of each shell.
- Return type:
np.array