diffpy.mpdf package

Submodules

diffpy.mpdf.magstructure module

functions and classes to create magnetic structures for mPDF calculations.

class diffpy.mpdf.magstructure.MagSpecies(struc=None, label='', magIdxs=[0], atoms=None, spins=None, rmaxAtoms=30.0, basisvecs=None, kvecs=None, gS=2.0, gL=0.0, ffparamkey=None, ffqgrid=None, ff=None, useDiffpyStruc=True, latVecs=None, atomBasis=None, spinBasis=None)[source]

Store information for a single species of magnetic atom.

This class takes a diffpy.Structure object and uses it to generate spins based on a set of propagation vectors and basis vectors. For more info about magnetic propagation vectors, visit e.g. http://andrewsteele.co.uk/physics/mmcalc/docs/magneticpropagationvector

Parameters:
  • struc (diffpy.Structure object) – provides lattice parameters and unit cell of desired structure.
  • label (string) – label for this particular magnetic species. Should be different from the labels for any other magnetic species you make.
  • magIdxs (python list) – list of integers giving indices of magnetic atoms in the unit cell
  • atoms (numpy array) – list of atomic coordinates of all the magnetic atoms in the structure; e.g. generated by generateAtomsXYZ()
  • spins (numpy array) – triplets giving the spin vectors of all the atoms, in the same order as the atoms array provided as input.
  • rmaxAtoms (float) – maximum distance from the origin of atomic positions generated by the makeAtoms method.
  • basisvecs (numpy array) – nested three-vector(s) giving the basis vectors to generate the spins. e.g. np.array([[0, 0, 1]]). Any phase factor should be included directly with the basisvecs.
  • kvecs (numpy array) – nested three-vector(s) giving the propagation vectors for the magnetic structure in r.l.u., e.g. np.array([[0.5, 0.5, 0.5]])
  • gS (float) – spin component of the Lande g-factor (g = gS+gL)
  • gL (float) – orbital component of the Lande g-factor
  • ffparamkey (string) – gives the appropriate key for getFFparams()
  • ffqgrid (numpy array) – grid of momentum transfer values used for calculating the magnetic form factor.
  • ff (numpy array) – magnetic form factor.
  • useDiffpyStruc (boolean) – True if atoms/spins to be generated from a diffpy structure object; False if a user-provided unit cell is to be used. Note that there may be some problems with user- provided unit cells with lattice angles strongly deviated from 90 degrees.
  • latVecs (numpy array) – Provides the unit cell lattice vectors as np.array([avec, bvec, cvec]). Only useful if useDiffpyStruc = False.
  • atomBasis (numpy array) – Provides positions of the magnetic atoms in fractional coordinates within the unit cell. Only useful if useDiffpyStruc = False. Example: np.array([[0, 0, 0], [0.5, 0.5, 0.5]])
  • spinBasis (numpy array) – Provides the orientations of the spins in the unit cell, in the same order as atomBasis. Only useful if useDiffpyStruc = False. Example: np.array([[0, 0, 1], [0, 0, -1]])
atomsFromSpins(spinvecs, fractional=True, returnIdxs=False)[source]

Return the atomic positions corresponding to specified spins.

This method calls the diffpy.mpdf.atomsFromSpins() method.

Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • spinvecs (list or array) – spin vectors for which the corresponding atoms should be returned.
  • fractional (boolean) – set as True if the atomic positions are to be returned as fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the atoms will also be returned.
Returns:

List of arrays of atoms corresponding to the spins.

copy()[source]

Return a deep copy of the MagSpecies object.

makeAtoms()[source]

Generate the Cartesian coordinates of the atoms for this species.

makeFF()[source]

Generate the magnetic form factor.

makeSpins()[source]

Generate the Cartesian coordinates of the spin vectors in the structure. Must provide propagation vector(s) and basis vector(s).

runChecks()[source]

Run some simple checks and raise a warning if a problem is found.

spinsFromAtoms(positions, fractional=True, returnIdxs=False)[source]
Return the spin vectors corresponding to specified atomic
positions.

