Spectra
These are generic objects for storing 1D/2D spectra e.g. density of states or S(Q,w) maps.
Metadata
All spectra objects have a metadata
attribute. By default it is an empty
dictionary, but it can contain keys/values to help describe the contained
spectra. The keys should be strings and values should be only strings or
integers. Note there are some special ‘functional’ keys, see the metadata
docstring below for each specific spectrum class for details.
Spectrum1D
Adding Spectrum1D
Two Spectrum1D
objects can be added together (provided they have the
same x_data
axes) with the +
operator. This will add their y_data
and return a single Spectrum1D
. Note that any metadata key/value pairs
that aren’t common to both spectra will be omitted from the new object. For
example:
from euphonic import Spectrum1D
pdos_si = Spectrum1D.from_json_file('si_pdos.json')
pdos_o = Spectrum1D.from_json_file('o_pdos.json')
total_dos = pdos_si + pdos_o
Broadening
A 1D spectrum can be broadened using
Spectrum1D.broaden
,
which broadens along the x-axis and returns a new Spectrum1D
object. It can broaden with either a Gaussian or Lorentzian and requires
a broadening FWHM in the same type of units as x_data
. For example:
from euphonic import ureg, Spectrum1D
dos = Spectrum1D.from_json_file('dos.json')
fwhm = 1.5*ureg('meV')
dos_broaden = dos.broaden(fwhm, shape='lorentz')
Variable-width broadening can be achieved by passing a Python “Callable” (i.e. function) that accepts a position value as input and returns a broadening width. For example, the energy-dependent resolution of the TOSCA spectrometer is typically treated as a Gaussian following the quadratic function \(\sigma = 2.5 + 0.005 \omega + 0.0000001 \omega\) where \(sigma\) and \(omega\) are the standard deviation and energy transfer in cm:math:^{-1}. To apply this resolution function to a Spectrum1D we package the polynomial into a suitable Callable:
from numpy.polynomial import Polynomial
def tosca_resolution(energy):
poly_in_wavenumber = Polynomial([2.5, 0.005, 1e-7])
return poly_in_wavenumber(energy.to('1/cm').magnitude) * ureg('1/cm')
dos_tosca = dos.broaden(tosca_resolution, shape='gauss', width_convention='std')
Plotting
See Plotting
Docstring
- class Spectrum1D(x_data, y_data, x_tick_labels=None, metadata=None)
For storing generic 1D spectra e.g. density of states
- Variables
x_data – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))
y_data – Shape (n_x_data,) float Quantity. The plot data in y
x_tick_labels – Sequence[Tuple[int, str]] or None. Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied to
metadata –
Dict[str, Union[int, str]]. Contains metadata about the spectrum. Keys should be strings and values should be strings or integers There are some functional keys:
’label’ : str. This is used label lines on a 1D plot
- __init__(x_data, y_data, x_tick_labels=None, metadata=None)
- Parameters
x_data (
Quantity
) – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))y_data (
Quantity
) – Shape (n_x_data,) float Quantity. The plot data in yx_tick_labels (
Optional
[Sequence
[Tuple
[int
,str
]]]) – Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied tometadata (
Optional
[Dict
[str
,Union
[int
,str
]]]) –Contains metadata about the spectrum. Keys should be strings and values should be strings or integers There are some functional keys:
’label’ : str. This is used label lines on a 1D plot
- copy()
Get an independent copy of spectrum
- Return type
TypeVar
(T
, bound= Spectrum1D)
- to_dict()
Convert to a dictionary. See Spectrum1D.from_dict for details on keys/values
- Return type
Dict
[str
,Any
]- Returns
dict
- to_text_file(filename, fmt=None)
Write to a text file. The header contains metadata and unit information, the first column contains x_data and the second column contains y_data. Note that text files written in this format cannot be read back in by Euphonic.
- Parameters
filename (
str
) – Name of the text file to write tofmt (
Union
[str
,Sequence
[str
],None
]) – A format specifier or sequence of specifiers (one for each column), to be passed to numpy.savetxt
- Return type
None
- classmethod from_dict(d)
Convert a dictionary to a Spectrum1D object
- Parameters
d (
Dict
[str
,Any
]) –A dictionary with the following keys/values:
’x_data’: (n_x_data,) or (n_x_data + 1,) float ndarray
’x_data_unit’: str
’y_data’: (n_x_data,) float ndarray
’y_data_unit’: str
There are also the following optional keys:
’x_tick_labels’: list of (int, string) tuples
’metadata’: dict
- Return type
TypeVar
(T
, bound= Spectrum1D)- Returns
spectrum
- classmethod from_castep_phonon_dos(filename, element=None)
Reads DOS from a CASTEP .phonon_dos file
- Parameters
filename (
str
) – The path and name of the .phonon_dos file to readelement (
Optional
[str
]) – Which element’s PDOS to read. If not supplied reads the total DOS.
- Return type
TypeVar
(T
, bound= Spectrum1D)
- broaden(x_width: Quantity, shape: Literal['gauss', 'lorentz'] = 'gauss', method: Optional[Literal['convolve']] = None) T
- broaden(x_width: Callable[[Quantity], Quantity], shape: Literal['gauss', 'lorentz'] = 'gauss', method: Optional[Literal['convolve']] = None, width_lower_limit: Optional[Quantity] = None, width_convention: Literal['fwhm', 'std'] = 'fwhm', width_interpolation_error: float = 0.01, width_fit: Literal['cheby-log', 'cubic'] = 'cheby-log') T
Broaden y_data and return a new broadened spectrum object
- Parameters
x_width – Scalar float Quantity, giving broadening FWHM for all y data, or a function handle accepting Quantity consistent with x-axis as input and returning FWHM corresponding to input array. This would typically be an energy-dependent resolution function.
shape – One of {‘gauss’, ‘lorentz’}. The broadening shape
method – Can be None or ‘convolve’. Currently the only broadening method available is convolution with a broadening kernel but this may not produce correct results for unequal bin widths. To use convolution anyway, explicitly set method=’convolve’
width_lower_limit – Set a lower bound to width obtained calling x_width function. By default, this is equal to bin width. To disable, set to -Inf.
width_convention – By default (‘fwhm’), x_width is interpreted as full-width half-maximum. Set to ‘std’ to instead define standard deviation.
width_interpolation_error – When x_width is a callable function, variable-width broadening is implemented by an approximate kernel-interpolation scheme. This parameter determines the target error of the kernel approximations.
width_fit – Select parametrisation of kernel width spacing to width_interpolation_error. ‘cheby-log’ is recommended: for shape ‘gauss’, ‘cubic’ is also available.
- Return type
TypeVar
(T
, bound= Spectrum1D)- Returns
broadened_spectrum – A new Spectrum1D object with broadened y_data
- Raises
ValueError – If shape is not one of the allowed strings
ValueError – If method is None and bins are not of equal size
- assert_regular_bins(message='', rtol=1e-05, atol=0.0)
Raise AssertionError if x-axis bins are not evenly spaced.
- Parameters
message (
str
) – Text appended to ValueError for more informative output.rtol (
float
) – Relative tolerance for ‘close enough’ valuesatol (
float
) – Absolute tolerance for ‘close enough’ values. Note this is a bare float and follows the stored units of the bins.
- Return type
None
- classmethod from_json_file(filename)
Read from a JSON file. See from_dict for required fields
- Parameters
filename (str) – The file to read from
- Return type
TypeVar
(T
, bound= Spectrum)
- get_bin_centres()
Get x-axis bin centres. If the size of x_data is the same size as y_data, x_data is assumed to contain bin centres, but if x_data is one element larger, x_data is assumed to contain bin edges and a conversion is made. In this conversion, the bin centres are assumed to be in the middle of each bin, which may not be an accurate assumption in the case of differently sized bins.
- Return type
Quantity
- get_bin_edges()
Get x-axis bin edges. If the size of x_data is one element larger than y_data, x_data is assumed to contain bin edges, but if x_data is the same size, x_data is assumed to contain bin centres and a conversion is made. In the conversion, the bin edges will not go outside the existing data bounds so the first and last bins may be half-size. In addition, each bin edge is assumed to be halfway between each bin centre, which may not be an accurate assumption in the case of differently sized bins.
- Return type
Quantity
- get_bin_widths()
Get x-axis bin widths
- Return type
Quantity
- split(indices=None, btol=None)
Split to multiple spectra
Data may be split by index. Alternatively, x-axis data may be split automatically by searching for unusually large gaps in energy values. These would usually correspond to disconnected special-points in a band-structure diagram.
- Parameters
indices (
Union
[Sequence
[int
],ndarray
,None
]) – positions in data of breakpointsbtol (
Optional
[float
]) – parameter used to identify breakpoints. This is a ratio between the gap in values and the median gap between neighbouring x-values. If neither indices nor btol is specified, this is set to 10.0.
- Return type
List
[TypeVar
(T
, bound= Spectrum)]- Returns
split_spectra – Separated spectrum regions. If passed to the appropriate functions in euphonic.plot this would be interpreted as a series of subplots.
- to_json_file(filename)
Write to a JSON file. JSON fields are equivalent to from_dict keys
- Parameters
filename (str) – Name of the JSON file to write to
- Return type
None
Spectrum1DCollection
This is an object for storing multiple 1D spectra which share the same x-axis, e.g. bands in a dispersion plot.
From Spectrum1D
If you have multiple Spectrum1D
objects with the same x_data
,
they can be grouped together into a Spectrum1DCollection
object.
For example:
from euphonic import Spectrum1D, Spectrum1DCollection
pdos_si = Spectrum1D.from_json_file('si_pdos.json')
pdos_o = Spectrum1D.from_json_file('o_pdos.json')
pdos_collection = Spectrum1DCollection.from_spectra([pdos_si, pdos_o])
Adding Spectrum1DCollection
Two Spectrum1DCollection
objects can be added together (provided they
have the same x_data
axes) with the +
operator. This will concatenate
their y_data
, returning a single Spectrum1DCollection
that contains
all the y_data
of both objects. For example:
from euphonic import Spectrum1DCollection
coh_pdos = Spectrum1DCollection.from_json_file('coherent_pdos.json')
incoh_pdos = Spectrum1DCollection.from_json_file('incoherent_pdos.json')
all_pdos = coh_pdos + incoh_pdos
Indexing
A Spectrum1DCollection
can be indexed just like a list to obtain
a specific spectrum as a Spectrum1D
, or a subset of spectra as
another Spectrum1DCollection
, for example, to plot only specific
spectra:
from euphonic import Spectrum1DCollection
from euphonic.plot import plot_1d
spec1d_col = Spectrum1DCollection.from_json_file('coherent_pdos.json')
# Plot the 1st spectrum
spec1d_0 = spec1d_col[0]
fig1 = plot_1d(spec1d_0)
# Plot the 2nd - 5th spectra
spec1d_col_1_5 = spec1d_col[1:5]
fig2 = plot_1d(spec1d_col_1_5)
Broadening
A collection of 1D spectra can also be broadened using
Spectrum1DCollection.broaden
,
which broadens each spectrum individually, giving the same result
as using Spectrum1D.broaden
on each contained spectrum.
Grouping By Metadata
You can group and sum specific spectra from a Spectrum1DCollection
based on their metadata using
Spectrum1DCollection.group_by
.
For example, if you have a collection spec1d_col
containing
8 spectra with the following metadata:
>>> spec1d_col.metadata
{'line_data': [{'index': 1, 'species': 'Si', 'weighting': 'coherent'}, {'index': 2, 'species': 'Si', 'weighting': 'coherent'}, {'index': 3, 'species': 'O', 'weighting': 'coherent'}, {'index': 4, 'species': 'O', 'weighting': 'coherent'}, {'index': 1, 'species': 'Si', 'weighting': 'incoherent'}, {'index': 2, 'species': 'Si', 'weighting': 'incoherent'}, {'index': 3, 'species': 'O', 'weighting': 'incoherent'}, {'index': 4, 'species': 'O', 'weighting': 'incoherent'}]}
If you want to group and sum spectra that have the same weighting
,
pass 'weighting'
to the group_by
method. This would produce a
collection containing 2 spectra with the following metadata (the species
and index
metadata are not common across all the grouped spectra, so have
been discarded):
>>> weighting_pdos = spec1d_col.group_by('weighting')
>>> weighting_pdos.metadata
{'line_data': [{'weighting': 'coherent'}, {'weighting': 'incoherent'}]}
You can also group by multiple keys, for example you can group and sum
spectra that have both the same weighting
and species
, producing
the following metadata:
>>> weighting_species_pdos = spec1d_col.group_by('weighting', 'species')
>>> weighting_species_pdos.metadata
{'line_data': [{'weighting': 'coherent', 'species': 'Si'}, {'species': 'O', 'weighting': 'coherent'}, {'weighting': 'incoherent', 'species': 'Si'}, {'weighting': 'incoherent', 'species': 'O'}]}
Selecting By Metadata
You can select specific spectra from a Spectrum1DCollection
based
on their metadata using
Spectrum1DCollection.select
.
For example, if you have a collection spec1d_col
containing
8 spectra with the following metadata:
>>> spec1d_col.metadata
{'line_data': [{'index': 1, 'species': 'Si', 'weighting': 'coherent'}, {'index': 2, 'species': 'Si', 'weighting': 'coherent'}, {'index': 3, 'species': 'O', 'weighting': 'coherent'}, {'index': 4, 'species': 'O', 'weighting': 'coherent'}, {'index': 1, 'species': 'Si', 'weighting': 'incoherent'}, {'index': 2, 'species': 'Si', 'weighting': 'incoherent'}, {'index': 3, 'species': 'O', 'weighting': 'incoherent'}, {'index': 4, 'species': 'O', 'weighting': 'incoherent'}]}
If you want to select only the spectra where weighting='coherent'
,
use the select
method, which would create a collection containing
4 spectra, with the following metadata:
>>> coh_pdos = spec1d_col.select(weighting='coherent')
>>> coh_pdos.metadata
{'weighting': 'coherent', 'line_data': [{'index': 1, 'species': 'Si'}, {'index': 2, 'species': 'Si'}, {'index': 3, 'species': 'O'}, {'index': 4, 'species': 'O'}]}
You can also select multiple values for a specific key. For example, to
select spectra where index=1
or index=2
:
>>> coh_or_incoh_pdos = spec1d_col.select(index=[1, 2])
>>> coh_or_incoh_pdos.metadata
{'species': 'Si', 'line_data': [{'index': 1, 'weighting': 'coherent'}, {'index': 1, 'weighting': 'incoherent'}, {'index': 2, 'weighting': 'coherent'}, {'index': 2, 'weighting': 'incoherent'}]}
You can also select by multiple key/values. To select only the spectra with
weighting='coherent'
and species='Si'
:
>>> coh_si_pdos = spec1d_col.select(weighting='coherent', species='Si')
>>> coh_si_pdos.metadata
{'weighting': 'coherent', 'species': 'Si', 'line_data': [{'index': 1}, {'index': 2}]}
Summing Spectra
All spectra in a Spectrum1DCollection
can be summed with
Spectrum1DCollection.sum
.
This produces a single Spectrum1D
object.
Plotting
See Plotting
Docstring
- class Spectrum1DCollection(x_data, y_data, x_tick_labels=None, metadata=None)
A collection of Spectrum1D with common x_data and x_tick_labels
Intended for convenient storage of band structures, projected DOS etc. This object can be indexed or iterated to obtain individual Spectrum1D.
- Variables
x_data – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))
y_data – Shape (n_entries, n_x_data) float Quantity. The plot data in y, in rows corresponding to separate 1D spectra
x_tick_labels – Sequence[Tuple[int, str]] or None. Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied to
metadata –
Dict[str, Union[int, str, LineData]] or None. Contains metadata about the spectra. Keys should be strings and values should be strings or integers. There are some functional keys:
- ’line_data’LineData
This is a Sequence[Dict[str, Union[int, str]], it contains metadata for each spectrum in the collection, and must be of length n_entries
- __init__(x_data, y_data, x_tick_labels=None, metadata=None)
- Parameters
x_data (
Quantity
) – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))y_data (
Quantity
) – Shape (n_entries, n_x_data) float Quantity. The plot data in y, in rows corresponding to separate 1D spectrax_tick_labels (
Optional
[Sequence
[Tuple
[int
,str
]]]) – Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied tometadata (
Optional
[Dict
[str
,Union
[str
,int
,Sequence
[Dict
[str
,Union
[int
,str
]]]]]]) –Contains metadata about the spectra. Keys should be strings and values should be strings or integers. There are some functional keys:
- ’line_data’LineData
This is a Sequence[Dict[str, Union[int, str]], it contains metadata for each spectrum in the collection, and must be of length n_entries
- copy()
Get an independent copy of spectrum
- Return type
TypeVar
(T
, bound= Spectrum1DCollection)
- to_dict()
Convert to a dictionary consistent with from_dict()
- Return type
Dict
[str
,Any
]- Returns
dict
- to_text_file(filename, fmt=None)
Write to a text file. The header contains metadata and unit information, the first column is x_data and each subsequent column is a y_data spectrum. Note that text files written in this format cannot be read back in by Euphonic.
- Parameters
filename (
str
) – Name of the text file to write tofmt (
Union
[str
,Sequence
[str
],None
]) – A format specifier or sequence of specifiers (one for each column), to be passed to numpy.savetxt
- Return type
None
- classmethod from_dict(d)
Convert a dictionary to a Spectrum1DCollection object
- Parameters
d (dict) –
A dictionary with the following keys/values:
’x_data’: (n_x_data,) or (n_x_data + 1,) float ndarray
’x_data_unit’: str
’y_data’: (n_x_data,) float ndarray
’y_data_unit’: str
There are also the following optional keys:
’x_tick_labels’: list of (int, string) tuples
’metadata’: dict
- Return type
TypeVar
(T
, bound= Spectrum1DCollection)- Returns
spectrum_collection
- classmethod from_castep_phonon_dos(filename)
Reads total DOS and per-element PDOS from a CASTEP .phonon_dos file
- Parameters
filename (
str
) – The path and name of the .phonon_dos file to read- Return type
TypeVar
(T
, bound= Spectrum1DCollection)
- broaden(x_width: Quantity, shape: Literal['gauss', 'lorentz'] = 'gauss', method: Optional[Literal['convolve']] = None) T
- broaden(x_width: Callable[[Quantity], Quantity], shape: Literal['gauss', 'lorentz'] = 'gauss', method: Optional[Literal['convolve']] = None, width_lower_limit: Optional[Quantity] = None, width_convention: Literal['fwhm', 'std'] = 'fwhm', width_interpolation_error: float = 0.01, width_fit: Literal['cheby-log', 'cubic'] = 'cheby-log') T
Individually broaden each line in y_data, returning a new Spectrum1DCollection
- Parameters
x_width – Scalar float Quantity, giving broadening FWHM for all y data, or a function handle accepting Quantity consistent with x-axis as input and returning FWHM corresponding to input array. This would typically be an energy-dependent resolution function.
shape – One of {‘gauss’, ‘lorentz’}. The broadening shape
method – Can be None or ‘convolve’. Currently the only broadening method available is convolution with a broadening kernel but this may not produce correct results for unequal bin widths. To use convolution anyway, explicitly set method=’convolve’
width_lower_limit – Set a lower bound to width obtained calling x_width function. By default, this is equal to bin width. To disable, set to -Inf.
width_convention – By default (‘fwhm’), x_width is interpreted as full-width half-maximum. Set to ‘std’ to instead define standard deviation.
width_interpolation_error – When x_width is a callable function, variable-width broadening is implemented by an approximate kernel-interpolation scheme. This parameter determines the target error of the kernel approximations.
width_fit – Select parametrisation of kernel width spacing to width_interpolation_error. ‘cheby-log’ is recommended: for shape ‘gauss’, ‘cubic’ is also available.
- Return type
TypeVar
(T
, bound= Spectrum1DCollection)- Returns
broadened_spectrum – A new Spectrum1DCollection object with broadened y_data
- Raises
ValueError – If shape is not one of the allowed strings
ValueError – If method is None and bins are not of equal size
- group_by(*line_data_keys)
Group and sum y_data for each spectrum according to the values mapped to the specified keys in metadata[‘line_data’]
- Parameters
line_data_keys (
str
) – The key(s) to group by. If only one line_data_key is supplied, if the value mapped to a key is the same for multiple spectra, they are placed in the same group and summed. If multiple line_data_keys are supplied, the values must be the same for all specified keys for them to be placed in the same group- Return type
TypeVar
(T
, bound= Spectrum1DCollection)- Returns
grouped_spectrum – A new Spectrum1DCollection with one line for each group. Any metadata in ‘line_data’ not common across all spectra in a group will be discarded
- sum()
Sum y_data over all spectra
- Return type
- Returns
summed_spectrum – A Spectrum1D created from the summed y_data. Any metadata in ‘line_data’ not common across all spectra will be discarded
- select(**select_key_values)
Select spectra by their keys and values in metadata[‘line_data’]
- Parameters
**select_key_values (
Union
[str
,int
,Sequence
[str
],Sequence
[int
]]) – Key-value/values pairs in metadata[‘line_data’] describing which spectra to extract. For example, to select all spectra where metadata[‘line_data’][‘species’] = ‘Na’ or ‘Cl’ use spectrum.select(species=[‘Na’, ‘Cl’]). To select ‘Na’ and ‘Cl’ spectra where weighting is also coherent, use spectrum.select(species=[‘Na’, ‘Cl’], weighting=’coherent’)- Return type
TypeVar
(T
, bound= Spectrum1DCollection)- Returns
selected_spectra – A Spectrum1DCollection containing the selected spectra
- Raises
ValueError – If no matching spectra are found
- assert_regular_bins(message='', rtol=1e-05, atol=0.0)
Raise AssertionError if x-axis bins are not evenly spaced.
- Parameters
message (
str
) – Text appended to ValueError for more informative output.rtol (
float
) – Relative tolerance for ‘close enough’ valuesatol (
float
) – Absolute tolerance for ‘close enough’ values. Note this is a bare float and follows the stored units of the bins.
- Return type
None
- count(value) integer -- return number of occurrences of value
- classmethod from_json_file(filename)
Read from a JSON file. See from_dict for required fields
- Parameters
filename (str) – The file to read from
- Return type
TypeVar
(T
, bound= Spectrum)
- get_bin_centres()
Get x-axis bin centres. If the size of x_data is the same size as y_data, x_data is assumed to contain bin centres, but if x_data is one element larger, x_data is assumed to contain bin edges and a conversion is made. In this conversion, the bin centres are assumed to be in the middle of each bin, which may not be an accurate assumption in the case of differently sized bins.
- Return type
Quantity
- get_bin_edges()
Get x-axis bin edges. If the size of x_data is one element larger than y_data, x_data is assumed to contain bin edges, but if x_data is the same size, x_data is assumed to contain bin centres and a conversion is made. In the conversion, the bin edges will not go outside the existing data bounds so the first and last bins may be half-size. In addition, each bin edge is assumed to be halfway between each bin centre, which may not be an accurate assumption in the case of differently sized bins.
- Return type
Quantity
- get_bin_widths()
Get x-axis bin widths
- Return type
Quantity
- index(value[, start[, stop]]) integer -- return first index of value.
Raises ValueError if the value is not present.
Supporting start and stop arguments is optional, but recommended.
- split(indices=None, btol=None)
Split to multiple spectra
Data may be split by index. Alternatively, x-axis data may be split automatically by searching for unusually large gaps in energy values. These would usually correspond to disconnected special-points in a band-structure diagram.
- Parameters
indices (
Union
[Sequence
[int
],ndarray
,None
]) – positions in data of breakpointsbtol (
Optional
[float
]) – parameter used to identify breakpoints. This is a ratio between the gap in values and the median gap between neighbouring x-values. If neither indices nor btol is specified, this is set to 10.0.
- Return type
List
[TypeVar
(T
, bound= Spectrum)]- Returns
split_spectra – Separated spectrum regions. If passed to the appropriate functions in euphonic.plot this would be interpreted as a series of subplots.
- to_json_file(filename)
Write to a JSON file. JSON fields are equivalent to from_dict keys
- Parameters
filename (str) – Name of the JSON file to write to
- Return type
None
Spectrum2D
Broadening
A 2D spectrum can be broadened using
Spectrum2D.broaden
, which
broadens along either or both of the x/y-axes and returns a new
Spectrum2D object. It can broaden with either a Gaussian or Lorentzian
and requires a broadening FWHM in the same type of units as
x_data
/y_data
for broadening along the x/y-axis respectively.
For example:
from euphonic import ureg, Spectrum2D
sqw = Spectrum2D.from_json_file('sqw.json')
x_fwhm = 0.05*ureg('1/angstrom')
y_fwhm = 1.5*ureg('meV')
sqw_broaden = sqw.broaden(x_width=x_fwhm, y_width=y_fwhm, shape='lorentz')
Plotting
See Plotting
Kinematic constraints
Inelastic neutron-scattering (INS) experiments are often performed on time-of-flight (ToF) spectrometers where a wide ToF range yields a wide range of energy transfers which are measured simultaneously. In “direct geometry” the incident energy is fixed (e.g. by a Fermi chopper) while in “indirect geometry” the scattered energy is fixed (e.g. by scattering from an “analyser” crystal). Conservation laws allow the overall energy and momentum transfer to be determined for a given scattering angle and crystal orientation. In powder measurements, the crystal orientation is not needed so the kinematic limits — the accessible \((q, \omega)\) range — determined by the conservation laws are given solely by these instrument parameters.
The function euphonic.spectra.apply_kinematic_constraints
applies these limits
to a powder-averaged Spectrum2D object with appropriate dimensions
(i.e. the x- and y-axes represent \(|q|\) and \(\omega\) respectively).
Inaccessible data values are set to NaN
; in Matplotlib colour maps
this will leave them unset.
Facility |
Instrument |
\(E_i\) / meV |
\(E_f\) / meV |
\(2\theta\) / \({}^\circ\) |
---|---|---|---|---|
ILL |
LAGRANGE |
4.5 |
10–90 |
|
ISIS |
LET |
1–25 |
5–140 |
|
ISIS |
MAPS |
15–2000 |
3–60 |
|
ISIS |
MARI |
7–1000 |
3–135 |
|
ISIS |
MERLIN |
7–2000 |
3–135 |
|
ILL |
PANTHER |
76, 112, 150 |
5–136 |
|
ISIS |
TOSCA |
3.97 |
45, 135 |
- apply_kinematic_constraints(spectrum, e_i=None, e_f=None, angle_range=(0, 180.0))
Set events to NaN which violate energy/momentum limits:
Energy transfer greater than e_i
q outside region accessible for given e_i and angle range
This requires x_data to be in wavevector units and y_data to be energy.
Either e_i or e_f should be set, according to direct/instrument geometry. The other values will be inferred, interpreting y_data as energy transfer.
- Parameters
spectrum (
Spectrum2D
) – input 2-D spectrum, with |q| x_data and energy y_datae_i (
Optional
[Quantity
]) – incident energy of direct-geometry spectrometere_f (
Optional
[Quantity
]) – final energy of indirect-geometry spectrometerangle_range (
Tuple
[float
]) – min and max scattering angles (2θ) of detector bank in degrees
- Return type
- Returns
Masked spectrum with inaccessible bins set to NaN in z_data.
Docstring
- class Spectrum2D(x_data, y_data, z_data, x_tick_labels=None, metadata=None)
For storing generic 2D spectra e.g. S(Q,w)
- Variables
x_data – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))
y_data – Shape (n_y_data,) or (n_y_data + 1,) float Quantity. The y_data bin points (if size == (n_y_data,)) or y_data bin edges (if size == (n_y_data + 1,))
z_data – Shape (n_x_data, n_y_data) float Quantity. The plot data in z
x_tick_labels – Sequence[Tuple[int, str]] or None. Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied to
metadata – Dict[str, Union[int, str]]. Contains metadata about the spectrum. Keys should be strings and values should be strings or integers
- __init__(x_data, y_data, z_data, x_tick_labels=None, metadata=None)
- Parameters
x_data (
Quantity
) – Shape (n_x_data,) or (n_x_data + 1,) float Quantity. The x_data points (if size == (n_x_data,)) or x_data bin edges (if size == (n_x_data + 1,))y_data (
Quantity
) – Shape (n_y_data,) or (n_y_data + 1,) float Quantity. The y_data bin points (if size == (n_y_data,)) or y_data bin edges (if size == (n_y_data + 1,))z_data (
Quantity
) – Shape (n_x_data, n_y_data) float Quantity. The plot data in zx_tick_labels (
Optional
[Sequence
[Tuple
[int
,str
]]]) – Special tick labels e.g. for high-symmetry points. The int refers to the index in x_data the label should be applied tometadata (
Optional
[Dict
[str
,Union
[int
,str
]]]) – Contains metadata about the spectrum. Keys should be strings and values should be strings or integers.
- broaden(x_width=None, y_width=None, shape='gauss', method=None, x_width_lower_limit=None, y_width_lower_limit=None, width_convention='fwhm', width_interpolation_error=0.01, width_fit='cheby-log')
Broaden z_data and return a new broadened Spectrum2D object
Callable functions can be used to access variable-width broadening. In this case, broadening is implemented with a kernel-interpolating approximate scheme.
- Parameters
x_width (
Union
[Quantity
,Callable
[[Quantity
],Quantity
],None
]) – Either a scalar float Quantity representing the broadening width, or a callable function that accepts and returns Quantity consistent with x-axis units.y_width (
Union
[Quantity
,Callable
[[Quantity
],Quantity
],None
]) – Either a scalar float Quantity representing the broadening width, or a callable function that accepts and returns Quantity consistent with y-axis units.shape (
Literal
[‘gauss’, ‘lorentz’]) – One of {‘gauss’, ‘lorentz’}. The broadening shapemethod (
Optional
[Literal
[‘convolve’]]) – Can be None or ‘convolve’. Currently the only broadening method available is convolution with a broadening kernel, but this may not produce correct results for unequal bin widths. To use convolution anyway, explicitly set method=’convolve’x_width_lower_limit (
Optional
[Quantity
]) – Set a lower bound to width obtained calling x_width function. By default, this is equal to x bin width. To disable, set to -Inf.y_width_lower_limit (
Optional
[Quantity
]) – Set a lower bound to width obtained calling y_width function. By default, this is equal to y bin width. To disable, set to -Inf.width_convention (
Literal
[‘fwhm’, ‘std’]) – By default (‘fwhm’), widths are interpreted as full-width half-maximum. Set to ‘std’ to instead define standard deviation.width_interpolation_error (
float
) – When x_width is a callable function, variable-width broadening is implemented by an approximate kernel-interpolation scheme. This parameter determines the target error of the kernel approximations.width_fit (
Literal
[‘cheby-log’, ‘cubic’]) – Select parametrisation of kernel width spacing to width_interpolation_error. ‘cheby-log’ is recommended: for shape ‘gauss’, ‘cubic’ is also available.
- Return type
TypeVar
(T
, bound= Spectrum2D)- Returns
broadened_spectrum – A new Spectrum2D object with broadened z_data
- Raises
ValueError – If shape is not one of the allowed strings
ValueError – If method is None and bins are not of equal size
- copy()
Get an independent copy of spectrum
- Return type
TypeVar
(T
, bound= Spectrum2D)
- classmethod from_json_file(filename)
Read from a JSON file. See from_dict for required fields
- Parameters
filename (str) – The file to read from
- Return type
TypeVar
(T
, bound= Spectrum)
- get_bin_edges(bin_ax='x')
Get bin edges for the axis specified by bin_ax. If the size of bin_ax is one element larger than z_data along that axis, bin_ax is assumed to contain bin edges, but if bin_ax is the same size, bin_ax is assumed to contain bin centres and a conversion is made. In the conversion, the bin edges will not go outside the existing data bounds so the first and last bins may be half-size. In addition, each bin edge is assumed to be halfway between each bin centre, which may not be an accurate assumption in the case of differently sized bins.
- Parameters
bin_ax (
Literal
[‘x’, ‘y’]) – The axis to get the bin edges for, ‘x’ or ‘y’- Return type
Quantity
- split(indices=None, btol=None)
Split to multiple spectra
Data may be split by index. Alternatively, x-axis data may be split automatically by searching for unusually large gaps in energy values. These would usually correspond to disconnected special-points in a band-structure diagram.
- Parameters
indices (
Union
[Sequence
[int
],ndarray
,None
]) – positions in data of breakpointsbtol (
Optional
[float
]) – parameter used to identify breakpoints. This is a ratio between the gap in values and the median gap between neighbouring x-values. If neither indices nor btol is specified, this is set to 10.0.
- Return type
List
[TypeVar
(T
, bound= Spectrum)]- Returns
split_spectra – Separated spectrum regions. If passed to the appropriate functions in euphonic.plot this would be interpreted as a series of subplots.
- to_json_file(filename)
Write to a JSON file. JSON fields are equivalent to from_dict keys
- Parameters
filename (str) – Name of the JSON file to write to
- Return type
None
- get_bin_centres(bin_ax='x')
Get bin centres for the axis specified by bin_ax. If the size of bin_ax is the same size as z_data along that axis, bin_ax is assumed to contain bin centres, but if bin_ax is one element larger, bin_ax is assumed to contain bin centres and a conversion is made. In this conversion, the bin centres are assumed to be in the middle of each bin, which may not be an accurate assumption in the case of differently sized bins.
- Parameters
bin_ax (
Literal
[‘x’, ‘y’]) – The axis to get the bin centres for, ‘x’ or ‘y’- Return type
Quantity
- get_bin_widths(bin_ax='x')
Get x-bin widths along specified axis
- Parameters
bin_ax (
Literal
[‘x’, ‘y’]) – Axis of interest, ‘x’ or ‘y’- Return type
Quantity
- assert_regular_bins(bin_ax, message='', rtol=1e-05, atol=0.0)
Raise AssertionError if x-axis bins are not evenly spaced.
- Parameters
bin_ax (
Literal
[‘x’, ‘y’]) – Axis of interest, ‘x’ or ‘y’message (
str
) – Text appended to ValueError for more informative output.rtol (
float
) – Relative tolerance for ‘close enough’ valuesatol (
float
) – Absolute tolerance for ‘close enough’ values. Note this is a bare float and follows the stored units of the bins.
- Return type
None
- to_dict()
Convert to a dictionary. See Spectrum2D.from_dict for details on keys/values
- Return type
Dict
[str
,Any
]- Returns
dict
- classmethod from_dict(d)
Convert a dictionary to a Spectrum2D object
- Parameters
d (dict) –
A dictionary with the following keys/values:
’x_data’: (n_x_data,) or (n_x_data + 1,) float ndarray
’x_data_unit’: str
’y_data’: (n_y_data,) or (n_y_data + 1,) float ndarray
’y_data_unit’: str
’z_data’: (n_x_data, n_y_data) float Quantity
’z_data_unit’: str
There are also the following optional keys:
’x_tick_labels’: list of (int, string) tuples
’metadata’: dict
- Return type
TypeVar
(T
, bound= Spectrum2D)- Returns
spectrum