diffpy.Structure package

Submodules

diffpy.Structure.SpaceGroups module

Space group classes and definitions from mmLib and sgtbx.

diffpy.Structure.SpaceGroups.GetSpaceGroup(sgid)

Returns the SpaceGroup instance for the given identifier.

sgid – space group symbol, either short_name or pdb_name,
whatever it means in mmlib. Can be also an integer.

Return space group instance. Raise ValueError when not found.

diffpy.Structure.SpaceGroups.IsSpaceGroupIdentifier(sgid)

Check if identifier can be used as an argument to GetSpaceGroup.

Return bool.

diffpy.Structure.StructureErrors module

Exceptions used in Structure package.

exception diffpy.Structure.StructureErrors.LatticeError

Bases: exceptions.Exception

Exception for impossible lattice parameters.

exception diffpy.Structure.StructureErrors.StructureFormatError

Bases: exceptions.Exception

Exception for failed IO from Structure file

exception diffpy.Structure.StructureErrors.SymmetryError

Bases: exceptions.Exception

Exception raised for invalid symmetry operations.

diffpy.Structure.SymmetryUtilities module

Symmetry utility functions such as expansion of asymmetric unit, and generation of positional constraints.

class diffpy.Structure.SymmetryUtilities.ExpandAsymmetricUnit(spacegroup, corepos, coreUijs=None, sgoffset=[0, 0, 0], eps=None)

Bases: object

Expand asymmetric unit and anisotropic thermal displacement

Data members:

spacegroup – instance of SpaceGroup corepos – list of positions in asymmetric unit,

it may contain duplicates

coreUijs – thermal factors for corepos (defaults to zeros) sgoffset – optional offset of space group origin [0, 0, 0] eps – cutoff for equivalent positions

Calculated data members:
multiplicity – multiplicity of each site in corepos Uisotropy – bool flags for isotropic sites in corepos expandedpos – list of equivalent positions per each site in corepos expandedUijs – list of thermal factors per each site in corepos
class diffpy.Structure.SymmetryUtilities.GeneratorSite(spacegroup, xyz, Uij=array([[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]]), sgoffset=[0, 0, 0], eps=None)

Bases: object

Storage of data related to a generator positions.

Data members:

xyz – fractional coordinates of generator site Uij – anisotropic thermal displacement at generator site sgoffset – offset of space group origin [0, 0, 0] eps – cutoff for equal positions eqxyz – list of equivalent positions eqUij – list of displacement matrices at equivalent positions symops – nested list of operations per each eqxyz multiplicity – generator site multiplicity Uisotropy – bool flag for isotropic thermal factors invariants – list of invariant operations for generator site null_space – null space of all possible differences of invariant

rotational matrices, this is a base of symmetry allowed shifts.

Uspace – 3D array of independent components of U matrices. pparameters – list of (xyz symbol, value) pairs Uparameters – list of (U symbol, value) pairs

UFormula(pos, Usymbols=['U11', 'U22', 'U33', 'U12', 'U13', 'U23'])

List of atom displacement formulas with custom parameter symbols.

pos – fractional coordinates of possibly equivalent site Usymbols – 6 symbols for possible U matrix parameters

Return 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.

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.]]])
eqIndex(pos)

Index of the nearest generator equivalent site

pos – fractional coordinates

Return integer.

idx2Usymbol = {0: 'U11', 1: 'U12', 2: 'U13', 3: 'U12', 4: 'U22', 5: 'U23', 6: 'U13', 7: 'U23', 8: 'U33'}
positionFormula(pos, xyzsymbols=('x', 'y', 'z'))

Formula of equivalent position with respect to generator site

pos – fractional coordinates of possibly equivalent site xyzsymbols – symbols for parametrized coordinates

Return 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”.

signedRatStr(x)

Convert floating point number to signed rational representation. Possible fractional are multiples of 1/3, 1/6, 1/7, 1/9, if these are not close, return “%+g” format.

Return string.

class diffpy.Structure.SymmetryUtilities.SymmetryConstraints(spacegroup, positions, Uijs=None, sgoffset=[0, 0, 0], eps=None)

Bases: object

Generate symmetry constraints for specified positions

Data members:
spacegroup – instance of SpaceGroup positions – all positions to be constrained Uijs – thermal factors for all positions (defaults to zeros) sgoffset – optional offset of space group origin [0, 0, 0] eps – cutoff for equivalent positions
Calculated data members:

corepos – list of of positions in the asymmetric unit coremap – dictionary mapping indices of asymmetric core positions

to indices of all symmetry related positions
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”.

pospars – list of (xyz symbol, value) pairs 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”