This method calls the diffpy.mpdf.spinsFromAtoms() method.

Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • positions (list or array) – atomic positions for which the corresponding spins should be returned.
  • fractional (boolean) – set as True if the atomic positions are in fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the spins will also be returned.
Returns:

Array consisting of the spins corresponding to the atomic positions.

class diffpy.mpdf.magstructure.MagStructure(struc=None, species=None, atoms=None, spins=None, gfactors=None, rmaxAtoms=30.0, ffqgrid=None, ff=None, label='')[source]

Build on the diffpy.Structure class to include magnetic attributes.

This class takes a diffpy.Structure object and packages additional info relating to magnetic structure, which can then be fed to an MPDFcalculator object.

Parameters:
  • struc (diffpy.Structure object) – provides lattice parameters and unit cell of desired structure.
  • species (python dictionary) – dictionary of magnetic species in the structure. The values are MagSpecies objects.
  • atoms (numpy array) – list of atomic coordinates of all the magnetic atoms in the structure; e.g. generated by generateAtomsXYZ()
  • spins (numpy array) – triplets giving the spin vectors of all the atoms, in the same order as the atoms array provided as input.
  • gfactors (numpy array) – Lande g-factors of the magnetic moments
  • rmaxAtoms (float) – maximum distance from the origin of atomic positions generated by the makeAtoms method.
  • ffqgrid (numpy array) – grid of momentum transfer values used for calculating the magnetic form factor.
  • ff (numpy array) – magnetic form factor. Should be same shape as ffqgrid.
  • label (string) – Optional descriptive string for the MagStructure.
atomsFromSpins(spinvecs, fractional=True, returnIdxs=False)[source]

Return the atomic positions corresponding to specified spins.

This method calls the diffpy.mpdf.atomsFromSpins() method.

Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • spinvecs (list or array) – spin vectors for which the corresponding atoms should be returned.
  • fractional (boolean) – set as True if the atomic positions are to be returned as fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the atoms will also be returned.
Returns:

List of arrays of atoms corresponding to the spins.

copy()[source]

Return a deep copy of the MagStructure object.

getCoordsFromSpecies()[source]

Read in atomic positions and spins from magnetic species.

This differs from makeSpins() and makeAtoms() because it simply loads the atoms and spins from the species without re-generating them from the structure.

getSpeciesIdxs()[source]

Return a dictionary with the starting index in the atoms and spins arrays corresponding to each magnetic species.

loadSpecies(magSpec)[source]

Load in an already-existing MagSpecies object

Parameters:magSpec (MagSpecies object) – The magnetic species to be imported into the structure.
makeAll()[source]

Shortcut method to generate atoms, spins, g-factors, and form factor for the magnetic structure all in one go.

makeAtoms()[source]

Generate the Cartesian coordinates of the atoms for this species.

Parameters:
  • fromUnitCell (boolean) – True if atoms/spins to be generated from a unit cell provided by the user; False if the diffpy structure object is to be used.
  • unitcell (numpy array) – Provides the unit cell lattice vectors as np.array((avec, bvec, cvec)).
  • atombasis (numpy array) – Provides positions of the magnetic atoms in fractional coordinates within the unit cell.
  • cell (spin) – Provides the orientations of the spins in the unit cell, in the same order as atombasis
makeFF()[source]

Generate the properly weighted average magnetic form factor of all the magnetic species in the structure.

makeGfactors()[source]

Generate an array of Lande g-factors in the same order as the spins in the MagStructure.

makeSpecies(label, magIdxs=None, atoms=None, spins=None, basisvecs=None, kvecs=None, gS=2.0, gL=0.0, ffparamkey=None, ffqgrid=None, ff=None)[source]

Create a MagSpecies object and add it to the species dictionary.

