HyperSpy API is changing in version 2.0, see the release notes!

components1D#

Components that can be used to define a 1D model for e.g. curve fitting.

Writing a new template is easy: see the user guide documentation on creating components.

For more details see each component docstring.#

Arctan

Arctan function component.

Bleasdale

Bleasdale function component.

Doniach

Doniach Sunjic lineshape component.

Erf

Error function component.

Exponential

Exponential function component.

Expression

Create a component from a string expression.

Gaussian

Normalized Gaussian function component.

GaussianHF

Normalized gaussian function component, with a fwhm parameter

HeavisideStep

The Heaviside step function.

Logistic

Logistic function (sigmoid or s-shaped curve) component.

Lorentzian

Cauchy-Lorentz distribution (a.k.a. Lorentzian function) component.

Offset

Component to add a constant value in the y-axis.

Polynomial

n-order polynomial component.

PowerLaw

Power law component.

RC

ScalableFixedPattern

Fixed pattern component with interpolation support.

SkewNormal

Skew normal distribution component.

SplitVoigt

Split pseudo-Voigt component.

Voigt

Voigt component.

class hyperspy.api.model.components1D.Arctan(A=1.0, k=1.0, x0=1.0, module=['numpy', 'scipy'], **kwargs)#

Bases: Expression

Arctan function component.

\[f(x) = A \cdot \arctan \left[ k \left( x-x_0 \right) \right]\]

Variable

Parameter

\(A\)

A

\(k\)

k

\(x_0\)

x0

Parameters:
Afloat

Amplitude parameter. \(\lim_{x\to -\infty}f(x)=-A\) and \(\lim_{x\to\infty}f(x)=A\)

kfloat

Slope (steepness of the step). The larger \(k\), the sharper the step.

x0float

Center parameter (position of zero crossing \(f(x_0)=0\)).

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

class hyperspy.api.model.components1D.Bleasdale(a=1.0, b=1.0, c=1.0, module=None, **kwargs)#

Bases: Expression

Bleasdale function component.

Also called the Bleasdale-Nelder function. Originates from the description of the yield-density relationship in crop growth.

\[f(x) = \left(a+b\cdot x\right)^{-1/c}\]
Parameters:
afloat, default=1.0

The value of Parameter a.

bfloat, default=1.0

The value of Parameter b.

cfloat, default=1.0

The value of Parameter c.

**kwargs

Extra keyword arguments are passed to Expression.

Notes

For \((a+b\cdot x)\leq0\), the component will be set to 0.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

grad_a(x)#

Returns d(function)/d(parameter_1)

grad_b(x)#

Returns d(function)/d(parameter_1)

grad_c(x)#

Returns d(function)/d(parameter_1)

class hyperspy.api.model.components1D.Doniach(centre=0.0, A=1.0, sigma=1.0, alpha=0.5, module=['numpy', 'scipy'], **kwargs)#

Bases: Expression

Doniach Sunjic lineshape component.

\[ f(x) = \frac{A \cos[ \frac{{\pi\alpha}}{2}+ (1-\alpha)\tan^{-1}(\frac{x-centre+dx}{\sigma})]} {(\sigma^2 + (x-centre+dx)^2)^{\frac{(1-\alpha)}{2}}} \] \[ dx = \frac{2.354820\sigma}{2 tan[\frac{\pi}{2-\alpha}]} \]

Variable

Parameter

\(A\)

A

\(\sigma\)

sigma

\(\alpha\)

alpha

\(centre\)

centre

Parameters:
Afloat

Height

sigmafloat

Variance parameter of the distribution

alphafloat

Tail or asymmetry parameter

centrefloat

Location of the maximum (peak position).

**kwargs

Extra keyword arguments are passed to the Expression component.

Notes

This is an asymmetric lineshape, originially design for xps but generally useful for fitting peaks with low side tails See Doniach S. and Sunjic M., J. Phys. 4C31, 285 (1970) or http://www.casaxps.com/help_manual/line_shapes.htm for a more detailed description

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the Donach by calculating the median (centre) and the variance parameter (sigma).

Note that an insufficient range will affect the accuracy of this method and that this method doesn’t estimate the asymmetry parameter (alpha).

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Returns True when the parameters estimation is successful

Examples

>>> g = hs.model.components1D.Lorentzian()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
class hyperspy.api.model.components1D.Erf(A=1.0, sigma=1.0, origin=0.0, module=['numpy', 'scipy'], **kwargs)#