Upars – list of (U symbol, value) pairs Uisotropy – list of bool flags for isotropic thermal displacements

UFormulas(Usymbols=None)

List of atom displacement formulas with custom parameter symbols.

Usymbols – list of custom symbols used in formula strings

Return 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”.

UFormulasPruned(Usymbols=None)

List of atom displacement formula dictionaries with constant items removed. See also UFormulas().

Usymbols – list of custom symbols used in formula strings

Return list of atom displacement formulas in tuples of (U11, U22, U33, U12, U13, U23).

UparSymbols()

Return list of standard atom displacement parameter symbols.

UparValues()

Return list of atom displacement parameters values.

positionFormulas(xyzsymbols=None)

List of position formulas with custom parameter symbols.

xyzsymbols – list of custom symbols used in formula strings

Return 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”.

positionFormulasPruned(xyzsymbols=None)

List of position formula dictionaries with constant items removed. See also positionFormulas().

xyzsymbols – list of custom symbols used in formula strings

Return list of coordinate formula dictionaries.

posparSymbols()

Return list of standard position parameter symbols.

posparValues()

Return list of position parameters values.

diffpy.Structure.SymmetryUtilities.equalPositions(xyz0, xyz1, eps)

Equality of two coordinates with optional tolerance.

xyz0, xyz1 – fractional coordinates eps – tolerance for equality of coordinates

Return bool.

diffpy.Structure.SymmetryUtilities.expandPosition(spacegroup, xyz, sgoffset=[0, 0, 0], eps=None)

Obtain unique equivalent positions and corresponding operations.

spacegroup – instance of SpaceGroup xyz – position to be expanded sgoffset – offset of space group origin [0, 0, 0] eps – cutoff for equal positions

Return a tuple with (list of unique equivalent positions, nested list of SpaceGroups.SymOp instances, site multiplicity).

diffpy.Structure.SymmetryUtilities.isSpaceGroupLatPar(spacegroup, a, b, c, alpha, beta, gamma)

Check if space group allows passed lattice parameters.

spacegroup – instance of SpaceGroup a, b, c, alpha, beta, gamma – lattice parameters

Return bool.

diffpy.Structure.SymmetryUtilities.isconstantFormula(s)

Check if formula string is constant. True when argument is a floating point number or a fraction of float with integer.

s – formula string

Return bool.

diffpy.Structure.SymmetryUtilities.nearestSiteIndex(sites, xyz)

Index of the nearest site to a specified position.

sites – list of coordinates or a 2-dimensional numpy.array xyz – single position

Return integer.

diffpy.Structure.SymmetryUtilities.nullSpace(A)

Null space of matrix A.

diffpy.Structure.SymmetryUtilities.positionDifference(xyz0, xyz1)

Smallest difference between two coordinates in periodic lattice.

xyz0, xyz1 – fractional coordinates

Return dxyz, a numpy.array dxyz with 0 <= dxyz <= 0.5.

diffpy.Structure.SymmetryUtilities.pruneFormulaDictionary(eqdict)

Remove constant items from formula dictionary.

eqdict – formula dictionary which maps standard variable symbols
(“x”, “U11”) to string formulas (“0”, “-x3”, “z7 +0.5”)

Return pruned formula dictionary.

diffpy.Structure.atom module

class Atom for storing properties of a single atom

class diffpy.Structure.atom.Atom(atype=None, xyz=None, label=None, occupancy=None, anisotropy=None, U=None, Uisoequiv=None, lattice=None)

Bases: object

Atom –> class for storing atom information

Data members:

element – type of the atom xyz – fractional coordinates, numpy array x, y, z – fractional coordiantes, float properties synced with xyz label – string atom label occupancy – fractional occupancy xyz_cartn – absolute Cartesian coordinates, property synced with xyz anisotropy – flag for anisotropic thermal displacements, property U – anisotropic thermal displacement tensor, property Uij – elements of U tensor, where i, j are from (1, 2, 3),

property
Uisoequiv – isotropic thermal displacement or equivalent value,
property
Bisoequiv – Debye-Waler isotropic temperature factor or equivalent
value, property
lattice – coordinate system for fractional coordinates,
an instance of Lattice or None for Cartesian system
Private data:
_U – storage of U property, 3x3 numpy matrix
B11

B11 element of Debye-Waler displacement tensor

B12

B12 element of Debye-Waler displacement tensor

B13

B13 element of Debye-Waler displacement tensor

B22

B22 element of Debye-Waler displacement tensor

B23

B23 element of Debye-Waler displacement tensor

B33

B33 element of Debye-Waler displacement tensor

Bisoequiv