Parameters:
  • label (string) – label for this particular magnetic species. Should be different from the labels for any other magnetic species you make.
  • magIdxs (python list) – list of integers giving indices of magnetic atoms in the unit cell
  • atoms (numpy array) – list of atomic coordinates of all the magnetic atoms in the structure; e.g. generated by generateAtomsXYZ()
  • spins (numpy array) – triplets giving the spin vectors of all the atoms, in the same order as the atoms array provided as input.
  • basisvecs (numpy array) – nested three-vector(s) giving the basis vectors to generate the spins. e.g. np.array([[0, 0, 1]]). Any phase factor should be included directly with the basisvecs.
  • kvecs (numpy array) – nested three-vector(s) giving the propagation vectors for the magnetic structure in r.l.u., e.g. np.array([[0.5, 0.5, 0.5]])
  • gS (float) – spin component of the Lande g-factor (g = gS+gL)
  • gL (float) – orbital component of the Lande g-factor
  • ffparamkey (string) – gives the appropriate key for getFFparams()
  • ffqgrid (numpy array) – grid of momentum transfer values used for calculating the magnetic form factor.
  • ff (numpy array) – magnetic form factor.
makeSpins()[source]

Generate the Cartesian coordinates of the spin vectors in the structure. Calls the makeSpins() method for each MagSpecies in the species dictionary and concatenates them together.

removeSpecies(label, update=True)[source]

Remove a magnetic species from the species dictionary.

Parameters:
  • label (string) – key for the dictionary entry to be removed.
  • update (boolean) – if True, the MagStructure will update its atoms and spins with the removed species now excluded.
runChecks()[source]

Run some simple checks and raise a warning if a problem is found.

spinsFromAtoms(positions, fractional=True, returnIdxs=False)[source]
Return the spin vectors corresponding to specified atomic
positions.

This method calls the diffpy.mpdf.spinsFromAtoms() method.

Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • positions (list or array) – atomic positions for which the corresponding spins should be returned.
  • fractional (boolean) – set as True if the atomic positions are in fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the spins will also be returned.
Returns:

Array consisting of the spins corresponding to the atomic positions.

visualize(atoms, spins, showcrystalaxes=False, axesorigin=array([0, 0, 0]))[source]

Generate a crude 3-d plot to visualize the selected spins.

Parameters:
  • atoms (numpy array) – array of atomic positions of spins to be visualized.
  • spins (numpy array) – array of spin vectors in same order as atoms.
  • showcrystalaxes (boolean) – if True, will display the crystal axes determined from the first magnetic species in the MagStructure
  • axesorigin (array) – position at which the crystal axes should be displayed
diffpy.mpdf.magstructure.atomsFromSpins(magstruc, spinvecs, fractional=True, returnIdxs=False)[source]

Return the atomic positions corresponding to specified spins.

Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • spinvecs (list or array) – spin vectors for which the corresponding atoms should be returned.
  • fractional (boolean) – set as True if the atomic positions are to be returned as fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the atoms will also be returned.
Returns:

List of arrays of atoms corresponding to the spins.

diffpy.mpdf.magstructure.generateAtomsXYZ(struc, rmax=30.0, magIdxs=[0], square=False)[source]

Generate array of atomic Cartesian coordinates from a given structure.

Parameters:
  • struc (diffpy.Structure object) – provides lattice parameters and unit cell of the desired structure
  • rmax (float) – largest distance from central atom that should be included
  • magIdxs (python list) – list of integers giving indices of magnetic atoms in the unit cell
  • square (boolean) – if not True, atoms within a given radius from the origin will be returned; if True, then the full grid will be returned rather than just a spherical selection.
Returns:

numpy array of triples giving the Cartesian coordinates of all the

magnetic atoms. Atom closest to the origin placed first in array.

Note: If square = True, this may have problems for structures that have
a several distorted unit cell (i.e. highly non-orthorhombic).
diffpy.mpdf.magstructure.generateFromUnitCell(unitcell, atombasis, spinbasis, rmax=30.0)[source]

Generate array of atomic Cartesian coordinates from a given structure.

