Spin Bath
BathArray
Documentation for the pycce.BathArray - central class, containing properties of the bath spins.
- class pycce.bath.array.BathArray(shape=None, array=None, names=None, hyperfines=None, quadrupoles=None, types=None, imap=None, ca=None, sn=None, hf=None, q=None, efg=None, state=None, center=1)
Subclass of
ndarraycontaining information about the bath spins.The subclass has fixed structured datatype:
_dtype_bath = np.dtype([('N', np.str_, 16), ('xyz', np.float64, (3,)), ('A', np.float64, (3, 3)), ('Q', np.float64, (3, 3))])
Accessing different fields results in the
ndarrayview.Each of the fields can be accessed as the attribute of the
BathArrayinstance and modified accordingly. In addition to the name fields, the information of the bath spin types is stored in thetypesattribute. All of the items intypescan be accessed as attributes of the BathArray itself.- Examples:
Generate empty
BathArrayinstance.>>> ba = BathArray((3,)) >>> print(ba) [('', [0., 0., 0.], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]]) ('', [0., 0., 0.], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]]) ('', [0., 0., 0.], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]], [[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]])]
Generate BathArray from the set of arrays:
>>> import numpy as np >>> ca = np.random.random((2, 3)) >>> sn = ['1H', '2H'] >>> hf = np.random.random((2, 3, 3)) >>> ba = BathArray(ca=ca, hf=hf, sn=sn) >>> print(ba.N, ba.types) ['1H' '2H'] SpinDict(1H: (1H, 0.5, 26.7519), 2H: (2H, 1, 4.1066, 0.00286))
Warning
Due to how structured arrays work, if one uses a boolean array to access an subarray, and then access the name field, the initial array will not change.
Example:
>>> ha = BathArray((10,), sn='1H') >>> print(ha.N) ['1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H'] >>> bool_mask = np.arange(10) % 2 == 0 >>> ha[bool_mask]['N'] = 'e' >>> print(ha.N) ['1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H' '1H']
To achieve the desired result, one should first access the name field and only then apply the boolean mask:
>>> ha['N'][bool_mask] = 'e' >>> print(ha.N) ['e' '1H' 'e' '1H' 'e' '1H' 'e' '1H' 'e' '1H']
Each bath spin can initiallized in some specific state accessing the
.stateattribute. It takes both state vectors and density matrices as values. See.stateattribute documentation for details.- Args:
shape (tuple): Shape of the array.
- array (array-like):
Either an unstructured array with shape (n, 3) containing coordinates of bath spins as rows OR structured ndarray with the same fields as the datatype of the bath.
- name (array-like):
Array of the bath spin name.
- hyperfines (array-like):
Array of the hyperfine tensors with shape (n, 3, 3).
- quadrupoles (array-like):
Array of the quadrupole tensors with shape (n, 3, 3).
- efg (array-like):
Array of the electric field gradients with shape (n, 3, 3) for each bath spin. Used to compute Quadrupole tensors for spins >= 1. Requires the spin types either be found in
common_isotopesor specified withtypesargument.- types (SpinDict):
SpinDict or input to create one. Contains either SpinTypes of the bath spins or tuples which will initialize those. See
pycce.bath.SpinDictdocumentation for details.- imap (InteractionMap):
Instance of InteractionMap containing user defined interaction tensors between bath spins stored in the array.
- ca (array-like):
Shorthand notation for
arrayargument.- sn (array-like):
Shorthand notation for
nameargument.- hf (array-like):
Shorthand notation for
hyperfinesargument.- q (array-like):
Shorthand notation for
quadrupolesargument.
- property A
ndarray: Array of hyperfine tensors for each spin in the array (
bath['A']).
- property N
ndarray: Array of name for each spin in the array (
bath['N']).
- property Q
ndarray: Array of quadrupole tensors for each spin in the array (
bath['Q']).
- add_interaction(i, j, tensor)
Add interactions tensor between bath spins with indexes
iandj.Note
If called from the subarray this method does not change the tensors of the total BathArray.
- Args:
- i (int or ndarray (n,) ):
Index of the first spin in the pair or array of the indexes of the first spins in n pairs.
- j (int or ndarray with shape (n,)):
Index of the second spin in the pair or array of the indexes of the second spins in n pairs.
- tensor (ndarray with shape (3,3) or (n, 3,3)):
Interaction tensor between the spins i and j or array of tensors.
- add_single_jump(operator, rate=1, units='rad', square_root=False, which=None)
Add single-spin jump operator for the given type of spins to be used in the Lindbladian master equation CCE.
- Args:
- operator (str or ndarray with shape (dim, dim)): Definition of the operator. Can be either of the following:
Pair of integers defining the Sven operator.
String where each symbol corresponds to the spin matrix or operation between them. Allowed symbols:
xyz+. If there is nothing between symbols, assume multiplication of the operators. If there is a+symbol, assume summation between terms. For example,xx+zwould correspond to the operator \(\hat S_x \hat S_x + \hat S_z\).String equal to
A. Then assumes that the correct matrix form of the operator has been provided by the user.
rate (float): Rate associated with the given jump operator. By default, is given in rad ms^-1. units (str): Units of the rate, can be either rad (for radial frequency units) or deg
(for rotational frequency).
- square_root (bool): True if the rate is given as a square root of the rate (to match how one sets up
collapse operators in Qutip). Default False.
- which (str): For which type of the spins add the jump operator. Default is None - if
there is only one spin type in the array then the jump operator is added, otherwise the exception is raised.
- add_type(*args, **kwargs)
Add spin type to the
typesdictionary.- Args:
For details and allowed inputs see
SpinDictdocumentation.- Returns:
SpinDict: A view of
self.typesinstance.
- property detuning
- ndarray: Array of the
detuning attribute for each spin in the array from
typesdictionary.
- ndarray: Array of the
- property dim
- ndarray: Array of the
dim(dimensions of the spin) attribute for each spin in the array from
typesdictionary.
- ndarray: Array of the
- dist(position=None)
Compute the distance of the bath spins from the given position.
- Args:
- position (ndarray with shape (3,)):
Cartesian coordinates of the position from which to compute the distance. Default is (0, 0, 0).
- Returns:
ndarray with shape (n,): Array of distances of each bath spin from the given position in angstrom.
- from_center(center, inplace=True, cube=None, which=0, **kwarg)
Generate hyperfine couplings using either the point dipole approximation or spin density in the .cube format, with the information from the CenterArray instance.
- Args:
center (CenterArray): Array, containing properties of the central spin
inplace (bool): True if changes parameters of the array in place. If False, returns copy of the array.
- cube (Cube or iterable of Cubes): An instance of
Cubeobject, which contains spatial distribution of spin density of central spins. For details see documentation of
Cubeclass.- which (int): If
cubeis a single Cube instance, this is an index of the central spin it corresponds to.
**kwarg: Additional arguments for .from_cube method.
- cube (Cube or iterable of Cubes): An instance of
- Returns:
BathArray: Updated BathArray instance.
- from_cube(cube, gyro_center=-17608.59705, inplace=True, which=0, **kwargs)
Generate hyperfine couplings, assuming that bath spins interaction with central spin can be approximated as a point dipole, interacting with given spin density distribution.
- Args:
- cube (Cube): An instance of Cube object, which contains spatial distribution of spin density.
For details see documentation of Cube class.
gyro_center (float): Gyromagnetic ratio of the central spin.
inplace (bool): True if changes parameters of the array in place. If False, returns copy of the array.
- Returns:
BathArray: Updated BathArray instance with changed hyperfine couplings.
- from_efg(efg, inplace=True)
Generate quadrupole splittings from electric field gradient tensors for spins >= 1.
- Args:
- efg (array-like): Array of the electric field gradients for each bath spin. The data for spins-1/2
should be included but can be any value.
inplace (bool): True if changes parameters of the array in place. If False, returns copy of the array.
- Returns:
BathArray: Updated BathArray instance with changed quadrupole tensors.
- from_func(func, *args, inplace=True, **kwargs)
Generate hyperfine couplings from user-defined function.
Args:
- Returns:
BathArray: Updated BathArray instance with changed hyperfine couplings.
- from_point_dipole(position, gyro_center=-17608.59705, inplace=True)
Generate hyperfine couplings, assuming that bath spins interaction with central spin is the same as the one between two magnetic point dipoles.
- Args:
position (ndarray with shape (3,)): position of the central spin
- gyro_center (float or ndarray with shape (3,3)):
gyromagnetic ratio of the central spin
OR
tensor corresponding to interaction between magnetic field and central spin.
inplace (bool): True if changes parameters of the array in place. If False, returns copy of the array.
- Returns:
BathArray: Updated BathArray instance with changed hyperfine couplings.
- property gyro
- ndarray: Array of the
gyro(gyromagnetic ratio) attribute for each spin in the array from
typesdictionary.
- ndarray: Array of the
- property h
dict: Dictionary with additional spin Hamiltonian parameters. Key denotes the product of spin operators as:
Either a string containing
x, y, z, +, -where each symbol is a corresponding spin operator:x== \(S_x\)y== \(S_y\)z== \(S_z\)p== \(S_+\)m== \(S_-\)
Several symbols is a product of those spin operators.
Or a tuple with indexes (k, q) for Stevens operators (see https://www.easyspin.org/documentation/stevensoperators.html).
The item is the coupling parameter in float.
Examples:
d['pm'] = 2000corresponds to the Hamiltonian term \(\hat H_{add} = A \hat S_+ \hat S_-\) with \(A = 2\) MHz.d[2, 0] = 1.5e6corresponds to Stevens operator \(B^q_k \hat O^q_k = 3 \hat S_z - s(s+1) \hat I\) with \(k = 2\), \(q = 0\), and \(B^q_k = 1.5\) GHz.
- property has_state
ndarray: Bool array. True if given spin was initialized with a state, False otherwise.
- property name
ndarray: Array of the
nameattribute for each spin in the array fromtypesdictionary.Note
While the value of this attribute should be the same as the
Nfield of the BathArray instance,.nameshould not be used for production as it creates a new array fromtypesdictionary.
- property nc
int: Number of central spins.
- property proj
ndarray: Array of \(S_z\) projections of the bath spin states.
- property q
- ndarray: Array of the
q(quadrupole moment) attribute for each spin in the array from
typesdictionary.
- ndarray: Array of the
- property s
ndarray: Array of the
spin(spin value) attribute for each spin in the array fromtypesdictionary.
- savetxt(filename, fmt='%18.8f', strip_isotopes=False, **kwargs)
Save name of the isotopes and their coordinates to the txt file of xyz format.
- Args:
filename (str or file): Filename or file handle. fmt (str): Format of the coordinate entry. strip_isotopes (bool): True if remove numbers from the name of bath spins. Default False. **kwargs: Additional keywords of the
numpy.savetxtfunction.
- sort(axis=-1, kind=None, order=None)
Sort array in-place. Is implemented only when imap is None. Otherwise use
np.sort.
- property state
BathState: Array of the bath spin states.
Can have three types of entries:
None. If entry is None, assumes fully random density matrix. Default value.
ndarray with shape (s,). If entry is vector, corresponds to the pure state of the spin.
ndarray with shape (s, s). If entry is a matrix, corresponds to the density matrix of the spin.
Examples:
>>> print(ba.state) [None None] >>> ba[0].state = np.array([0, 1]) >>> print(ba.state) [array([0, 1]) None]
- update(ext_bath, error_range=0.2, ignore_isotopes=True, inplace=True)
Update the properties of the spins in the array using data from other
BathArrayinstance. For each spin inext_bathcheck whether there is such spin in the array that has the same position within allowed error range given byerror_rangeand has the same name. If such spins is found in the array, then it’s coordinates, hyperfine tensor and quadrupole tensor are updated using the values of the spin in theext_bathobject.If
ignore_isotopesis true, then the name check ignores numbers in the name of the spins.- Args:
ext_bath (BathArray): Array of the new spins.
- error_range (float): +- distance in Angstrom within which two positions are considered to be the same.
Default is 0.2 A.
ignore_isotopes (bool): True if ignore numbers in the name of the spins. Default True.
inplace (bool): True if changes parameters of the array in place. If False, returns copy of the array.
- Returns:
BathArray: updated BathArray instance.
- property x
ndarray: Array of x coordinates for each spin in the array (
bath['xyz'][:, 0]).
- property xyz
ndarray: Array of coordinates for each spin in the array (
bath['xyz']).
- property y
ndarray: Array of y coordinates for each spin in the array (
bath['xyz'][:, 1]).
- property z
ndarray: Array of z coordinates for each spin in the array (
bath['xyz'][:, 2]).
- pycce.bath.array.argsort(a, *args, **kwargs)
Return
aindexes of a sorted array. Overridesnumpy.argsortfunction.
- pycce.bath.array.broadcast_array(array, root=0)
Using mpi4py broadcast
BathArrayorCenterArrayto all processes. Args:array (BathArray or CenterArray): Array to broadcast. root (int): Rank of the process to broadcast from.
- Returns:
BathArray or CenterArray: Broadcasted array.
- pycce.bath.array.check_gyro(gyro)
Check if gyro is matrix or scalar.
- Args:
gyro (ndarray or float): Gyromagnetic ratio matrix or float.
- Returns:
tuple: tuple containing:
ndarray or float: Gyromagnetic ratio.
bool: True if gyro is float, False otherwise.
- pycce.bath.array.point_dipole(pos, gyro_array, gyro_center)
Generate an array hyperfine couplings, assuming point dipole approximation.
- Args:
pos (ndarray with shape (n, 3)): Relative position of the bath spins. gyro_array (ndarray with shape (n,)): Array of the gyromagnetic ratios of the bath spins.
- gyro_center (float or ndarray with shape (3, 3)):
gyromagnetic ratio of the central spin
OR
tensor corresponding to interaction between magnetic field and central spin.
- Returns:
ndarray with shape (n, 3, 3): Array of hyperfine tensors.
- pycce.bath.array.same_bath_indexes(barray_1, barray_2, error_range=0.2, ignore_isotopes=True)
Find indexes of the same array elements in two
BathArrayinstances.- Args:
barray_1 (BathArray): First array. barray_2 (BathArray): Second array. error_range (float): If distance between positions in two arrays is smaller than
error_rangethey are assumed to be the same.
ignore_isotopes (bool): True if ignore numbers in the name of the spins. Default True.
- Returns:
tuple: tuple containing:
ndarray: Indexes of the elements in the first array found in the second.
ndarray: Indexes of the elements in the second array found in the first.
- pycce.bath.array.sort(a, axis=-1, kind=None, order=None)
Return a sorted copy of an array. Overrides numpy.sort function.
- utilities.rotmatrix(final_vector)
Generate 3D rotation matrix which applied on initial vector will produce vector, aligned with final vector.
Examples:
>>> R = rotmatrix([0,0,1], [1,1,1]) >>> R @ np.array([0,0,1]) array([0.577, 0.577, 0.577])
- Args:
initial_vector (ndarray with shape(3, )): Initial vector. final_vector (ndarray with shape (3, )): Final vector.
- Returns:
ndarray with shape (3, 3): Rotation matrix.
BathState
- class pycce.bath.state.BathState(size)
Class for storing the state of the bath spins. Usually is not generated directly, but accessed as an
BathArray.stateattribute.- Args:
size (int): Number of bath states to be stored.
- any(*args, **kawrgs)
Returns the output of
.has_state.anymethod. Args:- Returns:
bool: True if any entry is initialized. Otherwise False.
- gen_pure(rho, dim)
Generate pure states from the \(S_z\) projections to be stored in the given
BathStateobject.- Args:
rho (ndarray with shape (n, )): Array of the desired projections. dim (ndarray with shape (n,) ): Array of the dimensions of the spins.
- Returns:
BathState: View of the
BathStateobject.
- property has_state
ndarray: Bool property. True if given element was initialized as a state, False otherwise.
- property proj
ndarray: Projections of bath states on \(S_z\).
- project(rotation=None)
Generate projections of bath states on \(S_z\).
- Args:
- rotation (optional, ndarray with shape (3, 3)):
Matrix used to transform \(S_z\) matrix as \(S_z' = R^{\dagger} S_z R\).
- Returns:
ndarray with shape (n, ): Array with projections of the state.
- property pure
ndarray: Bool property. True if given entry is a pure state, False otherwise.
- property shape
tuple: Shape of the BathState underlying array.
- property size
int: Size of the BathState underlying array.
- property state
ndarray: Return an underlying object array.
Cube
- class pycce.bath.cube.Cube(filename)
Class to process the .cube datafiles with spin polarization.
- Args:
filename (str): Name of the .cube file.
- Attributes:
comments (str): First two lines of the .cube file. origin (ndarray with shape (3,)): Coordinates of the origin in angstrom. voxel (ndarray with shape (3,3)): Parameters of the voxel - unit of the 3D grid in angstrom. size (ndarray with shape (3,)): Size of the cube. atoms (BathArray with shape (n)): Array of atoms in the cube. data (ndarray with shape (size[0], size[1], size[2]): Data stored in cube. grid (ndarray with shape (size[0], size[1], size[2], 3): Coordinates of the points at which data is computed. integral (float): Data integrated over cube. spin (float): integral / 2 - total spin.
- integrate(position, gyro_n, gyro_e=-17608.59705, spin=None, parallel=False, root=0)
Integrate over polarization data, stored in Cube object, to obtain hyperfine dipolar-dipolar tensor.
- Args:
- position (ndarray with shape (3,) or (n, 3)): Position of the bath spin at which to compute
hyperfine tensor or array of positions.
gyro_n (float or ndarray with shape (n,) ): Gyromagnetic ratio of the bath spin or array of the ratios. gyro_e (float): Gyromagnetic ratio of central spin. spin (float): Total spin of the central spin. If not given, taken from the integral of the polarization.
- Returns:
ndarray with shape (3, 3) or (n, 3, 3): Hyperfine tensor or array of hyperfine tensors.
- transform(rotmatrix=None, shift=None)
Changes coordinates of the grid. DOES NOT ASSUME PERIODICITY.
- Args:
rotmatrix (ndarray with shape (3, 3)): Rotation matrix R:
\[\begin{split}R = &[n_1^{(1)} n_1^{(2)} n_1^{(3)}]\\ &[n_2^{(1)} n_2^{(2)} n_2^{(3)}]\\ &[n_3^{(1)} n_3^{(2)} n_3^{(3)}]\end{split}\]where \(n_i^{(j)}\) corresponds to the coefficient of initial basis vector \(i\) for \(j\) new basis vector:
\[e'_j = n_1^{(j)} \vec{e}_1 + n_2^{(j)} \vec{e}_2 + n_3^{(j)} \vec{e}_3\]In other words, columns of R are coordinates of the new basis in the old basis.
Given vector in initial basis v = [v1, v2, v3], vector in new basis is given as v’ = R.T @ v.
shift (ndarray with shape (3,)): Shift in the origin of coordinates (in the old basis).
SpinDict and SpinType
Documentation for the SpinDict - dict-like class which describes
the properties of the different types of the spins in the bath.
- class pycce.SpinDict(*args, **kwargs)
Wrapper class for dictionary tailored for containing properties of the spin types. Can take
np.voidorBathArrayinstances as keys. Every entry is instance of theSpinType.Each entry of the
SpinDictcan be initianlized as follows:As a Tuple containing name (optional), spin, gyromagnetic ratio, quadrupole constant (optional) and detuning (optional).
As a
SpinTypeinstance.
Examples:
>>> types = SpinDict() >>> types['1H'] = ('1H', 1 / 2, 26.7519) >>> types['2H'] = 1, 4.1066, 0.00286 >>> types['3H'] = SpinType('3H', 1 / 2, 28.535, 0) >>> print(types) SpinDict({'1H': (1H, 0.5, 26.7519, 0.0), '2H': (2H, 1, 4.1066, 0.00286), '3H': (3H, 0.5, 28.535, 0)})
If
SpinTypeof the given bath spin is not provided, when requestedSpinDictwill try to find information about the bath spins in thecommon_isotopes.If found, adds an entry to the given
SpinDictinstance and returns it. OtherwiseKeyErroris raised.To initiallize several
SpinTypeentries one can useadd_typesmethod.- Args:
*args: Any numbers of arguments which could initialize
SpinTypeinstances. **kwargs: Any numbers of keyword arguments which could initializeSpinTypeinstances.For details see
SpinDict.add_typemethod.
- add_type(*args, **kwargs)
Add one or several spin types to the spin dictionary.
- Args:
- *args:
Any numbers of arguments which could initialize
SpinTypeinstances. Accepted arguments:Tuple containing name, spin, gyromagnetic ratio, quadrupole constant (optional) and detuning (optional).
SpinTypeinstance.
Can also initialize one instance of
SpinTypeif each argument corresponds to each positional argument necessary to initiallize.- **kwargs: Any numbers of keyword arguments which could initialize
SpinTypeinstances. Usefull as an alternative for updating the dictionary. for each keyword argument adds an entry to the
SpinDictwith the same name as keyword.
- Examples:
>>> types = SpinDict() >>> types.add_type('1H', 1 / 2, 26.7519) >>> types.add_type(('1H_det', 1 / 2, 26.7519, 10), ('2H', 1, 4.1066, 0.00286), >>> SpinType('3H', 1 / 2, 28.535, 0), e=(1 / 2, 6.7283, 0)) >>> print(types) SpinDict(1H: (1H, 0.5, 26.7519), 1H_det: (1H_det, 0.5, 26.7519, 10), 2H: (2H, 1, 4.1066, 0.00286), 3H: (3H, 0.5, 28.535), e: (e, 0.5, 6.7283))
- class pycce.SpinType(name, s=0.0, gyro=0.0, q=0.0, detuning=0.0)
Class which contains properties of each spin type in the bath.
- Args:
name (str): Name of the bath spin. s (float): Total spin of the bath spin.
Default 0.
gyro (float): Gyromagnetic ratio in rad * kHz / G.
Default 0.
q (float): Quadrupole moment in barn (for s > 1/2).
Default 0.
- detuning (float): Energy detuning from the zeeman splitting in kHz,
included as an extra \(+\omega \hat S_z\) term in the Hamiltonian, where \(\omega\) is the detuning.
Default 0.
Attributes:
name (str): Name of the bath spin. s (float): Total spin of the bath spin. dim (int): Spin dimensionality = 2s + 1.
gyro (float): Gyromagnetic ratio in rad/(ms * G). q (float): Quadrupole moment in barn (for s > 1/2). detuning (float): Energy detuning from the zeeman splitting in kHz.
- property h
dict: Dictionary with additional spin Hamiltonian parameters. Key denotes the product of spin operators as:
Either a string containing
x, y, z, +, -where each symbol is a corresponding spin operator:x== \(S_x\)y== \(S_y\)z== \(S_z\)p== \(S_+\)m== \(S_-\)
Several symbols is a product of those spin operators.
Or a tuple with indexes (k, q) for Stevens operators (see https://www.easyspin.org/documentation/stevensoperators.html).
The item is the coupling parameter in float.
Examples:
d['pm'] = 2000corresponds to the Hamiltonian term \(\hat H_{add} = A \hat S_+ \hat S_-\) with \(A = 2\) MHz.d[2, 0] = 1.5e6corresponds to Stevens operator \(B^q_k \hat O^q_k = 3 \hat S_z - s(s+1) \hat I\) with \(k = 2\), \(q = 0\), and \(B^q_k = 1.5\) GHz.
- pycce.bath.array.common_isotopes = SpinDict(1H: (0.5, 26.7522), 2H: (1.0, 4.1066, 0.0029), 3He: (0.5, -20.3789), ...)
SpinDict: An instance of the
SpinDictdictionary, containing properties for the most of the common isotopes with nonzero spin. The isotope is considered common if it is stable and has nonzero concentration in nature.
- pycce.bath.array.common_concentrations = {element ('H', 'He',...) : { isotope ('1H', '2H', ..) : concentration}}
dict: Nested dict containing natural concentrations of the stable nuclear isotopes.
Random bath
Documentation for the pycce.random_bath function, used to generate random bath.
- pycce.bath.cell.random_bath(names, size, number=1000, density=None, types=None, density_units='cm-3', center=None, seed=None)
Generate random bath containing spins with names provided with argument
namein the box of sizesize. By default generates coordinates in range (-size/2; +size/2) but this behavior can be changed by providingcenterkeyword.Examples:
Generate 2000 \(^{13}\mathrm{C}\) nuclear spins in the cubic box with the side of 100 angstrom:
>>> atoms = random_bath('13C', 100, number=2000, seed=10) >>> print(atoms.size) 2000 >>> print(round(atoms.x.min()), round(atoms.x.max())) -50.0 50.0
Generate electron spin bath with density \(10^{17} \mathrm{cm}^{-3}\) in the cuboid box:
>>> electrons = random_bath('e', [1e3, 2e3, 3e3], density=1e17, >>> density_units='cm-3', seed=10) >>> print(electrons.size, round(electrons.x.min()), round(electrons.x.max())) 600 -494.0 500.0 >>> print(electrons.types) SpinDict(e: (e, 0.5, -17608.59705))
- Args:
names (str or array-like with length n): Name of the bath spin or array with the names of the bath spins, size (float or ndarray with shape (3,)): Size of the box. If float is given,
assumes 3D cube with the edge =
size. Otherwise the size specifies the dimensions of the box. Dimensionality is controlled by setting entries of the size array to 0.- number (int or array-like with length n): Number of the bath spins in the box
or array with the numbers of the bath spins. Has to have the same length as the
namearray.- density (float or array-like with length n): Concentration of the bath spin
or array with the concentrations. Has to have the same length as the
namearray.
types (SpinDict): Dictionary with SpinTypes or input to create one.
- density_units (str): If number of spins provided as density, defines units.
Values are accepted in the format
m, orm^xorm-xwhere m is the length unit, x is dimensionality of the bath (e.g. x = 1 for 1D, 2 for 2D etc). If onlymis provided the dimensions are inferred fromsizeargument. Accepted length units:mmeters;cmcentimeters;aangstroms.
- center (ndarray with shape (3,)): Coordinates of the (0, 0, 0) point of the final coordinate system
in the initial coordinates. Default is
size / 2- center is in the middle of the box.
seed (int): Seed for random number generator.
- Returns:
BathArray with shape (np.prod(number),): Array of the bath spins with random positions.
BathCell
Documentation for the pycce.BathCell - class for convenient generation of BathArray and the
necessary helper functions.
- class pycce.bath.cell.BathCell(a=None, b=None, c=None, alpha=None, beta=None, gamma=None, angle='rad', cell=None)
Generator of the bath spins positions from the unit cell of the material.
- Args:
a (float): a parameter of the primitive cell. b (float): b parameter of the primitive cell. c (float): c parameter of the primitive cell. alpha (float): \(\alpha\) angle of the primitive cell. beta (float): \(\beta\) angle of the primitive cell. gamma (float): \(\gamma\) angle of the primitive cell. angle (str): units of the \(\alpha\), \(\beta\), \(\gamma\) angles.
Can be either radians (
'rad'), or degrees ('deg').cell (ndarray with shape (3, 3)): Parameters of the cell.
cellis 3x3 matrix with columns of coordinates of crystallographic vectors in the cartesian reference frame. Seecellattribute.If provided, overrides a, b, and c.
Attributes:
- cell (ndarray with shape (3, 3)): Parameters of the cell.
cellis 3x3 matrix with entries:\[\begin{split}[&[a_x\ b_x\ c_x]\\ &[a_y\ b_y\ c_y]\\ &[a_z\ b_z\ c_z]]\end{split}\]where a, b, c are crystallographic vectors and x, y, z are their coordinates in the cartesian reference frame.
atoms (dict): Dictionary containing coordinates and occupancy of each lattice site:
{atom_1: [array([x1, y1, z1]), array([x2, y2, z2])], atom_2: [array([x3, y3, z3]), ...]}
- isotopes (dict):
Dictionary containing spin types and their concentration for each lattice site type:
{atom_1: {spin_1: concentration, spin_2: concentration}, atom_2: {spin_3: concentration ...}}
where
atom_iare lattice site types, andspin_iare spin types.
- add_atoms(*args, type='cell')
Add coordinates of the lattice sites to the unit cell.
- Args:
- args (tuple): List of tuples, each containing the type of atom N (*str),
and the xyz coordinates in the format (float, float, float):
(N, [x, y, z]).
type (str): Type of coordinates. Can take values of
['cell', 'angstrom'].If
type="cell", assumes crystallographic coordinates.If
type="angstrom"assumes that coordinates are given in the cartresian reference frame.
Returns:
- dict:
View of
cell.atomsdictionary, where each key is the type of lattice site, and each value is the list of coordinates in crystallographic frame.
Examples:
>>> cell = BathCell(10) >>> cell.add_atoms(('C', [0, 0, 0]), ('C', [5, 5, 5]), type='angstrom') >>> cell.add_atoms(('Si', [0, 0.5, 0.]), type='cell') >>> print(cell.atoms) {'C': [array([0., 0., 0.]), array([0.5, 0.5, 0.5])], 'Si': [array([0. , 0.5, 0. ])]}
- add_isotopes(*args)
Add spins that can populate each lattice site type.
Args:
*args (tuple or list of tuples): Each tuple can have any of the following formats:
Name of the lattice site N (str), name of the spin X (str), concentration c (float, in decimal):
(N, X, c).Isotope name X and concentration `c:
(X, c).In this case, the name of the isotope is given in the format
"{}{}".format(digits, atom_name)wheredigitsis any set of digits 0-9,atom_nameis the name of the corresponding lattice site. Convenient when generating nuclear spin bath.
- Returns:
- dict:
View of
cell.isotopesdictionary which contains information about lattice site types, spin types, and their concentrations:{atom_1: {spin_1: concentration, spin_2: concentration}, atom_2: {spin_3: concentration ...}}
Examples:
>>> cell = BathCell(10) >>> cell.add_atoms(('C', [0, 0, 0]), ('C', [5, 5, 5]), type='angstrom') >>> cell.add_isotopes(('C', 'X', 0.001), ('13C', 0.0107)) >>> print(cell.isotopes) {'C': {'X': 0.001, '13C': 0.0107}}
- classmethod from_ase(atoms_object)
Generate
BathCellinstance fromase.Atomsobject of Atomic Simulations Environment (ASE) package.- Args:
atoms_object (Atoms): Atoms object, used to generate new
BathCellinstance.- Returns:
BathCell: New instance of the
BathCellwith atoms read fromase.Atoms.
- gen_supercell(size, add=None, remove=None, seed=None, recenter=True)
Generate supercell populated with spins.
Note
If
isotopeswere not provided, assumes the natural concentration of nuclear spin isotopes for each lattice site type. However, if any isotope concentration is provided, then uses only user-defined ones.- Args:
- size (float): Approximate linear size of the supercell. The generated supercell will have
minimal distance between opposite sides larger than this parameter.
- add (tuple or list of tuples):
Tuple or list of tuples containing common_isotopes to add as a defect. Each tuple contains name of the new isotope and its coordinates in the cell basis:
(isotope_name, x_cell, y_cell, z_cell).- remove (tuple or list of tuples):
Tuple or list of tuples containing bath to remove in the defect. Each tuple contains name of the atom to remove and its coordinates in the cell basis:
(atom_name, x_cell, y_cell, z_cell).
seed (int): Seed for random number generator.
- recenter (bool):
True if place approximate center of the supercell at (0,0,0). False if start supercell at (0, 0, 0). Default True.
Note
While
addtakes the spin name as an argument,removetakes the lattice site name.- Returns:
BathArray: Array of the spins in the given supercell.
- rotate(rotation_matrix)
Rotate the BathCell using the rotation matrix provided.
- Args:
- rotation_matrix (ndarray with shape (3,)): Rotation matrix R which rotates the old basis of the
cartesian reference frame to the new basis.
- set_zdir(direction, type='cell')
Set z-direction of the cell.
- Args:
direction (ndarray with shape (3,)): Direction of the z axis.
- type (str): How coordinates in
directionare stored. Iftype="cell", assumes crystallographic coordinates. If
type="angstrom"assumes that z direction is given in the cartresian reference frame.
- type (str): How coordinates in
- to_cartesian(coord)
Transform coordinates from crystallographic basis to the cartesian reference frame.
- Args:
coord (ndarray with shape (3,) or (n, 3)): Coordinates in crystallographic basis or array of coordinates.
- Returns:
ndarray with shape (3,) or (n, 3): Cartesian coordinates in angstrom.
- to_cell(coord)
Transform coordinates from the cartesian coordinates of the reference frame to the cell coordinates.
- Args:
coord (ndarray with shape (3,) or (n, 3)): Cartesian coordinates in angstrom or array of coordinates.
- Returns:
ndarray with shape (3,) or (n, 3): Coordinates in the cell basis.
- property zdir
ndarray: z-direction of the reference cartesian coordinate frame in cell coordinates.
- pycce.bath.cell.defect(cell, atoms, add=None, remove=None)
Generate a defect in the given supercell.
The defect will be located in the unit cell, located roughly in the middle of the supercell, generated by
BathCell, such that (0, 0, 0) of cartesian reference frame is located at (0, 0, 0) position of this unit cell.- Args:
cell (ndarray with shape (3, 3)): parameters of the unit cell.
atoms (BathArray): Array of spins in the supercell.
- add (tuple or list of tuples):
Add spin type(s) to the supercell at specified positions to create point defect. Each tuple contains name of the new isotope and its coordinates in the cell basis:
(isotope_name, x_cell, y_cell, z_cell).- remove (tuple or list of tuples):
Remove lattice site from the supercell at specified position to create point defect. Each tuple contains name of the atom to remove and its coordinates in the cell basis:
(atom_name, x_cell, y_cell, z_cell).
- Returns:
BathArray: Array of spins with the defect added.
- pycce.bath.cell.read_ase(atoms_object)
Generate
BathCellinstance fromase.Atomsobject of Atomic Simulations Environment (ASE) package.- Args:
atoms_object (Atoms): Atoms object, used to generate new
BathCellinstance.- Returns:
BathCell: New instance of the
BathCellwith atoms read fromase.Atoms.