Debye-Waler isotropic thermal displacement or equivalent value

U

anisotropic thermal displacement tensor.

U11

U11 element of anisotropic displacement tensor

U12

U12 element of anisotropic displacement tensor

U13

U13 element of anisotropic displacement tensor

U22

U22 element of anisotropic displacement tensor

U23

U23 element of anisotropic displacement tensor

U33

U33 element of anisotropic displacement tensor

Uisoequiv

isotropic thermal displacement or equivalent value

anisotropy

flag for anisotropic thermal displacements.

element = ''
label = ''
lattice = None
msdCart(vc)

mean square displacement of an atom along cartesian vector

vc – vector in absolute cartesian coordinates

return mean square displacement

msdLat(vl)

mean square displacement of an atom along lattice vector

vl – vector in lattice coordinates

return mean square displacement

occupancy = 1.0
x

fractional coordinate x.

xyz_cartn

absolute Cartesian coordinates of an atom

y

fractional coordinate y.

z

fractional coordinate z.

diffpy.Structure.bratomsstructure module

definition of BRAtomsStructure class derived from Structure

class diffpy.Structure.bratomsstructure.BRAtomsStructure(*args, **kwargs)

Bases: diffpy.Structure.structure.Structure

BRAtomsStructure –> Structure with extra information for use with Bruce

Ravel’s atoms program.

Data members:
bratoms – dictionary for storing following extra parameters from
atoms .inp files. ‘space’, ‘output’, ‘rmax’, ‘core’, ‘edge’, ‘shift’, ‘nitrogen’, ‘argon’, ‘krypton’

see the following web site for descriptions.

http://leonardo.phys.washington.edu/~ravel/software/doc/Atoms/Atoms/node7.html

diffpy.Structure.lattice module

class Lattice stores properites and provides simple operations in lattice coordinate system.

Module variables:
cartesian – constant instance of Lattice, default Cartesian system
class diffpy.Structure.lattice.Lattice(a=None, b=None, c=None, alpha=None, beta=None, gamma=None, baserot=None, base=None)

Bases: object

Lattice –> general coordinate system

Data members:

a, b, c, alpha, beta, gamma – lattice parameters with cell angles
in degrees, read-write properties. Assignment to cell parameters calls the setLatPar method. The lattice parameters can be also changed by calling the setLatBase method.
ar, br, cr, alphar, betar, gammar – cell parameters for the
reciprocal lattice. Their values reflect the direct cell parameters and the calls to setLatPar or setLatBase methods. The reciprocal angles are in degrees. Read-only properties.
ca, cb, cg, sa, sb, sg – cosines and sines of the lattice angles,
read-only properties
car, cbr, cgr, sar, sbr, sgr – cosines and sines of the reciprocal
lattice angles, read-only properties.

metrics – metrics tensor base – matrix of row base vectors in Cartesian coordinates,

base = stdbase*baserot

stdbase – matrix of base vectors in standard orientation baserot – base rotation matrix recbase – inverse of base matrix, its columns are reciprocal

vectors in Cartesian coordinates

normbase – base with magnitudes of reciprocal vectors recnormbase – inverse of normbase isotropicunit – matrix for unit isotropic displacement parameters

in this coordinate system. Identity matrix when orthonormal.

Note: All data members except a, b, c, alpha, beta, gamma are read-only. Their values get updated by setting the lattice parameters or calling the setLatPar() or setLatBase() methods.

a

Unit cell length a

abcABG()

Return a tuple of 6 lattice parameters.

alpha

Cell angle alpha in degrees

alphar

Reciprocal lattice angle alpha in degrees

angle(u, v)

Return angle(u, v) –> angle of 2 lattice vectors in degrees.

ar

Cell length a of the reciprocal lattice

b

Unit cell length b

beta

Cell angle beta in degrees

betar

Reciprocal lattice angle beta in degrees

br

Cell length b of the reciprocal lattice

c

Unit cell length c

ca

Cosine of the cell angle alpha

car

Cosine of the reciprocal angle alpha

cartesian(u)

Return Cartesian coordinates of a lattice vector.

cb

Cosine of the cell angle beta

cbr

Cosine of the reciprocal angle beta

cg

Cosine of the cell angle gamma

cgr

Cosine of the reciprocal angle gamma

cr

Cell length c of the reciprocal lattice

dist(u, v)

Return distance between 2 points in lattice coordinates.

u – vector or N-by-3 matrix of fractional coordinates. v – vector or N-by-3 matrix of fractional coordinates.

Sizes of u, v must match when both of them are matrices.

Return float or an array of the same length as the matrix.

dot(u, v)

