diffpy.srfit.pdf package

Submodules

diffpy.srfit.pdf.basepdfgenerator module

PDF profile generator base class.

The BasePDFGenerator class interfaces with SrReal PDF calculators and is used as a base for the PDFGenerator and DebyePDFGenerator classes.

class diffpy.srfit.pdf.basepdfgenerator.BasePDFGenerator(name='pdf')

Bases: diffpy.srfit.fitbase.profilegenerator.ProfileGenerator

Base class for calculating PDF profiles using SrReal.

This works with diffpy.structure.Structure, pyobjcryst.crystal.Crystal and pyobjcryst.molecule.Molecule instances. Note that the managed Parameters are not created until the structure is added.

Attributes: _calc – PDFCalculator or DebyePDFCalculator instance for calculating

the PDF.

_phase – The structure ParameterSet used to calculate the profile. stru – The structure objected adapted by _phase. _lastr – The last value of r over which the PDF was calculated. This is

used to configure the calculator when r changes.

_pool – A multiprocessing.Pool for managing parallel computation.

Managed Parameters: scale – Scale factor delta1 – Linear peak broadening term delta2 – Quadratic peak broadening term qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

Managed ParameterSets: The structure ParameterSet (SrRealParSet instance) used to calculate the profile is named by the user.

Usable Metadata: stype – The scattering type “X” for x-ray, “N” for neutron (see

‘setScatteringType’).
qmax – The maximum scattering vector used to generate the PDF (see
setQmax).
qmin – The minimum scattering vector used to generate the PDF (see
setQmin).

scale – See Managed Parameters. delta1 – See Managed Parameters. delta2 – See Managed Parameters. qbroad – See Managed Parameters. qdamp – See Managed Parameters.

getQmax()

Get the qmax value.

getQmin()

Get the qmin value.

getScatteringType()

Get the scattering type. See ‘setScatteringType’.

parallel(ncpu, mapfunc=None)

Run calculation in parallel.

ncpu – Number of parallel processes. Revert to serial mode when 1. mapfunc – A mapping function to use. If this is None (default),

multiprocessing.Pool.imap_unordered will be used.

No return value.

processMetaData()

Process the metadata once it gets set.

setPhase(parset, periodic=True)

Set the phase that will be used to calculate the PDF.

Set the phase directly with a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts a structure object (from diffpy or pyobjcryst). The passed ParameterSet will be managed by this generator.

parset – A SrRealParSet that holds the structural information.
This can be used to share the phase between multiple BasePDFGenerators, and have the changes in one reflect in another.
periodic – The structure should be treated as periodic (default True).
Note that some structures do not support periodicity, in which case this will be ignored.
setQmax(qmax)

Set the qmax value.

setQmin(qmin)

Set the qmin value.

setScatteringType(stype='X')

Set the scattering type.

stype – “X” for x-ray, “N” for neutron, “E” for electrons,
or any registered type from diffpy.srreal from ScatteringFactorTable.getRegisteredTypes().

Raises ValueError for unknown scattering type.

setStructure(stru, name='phase', periodic=True)

Set the structure that will be used to calculate the PDF.

This creates a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts stru to a ParameterSet interface. See those classes (located in diffpy.srfit.structure) for how they are used. The resulting ParameterSet will be managed by this generator.

stru – diffpy.structure.Structure, pyobjcryst.crystal.Crystal or
pyobjcryst.molecule.Molecule instance. Default None.
name – A name to give to the managed ParameterSet that adapts stru
(default “phase”).
periodic – The structure should be treated as periodic (default
True). Note that some structures do not support periodicity, in which case this will have no effect on the PDF calculation.

diffpy.srfit.pdf.characteristicfunctions module

Form factors (characteristic functions) used in PDF nanoshape fitting.

These are used to calculate the attenuation of the PDF due to a finite size. For a crystal-like nanoparticle, one can calculate the PDF via Gnano(r) = f(r) Gcryst(r), where f(r) is the nanoparticle characteristic function and Gcryst(f) is the crystal PDF.

These functions are meant to be imported and added to a FitContribution using the ‘registerFunction’ method of that class.

diffpy.srfit.pdf.characteristicfunctions.sphericalCF(r, psize)