Bases: Expression

Error function component.

\[f(x) = \frac{A}{2}~\mathrm{erf}\left[\frac{(x - x_0)}{\sqrt{2} \sigma}\right]\]

Variable

Parameter

\(A\)

A

\(\sigma\)

sigma

\(x_0\)

origin

Parameters:
Afloat

The min/max values of the distribution are -A/2 and A/2.

sigmafloat

Width of the distribution.

originfloat

Position of the zero crossing.

**kwargs

Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

class hyperspy.api.model.components1D.Exponential(A=1.0, tau=1.0, module=None, **kwargs)#

Bases: Expression

Exponential function component.

\[f(x) = A\cdot\exp\left(-\frac{x}{\tau}\right)\]

Variable

Parameter

\(A\)

A

\(\tau\)

tau

Parameters:
A: float

Maximum intensity

tau: float

Scale parameter (time constant)

**kwargs

Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the parameters for the exponential component by splitting the signal window into two regions and using their geometric means

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool
class hyperspy.api.model.components1D.Expression(expression, name, position=None, module='numpy', autodoc=True, add_rotation=False, rotation_center=None, rename_pars={}, compute_gradients=True, linear_parameter_list=None, check_parameter_linearity=True, **kwargs)#

Bases: Component

Create a component from a string expression.

It automatically generates the partial derivatives and the class docstring.

Parameters:
expressionstr

Component function in SymPy text expression format with substitutions separated by ;. See examples and the SymPy documentation for details. In order to vary the components along the signal dimensions, the variables x and y must be included for 1D or 2D components. Also, if module is “numexpr” the functions are limited to those that numexpr support. See its documentation for details.

namestr

Name of the component.

positionstr, optional

The parameter name that defines the position of the component if applicable. It enables interative adjustment of the position of the component in the model. For 2D components, a tuple must be passed with the name of the two parameters e.g. (“x0”, “y0”).

moduleNone or str {"numpy" | "numexpr" | "scipy"}, default “numpy”

Module used to evaluate the function. numexpr is often faster but it supports fewer functions and requires installing numexpr. If None, the “numexpr” will be used if installed.

add_rotationbool, default False

This is only relevant for 2D components. If True it automatically adds rotation_angle parameter.

rotation_centerNone or tuple

If None, the rotation center is the center i.e. (0, 0) if position is not defined, otherwise the center is the coordinates specified by position. Alternatively a tuple with the (x, y) coordinates of the center can be provided.

rename_parsdict