Return dot product of 2 lattice vectors.

fractional(rc)

Return fractional coordinates of a Cartesian vector.

gamma

Cell angle gamma in degrees

gammar

Reciprocal lattice angle gamma in degrees

isanisotropic(umx)

True if ADP matrix is anisotropic in this coordinate system.

norm(xyz)

Calculate norm of a lattice vector.

xyz – vector or an N-by-3 array of fractional coordinates.

Return float or an array of the same length as xyz.

reciprocal()

Return the reciprocal lattice to self.

rnorm(hkl)

Calculate norm of a reciprocal vector.

hkl – vector or an N-by-3 array of reciprocal coordinates.

Return float or an array of the same length as hkl.

sa

Sine of the cell angle alpha

sar

Sine of the reciprocal angle alpha

sb

Sine of the cell angle beta

sbr

Sine of the reciprocal angle beta

setLatBase(base)

Set matrix of unit cell base vectors and calculate corresponding lattice parameters and stdbase, baserot and metrics tensors.

No return value.

setLatPar(a=None, b=None, c=None, alpha=None, beta=None, gamma=None, baserot=None)

set lattice parameters and all related tensors

a, b, c, alpha, beta, gamma – lattice parameters, unit cell angles
are in degrees.

baserot – unit cell rotation, base = stdbase*baserot

Note: parameters with value None will remain unchanged.

No return value.

sg

Sine of the cell angle gamma

sgr

Sine of the reciprocal angle gamma

unitvolume

Cell volume assuming a=b=c=1

volume

Lattice cell volume

diffpy.Structure.lattice.cosd(x)

Return the cosine of x (measured in degrees). Avoid round-off errors for exact cosine values.

diffpy.Structure.lattice.sind(x)

Return the sine of x (measured in degrees). Avoid round-off errors for exact sine values.

diffpy.Structure.mmlibspacegroups module

Space groups defined as a part of the pymmlib.

diffpy.Structure.pdffitstructure module

definition of PDFFitStructure class derived from Structure

class diffpy.Structure.pdffitstructure.PDFFitStructure(*args, **kwargs)

Bases: diffpy.Structure.structure.Structure

PDFFitStructure –> Structure with extra pdffit member

Data members:
pdffit – dictionary for storing following extra parameters from
PDFFit structure files:
‘scale’, ‘delta1’, ‘delta2’, ‘sratio’, ‘rcut’, ‘spcgr’, ‘dcell’, ‘ncell’
read(filename, format='auto')

Same as Structure.read, but update spcgr value in self.pdffit when parser can get spacegroup.

Return instance of StructureParser used to load the data. See Structure.read() for more info.

readStr(s, format='auto')

Same as Structure.readStr, but update spcgr value in self.pdffit when parser can get spacegroup.

Return instance of StructureParser used to load the data. See Structure.readStr() for more info.

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. Notable variables:

sgtbxSpaceGroupList – list of space group instances defined in this module

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, alt_name=None, point_group_name=None, crystal_system=None, pdb_name=None, symop_list=None)

Bases: object

Contains the various names and symmetry operations for one space group.

check_group_name(name)

Checks if the given name is a name for this space group, returns True or False. The space group name can be in several forms: the short name, the longer PDB-style name, or the space group number.

iter_equivalent_positions(vec)

Iterate the symmetry equivalent positions of the argument vector. The vector must already be in fractional coordinates, and the symmetry equivalent vectors are also in fractional coordinates.

iter_symops()

Iterates over all SymOps in the SpaceGroup.

class diffpy.Structure.spacegroupmod.SymOp(R, t)

Bases: object

A subclass of the tuple class for performing one symmetry operation.

is_identity()

Returns True if this SymOp is an identity symmetry operation (no rotation, no translation), otherwise returns False.

diffpy.Structure.structure module

This module defines class Structure.

class diffpy.Structure.structure.Structure(atoms=None, lattice=None, title=None, filename=None, format=None)

Bases: list

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.

Data members:
title – structure description lattice – coordinate system (instance of Lattice) pdffit – None or a dictionary of PDFFit-related metadata
B11

Array of B11 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

B12

Array of B12 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

B13

Array of B13 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

B22

Array of B22 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

B23

Array of B23 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

B33

Array of B33 elements of the Debye-Waller displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

Bisoequiv

Array of Debye-Waller isotropic thermal displacement or equivalent values. Assignment updates the U attribute of all atoms.

U

Array of anisotropic thermal displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U11

Array of U11 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U12

Array of U12 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U13

Array of U13 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U22

Array of U22 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U23

Array of U23 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

U33