Spherical nanoparticle characteristic function.

r – distance of interaction psize – The particle diameter

From Kodama et al., Acta Cryst. A, 62, 444-453 (converted from radius to diameter)

diffpy.srfit.pdf.characteristicfunctions.spheroidalCF(r, erad, prad)

Spheroidal characteristic function specified using radii.

Spheroid with radii (erad, erad, prad)

prad – polar radius erad – equatorial radius

erad < prad equates to a prolate spheroid erad > prad equates to a oblate spheroid erad == prad is a sphere

diffpy.srfit.pdf.characteristicfunctions.spheroidalCF2(r, psize, axrat)

Spheroidal nanoparticle characteristic function.

Form factor for ellipsoid with radii (psize/2, psize/2, axrat*psize/2)

r – distance of interaction psize – The equatorial diameter axrat – The ratio of axis lengths

From Lei et al., Phys. Rev. B, 80, 024118 (2009)

diffpy.srfit.pdf.characteristicfunctions.lognormalSphericalCF(r, psize, psig)

Spherical nanoparticle characteristic function with lognormal size distribution.

r – distance of interaction psize – The mean particle diameter psig – The log-normal width of the particle diameter

Here, r is the independent variable, mu is the mean of the distrubution (not of the particle size), and s is the width of the distribution. This is the characteristic function for the lognormal distribution of particle diameter:

F(r, mu, s) = 0.5*Erfc((-mu-3*s^2+Log(r))/(sqrt(2)*s))
  • 0.25*r^3*Erfc((-mu+Log(r))/(sqrt(2)*s))*exp(-3*mu-4.5*s^2)
  • 0.75*r*Erfc((-mu-2*s^2+Log(r))/(sqrt(2)*s))*exp(-mu-2.5*s^2)

The expectation value of the distribution gives the average particle diameter, psize. The variance of the distribution gives psig^2. mu and s can be expressed in terms of these as:

s^2 = log((psig/psize)^2 + 1) mu = log(psize) - s^2/2

Source unknown

diffpy.srfit.pdf.characteristicfunctions.sheetCF(r, sthick)

Nanosheet characteristic function.

r – distance of interaction sthick – Thickness of nanosheet

From Kodama et al., Acta Cryst. A, 62, 444-453

diffpy.srfit.pdf.characteristicfunctions.shellCF(r, radius, thickness)

Spherical shell characteristic function.

radius – Inner radius thickness – Thickness of shell

outer radius = radius + thickness

From Lei et al., Phys. Rev. B, 80, 024118 (2009)

diffpy.srfit.pdf.characteristicfunctions.shellCF2(r, a, delta)

Spherical shell characteristic function.

a – Central radius delta – Thickness of shell

outer radius = a + thickness/2

From Lei et al., Phys. Rev. B, 80, 024118 (2009)

class diffpy.srfit.pdf.characteristicfunctions.SASCF(name, model)

Bases: diffpy.srfit.fitbase.calculator.Calculator

Calculator class for characteristic functions from sas-models.

This class wraps a sas.models.BaseModel to calculate I(Q) related to nanoparticle shape. This I(Q) is inverted to f(r) according to: f(r) = 1 / (4 pi r) * SINFT(I(Q)), where “SINFT” represents the sine Fourier transform.

Attributes: _model – BaseModel object this adapts.

Managed Parameters: These depend on the parameters of the BaseModel object held by _model. They are created from the ‘params’ attribute of the BaseModel. If a dispersion is set for the BaseModel, the dispersion “width” will be accessible under “<parname>_width”, where <parname> is the name a parameter adjusted by dispersion.

diffpy.srfit.pdf.debyepdfgenerator module

PDF profile generator using the Debye equation.

The DebyePDFGenerator class can take a diffpy.structure, pyobjcryst.crystal.Crystal or pyobjcryst.molecule.Molecule object and calculate the PDF from it. This generator is especially appropriate for isolated scatterers, such as nanoparticles and molecules.

class diffpy.srfit.pdf.debyepdfgenerator.DebyePDFGenerator(name='pdf')

Bases: diffpy.srfit.pdf.basepdfgenerator.BasePDFGenerator

A class for calculating the PDF from an isolated scatterer.