Parameters:
  • unitcell (numpy array) – np.array([avec, bvec, cvec])
  • atombasis (numpy array) – gives positions of magnetic atoms in fractional coordinates; np.array([pos1, pos2, pos3, ...])
  • spinbasis (numpy array) – gives orientations of the magnetic moments in the unit cell, in the same order as atombasis
  • rmax (float) – largest distance from central atom that should be included
Returns:

atoms = numpy array of triples giving the Cartesian coordinates of all

the magnetic atoms. Atom closest to the origin placed first in array.

spins = numpy array of triples giving the Cartesian coordinates of all

the spins in the structure, in the same order as atoms.

Note: This will only work well for structures that can be expressed with a
unit cell that is close to orthorhombic or higher symmetry.
diffpy.mpdf.magstructure.generateSpinsXYZ(struc, atoms=array([], shape=(1, 0), dtype=float64), kvecs=array([[0, 0, 0]]), basisvecs=array([[0, 0, 1]]))[source]

Generate array of 3-vectors representing the spins in a structure.

Parameters:
  • struc (diffpy.Structure object) – provides lattice parameters and unit cell of the desired structure
  • atoms (numpy array) – list of atomic coordinates of all the magnetic atoms in the structure; e.g. generated by generateAtomsXYZ()
  • kvecs (numpy array) – list of three-vectors giving the propagation vector(s) of the magnetic structure
  • basisvecs (numpy array) – list of three-vectors describing the spin located at the spin origin.
Returns:

numpy array of triples giving the spin vectors of all the magnetic

atoms, in the same order as the atoms array provided as input.

diffpy.mpdf.magstructure.getFFparams(name, j2=False)[source]

Get list of parameters for approximation of magnetic form factor

Parameters:
  • name (str) – Name of magnetic ion in form ‘Mn2’ for Mn2+, etc.
  • j2 (boolean) – True of the j2 approximation should be calculated; otherwise, the j0 approximation is calculated.
Returns:

Python list of the 7 coefficients in the analytical approximation

given at e.g. https://www.ill.eu/sites/ccsl/ffacts/ffachtml.html

diffpy.mpdf.magstructure.spinsFromAtoms(magstruc, positions, fractional=True, returnIdxs=False)[source]
Return the spin vectors corresponding to specified atomic
positions.
Parameters:
  • magstruc – MagSpecies or MagStructure object containing atoms and spins
  • positions (list or array) – atomic positions for which the corresponding spins should be returned.
  • fractional (boolean) – set as True if the atomic positions are in fractional coordinates of the crystallographic lattice vectors.
  • returnIdxs (boolean) – if True, the indices of the spins will also be returned.
Returns:

Array consisting of the spins corresponding to the atomic positions.

diffpy.mpdf.magstructure.visualizeSpins(atoms, spins)[source]

Set up a 3d plot showing an arrangement of spins.

Parameters:
  • atoms (numpy array) – array of atomic positions of spins to be visualized.
  • spins (numpy array) – array of spin vectors in same order as atoms.
Returns:

matplotlib figure object with a quiver plot on 3d axes.

diffpy.mpdf.mpdfcalculator module

functions and classes to perform mPDF calculations

class diffpy.mpdf.mpdfcalculator.MPDFcalculator(magstruc=None, calcList=[0], maxextension=10.0, gaussPeakWidth=0.1, dampRate=0.0, dampPower=2.0, qmin=-1.0, qmax=-1.0, rmin=0.0, rmax=20.0, rstep=0.01, ordScale=0.3989422804014327, paraScale=1.0, rmintr=-5.0, rmaxtr=5.0, label='')[source]

Create an MPDFcalculator object to help calculate mPDF functions.

This class is loosely modelled after the PDFcalculator class in diffpy. At minimum, it requires a magnetic structure with atoms and spins, and it calculates the mPDF from that. Various other options can be specified for the calculated mPDF.