The desired name of a parameter may sometimes coincide with e.g. the name of a scientific function, what prevents using it in the expression. rename_parameters is a dictionary to map the name of the parameter in the expression` to the desired name of the parameter in the Component. For example: {“_gamma”: “gamma”}.

compute_gradientsbool, optional

If True, compute the gradient automatically using sympy. If sympy does not support the calculation of the partial derivatives, for example in case of expression containing a “where” condition, it can be disabled by using compute_gradients=False.

linear_parameter_listlist

A list of the components parameters that are known to be linear parameters.

check_parameter_linearitybool

If True, automatically check if each parameter is linear and set its corresponding attribute accordingly. If False, the default is to set all parameters, except for those who are specified in linear_parameter_list.

**kwargsdict

Keyword arguments can be used to initialise the value of the parameters.

Notes

As of version 1.4, Sympy’s lambdify function, that the Expression components uses internally, does not support the differentiation of some expressions, for example those containing a “where” condition. In such cases, the gradients can be set manually if required.

Examples

The following creates a Gaussian component and set the initial value of the parameters:

>>> hs.model.components1D.Expression(
... expression="height * exp(-(x - x0) ** 2 * 4 * log(2)/ fwhm ** 2)",
... name="Gaussian",
... height=1,
... fwhm=1,
... x0=0,
... position="x0",)
<Gaussian (Expression component)>

Substitutions for long or complicated expressions are separated by semicolumns:

>>> expr = 'A*B/(A+B) ; A = sin(x)+one; B = cos(y) - two; y = tan(x)'
>>> comp = hs.model.components1D.Expression(
...    expression=expr,
...    name='my function'
... )
>>> comp.parameters
(<Parameter one of my function component>,
 <Parameter two of my function component>)
Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

compile_function(module, position=False)#

Compile the function and calculate the gradient automatically when possible. Useful to recompile the function and gradient with a different module.

function_nd(*args)#

Returns a numpy array containing the value of the component for all indices. If enough memory is available, this is useful to quickly to obtain the fitted component without iterating over the navigation axes.

class hyperspy.api.model.components1D.Gaussian(A=1.0, sigma=1.0, centre=0.0, module=None, **kwargs)#

Bases: Expression

Normalized Gaussian function component.

\[f(x) = \frac{A}{\sigma \sqrt{2\pi}}\exp\left[ -\frac{\left(x-x_0\right)^{2}}{2\sigma^{2}}\right]\]

Variable

Parameter

\(A\)

A

\(\sigma\)

sigma

\(x_0\)

centre

Parameters:
Afloat

Area, equals height scaled by \(\sigma\sqrt{(2\pi)}\). GaussianHF implements the Gaussian function with a height parameter corresponding to the peak height.

sigmafloat

Scale parameter of the Gaussian distribution.

centrefloat

Location of the Gaussian maximum (peak position).

**kwargs

Extra keyword arguments are passed to the Expression component.

See also

GaussianHF
Attributes:
fwhmfloat

Convenience attribute to get and set the full width at half maximum.

heightfloat

Convenience attribute to get and set the height.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the Gaussian by calculating the momenta.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Notes

Adapted from https://scipy-cookbook.readthedocs.io/items/FittingData.html

Examples

>>> g = hs.model.components1D.Gaussian()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
class hyperspy.api.model.components1D.GaussianHF(height=1.0, fwhm=1.0, centre=0.0, module=None, **kwargs)#

Bases: Expression

Normalized gaussian function component, with a fwhm parameter instead of the sigma parameter, and a height parameter instead of the area parameter A (scaling difference of \(\sigma \sqrt{\left(2\pi\right)}\)). This makes the parameter vs. peak maximum independent of \(\sigma\), and thereby makes locking of the parameter more viable. As long as there is no binning, the height parameter corresponds directly to the peak maximum, if not, the value is scaled by a linear constant (signal_axis.scale).

\[f(x) = h\cdot\exp{\left[-\frac{4 \log{2} \left(x-c\right)^{2}}{W^{2}}\right]}\]

Variable

Parameter

\(h\)

height

\(W\)

fwhm

\(c\)

centre

Parameters:
height: float

The height of the peak. If there is no binning, this corresponds directly to the maximum, otherwise the maximum divided by signal_axis.scale

fwhm: float

The full width half maximum value, i.e. the width of the gaussian at half the value of gaussian peak (at centre).

centre: float

Location of the gaussian maximum, also the mean position.

**kwargs

Extra keyword arguments are passed to the Expression component.

See also

Gaussian
Attributes:
Afloat

Convenience attribute to get, set the area and defined for compatibility with Gaussian component.

sigmafloat

Convenience attribute to get, set the width and defined for compatibility with Gaussian component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the gaussian by calculating the momenta.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Notes

Adapted from https://scipy-cookbook.readthedocs.io/items/FittingData.html

Examples

>>> g = hs.model.components1D.GaussianHF()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
integral_as_signal()#

Utility function to get gaussian integral as Signal1D

class hyperspy.api.model.components1D.HeavisideStep(A=1.0, n=0.0, module='numpy', compute_gradients=True, **kwargs)#

Bases: Expression

The Heaviside step function.

Based on the corresponding numpy function using the half maximum definition for the central point:

\[\begin{split}f(x) = \begin{cases} 0 & x<n\\ A/2 & x=n\\ A & x>n \end{cases}\end{split}\]

Variable

Parameter

\(n\)

centre

\(A\)

height

Parameters:
nfloat

Location parameter defining the x position of the step.

Afloat

Height parameter for x>n.

**kwargs

Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

class hyperspy.api.model.components1D.Logistic(a=1.0, b=1.0, c=1.0, origin=0.0, module=None, **kwargs)#

Bases: Expression

Logistic function (sigmoid or s-shaped curve) component.

\[f(x) = \frac{a}{1 + b\cdot \mathrm{exp}\left[-c \left((x - x_0\right)\right]}\]

Variable

Parameter

\(A\)

a

\(b\)

b

\(c\)

c

\(x_0\)

origin

Parameters:
afloat

The curve’s maximum y-value, \(\mathrm{lim}_{x\to\infty}\left(y\right) = a\)

bfloat

Additional parameter: b>1 shifts origin to larger values; 0<b<1 shifts origin to smaller values; b<0 introduces an asymptote

cfloat

Logistic growth rate or steepness of the curve

originfloat

Position of the sigmoid’s midpoint

**kwargsdict

Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

class hyperspy.api.model.components1D.Lorentzian(A=1.0, gamma=1.0, centre=0.0, module=None, **kwargs)#

Bases: Expression

Cauchy-Lorentz distribution (a.k.a. Lorentzian function) component.

\[f(x)=\frac{A}{\pi}\left[\frac{\gamma}{\left(x-x_{0}\right)^{2} +\gamma^{2}}\right]\]

Variable

Parameter

\(A\)

A

\(\gamma\)

gamma

\(x_0\)

centre

Parameters:
Afloat

Area parameter, where \(A/(\gamma\pi)\) is the maximum (height) of peak.

gammafloat

Scale parameter corresponding to the half-width-at-half-maximum of the peak, which corresponds to the interquartile spread.

centrefloat

Location of the peak maximum.

**kwargs

Extra keyword arguments are passed to the Expression component.

For convenience the `fwhm` and `height` attributes can be used to get and set
the full-with-half-maximum and height of the distribution, respectively.
Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the Lorentzian by calculating the median (centre) and half the interquartile range (gamma).

Note that an insufficient range will affect the accuracy of this method.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Notes

Adapted from gaussian.py and https://en.wikipedia.org/wiki/Cauchy_distribution

Examples

>>> g = hs.model.components1D.Lorentzian()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
class hyperspy.api.model.components1D.Offset(offset=0.0)#

Bases: Component

Component to add a constant value in the y-axis.

\[f(x) = k\]

Variable

Parameter

\(k\)

offset

Parameters:
offsetfloat

The offset to be fitted

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the parameters by the two area method

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool
function_nd(axis)#

Returns a numpy array containing the value of the component for all indices. If enough memory is available, this is useful to quickly to obtain the fitted component without iterating over the navigation axes.

class hyperspy.api.model.components1D.Polynomial(order=2, module=None, **kwargs)#

Bases: Expression

n-order polynomial component.

Polynomial component consisting of order + 1 parameters. The parameters are named “a” followed by the corresponding order, i.e.

\[f(x) = a_{2} x^{2} + a_{1} x^{1} + a_{0}\]

Zero padding is used for polynomial of order > 10.

Parameters:
orderint

Order of the polynomial, must be different from 0.

**kwargs

Keyword arguments can be used to initialise the value of the parameters, i.e. a2=2, a1=3, a0=1. Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the parameters by the two area method

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool
class hyperspy.api.model.components1D.PowerLaw(A=1000000.0, r=3.0, origin=0.0, left_cutoff=0.0, module=None, compute_gradients=False, **kwargs)#

Bases: Expression

Power law component.

\[f(x) = A\cdot(x-x_0)^{-r}\]

Variable

Parameter

\(A\)

A

\(r\)

r

\(x_0\)

origin

Parameters:
Afloat

Height parameter.

rfloat

Power law coefficient.

originfloat

Location parameter.

**kwargs

Extra keyword arguments are passed to the Expression component.

Attributes:
left_cutofffloat

For x <= left_cutoff, the function returns 0. Default value is 0.0.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False, out=False)#

Estimate the parameters for the power law component by the two area method.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False, estimates the parameters for the full dataset.

outbool

If True, returns the result arrays directly without storing in the parameter maps/values. The returned order is (A, r).

Returns:
bool

Exit status required for the remove_background() function.

class hyperspy.api.model.components1D.RC(Vmax=1.0, V0=0.0, tau=1.0, module=None, **kwargs)#

Bases: Expression

RC function component (based on the time-domain capacitor voltage response of an RC-circuit)

\[f(x) = V_\mathrm{0} + V_\mathrm{max} \left[1 - \mathrm{exp}\left( -\frac{x}{\tau}\right)\right]\]

Variable

Parameter

\(V_\mathrm{max}\)

Vmax

\(V_\mathrm{0}\)

V0

\(\tau\)

tau

Parameters:
Vmaxfloat

maximum voltage, asymptote of the function for \(\mathrm{lim}_{x\to\infty}\)

V0float

vertical offset

taufloat

tau=RC is the RC circuit time constant (voltage rise time)

**kwargs

Extra keyword arguments are passed to the Expression component.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

class hyperspy.api.model.components1D.ScalableFixedPattern(signal1D, yscale=1.0, xscale=1.0, shift=0.0, interpolate=True)#

Bases: Component

Fixed pattern component with interpolation support.

\[f(x) = a \cdot s \left(b \cdot x - x_0\right) + c\]

Variable

Parameter

\(a\)

yscale

\(b\)

xscale

\(x_0\)

shift

Parameters:
yscalefloat

The scaling factor in y (intensity axis).

xscalefloat

The scaling factor in x.

shiftfloat

The shift of the component

interpolatebool

If False no interpolation is performed and only a y-scaled spectrum is returned.

Examples

The fixed pattern is defined by a Signal1D of navigation 0 which must be provided to the ScalableFixedPattern constructor, e.g.:

>>> s = hs.load('data.hspy') 
>>> my_fixed_pattern = hs.model.components1D.ScalableFixedPattern(s) 
Attributes:
yscaleParameter

The scaling factor in y (intensity axis).

xscaleParameter

The scaling factor in x.

shiftParameter

The shift of the component

interpolatebool

If False no interpolation is performed and only a y-scaled spectrum is returned.

Methods

prepare_interpolator(**kwargs)

Fine-tune the interpolation.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

function_nd(axis)#

Returns a numpy array containing the value of the component for all indices. If enough memory is available, this is useful to quickly to obtain the fitted component without iterating over the navigation axes.

gui(display=True, toolkit=None, **kwargs)#

Display or return interactive GUI element if available.

Parameters:
displaybool

If True, display the user interface widgets. If False, return the widgets container in a dictionary, usually for customisation or testing.

toolkitstr, iterable of str or None

If None (default), all available widgets are displayed or returned. If string, only the widgets of the selected toolkit are displayed if available. If an interable of toolkit strings, the widgets of all listed toolkits are displayed or returned.

prepare_interpolator(**kwargs)#

Fine-tune the interpolation.

Parameters:
xarray

The spectral axis of the fixed pattern

**kwargsdict

Keywords argument are passed to scipy.interpolate.make_interp_spline()

class hyperspy.api.model.components1D.SkewNormal(x0=0.0, A=1.0, scale=1.0, shape=0.0, module=['numpy', 'scipy'], **kwargs)#

Bases: Expression

Skew normal distribution component.

Asymmetric peak shape based on a normal distribution.

\[\begin{split}f(x) &= 2 A \phi(x) \Phi(x) \\ \phi(x) &= \frac{1}{\sqrt{2\pi}}\mathrm{exp}{\left[ -\frac{t(x)^2}{2}\right]} \\ \Phi(x) &= \frac{1}{2}\left[1 + \mathrm{erf}\left(\frac{ \alpha~t(x)}{\sqrt{2}}\right)\right] \\ t(x) &= \frac{x-x_0}{\omega}\end{split}\]

Variable

Parameter

\(x_0\)

x0

\(A\)

A

\(\omega\)

scale

\(\alpha\)

shape

Parameters:
x0float

Location of the peak position (not maximum, which is given by the mode property).

Afloat

Height parameter of the peak.

scalefloat

Width (sigma) parameter.

shape: float

Skewness (asymmetry) parameter. For shape=0, the normal distribution (Gaussian) is obtained. The distribution is right skewed (longer tail to the right) if shape>0 and is left skewed if shape<0.

**kwargs

Extra keyword arguments are passed to the Expression component.

Notes

The properties mean (position), variance, skewness and mode (position of maximum) are defined for convenience.

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the skew normal distribution by calculating the momenta.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Notes

Adapted from Lin, Lee and Yen, Statistica Sinica 17, 909-927 (2007) https://www.jstor.org/stable/24307705

Examples

>>> g = hs.model.components1D.SkewNormal()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager._axes[-1].offset = -10
>>> s.axes_manager._axes[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
property mean#

Mean (position) of the component.

property mode#

Mode (position of maximum) of the component.

property skewness#

Skewness of the component.

property variance#

Variance of the component.

class hyperspy.api.model.components1D.SplitVoigt(A=1.0, sigma1=1.0, sigma2=1.0, fraction=0.0, centre=0.0)#

Bases: Component

Split pseudo-Voigt component.

\[ pV(x,centre,\sigma) = (1 - \eta) G(x,centre,\sigma) + \eta L(x,centre,\sigma) \] \[ f(x) = \begin{cases} pV(x,centre,\sigma_1), & x \leq centre\\ pV(x,centre,\sigma_2), & x > centre \end{cases} \]

Variable

Parameter

\(A\)

A

\(\eta\)

fraction

\(\sigma_1\)

sigma1

\(\sigma_2\)

sigma2

\(centre\)

centre

Notes

This is a voigt function in which the upstream and downstream variance or sigma is allowed to vary to create an asymmetric profile In this case the voigt is a pseudo voigt consisting of a mixed gaussian and lorentzian sum

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#
Estimate the split voigt function by calculating the

momenta the gaussian.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Notes

Adapted from https://scipy-cookbook.readthedocs.io/items/FittingData.html

Examples

>>> g = hs.model.components1D.SplitVoigt()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True
function(x)#

Split pseudo voigt - a linear combination of gaussian and lorentzian

Parameters:
xarray

independent variable

Afloat

area of pvoigt peak

centerfloat

center position

sigma1float

standard deviation <= center position

sigma2float

standard deviation > center position

fractionfloat

weight for lorentzian peak in the linear combination, and (1-fraction) is the weight for gaussian peak.

function_nd(axis)#

Returns a numpy array containing the value of the component for all indices. If enough memory is available, this is useful to quickly to obtain the fitted component without iterating over the navigation axes.

class hyperspy.api.model.components1D.Voigt(centre=10.0, area=1.0, gamma=0.2, sigma=0.1, module=['numpy', 'scipy'], **kwargs)#

Bases: Expression

Voigt component.

Symmetric peak shape based on the convolution of a Lorentzian and Normal (Gaussian) distribution:

\[f(x) = G(x) \cdot L(x)\]

where \(G(x)\) is the Gaussian function and \(L(x)\) is the Lorentzian function. In this case using an approximate formula by David (see Notes). This approximation improves on the pseudo-Voigt function (linear combination instead of convolution of the distributions) and is, to a very good approximation, equivalent to a Voigt function:

\[\begin{split}z(x) &= \frac{x + i \gamma}{\sqrt{2} \sigma} \\ w(z) &= \frac{e^{-z^2} \text{erfc}(-i z)}{\sqrt{2 \pi} \sigma} \\ f(x) &= A \cdot \Re\left\{ w \left[ z(x - x_0) \right] \right\}\end{split}\]

Variable

Parameter

\(x_0\)

centre

\(A\)

area

\(\gamma\)

gamma

\(\sigma\)

sigma

Parameters:
centrefloat

Location of the maximum of the peak.

areafloat

Intensity below the peak.

gammafloat

\(\gamma\) = HWHM of the Lorentzian distribution.

sigma: float

\(2 \sigma \sqrt{(2 \log(2))}\) = FWHM of the Gaussian distribution.

**kwargs

Extra keyword arguments are passed to the Expression component.

Notes

For convenience the gwidth and lwidth attributes can also be used to set and get the FWHM of the Gaussian and Lorentzian parts of the distribution, respectively. For backwards compatability, FWHM is another alias for the Gaussian width.

W.I.F. David, J. Appl. Cryst. (1986). 19, 63-64, doi:10.1107/S0021889886089999

Parameters:
parameter_name_listlist

The list of parameter names.

linear_parameter_listlist, optional

The list of linear parameter. The default is None.

estimate_parameters(signal, x1, x2, only_current=False)#

Estimate the Voigt function by calculating the momenta of the Gaussian.

Parameters:
signalSignal1D
x1float

Defines the left limit of the spectral range to use for the estimation.

x2float

Defines the right limit of the spectral range to use for the estimation.

only_currentbool

If False estimates the parameters for the full dataset.

Returns:
bool

Exit status required for the remove_background() function.

Notes

Adapted from https://scipy-cookbook.readthedocs.io/items/FittingData.html

Examples

>>> g = hs.model.components1D.Voigt()
>>> x = np.arange(-10, 10, 0.01)
>>> data = np.zeros((32, 32, 2000))
>>> data[:] = g.function(x).reshape((1, 1, 2000))
>>> s = hs.signals.Signal1D(data)
>>> s.axes_manager[-1].offset = -10
>>> s.axes_manager[-1].scale = 0.01
>>> g.estimate_parameters(s, -10, 10, False)
True