This works with diffpy.structure.Structure, pyobjcryst.crystal.Crystal and pyobjcryst.molecule.Molecule instances. Note that the managed Parameters are not created until the structure is added.

Attributes: _calc – DebyePDFCalculator instance for calculating the PDF _phase – The structure ParameterSets used to calculate the profile. stru – The structure objected adapted by _phase. _lastr – The last value of r over which the PDF was calculated. This is

used to configure the calculator when r changes.

Managed Parameters: scale – Scale factor delta1 – Linear peak broadening term delta2 – Quadratic peak broadening term qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

Managed ParameterSets: The structure ParameterSet (SrRealStructure instance) used to calculate the profile is named by the user.

Usable Metadata: stype – The scattering type “X” for x-ray, “N” for neutron (see

‘setScatteringType’).
qmax – The maximum scattering vector used to generate the PDF (see
setQmax).
qmin – The minimum scattering vector used to generate the PDF (see
setQmin).

scale – See Managed Parameters. delta1 – See Managed Parameters. delta2 – See Managed Parameters. qbroad – See Managed Parameters. qdamp – See Managed Parameters.

setPhase(parset, periodic=False)

Set the phase that will be used to calculate the PDF.

Set the phase directly with a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts a structure object (from diffpy or pyobjcryst). The passed ParameterSet will be managed by this generator.

parset – A SrRealParSet that holds the structural information.
This can be used to share the phase between multiple BasePDFGenerators, and have the changes in one reflect in another.
periodic – The structure should be treated as periodic (default True).
Note that some structures do not support periodicity, in which case this will be ignored.
setStructure(stru, name='phase', periodic=False)

Set the structure that will be used to calculate the PDF.

This creates a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts stru to a ParameterSet interface. See those classes (located in diffpy.srfit.structure) for how they are used. The resulting ParameterSet will be managed by this generator.

stru – diffpy.structure.Structure, pyobjcryst.crystal.Crystal or
pyobjcryst.molecule.Molecule instance. Default None.
name – A name to give to the managed ParameterSet that adapts stru
(default “phase”).
periodic – The structure should be treated as periodic (default
False). Note that some structures do not support periodicity, in which case this will have no effect on the PDF calculation.

diffpy.srfit.pdf.pdfcontribution module

PDFContribution class.

This is a custom FitContribution that simplifies the creation of PDF fits.

class diffpy.srfit.pdf.pdfcontribution.PDFContribution(name)

Bases: diffpy.srfit.fitbase.fitcontribution.FitContribution

PDFContribution class.

PDFContribution is a FitContribution that is customized for PDF fits. Data and phases can be added directly to the PDFContribution. Setup of constraints and restraints requires direct interaction with the generator attributes (see setPhase).

Attributes name – A name for this FitContribution. profile – A Profile that holds the measured (and calculated)

signal.
_meta – Metadata dictionary. This is specific to this object,
and not shared with the profile. This is used to record configuration options, like qmax.

_calculators – A managed dictionary of Calculators, indexed by name. _constraints – A set of constrained Parameters. Constraints can be

added using the ‘constrain’ methods.

_generators – A managed dictionary of ProfileGenerators. _parameters – A managed OrderedDict of parameters. _restraints – A set of Restraints. Restraints can be added using the

‘restrain’ or ‘confine’ methods.

_parsets – A managed dictionary of ParameterSets. _eqfactory – A diffpy.srfit.equation.builder.EquationFactory

instance that is used to create constraints and restraints from string

_eq – The FitContribution equation that will be optimized. _reseq – The residual equation. _xname – Name of the x-variable _yname – Name of the y-variable _dyname – Name of the dy-variable

Managed Parameters: scale – Scale factor qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

addPhase(name, parset, periodic=True)

Add a phase that goes into the PDF calculation.

name – A name to give the generator that will manage the PDF
calculation from the passed parameter phase. The parset will be accessible via the name “phase” as an attribute of the generator, e.g., contribution.name.phase, where ‘contribution’ is this contribution and ‘name’ is passed name.
parset – A SrRealParSet that holds the structural information.
This can be used to share the phase between multiple BasePDFGenerators, and have the changes in one reflect in another.
periodic – The structure should be treated as periodic. If this is
True (default), then a PDFGenerator will be used to calculate the PDF from the phase. Otherwise, a DebyePDFGenerator will be used. Note that some structures do not support periodicity, in which case this may be ignored.