Array of U33 elements of the anisotropic displacement tensors. Assignment updates the U and anisotropy attributes of all atoms.

Uisoequiv

Array of isotropic thermal displacement or equivalent values. Assignment updates the U attribute of all atoms.

addNewAtom(*args, **kwargs)

Add new Atom instance to the end of this Structure.

All arguments are forwarded to Atom constructor.

No return value.

angle(aid0, aid1, aid2)

The bond angle at the second of three atoms in degrees.

aid0 – zero based index of the first atom or a string label
such as “Na1”
aid1 – index or string label for the second atom,
where the angle is formed

aid2 – index or string label for the third atom

Return float. Raise IndexError for invalid arguments.

anisotropy

Boolean array for anisotropic thermal displacement flags. Assignment updates the anisotropy attribute of all atoms.

append(a, copy=True)

Append atom to a structure and update its lattice attribute.

a – instance of Atom copy – flag for appending a copy of a.

When False, append a and update a.lattice.

No return value.

assignUniqueLabels()

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][+-]”.

No return value.

composition

Dictionary of chemical symbols and their total occupancies.

copy()

Return a deep copy of this Structure object.

distance(aid0, aid1)

Distance between 2 atoms, no periodic boundary conditions.

aid0 – zero based index of the first atom or a string label
such as “Na1”

aid1 – zero based index or string label of the second atom.

Return float. Raise IndexError for invalid arguments.

element

Character array of atom types. Assignment updates the element attribute of the respective atoms.

extend(atoms, copy=True)

Extend Structure by appending copies from a list of atoms.

atoms – list of Atom instances copy – flag for extending with copies of Atom instances.

When False extend with atoms and update their lattice attributes.

No return value.

getLastAtom()

Return Reference to the last Atom in this structure.

insert(idx, a, copy=True)

Insert atom a before position idx in this Structure.

idx – position in atom list a – instance of Atom copy – flag for inserting a copy of a.

When False, append a and update a.lattice.

No return value.

label

Character array of atom names. Assignment updates the label attribute of all atoms.

lattice

Coordinate system for this Structure.

occupancy

Array of atom occupancies. Assignment updates the occupancy attribute of all atoms.

pdffit = None
placeInLattice(new_lattice)

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

return self

read(filename, format='auto')

Load structure from a file, any original data become lost.

filename – file to be loaded format – all structure formats are defined in Parsers submodule,

when format == ‘auto’ all Parsers are tried one by one

Return instance of data Parser used to process file. This can be inspected for information related to particular format.

readStr(s, format='auto')

Load structure from a string, any original data become lost.

s – string with structure definition format – all structure formats are defined in Parsers submodule,

when format == ‘auto’ all Parsers are tried one by one

Return instance of data Parser used to process input string. This can be inspected for information related to particular format.

title = ''
tolist()

Return atoms in this Structure as a standard Python list.

write(filename, format)

Save structure to file in the specified format

No return value.

Note: available structure formats can be obtained by:
from Parsers import formats
writeStr(format)

return string representation of the structure in specified format

Note: available structure formats can be obtained by:
from Parsers import formats
x

Array of all fractional coordinates x. Assignment updates xyz attribute of all atoms.

xyz

Array of fractional coordinates of all atoms. Assignment updates xyz attribute of all atoms.

xyz_cartn

Array of absolute Cartesian coordinates of all atoms. Assignment updates the xyz attribute of all atoms.

y

Array of all fractional coordinates y. Assignment updates xyz attribute of all atoms.

z

Array of all fractional coordinates z. Assignment updates xyz attribute of all atoms.

diffpy.Structure.utils module

Small shared functions.

diffpy.Structure.utils.atomBareSymbol(smbl)

Remove atom type string stripped of isotope and ion charge symbols. This removes blank and leading [0-9-] or trailing [1-9][+-] characters.

smbl – atom type string such as “Cl-“, “Ca2+” or “12-C”.

Return bare element symbol.

diffpy.Structure.utils.isfloat(s)

True if argument can be converted to float

diffpy.Structure.version module

Definition of __version__, __date__, __gitsha__.

Module contents

classes related to structure of materials Classes:

Atom Lattice Structure PDFFitStructure
Exceptions:
StructureFormatError LatticeError SymmetryError
diffpy.Structure.loadStructure(filename, fmt='auto', **kw)

Load a new structure from a specified file.

filename – file to be loaded fmt – format of the structure file. Must be one of the formats

defined in the Parsers subpackage. When ‘auto’, all Parsers are tried in sequence.

kw – keyword arguments passed to the getParser factory function.

Return a new Structure object. Return PDFFitStructure object for ‘pdffit’ or ‘discus’ formats.