Parameters:
  • magstruc (MagStructure object) – provides information about the magnetic structure. Must have arrays of atoms and spins.
  • calcList (python list) – list giving the indices of the atoms array specifying the atoms to be used as the origin when calculating the mPDF.
  • maxextension (float) – extension of the r-grid on which the mPDF is calculated to properly account for contribution of pairs just outside the boundary.
  • gaussPeakWidth (float) – std deviation (in Angstroms) of Gaussian peak to be convoluted with the calculated mPDF to simulate thermal motion.
  • dampRate (float) – generalized (“stretched”) exponential damping rate of the mPDF.
  • dampPower (float) – power of the generalized exponential function.
  • qmin (float) – minimum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.
  • qmax (float) – maximum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.
  • rmin (float) – minimum value of r for which mPDF should be calculated.
  • rmax (float) – maximum value of r for which mPDF should be calculated.
  • rstep (float) – step size for r-grid of calculated mPDF.
  • ordScale (float) – overall scale factor for the mPDF function f(r).
  • paraScale (float) – scale factor for the paramagnetic part of the unnormalized mPDF function D(r).
  • rmintr (float) – minimum value of r for the Fourier transform of the magnetic form factor required for unnormalized mPDF.
  • rmaxtr (float) – maximum value of r for the Fourier transform of the magnetic form factor required for unnormalized mPDF.
  • drtr (float) – step size for r-grid used for calculating Fourier transform of magnetic form mactor.
  • label (string) – Optional descriptive string for the MPDFcalculator.
calc(normalized=True, both=False)[source]

Calculate the magnetic PDF.

Parameters:
  • normalized (boolean) – indicates whether or not the normalized mPDF should be returned.
  • both (boolean) – indicates whether or not both normalized and unnormalized mPDF quantities should be returned.
Returns: numpy array giving the r grid of the calculation, as well as
one or both the mPDF quantities.
copy()[source]

Return a deep copy of the MPDFcalculator object.

plot(normalized=True, both=False)[source]

Plot the magnetic PDF.

Parameters:
  • normalized (boolean) – indicates whether or not the normalized mPDF should be plotted.
  • both (boolean) – indicates whether or not both normalized and unnormalized mPDF quantities should be plotted.
rgrid()[source]

Return the current r grid of the mPDF calculator.

runChecks()[source]

Run some quick checks to help with troubleshooting.

diffpy.mpdf.mpdfcalculator.calculateDr(r, fr, q, ff, paraScale=1.0, rmintr=-5.0, rmaxtr=5.0, drtr=0.01, qmin=0, qmax=-1)[source]

Calculate the unnormalized mPDF quantity D(r).

This module requires a normalized mPDF as an input, as well as a magnetic form factor and associated q grid.

Parameters:
  • r (numpy array) – r grid for the properly normalized mPDF.
  • fr (numpy array) – the properly normalized mPDF.
  • q (numpy array) – grid of momentum transfer values used for calculating the magnetic form factor.
  • ff (numpy array) – magnetic form factor. Same shape as ffqgrid.
  • paraScale (float) – scale factor for the paramagnetic part of the unnormalized mPDF function D(r).
  • rmintr (float) – minimum value of r for the Fourier transform of the magnetic form factor required for unnormalized mPDF.
  • rmaxtr (float) – maximum value of r for the Fourier transform of the magnetic form factor required for unnormalized mPDF.
  • drtr (float) – step size for r-grid used for calculating Fourier transform of magnetic form mactor.
  • qmin (float) – minimum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.
  • qmax (float) – maximum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.

Returns: numpy array for the unnormalized mPDF Dr.

diffpy.mpdf.mpdfcalculator.calculatemPDF(xyz, sxyz, gfactors=array([ 2.]), calcList=array([0]), rstep=0.01, rmin=0.0, rmax=20.0, psigma=0.1, qmin=0, qmax=-1, dampRate=0.0, dampPower=2.0, maxextension=10.0, orderedScale=0.3989422804014327)[source]

Calculate the normalized mPDF.

At minimum, this module requires input lists of atomic positions and spins.