Returns the new phase (ParameterSet appropriate for what was passed in stru.)

addStructure(name, stru, periodic=True)

Add a phase that goes into the PDF calculation.

name – A name to give the generator that will manage the PDF
calculation from the passed structure. The adapted structure will be accessible via the name “phase” as an attribute of the generator, e.g. contribution.name.phase, where ‘contribution’ is this contribution and ‘name’ is passed name. (default), then the name will be set as “phase”.
stru – diffpy.structure.Structure, pyobjcryst.crystal.Crystal or
pyobjcryst.molecule.Molecule instance. Default None.
periodic – The structure should be treated as periodic. If this is
True (default), then a PDFGenerator will be used to calculate the PDF from the phase. Otherwise, a DebyePDFGenerator will be used. Note that some structures do not support periodicity, in which case this may be ignored.

Returns the new phase (ParameterSet appropriate for what was passed in stru.)

getQmax()

Get the qmax value.

getQmin()

Get the qmin value.

getScatteringType()

Get the scattering type. See ‘setScatteringType’.

loadData(data)

Load the data in various formats.

This uses the PDFParser to load the data and then passes it to the build-in profile with loadParsedData.

data – An open file-like object, name of a file that contains data
or a string containing the data.
savetxt(fname, **kwargs)

Call numpy.savetxt with x, ycalc, y, dy

This calls on the built-in Profile.

Arguments are passed to numpy.savetxt.

setCalculationRange(xmin=None, xmax=None, dx=None)

Set epsilon-inclusive calculation range.

Adhere to the observed xobs points when dx is the same as in the data. xmin and xmax are clipped at the bounds of the observed data.

Parameters:
  • xmin (float or "obs", optional) – The minimum value of the independent variable. Keep the current minimum when not specified. If specified as “obs” reset to the minimum observed value.
  • xmax (float or "obs", optional) – The maximum value of the independent variable. Keep the current maximum when not specified. If specified as “obs” reset to the maximum observed value.
  • dx (float or "obs", optional) – The sample spacing in the independent variable. When different from the data, resample the x as anchored at xmin.
  • that xmin is always inclusive (unless clipped) xmax is inclusive (Note) –
  • it is within the bounds of the observed data. (if) –
Raises:
  • AttributeError – If there is no observed data.
  • ValueError – When xmin > xmax or if dx <= 0. Also if dx > xmax - xmin.
setQmax(qmax)

Set the qmax value.

setQmin(qmin)

Set the qmin value.

setScatteringType(type='X')

Set the scattering type.

type – “X” for x-ray or “N” for neutron

Raises ValueError if type is not “X” or “N”

diffpy.srfit.pdf.pdfgenerator module

PDF profile generator.

The PDFGenerator class can take a diffpy.structure, pyobjcryst.crystal.Crystal or pyobjcryst.molecule.Molecule object and calculate the crystal PDF from it. The passed structure object is wrapped in a StructureParameter set, which makes its attributes refinable. See the class definition for more details and the examples for its use.

class diffpy.srfit.pdf.pdfgenerator.PDFGenerator(name='pdf')

Bases: diffpy.srfit.pdf.basepdfgenerator.BasePDFGenerator

A class for calculating the PDF from a single crystal structure.

This works with diffpy.structure.Structure, pyobjcryst.crystal.Crystal and pyobjcryst.molecule.Molecule instances. Note that the managed Parameters are not created until the structure is added.

Attributes: _calc – PDFCalculator instance for calculating the PDF _phase – The structure ParameterSet used to calculate the profile. _lastr – The last value of r over which the PDF was calculated. This is

used to configure the calculator when r changes.

Managed Parameters: scale – Scale factor delta1 – Linear peak broadening term delta2 – Quadratic peak broadening term qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

Managed ParameterSets: The structure ParameterSet (SrRealStructure instance) used to calculate the profile is named by the user.

Usable Metadata: stype – The scattering type “X” for x-ray, “N” for neutron (see

‘setScatteringType’).
qmax – The maximum scattering vector used to generate the PDF (see
setQmax).
qmin – The minimum scattering vector used to generate the PDF (see
setQmin).

