hyperspy.models.eelsmodel module
- class hyperspy.models.eelsmodel.EELSModel(signal1D, auto_background=True, auto_add_edges=True, ll=None, GOS=None, dictionary=None)
Bases:
Model1D
Build an EELS model
- Parameters:
spectrum (a Signal1D (or any Signal1D subclass) instance) –
auto_background (bool) – 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 (bool) – 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()
- _add_edges_from_subshells_names(e_shells=None)
Create the Edge instances and configure them appropriately
- Parameters:
e_shells (list of strings) –
- _classify_components()
Classify components between background and ionization edge components.
This method should be called every time that components are added and removed. An ionization edge becomes background when its onset falls to the left of the first non-masked energy channel. The ionization edges are stored in a list in the edges attribute. They are sorted by increasing onset_energy. The background components are stored in _background_components.
- _get_first_ionization_edge_energy(start_energy=None)
Calculate the first ionization edge energy.
- Returns:
iee – The first ionization edge energy or None if no edge is defined in the model.
- Return type:
float or None
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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 attribute 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- enable_background()
Enable the background components.
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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 freed.
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- fit(kind='std', **kwargs)
Fits the model to the experimental data.
Read more in the User Guide.
- Parameters:
kind ({"std", "smart"}, default "std") – If “std”, performs standard fit. If “smart”, performs a smart_fit - for more details see the User Guide.
optimizer (str or None, default None) –
The optimization algorithm used to perform the fitting.
”lm” performs least-squares optimization using the Levenberg-Marquardt algorithm, and supports bounds on parameters.
”trf” performs least-squares optimization using the Trust Region Reflective algorithm, and supports bounds on parameters.
”dogbox” performs least-squares optimization using the dogleg algorithm with rectangular trust regions, and supports bounds on parameters.
”odr” performs the optimization using the orthogonal distance regression (ODR) algorithm. It does not support bounds on parameters. See
scipy.odr
for more details.All of the available methods for
scipy.optimize.minimize()
can be used here. See the User Guide documentation for more details.”Differential Evolution” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.differential_evolution()
for more details on available options.”Dual Annealing” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.dual_annealing()
for more details on available options. Requiresscipy >= 1.2.0
.”SHGO” (simplicial homology global optimization” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.shgo()
for more details on available options. Requiresscipy >= 1.2.0
.
loss_function ({"ls", "ML-poisson", "huber", callable}, default "ls") –
The loss function to use for minimization. Only
"ls"
is available ifoptimizer
is one of["lm", "trf", "dogbox", "odr"]
.”ls” minimizes the least-squares loss function.
”ML-poisson” minimizes the negative log-likelihood for Poisson-distributed data. Also known as Poisson maximum likelihood estimation (MLE).
”huber” minimize the Huber loss function. The delta value of the Huber function is controlled by the
huber_delta
keyword argument (the default value is 1.0).callable supports passing your own minimization function.
grad ({"fd", "analytical", callable, None}, default "fd") –
Whether to use information about the gradient of the loss function as part of the optimization. This parameter has no effect if
optimizer
is a derivative-free or global optimization method.”fd” uses a finite difference scheme (if available) for numerical estimation of the gradient. The scheme can be further controlled with the
fd_scheme
keyword argument.”analytical” uses the analytical gradient (if available) to speed up the optimization, since the gradient does not need to be estimated.
callable should be a function that returns the gradient vector.
None means that no gradient information is used or estimated. Not available if
optimizer
is in["lm", "trf", ``"dogbox"
].
bounded (bool, default False) – If True, performs bounded parameter optimization if supported by
optimizer
.update_plot (bool, default False) – If True, the plot is updated during the optimization process. It slows down the optimization, but it enables visualization of the optimization progress.
print_info (bool, default False) – If True, print information about the fitting results, which are also stored in
model.fit_output
in the form of ascipy.optimize.OptimizeResult
object.return_info (bool, default True) – If True, returns the fitting results in the form of a
scipy.optimize.OptimizeResult
object.fd_scheme (str {"2-point", "3-point", "cs"}, default "2-point") – If
grad='fd'
, selects the finite difference scheme to use. Seescipy.optimize.minimize()
for details. Ignored ifoptimizer
is"lm"
,"trf"
or"dogbox"
.**kwargs (keyword arguments) – Any extra keyword argument will be passed to the chosen optimizer. For more information, read the docstring of the optimizer of your choice in
scipy.optimize
.
- Return type:
None
See also
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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- fix_fine_structure(edges_list=None)
Fixes the fine structure 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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 freed.
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- free_fine_structure(edges_list=None)
Frees the fine structure of the edges given in edges_list. If edges_list is None (default) all the edges will be freed.
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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.
See also
enable_edges
,disable_edges
,enable_background
,disable_background
,enable_fine_structure
,disable_fine_structure
,set_all_edges_intensities_positive
,unset_all_edges_intensities_positive
,enable_free_onset_energy
,disable_free_onset_energy
,fix_edges
,free_edges
,fix_fine_structure
,free_fine_structure
- 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 every time 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 adjustment of the core-loss edges fine structure width.
- Parameters:
update (bool, optional) – If True, also execute the automatic adjustment (default).
See also
- set_all_edges_intensities_positive()
Set all edges intensities positive by setting
ext_force_positive
andext_bounded
toTrue
.- Return type:
None.
- smart_fit(start_energy=None, **kwargs)
Fits EELS edges in a cascade style.
The fitting procedure acts in iterative manner along the energy-loss-axis. First it fits only the background up to the first edge. It continues by deactivating all edges except the first one, then performs the fit. Then it only activates the first two, fits, and repeats this until all edges are fitted simultaneously.
Other, non-EELSCLEdge components, are never deactivated, and fitted on every iteration.
- Parameters:
start_energy ({float, None}) – If float, limit the range of energies from the left to the given value.
optimizer (str or None, default None) –
The optimization algorithm used to perform the fitting.
”lm” performs least-squares optimization using the Levenberg-Marquardt algorithm, and supports bounds on parameters.
”trf” performs least-squares optimization using the Trust Region Reflective algorithm, and supports bounds on parameters.
”dogbox” performs least-squares optimization using the dogleg algorithm with rectangular trust regions, and supports bounds on parameters.
”odr” performs the optimization using the orthogonal distance regression (ODR) algorithm. It does not support bounds on parameters. See
scipy.odr
for more details.All of the available methods for
scipy.optimize.minimize()
can be used here. See the User Guide documentation for more details.”Differential Evolution” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.differential_evolution()
for more details on available options.”Dual Annealing” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.dual_annealing()
for more details on available options. Requiresscipy >= 1.2.0
.”SHGO” (simplicial homology global optimization” is a global optimization method. It does support bounds on parameters. See
scipy.optimize.shgo()
for more details on available options. Requiresscipy >= 1.2.0
.
loss_function ({"ls", "ML-poisson", "huber", callable}, default "ls") –
The loss function to use for minimization. Only
"ls"
is available ifoptimizer
is one of["lm", "trf", "dogbox", "odr"]
.”ls” minimizes the least-squares loss function.
”ML-poisson” minimizes the negative log-likelihood for Poisson-distributed data. Also known as Poisson maximum likelihood estimation (MLE).
”huber” minimize the Huber loss function. The delta value of the Huber function is controlled by the
huber_delta
keyword argument (the default value is 1.0).callable supports passing your own minimization function.
grad ({"fd", "analytical", callable, None}, default "fd") –
Whether to use information about the gradient of the loss function as part of the optimization. This parameter has no effect if
optimizer
is a derivative-free or global optimization method.”fd” uses a finite difference scheme (if available) for numerical estimation of the gradient. The scheme can be further controlled with the
fd_scheme
keyword argument.”analytical” uses the analytical gradient (if available) to speed up the optimization, since the gradient does not need to be estimated.
callable should be a function that returns the gradient vector.
None means that no gradient information is used or estimated. Not available if
optimizer
is in["lm", "trf", ``"dogbox"
].
bounded (bool, default False) – If True, performs bounded parameter optimization if supported by
optimizer
.update_plot (bool, default False) – If True, the plot is updated during the optimization process. It slows down the optimization, but it enables visualization of the optimization progress.
print_info (bool, default False) – If True, print information about the fitting results, which are also stored in
model.fit_output
in the form of ascipy.optimize.OptimizeResult
object.return_info (bool, default True) – If True, returns the fitting results in the form of a
scipy.optimize.OptimizeResult
object.fd_scheme (str {"2-point", "3-point", "cs"}, default "2-point") – If
grad='fd'
, selects the finite difference scheme to use. Seescipy.optimize.minimize()
for details. Ignored ifoptimizer
is"lm"
,"trf"
or"dogbox"
.**kwargs (keyword arguments) – Any extra keyword argument will be passed to the chosen optimizer. For more information, read the docstring of the optimizer of your choice in
scipy.optimize
.
See also
fit()
- suspend_auto_fine_structure_width()
Disable the automatic adjustment of the core-loss edges fine structure width.
See also
- two_area_background_estimation(E1=None, E2=None, powerlaw=None)
Estimates the parameters of a power law background with the two area method.
- unset_all_edges_intensities_positive()
Unset all edges intensities positive by setting
ext_force_positive
andext_bounded
toFalse
.- Return type:
None.