API

The multiphonon package provides tools to convert powder spectrum meausred at direct-geometry inelastic neutron spectrometers to phonon density of states (DOS).

Most users will only need functions in Convenient functions, which reduces raw data in the form of NeXus files to I(Q,E) using Mantid, and then obtain phonon DOS.

Some users may obtain I(Q,E) spectrum through other routes, and only need methods to convert I(Q,E) to phonon DOS. These methods are in the backward transformation section.

Sometimes it may be useful to “simulate” the neutron powder spectrum from a phonon DOS, probably obtained from ab initio calculation or other modeling techniques. These methods are in forward transformation.

Some helper functions users may found useful are in helper functions.

Convenient functions

Calculate DOS from NeXus files

multiphonon.getdos.getDOS(sample_nxs, mt_nxs=None, mt_fraction=0.9, const_bg_fraction=0.0, Emin=-100, Emax=100, dE=1.0, Qmin=0, Qmax=15.0, dQ=0.1, T=300, Ecutoff=50.0, elastic_E_cutoff=(-20.0, 7), M=50.94, C_ms=0.3, Ei=116.446, initdos=None, update_strategy_weights=None, workdir='work', iqe_nxs='iqe.nxs', iqe_h5='iqe.h5', maxiter=10)

Compute DOS from direct-geometry powder neutron scattering spectrum by performing multiphonon and multiple-scattering corrections. This is an iterator. Please call it with an evaluation of the iteration. For example:

>>> output = list(getDOS(...))
Parameters:

sample_nxs : str

Sample Nexus file

mt_nxs : str

Empty can Nexus file

mt_fraction : float

0<=mt_fraction<=1. Amount of empty can data to be subtracted from sample data

const_bg_fraction : float

Constant background fraction

Emin, Emax, dE : floats

Energy transfer axis setting

Qmin, Qmax, dQ : floats

Momentum transfer axis setting

T : float

Temperature (Kelvin)

Ecutoff : float

Maximum phonon energy

elastic_E_cutoff: 2-tuple of floats

cutoff for elastic peak (meV)

M : float

Average atomic mass (u)

C_ms: float

MS = C_ms * MP

Ei : float

Incident energy (meV)

initdos : histogram

initial guess of DOS

update_strategy_weights : floats

Weights for the update strategies (force continuity, area conservation). Useful only if multiple Ei.

work : str

Work directory

maxiter: int

Max iteration

backward transformation

S(Q,E) -> DOS

This is the core functionality of this package.

multiphonon.backward.sqe2dos.sqe2dos(sqe, T, Ecutoff, elastic_E_cutoff, M, C_ms=None, Ei=None, workdir='work', MAX_ITERATION=20, TOLERATION=0.0001, initdos=None, update_strategy_weights=None)

Given a SQE, compute DOS

This is an iterator.

  • Start with an initial guess of DOS and a SQE
  • Calculate SQE of multiphonon scattering
  • Calculate SQE of multiple scattering using C_ms and multiphonon scattering SQE
  • Subtract MS and MP SQE from the experimental SQE to obtain single-phonon SQE
  • Compute a new DOS from the single-phonon SQE
  • Compare the new DOS to the previous one and calculate the difference
  • If difference is large, continue the iteration. Otherwise the new DOS is what we want
Parameters:

sqe : histogram

S(Q, E)

T : float

Temperature (Kelvin)

Ecutoff : float

Maximum phonon energy

elastic_E_cutoff: 2-tuple of floats

cutoff for elastic peak (meV)

M : float

Average atomic mass (u)

C_ms: float

MS = C_ms * MP

Ei : float

Incident energy (meV)

workdir : str

Work directory

initdos : histogram

initial guess of DOS

update_strategy_weights : 2-tuple of floats

Weights for the update strategies (force continuity, area conservation). Useful only if multiple Ei.

MAX_ITERATION: int

Max iteration

TOLERATION: float

Toleration for convergence test

multiphonon.backward.singlephonon_sqe2dos.sqe2dos(sqe, T, Ecutoff, elastic_E_cutoff, M, initdos=None, update_weights=None)

Given a single-phonon SQE, compute DOS

The basic procedure is
  • construct an initial guess of DOS
  • use this DOS to compute 1-phonon SQE
  • for both exp and sim SQE, integrate along Q to obtain S(E)
  • scale the initial guess DOS by the S(E) ratio
  • optionally we can do this again
Parameters:

sqe:histogram

S(Q,E)

T:float

Temperature (kelvin)

Ecutoff:float

Cutoff energy beyond which DOS must be zero

Elastic_E_cutoff: 2-tuple of floats

Cutoff energy bracket for removing the elastic line (unit: meV)

M:float

Atomic mass

initdos:histogram

initial guess of DOS

update_weights:2-tuple of floats

weights for DOS update strategies (continuity, area conservation)

forward transformation

DOS -> S(Q,E)