scale – See Managed Parameters. delta1 – See Managed Parameters. delta2 – See Managed Parameters. qbroad – See Managed Parameters. qdamp – See Managed Parameters.

diffpy.srfit.pdf.pdfparser module

This module contains parsers for PDF data.

PDFParser is suitable for parsing data generated from PDFGetN and PDFGetX.

See the class documentation for more information.

class diffpy.srfit.pdf.pdfparser.PDFParser

Bases: diffpy.srfit.fitbase.profileparser.ProfileParser

Class for holding a diffraction pattern.

Attributes

_format – Name of the data format that this parses (string, default
“”). The format string is a unique identifier for the data format handled by the parser.
_banks – The data from each bank. Each bank contains a (x, y, dx, dy)

tuple: x – A numpy array containing the independent

variable read from the file.
y – A numpy array containing the profile
from the file.
dx – A numpy array containing the uncertainty in x
read from the file. This is 0 if the uncertainty cannot be read.
dy – A numpy array containing the uncertainty read
from the file. This is 0 if the uncertainty cannot be read.

_x – Indpendent variable from the chosen bank _y – Profile from the chosen bank _dx – Uncertainty in independent variable from the chosen bank _dy – Uncertainty in profile from the chosen bank _meta – A dictionary containing metadata read from the file.

General Metadata

filename – The name of the file from which data was parsed. This key
will not exist if data was not read from file.

nbanks – The number of banks parsed. bank – The chosen bank number.

Metadata - These may appear in the metadata dictionary

stype – The scattering type (“X”, “N”) qmin – Minimum scattering vector (float) qmax – Maximum scattering vector (float) qdamp – Resolution damping factor (float) qbroad – Resolution broadening factor (float) spdiameter – Nanoparticle diameter (float) scale – Data scale (float) temperature – Temperature (float) doping – Doping (float)

parseString(patstring)

Parse a string and set the _x, _y, _dx, _dy and _meta variables.

When _dx or _dy cannot be obtained in the data format it is set to 0.

This wipes out the currently loaded data and selected bank number.

Arguments patstring – A string containing the pattern

Raises ParseError if the string cannot be parsed

Module contents

PDF calculation tools.

class diffpy.srfit.pdf.PDFGenerator(name='pdf')

Bases: diffpy.srfit.pdf.basepdfgenerator.BasePDFGenerator

A class for calculating the PDF from a single crystal structure.

This works with diffpy.structure.Structure, pyobjcryst.crystal.Crystal and pyobjcryst.molecule.Molecule instances. Note that the managed Parameters are not created until the structure is added.

Attributes: _calc – PDFCalculator instance for calculating the PDF _phase – The structure ParameterSet used to calculate the profile. _lastr – The last value of r over which the PDF was calculated. This is

used to configure the calculator when r changes.

Managed Parameters: scale – Scale factor delta1 – Linear peak broadening term delta2 – Quadratic peak broadening term qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

Managed ParameterSets: The structure ParameterSet (SrRealStructure instance) used to calculate the profile is named by the user.

Usable Metadata: stype – The scattering type “X” for x-ray, “N” for neutron (see

‘setScatteringType’).
qmax – The maximum scattering vector used to generate the PDF (see
setQmax).
qmin – The minimum scattering vector used to generate the PDF (see
setQmin).

scale – See Managed Parameters. delta1 – See Managed Parameters. delta2 – See Managed Parameters. qbroad – See Managed Parameters. qdamp – See Managed Parameters.

class diffpy.srfit.pdf.DebyePDFGenerator(name='pdf')

Bases: diffpy.srfit.pdf.basepdfgenerator.BasePDFGenerator

A class for calculating the PDF from an isolated scatterer.

This works with diffpy.structure.Structure, pyobjcryst.crystal.Crystal and pyobjcryst.molecule.Molecule instances. Note that the managed Parameters are not created until the structure is added.

Attributes: _calc – DebyePDFCalculator instance for calculating the PDF _phase – The structure ParameterSets used to calculate the profile. stru – The structure objected adapted by _phase. _lastr – The last value of r over which the PDF was calculated. This is

