hyperspy.models.eelsmodel module

class hyperspy.models.eelsmodel.EELSModel(signal1D, auto_background=True, auto_add_edges=True, ll=None, GOS=None, dictionary=None)

Bases: hyperspy.models.model1d.Model1D

Build an EELS model

Parameters:
  • spectrum (a Signal1D (or any Signal1D subclass) instance) –
  • auto_background (boolean) – If True, and if spectrum is an EELS instance adds automatically a powerlaw to the model and estimate the parameters by the two-area method.
  • auto_add_edges (boolean) – If True, and if spectrum is an EELS instance, it will automatically add the ionization edges as defined in the Signal1D instance. Adding a new element to the spectrum using the components.EELSSpectrum.add_elements method automatically add the corresponding ionisation edges to the model.
  • ll ({None, EELSSpectrum}) – If an EELSSPectrum is provided, it will be assumed that it is a low-loss EELS spectrum, and it will be used to simulate the effect of multiple scattering by convolving it with the EELS spectrum.
  • GOS ({'hydrogenic', 'Hartree-Slater', None}) – The GOS to use when auto adding core-loss EELS edges. If None it will use the Hartree-Slater GOS if they are available, otherwise it will use the hydrogenic GOS.
  • dictionary ({dict, None}) – A dictionary to be used to recreate a model. Usually generated using hyperspy.model.as_dictionary()
append(component)

Add component to Model.

Parameters:thing (Component instance.) –
disable_background()

Disable the background components.

disable_edges(edges_list=None)

Disable the edges listed in edges_list. If edges_list is None (default) all the edges with onset in the spectrum energy region will be disabled.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
disable_fine_structure(edges_list=None)

Disable the fine structure of the edges listed in edges_list. If edges_list is None (default) the fine structure of all the edges with onset in the spectrum energy region will be disabled.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
disable_free_onset_energy(edges_list=None)

Disable the automatic freeing of the onset_energy parameter during a smart fit for the edges listed in edges_list. If edges_list is None (default) the onset_energy of all the edges with onset in the spectrum energy region will not be freed. Note that if their atribute edge.onset_energy.free is True, the parameter will be free during the smart fit.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
enable_background()

Enable the background componets.

enable_edges(edges_list=None)

Enable the edges listed in edges_list. If edges_list is None (default) all the edges with onset in the spectrum energy region will be enabled.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
enable_fine_structure(edges_list=None)

Enable the fine structure of the edges listed in edges_list. If edges_list is None (default) the fine structure of all the edges with onset in the spectrum energy region will be enabled.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
enable_free_onset_energy(edges_list=None)

Enable the automatic freeing of the onset_energy parameter during a smart fit for the edges listed in edges_list. If edges_list is None (default) the onset_energy of all the edges with onset in the spectrum energy region will be freeed.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
fit(fitter=None, method='ls', grad=False, bounded=False, ext_bounding=False, update_plot=False, kind='std', **kwargs)

Fits the model to the experimental data

Parameters:
  • fitter ({None, "leastsq", "odr", "mpfit", "fmin"}) – The optimizer to perform the fitting. If None the fitter defined in the Preferences is used. leastsq is the most stable but it does not support bounding. mpfit supports bounding. fmin is the only one that supports maximum likelihood estimation, but it is less robust than the Levenberg–Marquardt based leastsq and mpfit, and it is better to use it after one of them to refine the estimation.
  • method ({'ls', 'ml'}) – Choose ‘ls’ (default) for least squares and ‘ml’ for maximum-likelihood estimation. The latter only works with fitter = ‘fmin’.
  • grad (bool) – If True, the analytical gradient is used if defined to speed up the estimation.
  • ext_bounding (bool) – If True, enforce bounding by keeping the value of the parameters constant out of the defined bounding area.
  • bounded (bool) – If True performs bounded optimization if the fitter supports it. Currently only mpfit support bounding.
  • update_plot (bool) – If True, the plot is updated during the optimization process. It slows down the optimization but it permits to visualize the optimization evolution.
  • kind ({'std', 'smart'}) – If ‘std’ (default) performs standard fit. If ‘smart’ performs smart_fit
  • **kwargs (key word arguments) – Any extra key word argument will be passed to the chosen fitter

See also

multifit(), smart_fit()

fit_background(start_energy=None, only_current=True, **kwargs)

Fit the background to the first active ionization edge in the energy range.

Parameters:
  • start_energy ({float, None}, optional) – If float, limit the range of energies from the left to the given value. Default None.
  • only_current (bool, optional) – If True, only fit the background at the current coordinates. Default True.
  • **kwargs (extra key word arguments) – All extra key word arguments are passed to fit or multifit.
fix_edges(edges_list=None)

Fixes all the parameters of the edges given in edges_list. If edges_list is None (default) all the edges will be fixed.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
fix_fine_structure(edges_list=None)

Fixes all the parameters of the edges given in edges_list. If edges_list is None (default) all the edges will be fixed.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
free_edges(edges_list=None)

Frees all the parameters of the edges given in edges_list. If edges_list is None (default) all the edges will be freeed.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
free_fine_structure(edges_list=None)

Frees all the parameters of the edges given in edges_list. If edges_list is None (default) all the edges will be freeed.

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
quantify()

Prints the value of the intensity of all the independent active EELS core loss edges defined in the model

remove(component)

Remove component from model.

Examples

>>> s = hs.signals.Signal1D(np.empty(1))
>>> m = s.create_model()
>>> g = hs.model.components1D.Gaussian()
>>> m.append(g)

You could remove g like this

>>> m.remove(g)

Like this:

>>> m.remove("Gaussian")

Or like this:

>>> m.remove(0)
remove_fine_structure_data(edges_list=None)

Remove the fine structure data from the fitting routine as defined in the fine_structure_width parameter of the component.EELSCLEdge

Parameters:edges_list (None or list of EELSCLEdge or list of edge names) – If None, the operation is performed on all the edges in the model. Otherwise, it will be performed only on the listed components.
resolve_fine_structure(preedge_safe_window_width=2, i1=0)

Adjust the fine structure of all edges to avoid overlapping

This function is called automatically everytime the position of an edge changes

Parameters:preedge_safe_window_width (float) – minimum distance between the fine structure of an ionization edge and that of the following one. Default 2 (eV).
resume_auto_fine_structure_width(update=True)

Enable the automatic adjustament of the core-loss edges fine structure width.

Parameters:update (bool, optional) – If True, also execute the automatic adjustment (default).
set_all_edges_intensities_positive()
signal1D
smart_fit(start_energy=None, **kwargs)

Fits everything in a cascade style.

Parameters:
  • start_energy ({float, None}) – If float, limit the range of energies from the left to the given value.
  • **kwargs (key word arguments) – Any extra key word argument will be passed to the fit method. See the fit method documentation for a list of valid arguments.

See also

fit(), multifit()

suspend_auto_fine_structure_width()

Disable the automatic adjustament of the core-loss edges fine structure width.

two_area_background_estimation(E1=None, E2=None, powerlaw=None)

Estimates the parameters of a power law background with the two area method.

Parameters:
  • E1 (float) –
  • E2 (float) –
  • powerlaw (PowerLaw component or None) – If None, it will try to guess the right component from the background components of the model
unset_all_edges_intensities_positive()