| |
- builtins.object
-
- Atom
- Atoms
- Cell
- Crystal
-
- Superstructure
- Symmetry
class Atom(builtins.object) |
|
Atom(u, v, w, atom_type, label, occupancy=1.0, uiso=0.001, mxmymz=None)
Atom class
Contains site information |
|
Methods defined here:
- __init__(self, u, v, w, atom_type, label, occupancy=1.0, uiso=0.001, mxmymz=None)
- Initialize self. See help(type(self)) for accurate signature.
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- info(self)
- Display atomic properties
- total_moment(self)
- Return the total moment along a, b, c directions
- uvw(self)
- Returns a [1x3] array of current positions
:return: np.array([1x3])
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class Atoms(builtins.object) |
|
Atoms(u=[0], v=[0], w=[0], type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
Contains properties of atoms within the crystal
Each atom has properties:
u,v,w >> atomic coordinates, in the basis of the unit cell
type >> element species, given as element name, e.g. 'Fe'
label >> Name of atomic position, e.g. 'Fe1'
occupancy >> Occupancy of this atom at this atomic position
uiso >> atomic displacement factor (ADP) <u^2>
mxmymz >> magnetic moment direction [x,y,z] |
|
Methods defined here:
- __call__(self, u=[0], v=[0], w=[0], type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Re-initialises the class, generating new atomic positions
- __getitem__(self, idx)
- __init__(self, u=[0], v=[0], w=[0], type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Initialisation, defines Atoms defaults
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- addatom(self, u=0, v=0, w=0, type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Adds a new atom
:param u:
:param v:
:param w:
:param type:
:param label:
:param occupancy:
:param uiso:
:param mxmymz:
:return:
- atom(self, idx)
- Create Atom object for atom site
- changeatom(self, idx=None, u=None, v=None, w=None, type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Change an atoms properties
:param idx:
:param u:
:param v:
:param w:
:param type:
:param label:
:param occupancy:
:param uiso:
:param mxmymz:
:return: None
- check(self)
- Checks the validity of the contained attributes
:return: None
- findatom(self, u=None, v=None, w=None, type=None, label=None, occupancy=None, uiso=None, mxmymz=None, tol=0.01)
- Find atom using parameters, return idx
:param u: float
:param v: float
:param w: float
:param type: str
:param label: str
:param occupancy: float
:param uiso: float
:param mxmymz: [mx,my,mz]
:param tol: float, tolerance to match value
:return: array of indexes
- fitincell(self)
- Adjust all atom positions to fit within unit cell
- fromcif(self, cifvals)
- Import atom parameters from a cif dictionary
Required cif keys:
_atom_site_label
_atom_site_fract_x
_atom_site_fract_y
_atom_site_fract_z
Optional cif keys:
_atom_site_type_symbol
_atom_site_U_iso_or_equiv
_atom_site_B_iso_or_equiv
_atom_site_occupancy
_atom_site_moment_label
_atom_site_moment_crystalaxis_x
_atom_site_moment_crystalaxis_y
_atom_site_moment_crystalaxis_z
:param cifvals: dict
:return: none
- generate_lattice(self, U=1, V=1, W=0, centred=True)
- Expand the atomic positions beyond the unit cell, creating a lattice
uvw,type,label,occ,uiso,mxmymz = self.generate_lattice(U,V,W,centred)
U,V,W = maximum lattice index to loop to
centred = if True, positions will loop from e.g. -U to U,
otherwise, will loop from e.g. 0 to U
uvw,type,label,occ,uiso,mxmymz = standard array outputs of Atoms
- get(self)
- Returns the structure arrays
uvw, type, label, occupancy, uiso, mxmymz = Atoms.get()
- info(self, idx=None, type=None)
- Prints properties of all atoms
:param idx: None or array of atoms to display
:param type: None or str type of atom to dispaly
:return: str
- ismagnetic(self)
- Returns True if any ions have magnetic moments assigned
- mass_fraction(self)
- Return the mass fraction per element
:return: float
- mxmymz(self)
- Returns a [nx3] array of magnetic vectors
:return: np.array([nx3])
- remove_duplicates(self, min_distance=0.01, all_types=False)
- Remove atoms of the same type that are too close to each other
:param min_distance: remove atoms within this distance, in fractional units
:param all_types: if True, also remove atoms of different types
:return: None
- removeatom(self, idx)
- Removes atom number idx from the list
:param idx: int, atom index
:return: None
- total_moment(self)
- Return the total moment along a, b, c directions
- update_cif(self, cifvals)
- Update cif dict with stored values
:param cifvals: cif dict from readcif
:return: cifvals
- uvw(self)
- Returns a [nx3] array of current positions
:return: np.array([nx3])
- weight(self)
- Calculate the molecular weight in g/mol of all the atoms
:return: float
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class Cell(builtins.object) |
|
Cell(a=1.0, b=1.0, c=1.0, alpha=90.0, beta=90.0, gamma=90.0)
Contains lattice parameters and unit cell
Provides tools to convert between orthogonal and lattice bases in real and reciprocal space.
E.G.
UC = Cell() # instantiate the Cell object
UC.latt([2.85,2.85,10.8,90,90,120]) # Define the lattice parameters from a list
UC.tth([0,0,12],energy_kev=8.0) # Calculate the two-theta of a reflection
UC.lp() # Returns the current lattice parameters
UC.orientation # class to chanage cell orientation in space |
|
Methods defined here:
- Bmatrix(self)
- Calculate the Busing and Levy B matrix from a real space UV
"choose the x-axis parallel to a*, the y-axis in the plane of a* and b*, and the z-axis perpendicular to
that plane"
W. R. Busing & H. A. Levy, Acta Cryst. (1967). 22, 457
- Qmag(self, HKL)
- Returns the magnitude of wave-vector transfer of [h,k,l], in A-1
:param HKL: list of hkl reflections
:return: list of Q values
- UV(self)
- Returns the unit cell as a [3x3] array, [A,B,C]
The vector A is directed along the x-axis
- UVstar(self)
- Returns the reciprocal unit cell as a [3x3] array, [A*,B*,C*]
:return: [a*;b*;c*]
- __init__(self, a=1.0, b=1.0, c=1.0, alpha=90.0, beta=90.0, gamma=90.0)
- Initialize self. See help(type(self)) for accurate signature.
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- all_hkl(self, energy_kev=8.048, max_angle=180.0, wavelength_a=None, maxq=None)
- Returns an array of all (h,k,l) reflections at this energy
:param energy_kev: energy in keV
:param max_angle: max two-theta angle
:param wavelength_a: wavelength in A
:param maxq: maximum wavevetor transfere in A-1 (suplants all above)
:return: array hkl[:,3]
- angle(self, hkl1, hkl2)
- Return the angle between two reflections
:param hkl1: [h,k,l] reflection 1
:param hkl2: [h,k,l] reflection 2
:return: angle in degrees
- calculateQ(self, HKL)
- Convert coordinates [h,k,l], in the basis of the reciprocal lattice, to
coordinates [x,y,z], in an orthogonal basis, in units of A-1
Q(x,y,z) = hA* + kB* + lC*
E.G.
Q = Cell.calculateQ([1,0,0]) # for a hexagonal system, a = 2.85
> Q = array([[2.2046264, 1.2728417, 0.0000000]])
- calculateR(self, UVW)
- Convert coordinates [u,v,w], in the basis of the unit cell, to
coordinates [x,y,z], in an orthogonal basis, in units of A
R(x,y,z) = uA + vB + wC
E.G.
R = Cell.calculateR([0.1,0,0]) # for a hexagonal system, a = 2.85
> R = array([[0.285, 0, 0]])
- diff6circle(self, delta=0, gamma=0, energy_kev=None, wavelength=1.0)
- Calcualte wavevector in diffractometer axis using detector angles
:param delta: float angle in degrees in vertical direction (about diff-z)
:param gamma: float angle in degrees in horizontal direction (about diff-x)
:param energy_kev: float energy in KeV
:param wavelength: float wavelength in A
:return: q[1*3], ki[1*3], kf[1*3]
- diff6circle2hkl(self, phi=0, chi=0, eta=0, mu=0, delta=0, gamma=0, energy_kev=None, wavelength=1.0)
- Return [h,k,l] position of diffractometer axes at given energy
:param phi: float sample angle in degrees
:param chi: float sample angle in degrees
:param eta: float sample angle in degrees
:param mu: float sample angle in degrees
:param delta: float detector angle in degrees
:param gamma: float detector angle in degrees
:param energy_kev: float energy in KeV
:param wavelength: float wavelength in A
:return: [h,k,l]
- diff6circle_match(self, phi=0, chi=0, eta=0, mu=0, delta=0, gamma=0, energy_kev=None, wavelength=1.0, fwhm=0.5)
- Return the closest hkl and intensity factor
:param phi: float sample angle in degrees
:param chi: float sample angle in degrees
:param eta: float sample angle in degrees
:param mu: float sample angle in degrees
:param delta: float detector angle in degrees
:param gamma: float detector angle in degrees
:param energy_kev: float energy in KeV
:param wavelength: float wavelength in A
:param fwhm: float peak width in A-1
:return: [h,k,l], If
- dspace(self, hkl)
- Calculate the d-spacing in A
:param hkl: array : list of reflections
:return: d-spacing
- find_close_reflections(self, hkl, energy_kev, max_twotheta=2, max_angle=10)
- Find reflections near to given HKL for a given two-theta or reflection angle
:param hkl: [h,k,l] indices of reflection to start from
:param energy_kev: energy in keV
:param max_twotheta: matches reflections within two-theta of hkl
:param max_angle: matches reflections within max_angle of hkl
:return: list of matching [[h,k,l]] reflections
- fromcif(self, cifvals)
- Import lattice parameters from a cif dictionary
Required CIF keys:
_cell_length_a
_cell_length_b
_cell_length_c
_cell_angle_alpha
_cell_angle_beta
_cell_angle_gamma
:param cifvals: dict from readcif
:return: None
- generate_lattice(self, U, V, W)
- Returns the lattice parameters of a larger lattice
- indexQ(self, Q)
- Convert coordinates [x,y,z], in an orthogonal basis, to
coordinates [h,k,l], in the basis of the reciprocal lattice
H(h,k,l) = Q(x,y,z) / [A*,B*,C*]
E.G.
HKL = indexQ([2.2046264, 1.2728417, 0.0000000]) # for a hexagonal system, a = 2.85
> HKL = [1,0,0]
- indexR(self, R)
- Convert coordinates [x,y,z], in an orthogonal basis, to
coordinates [u,v,w], in the basis of the unit cell
U(u,v,w) = R(x,y,z) / [A,B,C]
E.G.
UVW = indexR([0.285, 0, 0]) # for a hexagonal system, a = 2.85
> UVW = [0.1,0,0]
- info(self)
- Prints the lattice parameters and cell volume"
:return: str
- labwavevector(self, hkl)
- Calculate the lab wavevector using the unit-vector, oritenation matrix and rotation matrix
Returns vectors in the lab coordinate system, by default defined like Diamond Light Source:
x-axis : away from synchrotron ring, towards wall
y-axis : towards ceiling
z-axis : along beam direction
:param hkl: [3xn] array of (h, k, l) reciprocal lattice vectors
:return: [3xn] array of Q vectors in the lab coordinate system
- latt(self, lattice_parameters=(), *args, **kwargs)
- Generate lattice parameters with list
latt(1) -> a=b=c=1,alpha=beta=gamma=90
latt([1,2,3]) -> a=1,b=2,c=3,alpha=beta=gamma=90
latt([1,2,3,120]) -> a=1,b=2,c=3,alpha=beta=90,gamma=120
latt([1,2,3,10,20,30]) -> a=1,b=2,c=3,alpha=10,beta=20,gamma=30
latt(1,2,3,10,20,30) -> a=1,b=2,c=3,alpha=10,beta=20,gamma=30
latt(a=1,b=2,c=3,alpha=10,beta=20,gamma=30]) -> a=1,b=2,c=3,alpha=10,beta=20,gamma=30
- lp(self)
- Returns the lattice parameters
:return: a,b,c,alpha,beta,gamma
- max_hkl(self, energy_kev=8.048, max_angle=180.0, wavelength_a=None, maxq=None)
- Returns the maximum index of h, k and l for a given energy
:param energy_kev: energy in keV
:param max_angle: maximum two-theta at this energy
:param wavelength_a: wavelength in A
:param maxq: maximum wavevetor transfere in A-1 (suplants all above)
:return: maxh, maxk, maxl
- moment(self, mxmymz)
- Calcualte moment from value stored in cif
- powder_average(self, hkl)
- Returns the powder average correction for the given hkl
:param hkl: array : list of reflections
:return: correction
- reciprocal_space_plane(self, x_axis=[1, 0, 0], y_axis=[0, 1, 0], centre=[0, 0, 0], q_max=4.0, cut_width=0.05)
- Returns positions within a reciprocal space plane
x_axis = direction along x, in units of the reciprocal lattice (hkl)
y_axis = direction along y, in units of the reciprocal lattice (hkl)
centre = centre of the plot, in units of the reciprocal lattice (hkl)
q_max = maximum distance to plot to - in A-1
Returns: X,Y,HKL
Qx = [nx1] array of x positions for each HKL in the plane
Qy = [nx1] array of y positions for each HKL in the plane
HKL= [nx1] array of each HKL in the plane
- reflection_hkl(self, energy_kev=8.048, max_angle=180.0, specular=(0, 0, 1), theta_offset=0, min_theta=0, max_theta=180.0)
- Returns an array of all (h,k,l) reflections in reflection geometry
:param energy_kev: energy in keV
:param max_angle: max two-theta angle
:param specular: (h,k,l) of direction normal to surface and the incident beam
:param theta_offset: float : angle (deg) of surface relative to specular normal
:param min_theta: float : cut hkl reflections with reflection-theta lower than min_theta
:param max_theta: flaot : cut hkl reflections with reflection-theta greater than max_theta
:return: array of hkl
- sort_hkl(self, hkl, ascend=True)
- Returns array of (h,k,l) sorted by two-theta
:param hkl: array : list of [h,k,l] values
:param ascend: True*/False : if False, lowest two-theta
:return: HKL[sorted,:]
- theta_reflection(self, HKL, energy_kev=8.048, specular=[0, 0, 1], theta_offset=0)
- Calculate the sample angle for diffraction in reflection geometry given a particular specular direction
- theta_transmission(self, HKL, energy_kev=8.048, parallel=[0, 0, 1], theta_offset=0)
- Calculate the sample angle for diffraction in transmission geometry given
a particular direction parallel to the beam
- transmission_hkl(self, energy_kev=8.048, max_angle=180.0, parallel=(0, 0, 1), theta_offset=0, min_theta=0, max_theta=180.0)
- Returns an array of all (h,k,l) reflections in reflection geometry
:param energy_kev: energy in keV
:param max_angle: max two-theta angle
:param parallel: (h,k,l) of direction normal to surface, parallel to the incident beam
:param theta_offset: float : angle (deg) of surface relative to specular normal
:param min_theta: float : cut hkl reflections with reflection-theta lower than min_theta
:param max_theta: flaot : cut hkl reflections with reflection-theta greater than max_theta
:return: array of hkl
- tth(self, HKL, energy_kev=8.048, wavelength_a=None)
- Returns the two-theta angle, in deg, of [h,k,l] at specified energy in keV
:param HKL: list of hkl reflections
:param energy_kev: photon energy in keV
:param wavelength_a: wavelength in Angstroms
:return: two-theta angles
- ubmatrix(self)
- Return UB matrix from Busing & Levy in the diffractometer frame
- update_cif(self, cifvals)
- Update cif dict with current values
:param cifvals: dict from readcif
:return: cifvals
- volume(self)
- Returns the volume of the unit cell, in A^3
:return: volume
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class Crystal(builtins.object) |
|
Crystal(filename=None)
Reads the structure information from a cif file and generates the full structure.
Allows the adjustment of the structure through the lattice parameters, symmetry
or atomic displacement.
Can calculate reflection intensities and two-theta values.
E.G.
xtl = Crystal('Diamond.cif')
xtl.Cell.lp() >> give the lattice parameters
xtl.Atoms.uvw() >> give the symmetric atomic positions
xtl.Symmetry.symmetry_operations >> give symmetry operations
xtl.Structure.uvw() >> give the full, unsymmetrised structure
xtl.Scatter.hkl([1,0,0],8.00) >> prints the intensity and two-theta of this reflection at this energy
xtl.Scatter.print_all_reflections(8.00) >> print all allowed reflections, with intensities, at this energy
xtl.write_cif('Diamond2.cif') >> write updated structure to file
To create your own crystal (BCC example):
xtl = Crystal()
xtl.new_latt([2.866])
xtl.new_atoms(u=[0,0.5],
v=[0,0.5],
w=[0,0.5],
type=['Fe','Fe'])
xtl.hkl([[1,0,0],[1,1,0],[1,1,1],[2,0,0]])
Also, see:
help(xtl.Cell)
help(xtl.Atoms)
help(xtl.Symmetry)
help(xtl.Scatter) |
|
Methods defined here:
- __add__(self, other)
- __init__(self, filename=None)
- Initialize self. See help(type(self)) for accurate signature.
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- add_parent(self, parent, P)
- Add parent structure, returning Crystal as superstructure
parent = Crystal(cif_parent)
xtl = Crystal(cif_superstructure)
su = xtl.add_parent(parent, [[1,1,0],[0,2,0],[0,0,1]])
- generate_lattice(self, U=1, V=1, W=0)
- Generate a repeated lattice of the crystal structure
latt = xtl.generate_lattice(2,0,0)
:param U: Repeat of the cell along the a axis
:param V: Repeat of the cell along the b axis
:param W: Repeat of the cell along the c axis
:return: Crystal object
- generate_structure(self)
- Combines the atomic positions with symmetry operations, returning the full structure as an Atoms class
:return: None
- generate_superstructure(self, P)
- Generate a superstructure of the current cell
a' = n1a + n2b + n3c
b' = m1a + m2b + m3c
c' = o1a + o2b + o3c
OR
[a',b',c'] = P[a,b,c]
Where
P = [[n1,n2,n3],
[m1,m2,m3],
[o1,o2,o3]]
Returns a superstructure Crystal class:
su = xtl.generate_superstructrue([[2,0,0],[0,2,0],[0,0,1]])
Superstructure Crystal classes have additional attributes:
su.P = P as given
su.Parent = the parent Crystal Class
Use >>hasattr(su,'Parent') to check if the current object is a
superstructure Crystal class
- info(self)
- Returns information about the crystal structure
:return: str
- invert_structure(self)
- Convert handedness of structure, transform from left-handed to right handed, or visa-versa
Equivlent to xtl.transform([[-1,0,0], [0,-1,0], [0,0,-1]])
:return: Superstructure Crystal class
- new_atoms(self, u=[0], v=[0], w=[0], type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Replace current atomic positions with new ones and regenerate structure
:param u: array : atomic positions u
:param v: array : atomic positions v
:param w: array : atomic positions w
:param type: array : atomic types
:param label: array : atomic labels
:param occupancy: array : atomic occupation
:param uiso: array : atomic isotropic thermal parameters
:param mxmymz: array : atomic magnetic vectors [mu,mv,mw]
:return: None
- new_cell(self, lattice_parameters=(), *args, **kwargs)
- Replace the lattice parameters
:param lattice_parameters: [a,b,c,alpha,beta,gamma]
:return: None
- start_gui(self)
- Start Crystal GUI
:return: None
- transform(self, P)
- Transform the current cell
a' = n1a + n2b + n3c
b' = m1a + m2b + m3c
c' = o1a + o2b + o3c
OR
[a',b',c'] = P[a,b,c]
Where
P = [[n1,n2,n3],
[m1,m2,m3],
[o1,o2,o3]]
Returns a superstructure Crystal class:
su = xtl.transform([[0,1,0],[1,0,0],[0,0,1]])
Superstructure Crystal classes have additional attributes:
su.P = P as given
su.Parent = the parent Crystal Class
Use >>hasattr(su,'Parent') to check if the current object is a
superstructure Crystal class
- update_cif(self, cifvals=None)
- Update self.cif dict with current values
:param cifvals: cif dict from readcif (None to use self.cif)
:return: cifvals
- write_cif(self, filename=None, comments=None)
- Write crystal structure to CIF (Crystallographic Information File)
Only basic information is saved to the file, but enough to open in VESTA etc.
If magnetic ions are defined, a magnetic cif (*.mcif) will be produce
:param filename: name to write too, if None, use writes to self.name (.cif/.mcif)
:param comments: str comments to add to file header
:return: None
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
Data and other attributes defined here:
- cif = {}
- filename = ''
- name = 'Crystal'
- scale = 1.0
|
class Superstructure(Crystal) |
|
Superstructure(Parent, P)
Generate a superstructure of the current cell
a' = n1a + n2b + n3c
b' = m1a + m2b + m3c
c' = o1a + o2b + o3c
OR
[a',b',c'] = P[a,b,c]
Where
P = [[n1,n2,n3],
[m1,m2,m3],
[o1,o2,o3]]
Returns a superstructure Crystal class:
xtl = Crystal()
su = Superstructure(xtl,[[2,0,0],[0,2,0],[0,0,1]])
Superstructure Crystal classes have additional attributes compared with Crystal classes:
su.P = P as given
su.Parent = the parent Crystal Class
And additional functions:
su.calculateQ_parent >> indexes (h,k,l) coordinates in the frame same cartesian frame as the Parent structure
su.superhkl2parent >> indexes (h,k,l) coordinates in the frame of the Parent structure
su.parenthkl2super >> indexes parent (h,k,l) coordinates in the frame of superstructure
Use >>hasattr(su,'Parent') to check is the current object is a
superstructure Crystal class |
|
- Method resolution order:
- Superstructure
- Crystal
- builtins.object
Methods defined here:
- __init__(self, Parent, P)
- Initialise
- calculateQ_parent(self, super_hkl)
- Indexes (h,k,l) coordinates in the frame same cartesian frame as the Parent structure
Q = h'*a'* + k'*b'* + l'*c'*
Where a'*, b'*, c'* are defined relative to the parent lattice, a*,b*,c*
[qx,qy,qz] = calculateQ_parent([h',k',l'])
- generate_super_positions(self)
- Generate the supercell and superstructure based on P and parent structure
:return: None, set new atom positions
- parentUV(self)
- Returns the parent unit vectors defined relative to the supercell
- parentUVstar(self)
- Returns the parent reciprocal cell unit vectors defined relative to the supercell
- parenthkl2super(self, parent_HKL)
- Indexes (h,k,l) coordinates in the frame of the Parent structure
Q = h*a* + k*b* + l*c* = h'*a'* + k'*b'* + l'*c'*
[h',k',l'] = Q/[a'*,b'*,c'*]
[h',k',l'] = parenthkl2super([h,k,l])
- set_scale(self)
- Set scale parameter automatically
Based on ratio of parent - to superstructure volume
- superUV(self)
- Returns the supercell unit vectors defined relative to the Parent cell
- superUVstar(self)
- Returns the reciprocal supercell unit vectors defined relative to the Parent cell
- superhkl2parent(self, super_HKL)
- Indexes (h,k,l) coordinates in the frame of the Parent structure
Q = h*a* + k*b* + l*c* = h'*a'* + k'*b'* + l'*c'*
[h',k',l'] = Q/[a'*,b'*,c'*]
[h,k,l] = superhkl2parent([h',k',l'])
Methods inherited from Crystal:
- __add__(self, other)
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- add_parent(self, parent, P)
- Add parent structure, returning Crystal as superstructure
parent = Crystal(cif_parent)
xtl = Crystal(cif_superstructure)
su = xtl.add_parent(parent, [[1,1,0],[0,2,0],[0,0,1]])
- generate_lattice(self, U=1, V=1, W=0)
- Generate a repeated lattice of the crystal structure
latt = xtl.generate_lattice(2,0,0)
:param U: Repeat of the cell along the a axis
:param V: Repeat of the cell along the b axis
:param W: Repeat of the cell along the c axis
:return: Crystal object
- generate_structure(self)
- Combines the atomic positions with symmetry operations, returning the full structure as an Atoms class
:return: None
- generate_superstructure(self, P)
- Generate a superstructure of the current cell
a' = n1a + n2b + n3c
b' = m1a + m2b + m3c
c' = o1a + o2b + o3c
OR
[a',b',c'] = P[a,b,c]
Where
P = [[n1,n2,n3],
[m1,m2,m3],
[o1,o2,o3]]
Returns a superstructure Crystal class:
su = xtl.generate_superstructrue([[2,0,0],[0,2,0],[0,0,1]])
Superstructure Crystal classes have additional attributes:
su.P = P as given
su.Parent = the parent Crystal Class
Use >>hasattr(su,'Parent') to check if the current object is a
superstructure Crystal class
- info(self)
- Returns information about the crystal structure
:return: str
- invert_structure(self)
- Convert handedness of structure, transform from left-handed to right handed, or visa-versa
Equivlent to xtl.transform([[-1,0,0], [0,-1,0], [0,0,-1]])
:return: Superstructure Crystal class
- new_atoms(self, u=[0], v=[0], w=[0], type=None, label=None, occupancy=None, uiso=None, mxmymz=None)
- Replace current atomic positions with new ones and regenerate structure
:param u: array : atomic positions u
:param v: array : atomic positions v
:param w: array : atomic positions w
:param type: array : atomic types
:param label: array : atomic labels
:param occupancy: array : atomic occupation
:param uiso: array : atomic isotropic thermal parameters
:param mxmymz: array : atomic magnetic vectors [mu,mv,mw]
:return: None
- new_cell(self, lattice_parameters=(), *args, **kwargs)
- Replace the lattice parameters
:param lattice_parameters: [a,b,c,alpha,beta,gamma]
:return: None
- start_gui(self)
- Start Crystal GUI
:return: None
- transform(self, P)
- Transform the current cell
a' = n1a + n2b + n3c
b' = m1a + m2b + m3c
c' = o1a + o2b + o3c
OR
[a',b',c'] = P[a,b,c]
Where
P = [[n1,n2,n3],
[m1,m2,m3],
[o1,o2,o3]]
Returns a superstructure Crystal class:
su = xtl.transform([[0,1,0],[1,0,0],[0,0,1]])
Superstructure Crystal classes have additional attributes:
su.P = P as given
su.Parent = the parent Crystal Class
Use >>hasattr(su,'Parent') to check if the current object is a
superstructure Crystal class
- update_cif(self, cifvals=None)
- Update self.cif dict with current values
:param cifvals: cif dict from readcif (None to use self.cif)
:return: cifvals
- write_cif(self, filename=None, comments=None)
- Write crystal structure to CIF (Crystallographic Information File)
Only basic information is saved to the file, but enough to open in VESTA etc.
If magnetic ions are defined, a magnetic cif (*.mcif) will be produce
:param filename: name to write too, if None, use writes to self.name (.cif/.mcif)
:param comments: str comments to add to file header
:return: None
Data descriptors inherited from Crystal:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
Data and other attributes inherited from Crystal:
- cif = {}
- filename = ''
- name = 'Crystal'
- scale = 1.0
|
class Symmetry(builtins.object) |
|
Symmetry(symmetry_operations=None, symmetry_operations_magnetic=None)
Contains symmetry information about the crystal, including the symmetry operations. |
|
Methods defined here:
- __init__(self, symmetry_operations=None, symmetry_operations_magnetic=None)
- Initialises the symmetry group
- __repr__(self)
- Return repr(self).
- __str__(self)
- Return str(self).
- addcen(self, operations, mag_operations=None)
- Apply centring operations to current symmetry operations
- addsym(self, operations, mag_operations=None)
- Add symmetry operations
Symmetry.addsym('x,y,z+1/2') >> adds single symmetry operation, magnetic operation is infered
Symmetry.addsym(['x,y,z+1/2','z,x,y']) >> adds multiple symmetry operation, magnetic operations are infered
Symmetry.addsym('x,y,z+1/2','x,y,-z') >> adds single symmetry operation + magnetic operation
- average_symmetric_intensity(self, hkl_list, intensity_list, tolerance=0.01)
- Return a list of reflections with symmetric reflections removed, matching reflections will be averaged
:param hkl_list: list of [h,k,l] reflections
:param intensity_list: list of intensities
:param tolerance: tolerance for matching reflections
:return: array of [h,k,l]
- axial_vector(self, uvw, remove_identical=True)
- Perform symmetry operations on an axial vector uvw
:param uvw: 3 element array/ list in cell coordinates
:param remove_identical: True/ False, if True, identical operations are removed
:return: [S*3] array of transformed coordinates
- changesym(self, idx, operation)
- Change a symmetry operation
:param idx: symmetry index
:param operation: str e.g. 'x,-y,z'
- fromcif(self, cifvals)
- Import symmetry information from a cif dictionary
Required cif keys:
None
Optional cif keys:
_symmetry_equiv_pos_as_xyz
_space_group_symop_operation_xyz
_space_group_symop_magn_operation_xyz
_space_group_symop_magn_operation_mxmymz
_space_group_symop_magn_centering_xyz
_space_group_symop_magn_centering_mxmymz
_symmetry_space_group_name_H-M
_space_group_name_H-M_alt
_space_group_magn_name_BNS
_symmetry_Int_Tables_number
_space_group_IT_number
_space_group_magn_number_BNS
:param cifvals: dict of values from cif
:return:
- generate_matrices(self)
- Generates the symmetry matrices from string symmetry operations
- info(self)
- Prints the symmetry information
:return: str
- invert_magsym(self, idx)
- Invert the time symmetry of a magnetic symmetry
:param idx: symmetry index, 0:Nsym
:return:
- is_symmetric_reflection(self, hkl1, hkl2, tolerance=0.01)
- Check if reflection 1 is a symmetric equivalent of reflection 2
:param hkl1: [h,k,l] reflection 1
:param hkl2: [h,k,l] reflection 2
:param tolerance: tolerance for matching reflections
:return: True/ False
- load_magnetic_spacegroup(self, msg_number=None, sg_dict=None)
- Load symmetry operations from a magnetic spacegroup from Bilbao crystallographic server
Replaces the current symmetry operators and the magnetic symmetry operators.
See functions_crystallography.spacegroup_magnetic for more details
:param msg_number: magnetic space group number e.g. 61.433
:param sg_dict: alternative inuput: spacegroup dict from fc.spacegroup_magnetic
:return: None
- load_spacegroup(self, sg_number=None, sg_dict=None)
- Load symmetry operations from a spacegroup from the International Tables of Crystallogrphy
See functions_crystallography.spacegroup for more details
:param sg_number: space group number (1-230)
:param sg_dict: alternative input: spacegroup dict from fc.spacegroup
:return: None
- parity_time_info(self)
- Returns string of parity and time operations for symmetry operations
:return: str
- print_magnetic_spacegroups(self)
- Return str of available magnetic spacegroups for this spacegroup
:return: str
- print_subgroups(self)
- Return str of subgroups of this spacegroup
:return: str
- print_symmetric_coordinate_operations(self, UVW, remove_identical=True)
- Returns array of symmetric operations for given position
Uses fc.gen_sym_pos
Returns list of identical symmetry operations, with identical positions removed
All positions returned if remove_identical=False
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.print_symmetric_coordinate_operations([0.1,0,0])
> n u v w Symmetry Magnetic Symmetry
0 0.1000 0.0000 0.0000 x,y,z x,y,z
1 0.9000 0.0000 0.0000 -x,-y,-z -x,-y,-z
2 0.0000 0.1000 0.0000 y,x,z y,x,z
- print_symmetric_vectors(self, HKL)
- Print symmetric vectors
:param HKL: [h,k,l] reflection
:return: str
- print_wyckoff_sites(self)
- Return info str about Wycoff positions for this spacegroup
- reflection_multiplyer(self, HKL)
- Returns the number of symmetric reflections for each hkl
:param HKL: [nx3] array of [h,k,l]
:return: [n] array of multiplyers
- remove_symmetric_reflections(self, hkl_list, tolerance=0.01)
- Return a list of reflections with symmetric reflections removed
:param hkl_list: list of [h,k,l] reflections
:param tolerance: tolerance for matching reflections
:return: array of [h,k,l]
- spacegroup_name(self)
- Return the spacegroup name and number as str
- symmetric_coordinate_operations(self, UVW, MXMYMZ=None)
- Returns array of symmetric operations for given position
Uses fc.gen_sym_pos
Returns list of identical symmetry operations, with identical positions removed.
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.symmetric_coordinate_operations([0.1,0,0])
> array(['x,y,z', '-x,-y,-z', 'y,x,z'])
- symmetric_coordinates(self, UVW, MXMYMZ=None, remove_identical=True)
- Returns array of symmetric coordinates
Uses fc.gen_sym_pos
Returns coordinates wrapped within the unit cell, with identical positions removed.
All positions returned if remove_identical=False
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.symmetric_coordinates([0.1,0,0])
> array([[0.1, 0.0, 0.0],
[0.9, 0.0, 0.0],
[0.0, 0.1, 0.0]])
- symmetric_intensity(self, HKL, I, dI=None)
- Returns symmetric reflections with summed intensities of repeated reflections
Assumes HKL reflections are unique, repeated reflections will be incorrectly added together.
Uses fc.gen_sym_mat
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.generate_matrixes()
HKL, I = Symmetry.symmetric_intensity([1,1,0],10)
> HKL = array([[1, 1, 0],
[-1,-1, 0]])
> I = array([20,10])
OR
HKL, I, dI = Symmetry.symmetric_intensity([1,1,0],10,1)
> HKL = array([[1, 1, 0],
[-1,-1, 0]])
> I = array([20,10])
> dI = array([2,1])
- symmetric_magnetic_vectors(self, MXMYMZ)
- NOT COMPLETE
- symmetric_reflections(self, HKL)
- Returns array of symmetric reflection indices
Uses fc.gen_sym_mat
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.generate_matrixes()
Symmetry.symmetric_reflections([1,0,0])
> array([[1, 0, 0],
[-1, 0, 0],
[0, -1, 0]])
- symmetric_reflections_count(self, HKL)
- Returns array of symmetric reflection indices,
identical reflections are removed, plus the counted sum of each reflection
Uses fc.gen_sym_mat
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.generate_matrixes()
HKL, count = Symmetry.symmetric_reflections([1,1,0])
> HKL = array([[1, 1, 0],
[-1,-1, 0]])
> count = array([2,1])
- symmetric_reflections_unique(self, HKL)
- Returns array of symmetric reflection indices, with identical reflections removed
Uses fc.gen_sym_mat
E.G.
Symmetry.symmetry_operations = ['x,y,z','-x,-y,-z','y,x,z']
Symmetry.generate_matrixes()
Symmetry.symmetric_reflections([1,1,0])
> array([[1, 1, 0],
[-1,-1, 0]])
- update_cif(self, cifvals)
- Update cifvals dict with current symmetry operations
:param cifvals: cif dict from functions_crystallography.readcif
:return: cifvals
- wyckoff_labels(self, UVW)
- Return Wyckoff site for position
:param UVW: (u,v,w) or None to use xtl.Atoms.uvw()
:return: list of Wyckoff site letters
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
Data and other attributes defined here:
- spacegroup = 'P1'
- spacegroup_dict = {'general positions': ['x,y,z'], 'magnetic space groups': ['1.1', '1.2', '1.3', '1.3.3'], 'positions centring': ['x,y,z'], 'positions coordinates': [['x,y,z']], 'positions multiplicity': [1], 'positions symmetry': ['1'], 'positions wyckoff letter': ['a'], 'space group name': 'P1', 'space group name html': '<i>P</i>1 (No. 1)', 'space group number': 1, ...}
- spacegroup_number = 1
- symmetry_matrices = array([[1., 0., 0., 0.],
[0., 1., 0., 0.],
[0., 0., 1., 0.],
[0., 0., 0., 1.]])
- symmetry_operations = ['x,y,z']
- symmetry_operations_magnetic = ['x,y,z']
- symmetry_operations_time = [1]
| |