used to configure the calculator when r changes.

Managed Parameters: scale – Scale factor delta1 – Linear peak broadening term delta2 – Quadratic peak broadening term qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

Managed ParameterSets: The structure ParameterSet (SrRealStructure instance) used to calculate the profile is named by the user.

Usable Metadata: stype – The scattering type “X” for x-ray, “N” for neutron (see

‘setScatteringType’).
qmax – The maximum scattering vector used to generate the PDF (see
setQmax).
qmin – The minimum scattering vector used to generate the PDF (see
setQmin).

scale – See Managed Parameters. delta1 – See Managed Parameters. delta2 – See Managed Parameters. qbroad – See Managed Parameters. qdamp – See Managed Parameters.

setPhase(parset, periodic=False)

Set the phase that will be used to calculate the PDF.

Set the phase directly with a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts a structure object (from diffpy or pyobjcryst). The passed ParameterSet will be managed by this generator.

parset – A SrRealParSet that holds the structural information.
This can be used to share the phase between multiple BasePDFGenerators, and have the changes in one reflect in another.
periodic – The structure should be treated as periodic (default True).
Note that some structures do not support periodicity, in which case this will be ignored.
setStructure(stru, name='phase', periodic=False)

Set the structure that will be used to calculate the PDF.

This creates a DiffpyStructureParSet, ObjCrystCrystalParSet or ObjCrystMoleculeParSet that adapts stru to a ParameterSet interface. See those classes (located in diffpy.srfit.structure) for how they are used. The resulting ParameterSet will be managed by this generator.

stru – diffpy.structure.Structure, pyobjcryst.crystal.Crystal or
pyobjcryst.molecule.Molecule instance. Default None.
name – A name to give to the managed ParameterSet that adapts stru
(default “phase”).
periodic – The structure should be treated as periodic (default
False). Note that some structures do not support periodicity, in which case this will have no effect on the PDF calculation.
class diffpy.srfit.pdf.PDFContribution(name)

Bases: diffpy.srfit.fitbase.fitcontribution.FitContribution

PDFContribution class.

PDFContribution is a FitContribution that is customized for PDF fits. Data and phases can be added directly to the PDFContribution. Setup of constraints and restraints requires direct interaction with the generator attributes (see setPhase).

Attributes name – A name for this FitContribution. profile – A Profile that holds the measured (and calculated)

signal.
_meta – Metadata dictionary. This is specific to this object,
and not shared with the profile. This is used to record configuration options, like qmax.

_calculators – A managed dictionary of Calculators, indexed by name. _constraints – A set of constrained Parameters. Constraints can be

added using the ‘constrain’ methods.

_generators – A managed dictionary of ProfileGenerators. _parameters – A managed OrderedDict of parameters. _restraints – A set of Restraints. Restraints can be added using the

‘restrain’ or ‘confine’ methods.

_parsets – A managed dictionary of ParameterSets. _eqfactory – A diffpy.srfit.equation.builder.EquationFactory

instance that is used to create constraints and restraints from string

_eq – The FitContribution equation that will be optimized. _reseq – The residual equation. _xname – Name of the x-variable _yname – Name of the y-variable _dyname – Name of the dy-variable

Managed Parameters: scale – Scale factor qbroad – Resolution peak broadening term qdamp – Resolution peak dampening term

addPhase(name, parset, periodic=True)

Add a phase that goes into the PDF calculation.

name – A name to give the generator that will manage the PDF
calculation from the passed parameter phase. The parset will be accessible via the name “phase” as an attribute of the generator, e.g., contribution.name.phase, where ‘contribution’ is this contribution and ‘name’ is passed name.
parset – A SrRealParSet that holds the structural information.
This can be used to share the phase between multiple BasePDFGenerators, and have the changes in one reflect in another.
periodic – The structure should be treated as periodic. If this is
True (default), then a PDFGenerator will be used to calculate the PDF from the phase. Otherwise, a DebyePDFGenerator will be used. Note that some structures do not support periodicity, in which case this may be ignored.

Returns the new phase (ParameterSet appropriate for what was passed in stru.)

addStructure(name, stru, periodic=True)

Add a phase that goes into the PDF calculation.