Parameters:
  • xyz (numpy array) – list of atomic coordinates of all the magnetic atoms in the structure.
  • sxyz (numpy array) – triplets giving the spin vectors of all the atoms, in the same order as the xyz array provided as input.
  • gfactors (numpy array) – Lande g-factors of spins in same order as spin array.
  • calcList (python list) – list giving the indices of the atoms array specifying the atoms to be used as the origin when calculating the mPDF.
  • rstep (float) – step size for r-grid of calculated mPDF.
  • rmin (float) – minimum value of r for which mPDF should be calculated.
  • rmax (float) – maximum value of r for which mPDF should be calculated.
  • psigma (float) – std deviation (in Angstroms) of Gaussian peak to be convoluted with the calculated mPDF to simulate thermal motion.
  • qmin (float) – minimum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.
  • qmax (float) – maximum experimentally accessible q-value (to be used for simulating termination ripples). If <0, no termination effects are included.
  • dampRate (float) – generalized (“stretched”) exponential damping rate of the mPDF.
  • dampPower (float) – power of the generalized exponential function.
  • maxextension (float) – extension of the r-grid on which the mPDF is calculated to properly account for contribution of pairs just outside the boundary.
  • ordScale (float) – overall scale factor for the mPDF function f(r).

Returns: numpy arrays for r and the mPDF fr.

diffpy.mpdf.mpdfcalculator.cosTransform(q, fq, rmin=0.0, rmax=50.0, rstep=0.1)[source]

Compute the cosine Fourier transform of a function.

This method uses direct integration rather than an FFT and doesn’t require an even grid. The grid for the Fourier transform is even and specifiable.

Parameters:
  • q (numpy array) – independent variable for function to be transformed
  • fq (numpy array) – dependent variable for function to be transformed
  • rmin (float, default = 0.0) – min value of conjugate independent variable grid
  • rmax (float, default = 50.0) – maximum value of conjugate independent variable grid
  • rstep (float, default = 0.1) – grid spacing for conjugate independent variable
Returns:

r – independent variable grid for transformed quantity

fr (numpy array): cosine Fourier transform of fq

Return type:

numpy array

diffpy.mpdf.mpdfcalculator.cv(x1, y1, x2, y2)[source]

Perform the convolution of two functions and give the correct output.

Parameters:
  • x1 (numpy array) – independent variable of first function; must be in ascending order
  • y1 (numpy array) – dependent variable of first function
  • x2 (numpy array) – independent variable of second function; must have same grid spacing as x1
  • y2 (numpy array) – dependent variable of second function
Returns:

xcv

independent variable of convoluted function, has

dimension len(x1) + len(x2) - 1

ycv (numpy array): convolution of y1 and y2, same shape as xcv

Return type:

numpy array

diffpy.mpdf.mpdfcalculator.getDiffData(fileNames=[], fmt='pdfgui', writedata=False)[source]

Extract the fit residual from a structural PDF fit.

Parameters:
  • fileNames (python list) – list of paths to the files containing the fit information (e.g. calculated and experimental PDF, as in the .fgr files from PDFgui exported fits)
  • fmt (string) – string identifying the format of the file(s). Options are currently just ‘pdfgui’.
  • writedata (boolean) – whether or not the output should be saved to file
Returns:

r – same r-grid as contained in the fit file

diff (numpy array): the structural PDF fit residual (i.e. the mPDF)

Return type:

numpy array

diffpy.mpdf.mpdfcalculator.jCalc(q, params=[0.2394, 26.038, 0.4727, 12.1375, 0.3065, 3.0939, -0.01906], j2=False)[source]

Calculate the magnetic form factor j0.

This method for calculating the magnetic form factor is based on the approximate empirical forms based on the tables by Brown, consisting of the sum of 3 Gaussians and a constant.

Parameters:
  • q (numpy array) – 1-d grid of momentum transfer on which the form factor is to be computed
  • params (python list) – provides the 7 numerical coefficients. The default is an average form factor of 3d j0 approximations.
  • j2 (boolean) – if True, calculate the j2 approximation for orbital angular momentum contributions
Returns:

numpy array with same shape as q giving the magnetic form factor j0 or j2.

Module contents

Tools for magnetic pair distribution function analysis