diffpy.structure package
Crystal structure container and parsers for structure formats.
- Classes related to the structure of materials:
Atom
Lattice
Structure
PDFFitStructure
- Other classes:
SpaceGroup
SymOp
ExpandAsymmetricUnit
GeneratorSite
SymmetryConstraints
- Exceptions:
StructureFormatError
LatticeError
SymmetryError
- diffpy.structure.loadStructure(filename, fmt='auto', **kw)[source]
Load new structure object from the specified file.
- Parameters:
filename (str) – Path to the file to be loaded.
fmt (str, Optional) – Format of the structure file such as ‘cif’ or ‘xyz’. Must be one of the formats listed by the parsers.inputFormats function. When ‘auto’, all supported formats are tried in a sequence.
kw (Optional) – Extra keyword arguments that are passed to parsers.getParser function. These configure the dedicated Parser object that is used to read content in filename.
- Returns:
stru – The new Structure object loaded from the specified file. Return a more specific PDFFitStructure type for ‘pdffit’ and ‘discus’ formats.
- Return type:
Structure, PDFFitStructure
Subpackages
Submodules
diffpy.structure.spacegroups module
Space group classes and definitions from mmLib and sgtbx.
- diffpy.structure.spacegroups.SpaceGroupList
List of all spacegroup definitions.
- Type:
list
- diffpy.structure.spacegroups.FindSpaceGroup(symops, shuffle=False)[source]
Lookup SpaceGroup from a given list of symmetry operations.
- Parameters:
symops (list) – The list of SymOp objects for which to find SpaceGroup.
shuffle (bool, Optional) – Flag for allowing different order of symops in the returned SpaceGroup. The default is
False
.
- Returns:
The SpaceGroup object with equivalent list of symmetry operations. Return predefined SpaceGroup instance when symmetry operations have the same order or when the shuffle flag is set.
- Return type:
- Raises:
ValueError – When symops do not match any known SpaceGroup.
- diffpy.structure.spacegroups.GetSpaceGroup(sgid)[source]
Returns the SpaceGroup instance for the given identifier.
- Parameters:
sgid (str, int) – space group symbol, either short_name or pdb_name, whatever it means in mmlib. Can be also an integer.
- Returns:
The SpaceGroup object for the given identifier.
- Return type:
- Raises:
ValueError – When the identifier is not found.
diffpy.structure._legacy_importer module
Support import of old camel-case module names with DeprecationWarning.
The imported camel-case modules are aliases for the current module instances. Their __name__ attributes are thus all in lower-case.
Note
this module must be only imported from diffpy.Structure.
Warning
This module is deprecated and will be removed in the future.
- class diffpy.structure._legacy_importer.FindRenamedStructureModule[source]
Bases:
MetaPathFinder
- prefix = 'diffpy.Structure.'
- class diffpy.structure._legacy_importer.MapRenamedStructureModule[source]
Bases:
Loader
Loader for old camel-case module names. Import the current module and alias it under the old name.
diffpy.structure.structureerrors module
Exceptions used in Structure package.
- exception diffpy.structure.structureerrors.LatticeError[source]
Bases:
Exception
Exception for impossible lattice parameters.
diffpy.structure.spacegroupmod module
Symmetry operations as functions on vectors or arrays.
- class diffpy.structure.spacegroupmod.SpaceGroup(number=None, num_sym_equiv=None, num_primitive_sym_equiv=None, short_name=None, point_group_name=None, crystal_system=None, pdb_name=None, symop_list=None)[source]
Bases:
object
Definition and basic operations for a specific space group.
Provide standard names and all symmetry operations contained in one space group.
- Parameters:
number (int) – The space group number.
num_sym_equiv (int) – The number of symmetry equivalent sites for a general position.
num_primitive_sym_equiv (int) – The number of symmetry equivalent sites in a primitive unit cell.
short_name (str) – The short Hermann-Mauguin symbol of the space group.
point_group_name (str) – The point group of this space group.
crystal_system (str) – The crystal system of this space group.
pdb_name (str) – The full Hermann-Mauguin symbol of the space group.
symop_list (list of SymOp) – The symmetry operations contained in this space group.
- number
A unique space group number. This may be incremented by several thousands to facilitate unique values for multiple settings of the same space group. Use
number % 1000
to get the standard space group number from International Tables.- Type:
int
- num_sym_equiv
The number of symmetry equivalent sites for a general position.
- Type:
int
- num_primitive_sym_equiv
The number of symmetry equivalent sites in a primitive unit cell.
- Type:
int
- short_name
The short Hermann-Mauguin symbol of the space group.
- Type:
str
- point_group_name
The point group to which this space group belongs to.
- Type:
str
- crystal_system
The crystal system of this space group. The possible values are
"TRICLINIC", "MONOCLINIC", "ORTHORHOMBIC", "TETRAGONAL", "TRIGONAL" "HEXAGONAL", "CUBIC"
.- Type:
str
- pdb_name
The full Hermann-Mauguin symbol of the space group.
- Type:
str
- symop_list
A list of SymOp objects for all symmetry operations in this space group.
- Type:
list of SymOp
- check_group_name(name)[source]
Check if given name matches this space group.
- Parameters:
name (str or int) – The space group identifier, a string name or number.
- Returns:
True
if the specified name matches one of the recognized names of this space group or if it equals its number. ReturnFalse
otherwise.- Return type:
bool
- iter_equivalent_positions(vec)[source]
Generate symmetry equivalent positions for the specified position.
The initial position must be in fractional coordinates and so are the symmetry equivalent positions yielded by iteration. This generates num_sym_equiv positions regardless of initial coordinates being a special symmetry position or not.
- Parameters:
vec (numpy.ndarray) – The initial position in fractional coordinates.
- Yields:
numpy.ndarray – The symmetry equivalent positions in fractional coordinates. The positions may be duplicate or outside of the
0 <= x < 1
unit cell bounds.
- class diffpy.structure.spacegroupmod.SymOp(R, t)[source]
Bases:
object
The transformation of coordinates to a symmetry-related position.
The SymOp operation involves rotation and translation in cell coordinates.
- Parameters:
R (numpy.ndarray) – The 3x3 matrix of rotation for this symmetry operation.
t (numpy.ndarray) – The vector of translation in this symmetry operation.
- R
The 3x3 matrix of rotation pertaining to unit cell coordinates. This may be identity, simple rotation, improper rotation, mirror or inversion. The determinant of R is either +1 or -1.
- Type:
numpy.ndarray
- t
The translation of cell coordinates applied after rotation R.
- Type:
numpy.ndarray
diffpy.structure.utils module
Small shared functions.
- diffpy.structure.utils.atomBareSymbol(smbl)[source]
Remove atom type string stripped of isotope and ion charge symbols.
This function removes any blank, isotope numbers (0-9), leading hyphens (-), and ion charge symbols (1-9)(+-) from the given atom type string, returning only the bare element symbol.
- Parameters:
smbl (str) – Atom type string that may include isotope numbers, ion charges, or hyphens.
- Returns:
The bare element symbol.
- Return type:
str
Examples
>>> atomBareSymbol("Cl-") 'Cl' >>> atomBareSymbol("Ca2+") 'Ca' >>> atomBareSymbol("12-C") 'C'
diffpy.structure.lattice module
Class Lattice stores properties and provides simple operations in lattice coordinate system.
- diffpy.structure.lattice.cartesian
Constant instance of Lattice, default Cartesian system.
- Type:
- class diffpy.structure.lattice.Lattice(a=None, b=None, c=None, alpha=None, beta=None, gamma=None, baserot=None, base=None)[source]
Bases:
object
General coordinate system and associated operations.
- Parameters:
a (float or Lattice, Optional) – The cell length a. When present, other cell parameters must be also specified. When of the Lattice type, create a duplicate Lattice.
b (float) – The cell length b.
c (float) – The cell length c.
alpha (float) – The angle between the b and c axes in degrees.
beta (float) – The angle between the b and c axes in degrees.
gamma (float) – The angle between the a and b axes in degrees.
baserot (array_like, Optional) – The 3x3 rotation matrix of the base vectors with respect to their standard setting.
base (array_like, Optional) – The 3x3 array of row base vectors. This must be the only argument when present.
- metrics
The metrics tensor.
- Type:
numpy.ndarray
- base
The 3x3 matrix of row base vectors in Cartesian coordinates, which may be rotated, i.e.,
base = stdbase @ baserot
.- Type:
numpy.ndarray
- stdbase
The 3x3 matrix of row base vectors in standard orientation.
- Type:
numpy.ndarray
- baserot
The rotation matrix for the base.
- Type:
numpy.ndarray
- recbase
The inverse of the base matrix, where the columns give reciprocal vectors in Cartesian coordinates.
- Type:
numpy.ndarray
- normbase
The base vectors scaled by magnitudes of reciprocal cell lengths.
- Type:
numpy.ndarray
- recnormbase
The inverse of the normbase matrix.
- Type:
numpy.ndarray
- isotropicunit
The 3x3 tensor for a unit isotropic displacement parameters in this coordinate system. This is an identity matrix when this Lattice is orthonormal.
- Type:
numpy.ndarray
Note
The array attributes are read-only. They get updated by changing some lattice parameters or by calling the setLatPar() or setLatBase() methods.
Examples
Create a Cartesian coordinate system:
>>> Lattice()
Create coordinate system with the cell lengths a, b, c and cell angles alpha, beta, gamma in degrees:
>>> Lattice(a, b, c, alpha, beta, gamma)
Create a duplicate of an existing Lattice lat:
>>> Lattice(lat)
Create coordinate system with the base vectors given by rows of the abc matrix:
>>> Lattice(base=abc)
- property a
The unit cell length a.
- abcABG()[source]
Return the cell parameters in the standard setting. :returns: A tuple of
(a, b, c, alpha, beta, gamma)
. :rtype: tuple
- property alpha
The cell angle alpha in degrees.
- property alphar
The reciprocal cell angle alpha in degrees.
- angle(u, v)[source]
Calculate angle between 2 lattice vectors in degrees.
- Parameters:
u (array_like) – The first lattice vector.
v (array_like) – The second lattice vector.
- Returns:
The angle between lattice vectors u and v in degrees.
- Return type:
float
- property ar
The cell length a of the reciprocal lattice.
- property b
The unit cell length b.
- property beta
The cell angle beta in degrees.
- property betar
The reciprocal cell angle beta in degrees
- property br
The cell length b of the reciprocal lattice.
- property c
The unit cell length c.
- property ca
The cosine of the cell angle alpha.
- property car
The cosine of the reciprocal angle alpha.
- cartesian(u)[source]
Transform lattice vector to Cartesian coordinates.
- Parameters:
u (array_like) – Vector of lattice coordinates or an Nx3 array of lattice vectors.
- Returns:
rc – Cartesian coordinates of the u vector.
- Return type:
numpy.ndarray
- property cb
The cosine of the cell angle beta.
- property cbr
The cosine of the reciprocal angle beta.
- property cg
The cosine of the cell angle gamma.
- property cgr
The cosine of the reciprocal angle gamma.
- property cr
The cell length c of the reciprocal lattice.
- dist(u, v)[source]
Calculate distance between 2 points in lattice coordinates.
- Parameters:
u (array_like) – A vector or an Nx3 matrix of fractional coordinates.
v (numpy.ndarray) – A vector or an Nx3 matrix of fractional coordinates.
Note
u and v must be of the same shape when matrices.
- Returns:
The distance between lattice points u and v.
- Return type:
float or numpy.ndarray
- dot(u, v)[source]
Calculate dot product of 2 lattice vectors.
- Parameters:
u (array_like) – The first lattice vector or an Nx3 array.
v (array_like) – The second lattice vector or an array of the same shape as u.
- Returns:
The dot product of lattice vectors u, v.
- Return type:
float or numpy.ndarray
- fractional(rc)[source]
Transform Cartesian vector to fractional lattice coordinates.
- Parameters:
rc (array_like) – A vector of Cartesian coordinates or an Nx3 array of Cartesian vectors.
- Returns:
u – Fractional coordinates of the Cartesian vector rc.
- Return type:
numpy.ndarray
- property gamma
The cell angle gamma in degrees.
- property gammar
The reciprocal cell angle gamma in degrees
- isanisotropic(umx)[source]
True if displacement parameter matrix is anisotropic.
This checks if the specified matrix of anisotropic displacement parameters (ADP) differs from isotropic values for this lattice by more than a small round-off error.
- Parameters:
umx (array_like) – The 3x3 matrix of displacement parameters.
- Returns:
True when umx is anisotropic by more than a round-off error.
- Return type:
bool
- norm(xyz)[source]
Calculate norm of a lattice vector.
- Parameters:
xyz (array_like) – A vector or an Nx3 array of fractional coordinates.
- Returns:
The magnitude of the lattice vector xyz.
- Return type:
float or numpy.ndarray
- reciprocal()[source]
Return the reciprocal lattice of the current lattice. :returns: The reciprocal lattice of the current lattice. :rtype: Lattice
- rnorm(hkl)[source]
Calculate norm of a reciprocal vector.
- Parameters:
hkl (array_like) – A vector or an Nx3 array of reciprocal coordinates.
- Returns:
The magnitude of the reciprocal vector hkl.
- Return type:
float or numpy.ndarray
- property sa
The sine of the cell angle alpha.
- property sar
The sine of the reciprocal angle alpha.
- property sb
The sine of the cell angle beta.
- property sbr
The sine of the reciprocal angle beta.
- setLatBase(base)[source]
Set new base vectors for this lattice.
This updates the cell lengths and cell angles according to the new base. The stdbase, baserot, and metrics attributes are also updated.
- Parameters:
base (array_like) – The 3x3 matrix of row base vectors expressed in Cartesian coordinates.
- setLatPar(a=None, b=None, c=None, alpha=None, beta=None, gamma=None, baserot=None)[source]
Set one or more lattice parameters.
This updates all attributes that depend on the lattice parameters.
- Parameters:
a (float, Optional) – The new value of the cell length a.
b (float, Optional) – The new value of the cell length b.
c (float, Optional) – The new value of the cell length c.
alpha (float, Optional) – The new value of the cell angle alpha in degrees.
beta (float, Optional) – The new value of the cell angle beta in degrees.
gamma (float, Optional) – The new value of the cell angle gamma in degrees.
baserot (array_like, Optional) – The new 3x3 rotation matrix of the base vectors with respect to their standard setting in Cartesian coordinates.
Note
Parameters that are not specified will keep their initial values.
- property sg
The sine of the cell angle gamma.
- property sgr
The sine of the reciprocal angle gamma.
- property unitvolume
The unit cell volume when a = b = c = 1.
- property volume
The unit cell volume.
diffpy.structure.structure module
This module defines class Structure.
- class diffpy.structure.structure.Structure(atoms=None, lattice=None, title=None, filename=None, format=None)[source]
Bases:
list
Define group of atoms in a specified lattice. Structure –> group of atoms.
Structure class is inherited from Python list. It contains a list of Atom instances. Structure overloads setitem and setslice methods so that the lattice attribute of atoms get set to lattice.
- Parameters:
atoms (list of Atom or Structure, Optional) – List of Atom instances to be included in this Structure. When atoms argument is an existing Structure instance, the new structure is its copy.
lattice (Lattice, Optional) – Instance of Lattice defining coordinate systems, property.
title (str, Optional) – String description of the structure.
filename (str, Optional) – Name of a file to load the structure from.
format (str, Optional) – Structure format of the loaded filename. By default all structure formats are tried one by one. Ignored when filename has not been specified.
Note
Cannot use filename and atoms arguments together. Overrides atoms argument when filename is specified.
- title
String description of the structure.
- Type:
str
- pdffit
Dictionary of PDFFit-related metadata.
- Type:
None or dict
Examples
Structure(stru)
create a copy of Structure instance stru.>>> stru = Structure() >>> copystru = Structure(stru)
Structure is inherited from a list it can use list expansions.
>>> oxygen_atoms = [a for a in stru if a.element == "O" ] >>> oxygen_stru = Structure(oxygen_atoms, lattice=stru.lattice)
- property B11
Array of B11 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property B12
Array of B12 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property B13
Array of B13 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property B22
Array of B22 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property B23
Array of B23 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property B33
Array of B33 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property Bisoequiv
Array of Debye-Waller isotropic thermal displacement or equivalent values. Assignment updates the U attribute of all Atoms.
- property U
Array of anisotropic thermal displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U11
Array of U11 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U12
Array of U12 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U13
Array of U13 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U22
Array of U22 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U23
Array of U23 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property U33
Array of U33 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all Atoms.
- property Uisoequiv
Array of isotropic thermal displacement or equivalent values. Assignment updates the U attribute of all Atoms.
- addNewAtom(*args, **kwargs)[source]
Add new Atom instance to the end of this Structure.
- Parameters:
*args – See Atom class constructor.
**kwargs – See Atom class constructor.
- angle(aid0, aid1, aid2)[source]
The bond angle at the second of three Atoms in degrees.
- Parameters:
aid0 (int or str) – Zero based index of the first Atom or a string label.
aid1 (int or str) – Index or string label for the second atom, where the angle is formed.
aid2 (int or str) – Index or string label for the third atom.
- Returns:
The bond angle in degrees.
- Return type:
float
- Raises:
IndexError – If any of the arguments are invalid.
- property anisotropy
Boolean array for anisotropic thermal displacement flags. Assignment updates the anisotropy attribute of all Atoms.
- append(a, copy=True)[source]
Append Atom to a structure and update its lattice attribute.
- Parameters:
a (Atom) – Instance of Atom to be appended.
copy (bool, Optional) – Flag for appending a copy of a. When
False
, append a and update a.lattice.
- assignUniqueLabels()[source]
Set a unique label string for each Atom in this structure.
The label strings are formatted as “%(baresymbol)s%(index)i”, where baresymbol is the element right-stripped of “[0-9][+-]”.
- property composition
Dictionary of chemical symbols and their total occupancies.
- distance(aid0, aid1)[source]
Calculate distance between 2 Atoms, no periodic boundary conditions.
- Parameters:
aid0 (int or str) – Zero based index of the first Atom or a string label.
aid1 (int or str) – Zero based index or string label of the second atom.
- Returns:
Distance between the two Atoms in Angstroms.
- Return type:
float
- Raises:
IndexError – If any of the Atom indices or labels are invalid.
- property element
Character array of Atom types. Assignment updates the element attribute of the respective Atoms. Set the maximum length of the element string to 5 characters.
- extend(atoms, copy=None)[source]
Extend Structure with an iterable of atoms.
Update the lattice attribute of all added atoms.
- Parameters:
atoms (Iterable) – The Atom objects to be appended to this Structure.
copy (bool, Optional) – Flag for adding copies of Atom objects. Make copies when
True
, append atoms unchanged whenFalse
. The default behavior is to make copies when atoms are of Structure type or if new atoms introduce repeated objects.
- insert(idx, a, copy=True)[source]
Insert Atom a before position idx in this Structure.
- Parameters:
idx (int) – Position in Atom list.
a (Atom) – Instance of Atom to be inserted.
copy (bool, Optional) – Flag for inserting a copy of a. When
False
, append a and update a.lattice.
- property label
Character array of Atom names. Assignment updates the label attribute of all Atoms. Set the maximum length of the label string to 5 characters.
- property lattice
Coordinate system for this Structure.
- property occupancy
Array of Atom occupancies. Assignment updates the occupancy attribute of all Atoms.
- pdffit = None
default values for pdffit.
- Type:
None
- placeInLattice(new_lattice)[source]
place structure into new_lattice coordinate system.
Sets lattice to new_lattice and recalculate fractional coordinates of all Atoms so their absolute positions remain the same.
- read(filename, format='auto')[source]
Load structure from a file, any original data become lost.
- Parameters:
filename (str) – File to be loaded.
format (str, Optional) – All structure formats are defined in parsers submodule, when
format == 'auto'
all parsers are tried one by one.
- Returns:
Return instance of data Parser used to process input string. This can be inspected for information related to particular format.
- Return type:
Parser
- readStr(s, format='auto')[source]
Read structure from a string.
- Parameters:
s (str) – String with structure definition.
format (str, Optional) – All structure formats are defined in parsers submodule. When
format == 'auto'
, all parsers are tried one by one.
- Returns:
Return instance of data Parser used to process input string. This can be inspected for information related to particular format.
- Return type:
Parser
- title = ''
default values for title.
- Type:
str
- write(filename, format)[source]
Save structure to file in the specified format.
- Parameters:
filename (str) – File to save the structure to.
format (str) – Structure format to use for saving.
Note
Available structure formats can be obtained by:
from parsers import formats
- writeStr(format)[source]
return string representation of the structure in specified format.
Note
Available structure formats can be obtained by:
from parsers import formats
- property x
Array of all fractional coordinates x. Assignment updates xyz attribute of all Atoms.
- property xyz
Array of fractional coordinates of all Atoms. Assignment updates xyz attribute of all Atoms.
- property xyz_cartn
Array of absolute Cartesian coordinates of all Atoms. Assignment updates the xyz attribute of all Atoms.
- property y
Array of all fractional coordinates y. Assignment updates xyz attribute of all Atoms.
- property z
Array of all fractional coordinates z. Assignment updates xyz attribute of all Atoms.
diffpy.structure.mmlibspacegroups module
Space groups defined as a part of the pymmlib.
diffpy.structure.symmetryutilities module
Symmetry utility functions such as expansion of asymmetric unit, and generation of positional constraints.
- diffpy.structure.symmetryutilities.epsilon
Default tolerance for equality of 2 positions, also used for identification of special positions.
- Type:
float
- diffpy.structure.symmetryutilities.stdUsymbols
Standard symbols denoting elements of anisotropic thermal displacement tensor.
- Type:
list
- class diffpy.structure.symmetryutilities.ExpandAsymmetricUnit(spacegroup, corepos, coreUijs=None, sgoffset=[0, 0, 0], eps=None)[source]
Bases:
object
Expand asymmetric unit and anisotropic thermal displacement.
- Parameters:
spacegroup (SpaceGroup) – Instance of SpaceGroup.
corepos (array_like) – List of positions in asymmetric unit, it may contain duplicates.
coreUijs (numpy.ndarray, Optional) – Thermal factors for corepos.
sgoffset (list, Optional) – Offset of space group origin
[0, 0, 0]
. Default is[0, 0, 0]
.eps (float, Optional) – Cutoff for duplicate positions. Default is
1.0e-5
.
- spacegroup
Instance of SpaceGroup.
- Type:
- corepos
List of positions in asymmetric unit, it may contain duplicates.
- Type:
array_like
- coreUijs
Thermal factors for corepos. Defaults to zeros.
- Type:
numpy.ndarray
- sgoffset
Offset of space group origin
[0, 0, 0]
. Default to zeros.- Type:
numpy.ndarray
- eps
Cutoff for equivalent positions. Default is
1.0e-5
.- Type:
float
- multiplicity
Multiplicity of each site in corepos.
- Type:
list
- Uisotropy
Bool flags for isotropic sites in corepos.
- Type:
list
- expandedpos
List of equivalent positions per each site in corepos.
- Type:
list
- expandedUijs
List of thermal factors per each site in corepos.
- Type:
list
- class diffpy.structure.symmetryutilities.GeneratorSite(spacegroup, xyz, Uij=array([[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]]), sgoffset=[0, 0, 0], eps=None)[source]
Bases:
object
Storage of data related to a generator positions.
- Parameters:
spacegroup (SpaceGroup) – Instance of SpaceGroup.
xyz (array_like) – Generating site. When xyz is close to special position self.xyz will be adjusted.
Uij (array_like, Optional) – Thermal factors at generator site. Yields self.Uij after adjusting to spacegroup symmetry. Default is zeros.
sgoffset (list, Optional) – Offset of space group origin
[0, 0, 0]
. Default is[0, 0, 0]
.eps (float, Optional) – Cutoff for equal positions. Default is
1.0e-5
.
- xyz
Fractional coordinates of generator site.
- Type:
numpy.ndarray
- Uij
Anisotropic thermal displacement at generator site.
- Type:
numpy.ndarray
- sgoffset
Offset of space group origin
[0, 0, 0]
.- Type:
numpy.ndarray
- eps
Cutoff for equal positions.
- Type:
float
- eqxyz
List of equivalent positions.
- Type:
list
- eqUij
List of displacement matrices at equivalent positions.
- Type:
list
- symops
Nested list of operations per each eqxyz.
- Type:
list
- multiplicity
Generator site multiplicity.
- Type:
int
- Uisotropy
Bool flag for isotropic thermal factors.
- Type:
bool
- invariants
List of invariant operations for generator site.
- Type:
list
- null_space
Null space of all possible differences of invariant rotational matrices, this is a base of symmetry allowed shifts.
- Type:
numpy.ndarray
- Uspace
3D array of independent components of U matrices.
- Type:
numpy.ndarray
- pparameters
List of
(xyz symbol, value)
pairs.- Type:
list
- Uparameters
List of
(U symbol, value)
pairs.- Type:
list
- UFormula(pos, Usymbols=['U11', 'U22', 'U33', 'U12', 'U13', 'U23'])[source]
List of atom displacement formulas with custom parameter symbols.
- Parameters:
pos (array_like) – Fractional coordinates of possibly equivalent site.
Usymbols (list, Optional) – 6 symbols for possible U matrix parameters, default is
["U11", "U22", "U33", "U12", "U13", "U23"]
.
- Returns:
Uformula – U element formulas in a dictionary where keys are from
('U11','U22','U33','U12','U13','U23')
or empty dictionary when pos is not equivalent to generator.- Return type:
dict
- Ucomponents = array([[[1., 0., 0.], [0., 0., 0.], [0., 0., 0.]], [[0., 0., 0.], [0., 1., 0.], [0., 0., 0.]], [[0., 0., 0.], [0., 0., 0.], [0., 0., 1.]], [[0., 1., 0.], [1., 0., 0.], [0., 0., 0.]], [[0., 0., 1.], [0., 0., 0.], [1., 0., 0.]], [[0., 0., 0.], [0., 0., 1.], [0., 1., 0.]]])
6x3x3 array of independent components of U matrices.
- Type:
numpy.ndarray
- eqIndex(pos)[source]
Index of the nearest generator equivalent site.
- Parameters:
pos (array_like) – Fractional coordinates.
- Returns:
Index of the nearest generator equivalent site.
- Return type:
int
- idx2Usymbol = {0: 'U11', 1: 'U12', 2: 'U13', 3: 'U12', 4: 'U22', 5: 'U23', 6: 'U13', 7: 'U23', 8: 'U33'}
Mapping of index to standard U symbol.
- Type:
dict
- positionFormula(pos, xyzsymbols=('x', 'y', 'z'))[source]
Formula of equivalent position with respect to generator site.
- Parameters:
pos (array_like) – Fractional coordinates of possibly equivalent site.
xyzsymbols (tuple, Optional) – Symbols for parametrized coordinates.
- Returns:
Position formulas in a dictionary with keys equal
("x", "y", "z")
or an empty dictionary when pos is not equivalent to generator. Formulas are formatted as[[-][%g*]{x|y|z}] [{+|-}%g]
, for example-x
,z +0.5
,0.25
.- Return type:
dict
- class diffpy.structure.symmetryutilities.SymmetryConstraints(spacegroup, positions, Uijs=None, sgoffset=[0, 0, 0], eps=None)[source]
Bases:
object
Generate symmetry constraints for specified positions.
- Parameters:
spacegroup (SpaceGroup) – Instance of SpaceGroup.
positions (array_like) – List of all positions to be constrained.
Uijs (array_like, Optional) – List of U matrices for all constrained positions.
sgoffset (list, Optional) – Offset of space group origin
[0, 0, 0]
. Default is[0, 0, 0]
.eps (float, Optional) – Cutoff for duplicate positions. Default is
1.0e-5
.
- spacegroup
Instance of SpaceGroup.
- Type:
- positions
All positions to be constrained.
- Type:
numpy.ndarray
- Uijs
Thermal factors for all positions. Defaults to zeros.
- Type:
numpy.ndarray
- sgoffset
Optional offset of space group origin
[0, 0, 0]
.- Type:
numpy.ndarray
- eps
Cutoff for equivalent positions. Default is
1.0e-5
.- Type:
float
- corepos
List of of positions in the asymmetric unit.
- Type:
list
- coremap
Dictionary mapping indices of asymmetric core positions to indices of all symmetry related positions.
- Type:
dict
- poseqns
List of coordinate formula dictionaries per each site. Formula dictionary keys are from
("x", "y", "z")
and the values are formatted as[[-]{x|y|z}%i] [{+|-}%g]
, for example:x0
,-x3
,z7 +0.5
,0.25
.- Type:
list
- pospars
List of
(xyz symbol, value)
pairs.- Type:
list
- Ueqns
List of anisotropic atomic displacement formula dictionaries per each position. Formula dictionary keys are from
('U11','U22','U33','U12','U13','U23')
and the values are formatted as{[%g*][Uij%i]|0}
, for example:U110
,0.5*U2213
,0
.- Type:
list
- Upars
List of
(U symbol, value)
pairs.- Type:
list
- Uisotropy
List of bool flags for isotropic thermal displacements.
- Type:
list
- UFormulas(Usymbols=None)[source]
List of atom displacement formulas with custom parameter symbols.
- Parameters:
Usymbols (list, Optional) – List of custom symbols used in formula strings.
- Returns:
List of atom displacement formula dictionaries per each site. Formula dictionary keys are from
('U11','U22','U33','U12','U13','U23')
and the values are formatted as{[%g*][Usymbol]|0}
, for example:U11
,0.5*@37
,0
.- Return type:
list
- UFormulasPruned(Usymbols=None)[source]
List of atom displacement formula dictionaries with constant items removed.
See also
- Parameters:
Usymbols (list, Optional) – List of custom symbols used in formula strings.
- Returns:
List of atom displacement formulas in tuples of
(U11, U22, U33, U12, U13, U23)
.- Return type:
list
- positionFormulas(xyzsymbols=None)[source]
List of position formulas with custom parameter symbols.
- Parameters:
xyzsymbols (list, Optional) – List of custom symbols used in formula strings.
- Returns:
List of coordinate formulas dictionaries. Formulas dictionary keys are from
("x", "y", "z")
and the values are formatted as[[-]{symbol}] [{+|-}%g]
, for example:x0
,-sym
,@7 +0.5
,0.25
.- Return type:
list
- diffpy.structure.symmetryutilities.equalPositions(xyz0, xyz1, eps)[source]
Equality of two coordinates with optional tolerance.
- Parameters:
xyz0 (array_like) – Fractional coordinates.
xyz1 (array_like) – Fractional coordinates.
eps (float) – Tolerance for equality of coordinates.
- Returns:
True
when two coordinates are closer than eps.- Return type:
bool
- diffpy.structure.symmetryutilities.expandPosition(spacegroup, xyz, sgoffset=[0, 0, 0], eps=None)[source]
Obtain unique equivalent positions and corresponding operations.
- Parameters:
spacegroup (SpaceGroup) – Instance of SpaceGroup.
xyz (list) – Position to be expanded.
sgoffset (list, Optional) – Offset of space group origin
[0, 0, 0]
. Default is[0, 0, 0]
.eps (float, Optional) – Cutoff for equal positions, default is
1.0e-5
.
- Returns:
A tuple with
(list of unique equivalent positions, nested list of `SpaceGroups.SymOp` instances, site multiplicity)
.- Return type:
tuple
- diffpy.structure.symmetryutilities.isSpaceGroupLatPar(spacegroup, a, b, c, alpha, beta, gamma)[source]
Check if space group allows passed lattice parameters.
- Parameters:
spacegroup (SpaceGroup) – Instance of SpaceGroup.
a (float) – Lattice parameters.
b (float) – Lattice parameters.
c (float) – Lattice parameters.
alpha (float) – Lattice parameters.
beta (float) – Lattice parameters.
gamma (float) – Lattice parameters.
- Returns:
True
when lattice parameters are allowed by space group.- Return type:
bool
Note
Crystal system rules:
Benjamin, W. A., Introduction to crystallography, New York (1969), p.60.
- diffpy.structure.symmetryutilities.isconstantFormula(s)[source]
Check if formula string is constant.
- Parameters:
s (str) – Formula string.
- Returns:
True
when argument is a floating point number or a fraction of float with integer.- Return type:
bool
- diffpy.structure.symmetryutilities.nearestSiteIndex(sites, xyz)[source]
Index of the nearest site to a specified position.
- Parameters:
sites (array_like) – List of coordinates.
xyz (array_like) – Single position.
- Returns:
Index of the nearest site.
- Return type:
int
- diffpy.structure.symmetryutilities.positionDifference(xyz0, xyz1)[source]
Smallest difference between two coordinates in periodic lattice.
- Parameters:
xyz0 (array_like) – Fractional coordinates.
xyz1 (array_like) – Fractional coordinates.
- Returns:
dxyz – Smallest difference between two coordinates in periodic lattice with
0 <= dxyz <= 0.5
.- Return type:
numpy.ndarray
- diffpy.structure.symmetryutilities.pruneFormulaDictionary(eqdict)[source]
Remove constant items from formula dictionary.
- Parameters:
eqdict (dict) – Formula dictionary which maps standard variable symbols
("x", "U11")
to string formulas("0", "-x3", "z7 +0.5")
.- Returns:
Pruned formula dictionary.
- Return type:
dict
diffpy.structure.atom module
Provide class Atom for managing properties of an atom in structure model.
- class diffpy.structure.atom.Atom(atype=None, xyz=None, label=None, occupancy=None, anisotropy=None, U=None, Uisoequiv=None, lattice=None)[source]
Bases:
object
Storage of structure information relevant for a single atom.
This class manages atom information such as element symbol, position in fractional and Cartesian coordinates, atomic displacement parameters and so forth.
- Parameters:
atype (str or Atom, Optional) – The string atom type to be set as the element attribute. By default an empty string. When of the Atom type, create a copy of atype and adjust it per other arguments.
xyz (numpy.ndarray, Optional) – Fractional coordinates within the associated lattice. By default
[0, 0, 0]
.label (str, Optional) – A unique string label for referring to this Atom. By default an empty string.
occupancy (float, Optional) – The initial occupancy of this atom, by default
1
.anisotropy (bool, Optional) – The flag for anisotropic thermal displacements parameters. This overrides anisotropy implied by presence of the U or Uisoequiv arguments. Defaults to
False
when not set in any other way.U (numpy.ndarray, Optional) – The 3x3 matrix of anisotropic thermal displacement parameters. When present anisotropy defaults to
True
.Uisoequiv (float, Optional) – The isotropic atomic displacement parameter. The anisotropy defaults to
False
when present. Only one of the U and Uisoequiv arguments may be provided at the same time. Assume zero atomic displacements when U and Uisoequiv are unset.lattice (Lattice, Optional) – Coordinate system for the fractional coordinates xyz. Use the absolute Cartesian system when
None
.
- element
The string type of the atom. An element or ion symbol.
- Type:
str
- xyz
The fractional coordinates in the associated lattice.
- Type:
numpy.ndarray
- label
A unique string label referring to this atom, for example, “C_1”. The label can be used to reference this atom when contained in a Structure object.
- Type:
str
- occupancy
The fractional occupancy of this atom.
- Type:
float
- lattice
Coordinate system for the fractional coordinates xyz and the tensor of atomic displacement parameters U. Use the absolute Cartesian coordinates when
None
.- Type:
Note
Cannot use both U and Uisoequiv arguments at the same time.
- property B11
The
B11
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U11
. When anisotropy isFalse
setting a new value updates entire tensor U.- Type:
float
- property B12
The
B12
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U12
. Setting a new value updates U in a symmetric way. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property B13
The
B13
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U13
. Setting a new value updates U in a symmetric way. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property B22
The
B22
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U22
. When anisotropy isFalse
setting a new value updates entire tensor U.- Type:
float
- property B23
The
B23
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U23
. Setting a new value updates U in a symmetric way. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property B33
The
B33
element of the Debye-Waller matrix.This is equivalent to
8 * pi**2 * U33
. When anisotropy isFalse
setting a new value updates entire tensor U.- Type:
float
- property Bisoequiv
The Debye-Waller isotropic displacement or an equivalent value.
This equals
8 * pi**2 * Uisoequiv
. Setting a new value rescales U tensor to yield equivalent direction-average of Debye-Waller displacements.- Type:
float
- property U
The 3x3 matrix of anisotropic atomic displacements.
For isotropic displacements (when anisotropy is
False
) assignment to U uses only the firstUnew[0, 0]
element and the remaining components of U are adjusted to obtain isotropic tensor in the active lattice.Note
Elements of the U tensor such as
U[0, 1]
should be considered read-only as setting them directly leads to undefined behavior. Use the U11, U22, …, or B11, B22, …, descriptors to set only some U components.- Type:
numpy.ndarray
- property U11
The
U[0, 0]
component of the displacement tensor U.When anisotropy is
False
setting a new value updates entire tensor U.- Type:
float
- property U12
The
U[0, 1]
element of the displacement tensor U.Sets
U[1, 0]
together withU[0, 1]
. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property U13
The
U[0, 2]
element of the displacement tensor U.Sets
U[2, 0]
together withU[0, 2]
. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property U22
The
U[1, 1]
component of the displacement tensor U.When anisotropy is
False
setting a new value updates entire tensor U.- Type:
float
- property U23
The
U[1, 2]
element of the displacement tensor U.Sets
U[2, 1]
together withU[1, 2]
. Assignment has no effect when anisotropy isFalse
.- Type:
float
- property U33
The
U[2, 2]
component of the displacement tensor U.When anisotropy is
False
setting a new value updates entire tensor U.- Type:
float
- property Uisoequiv
The isotropic displacement parameter or an equivalent value.
Setting a new value rescales tensor U so it yields equivalent direction-averaged displacements.
- Type:
float
- property anisotropy
Flag for allowing anisotropic displacement parameters.
When
False
the tensor of thermal displacement parameters U must be isotropic and only its diagonal elements are taken into account.- Type:
bool
- element = ''
Default values of element.
- Type:
str
- label = ''
Default values of label.
- Type:
str
- lattice = None
Default values of lattice.
- Type:
None
- msdCart(vc)[source]
Calculate mean square displacement along the Cartesian vector.
- Parameters:
vc (array_like) – Vector in Cartesian coordinates.
- Returns:
The mean square displacement along vc.
- Return type:
float
- msdLat(vl)[source]
Calculate mean square displacement along the lattice vector.
- Parameters:
vl (array_like) – The vector in lattice coordinates.
- Returns:
The mean square displacement along vl.
- Return type:
float
- occupancy = 1.0
Default values of occupancy.
- Type:
float
- property x
fractional coordinate x, same as
xyz[0]
.- Type:
float
- property xyz_cartn
Atom position in absolute Cartesian coordinates.
This is computed from fractional coordinates xyz and the current lattice setup. Assignment to xyz_cartn or its components is applied on fractional coordinates xyz.
- Type:
numpy.ndarray
- property y
fractional coordinate y, same as
xyz[1]
.- Type:
float
- property z
fractional coordinate z, same as
xyz[2]
.- Type:
float
diffpy.structure.pdffitstructure module
Definition of PDFFitStructure class derived from Structure
- class diffpy.structure.pdffitstructure.PDFFitStructure(*args, **kwargs)[source]
Bases:
Structure
PDFFitStructure –> Structure with extra pdffit member.
- Parameters:
*args – See Structure class constructor.
**kwargs – See Structure class constructor.
- pdffit
Dictionary for storing following extra parameters from PDFFit structure files:
‘scale’, ‘delta1’, ‘delta2’, ‘sratio’, ‘rcut’, ‘spcgr’, ‘dcell’, ‘ncell’
- Type:
dict
- read(filename, format='auto')[source]
Same as Structure.read, but update spcgr value in self.pdffit when parser can get spacegroup.
See Structure.read() for more info.
- Parameters:
filename (str) – File to be loaded.
format (str, Optional) – All structure formats are defined in parsers submodule, when
format == 'auto'
all parsers are tried one by one.
- Returns:
Instance of StructureParser used to load the data.
- Return type:
- readStr(s, format='auto')[source]
Same as Structure.readStr, but update spcgr value in self.pdffit when parser can get spacegroup.
See Structure.readStr() for more info.
- Parameters:
s (str) – String with structure definition.
format (str, Optional) – All structure formats are defined in parsers submodule. When
format == 'auto'
, all parsers are tried one by one.
- Returns:
Instance of StructureParser used to load the data.
- Return type:
diffpy.structure.sgtbxspacegroups module
Extra space group representations generated using sgtbx module from cctbx. Import of this module extends the SpaceGroupList in the SpaceGroups module.
- diffpy.structure.sgtbxspacegroups.sgtbxSpaceGroupList
List of space group instances defined in this module.
- Type:
list