name – A name to give the generator that will manage the PDF
calculation from the passed structure. The adapted structure will be accessible via the name “phase” as an attribute of the generator, e.g. contribution.name.phase, where ‘contribution’ is this contribution and ‘name’ is passed name. (default), then the name will be set as “phase”.
stru – diffpy.structure.Structure, pyobjcryst.crystal.Crystal or
pyobjcryst.molecule.Molecule instance. Default None.
periodic – The structure should be treated as periodic. If this is
True (default), then a PDFGenerator will be used to calculate the PDF from the phase. Otherwise, a DebyePDFGenerator will be used. Note that some structures do not support periodicity, in which case this may be ignored.

Returns the new phase (ParameterSet appropriate for what was passed in stru.)

getQmax()

Get the qmax value.

getQmin()

Get the qmin value.

getScatteringType()

Get the scattering type. See ‘setScatteringType’.

loadData(data)

Load the data in various formats.

This uses the PDFParser to load the data and then passes it to the build-in profile with loadParsedData.

data – An open file-like object, name of a file that contains data
or a string containing the data.
savetxt(fname, **kwargs)

Call numpy.savetxt with x, ycalc, y, dy

This calls on the built-in Profile.

Arguments are passed to numpy.savetxt.

setCalculationRange(xmin=None, xmax=None, dx=None)

Set epsilon-inclusive calculation range.

Adhere to the observed xobs points when dx is the same as in the data. xmin and xmax are clipped at the bounds of the observed data.

Parameters:
  • xmin (float or "obs", optional) – The minimum value of the independent variable. Keep the current minimum when not specified. If specified as “obs” reset to the minimum observed value.
  • xmax (float or "obs", optional) – The maximum value of the independent variable. Keep the current maximum when not specified. If specified as “obs” reset to the maximum observed value.
  • dx (float or "obs", optional) – The sample spacing in the independent variable. When different from the data, resample the x as anchored at xmin.
  • that xmin is always inclusive (unless clipped) xmax is inclusive (Note) –
  • it is within the bounds of the observed data. (if) –
Raises:
  • AttributeError – If there is no observed data.
  • ValueError – When xmin > xmax or if dx <= 0. Also if dx > xmax - xmin.
setQmax(qmax)

Set the qmax value.

setQmin(qmin)

Set the qmin value.

setScatteringType(type='X')

Set the scattering type.

type – “X” for x-ray or “N” for neutron

Raises ValueError if type is not “X” or “N”

class diffpy.srfit.pdf.PDFParser

Bases: diffpy.srfit.fitbase.profileparser.ProfileParser

Class for holding a diffraction pattern.

Attributes

_format – Name of the data format that this parses (string, default
“”). The format string is a unique identifier for the data format handled by the parser.
_banks – The data from each bank. Each bank contains a (x, y, dx, dy)

tuple: x – A numpy array containing the independent

variable read from the file.
y – A numpy array containing the profile
from the file.
dx – A numpy array containing the uncertainty in x
read from the file. This is 0 if the uncertainty cannot be read.
dy – A numpy array containing the uncertainty read
from the file. This is 0 if the uncertainty cannot be read.

_x – Indpendent variable from the chosen bank _y – Profile from the chosen bank _dx – Uncertainty in independent variable from the chosen bank _dy – Uncertainty in profile from the chosen bank _meta – A dictionary containing metadata read from the file.

General Metadata

filename – The name of the file from which data was parsed. This key
will not exist if data was not read from file.

nbanks – The number of banks parsed. bank – The chosen bank number.

Metadata - These may appear in the metadata dictionary

stype – The scattering type (“X”, “N”) qmin – Minimum scattering vector (float) qmax – Maximum scattering vector (float) qdamp – Resolution damping factor (float) qbroad – Resolution broadening factor (float) spdiameter – Nanoparticle diameter (float) scale – Data scale (float) temperature – Temperature (float) doping – Doping (float)

parseString(patstring)

Parse a string and set the _x, _y, _dx, _dy and _meta variables.

When _dx or _dy cannot be obtained in the data format it is set to 0.

This wipes out the currently loaded data and selected bank number.

Arguments patstring – A string containing the pattern

Raises ParseError if the string cannot be parsed