Simulate S(Q,E) spectrum from phonon DOS using incoherent approximation

multiphonon.forward.dos2sqe(dos, C_ms, sqe, T, M, Ei)

calculate SQE from DOS.

The computed SQE has similar props as the given (experimental) SQE.

Parameters:

dos: histogram

Phonon density of states

C_ms:float

MS = C_ms * MP

sqe:histogram

S(Q,E)

T: float

Temperature (Kelvin)

M: float

Atomic mass

Ei: float

Incident energy (meV)

multiphonon.forward.phonon.sqe(E, g, Qmax=None, Qmin=0, dQ=None, T=300, M=50, N=5, starting_order=2, Emax=None)

compute sum of multiphonon SQE from dos

S = sum_{i=2,N} S_i(Q,E)

Note: single phonon scattering is not included. only 2-phonons and up

Parameters:

E:numpy array of floats

energies in meV

g:numpy array of floats

density of states at the given energies

Qmax:float

maximum value for momentum transfer axis in inverse angstrom

Qmin:float

minimum value for momentum transfer axis in inverse angstrom

dQ:float

the step size for momentum transfer axis in inverse angstrom

T:float

temperature (Kelvin)

M:float

atomic mass

N:integer

maximum number of order for multi-phonon scattering

starting_order:integer

starting number for phonon scattering order

multiphonon.forward.phonon.sqehist(E, g, **kwds)

a simple wrapper of method sqe to return a histogram

Please see method sqe for details of all keyword parameters

Parameters:

E:numpy array of floats

energies in meV

g:numpy array of floats

density of states at the specified energies

helper functions

multiphonon.sqe.interp(iqehist, newE)

compute a new IQE histogram from the given IQE using the new energy array by interpolation

Parameters:

iqehist: histogram

input IQE

newE:numpy array of floats

new energy centers in meV

multiphonon.getdos.reduce2iqe(sample_nxs, Emin, Emax, dE, Qmin, Qmax, dQ, mt_nxs=None, iqe_nxs='iqe.nxs', iqe_h5='iqe.h5', workdir='work')

Reduce sample and (optionally) empty can nxs files and generate I(Q,E) histograms

This is an iterator of processing messages. Please call in this form:

>>> for msg in reduce2iqe(...): print msg
Parameters:

sample_nxs : str

Sample Nexus file

Emin : float

Energy tranfer axis minimum

Emax : float

Energy tranfer axis maximum

dE : float

Energy tranfer axis step size

Qmin : float

Momentum tranfer axis minimum

Qmax : float

Momentum tranfer axis maximum

dQ : float

Momentum tranfer axis step size

mt_nxs : str

Empty can Nexus file

iqe_nxs : str

nexus filename of reduced I(Q,E). temporary file.

iqe_h5: str

output histogram filename of reduced I(Q,E)

workdir: str

path to working directory

multiphonon.redutils.reduce(nxsfile, qaxis, outfile, use_ei_guess=False, ei_guess=None, eaxis=None, tof2E=True, ibnorm='ByCurrent')

reduce a NeXus file to a I(Q,E) histogram using Mantid

This is a wrapper of Mantid algorithms to reduce a NeXus file to IQE histogram.

Parameters:

nxsfile: str

path to nxs file

qaxis: 3-tuple of floats

Momentum transfer axis. (Qmin, dQ, Qmax). unit: inverse angstrom

outfile: str

path to save nxs data

use_ei_guess: boolean

Use incident energy guess

ei_guess: float

Initial guess of incident energy (meV)

eaxis: 3-tuple of floats

Energy transfer axis. (Emin, dE, Emax). unit: meV

tof2E: boolean

Conversion from time of flight axis to energy axis or not. If the NeXus file is in time of flight, tof2E=True If the NeXus file is processed and in energy transfer, tof2E=False

ibnorm: str

Incident beam normalization choice. Allowed values: None, ByCurrent, ToMonitor For more details, see http://docs.mantidproject.org/nightly/algorithms/DgsReduction-v1.html

multiphonon.backward.plotutils.plot_dos_iteration(curdir, total_rounds=None)

plot the DOS for each iteration

Parameters:

curdir : str

path to the root of the working directory for SQE->DOS calculation

total_rounds : integer

number of iterations

multiphonon.backward.plotutils.plot_intermediate_result_se(curdir)

plot the intermediate S(E)

Parameters:

curdir: str

path to one of the iteration working directory for SQE->DOS calculation, for example, work/round-5

multiphonon.backward.plotutils.plot_intermediate_result_sqe(curdir)

plot the intermediate S(Q,E)

Parameters:

curdir: str

path to one of the iteration working directory for SQE->DOS calculation, for example, work/round-5

multiphonon.backward.plotutils.plot_residual(curdir)

plot the residual DOS

Parameters:

curdir : str

path to the root of the working directory for SQE->DOS calculation