Tutorial

The following pages show an example analysis with emgfit broken down into the essential steps. Many of the more advanced features of emgfit are left out or only briefly mentioned in passing, so feel free to explore the documentation further!

This tutorial was created in the Jupyter Notebook emgfit_tutorial.ipynb which can be found in the emgfit/emgfit/examples/tutorial/ directory of the emgfit distribution. Feel free to copy the tutorial folder to a different directory (outside of the emgfit/ directory!) and follow along with the tutorial by actually running the code. You can also use this notebook as a template for your own analyses (consider removing some of the explanations). It is recommended to use a separate notebook for each spectrum to be fitted. This enables you to go back to the notebook at any time and check on all the details of how the fits were performed.

emgfit is optimized to be run within Jupyter Notebooks. There is dozens of decent introductions to using Jupyter Notebooks, a nice overview can e.g. be found at https://realpython.com/jupyter-notebook-introduction/. Naturally, the first step of an analysis with emgfit is starting up your notebook server by running jupyter notebook in your command-line interface. This should automatically make the Jupyter interface pop up in a browser window. From there you can navigate to different directories and create a new notebook (new panel on the top right) or open an existing notebook (.ipynb files). If you installed emgfit into a virtual environment be sure that the correct kernel is selected before executing your notebook cells.

Import the package

Assuming you have setup emgfit following the installation instructions, the first step after launching your Jupyter Notebook will be importing the emgfit package:

[1]:
### Import fit package
import emgfit as emg
print("emgfit version:",emg.__version__)
emgfit version: 0.5.0

How to access the documentation

Before we actually start processing a spectrum it is important to know how to get access to emgfit’s documentation. There is multiple options for this:

  1. The html documentation can be viewed in any browser. It contains usage examples, detailed explanations of the crucial components and API docs with the different modules and all their methods. The search option and cross references enable quick and easy browsing for help.

  2. Once you have imported emgfit you can access the docs directly from the Jupyter Notebook:

    • print all available methods of e.g. the spectrum class by running dir(emg.spectrum)

    • print documentation of a method using help(), e.g. the docs of the add_peak method are printed by running help(emg.spectrum.add_peak) in a code cell

    • keyboard shortcuts can be even more convenient:

      • Use TAB to get suggestions for auto-completion of method and variable names

      • Place the cursor inside the brackets of a function/method and press SHIFT + TAB to have a window with the function/method documention pop up

      • Pressing the H key inside a Jupyter Notebook shows you all available keyboard shortcuts)

Import data

The following code imports the mass data and creates an emgfit spectrum object called spec. The input file must be a TXT or CSV-file with the bin centers and counts per bin as the respective columns (this is in line e.g. with the format of the hist export mode in the mass acquisition software MAc). From here on the analysis of the spectrum proceeds by calling the various methods on our spectrum object spec.

[2]:
### Import mass data, plot full spectrum and indicate chosen fit range
filename = "2019-09-13_004-_006 SUMMED High stats 62Ga"
skiprows = 38 # number of header rows to skip upon data import
m_start = 61.9243 # low-mass cut off
m_stop = 61.962 # high-mass cut off
resolving_power = 300000 # resolving power to use for shape parameter initialization

spec = emg.spectrum(filename+'.txt', m_start, m_stop, skiprows=skiprows,
                    resolving_power=resolving_power)
_images/emgfit_tutorial_5_0.png

Add peaks to the spectrum

This can be done with the automatic peak detection spectrum method) and/or by manually adding peaks (add_peak() spectrum method). The plots shown below are (optional) outputs of the detect_peaks() method and depicts the different stages of the automatic peak detection.

All information about the peaks associated with the spectrum are compiled in the peak properties table. The table’s left-most column shows the respective peak indeces. In all fits, the peaks’ x_pos will be used as the initial values for the peak position parameters mu (to be exact: mu marks the centroid of the underlying Gaussians).

[3]:
### Detect peaks and add them to spectrum object 'spec'
spec.detect_peaks() # automatic peak detection
#spec.add_peak(61.925,species='?') # manually add a peak at x_pos = 61.925u
#spec.remove_peak(peak_index=0) # manually remove the peak with index 0
_images/emgfit_tutorial_7_0.png
_images/emgfit_tutorial_7_1.png
_images/emgfit_tutorial_7_2.png
_images/emgfit_tutorial_7_3.png
Peak properties table after peak detection:
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 ? - None None None None False None None None None None None None None None None None None None
1 61.932021 ? - None None None None False None None None None None None None None None None None None None
2 61.934369 ? - None None None None False None None None None None None None None None None None None None
3 61.943618 ? - None None None None False None None None None None None None None None None None None None
4 61.946994 ? - None None None None False None None None None None None None None None None None None None
5 61.949527 ? - None None None None False None None None None None None None None None None None None None
6 61.956611 ? - None None None None False None None None None None None None None None None None None None
7 61.958997 ? - None None None None False None None None None None None None None None None None None None

Assign species to the peaks (optional)

Although this step is optional, it is highly recommended that it is not skipped. By assigning species labels to your peaks you do not only gain more overview over your spectrum, but also allow for literature values to be automatically fetched from the AME database and entered into the peak properties table. Once a species label has been assigned, you can refer to this peak not only via its index but also via the label.

The assign_species() method allows to assign species identifications either to a single selected peak or to all peaks at once. Here the second option was used by passing a list of species labels to assign_species(). The list must have the same length as the number of peaks associated with the spectrum object. If there are peaks whose labels should not be changed (e.g. unidentified peaks), simply insert None as a placeholder at the corresponding spots (as done for peaks 2 and 7 below). The syntax for species labels follows the :-notation. It is important not to forget to subtract the number of electrons corresponding to the ion’s charge state! Otherwise the analysis would mistakenly proceed with the atomic instead of the ionic mass. Note that currently only singly charged species are supported by emgfit. Tentative peak identifications can be indicated by adding a '?' to the end of the species string. In this case the literature values are not fetched. The user can also define custom literature values (e.g. to handle isomers or if there are recent measurements that have not entered the AME yet). For more details see the documentation of assign_species().

This is also a good point in time to add any comments to the peaks using the add_peak_comment() method. These comments can be particularly helpful for post-processing in Excel since they are also written into the output file with the fit results (as is the entire peak properties table). More general comments that concern the entire spectrum can instead be added with the add_spectrum_comment() method.

[4]:
### Assign species and add peak comments
spec.assign_species(['Ni62:-1e', 'Cu62:-1e?', None, 'Ga62:-1e', 'Ti46:O16:-1e',
                     'Sc46:O16:-1e','Ca43:F19:-1e', None])
spec.add_peak_comment('Non-isobaric', peak_index=2)
spec.show_peak_properties() # check the changes by printing the peak properties table
Species of peak 0 assigned as Ni62:-1e
Species of peak 1 assigned as Cu62:-1e?
Species of peak 3 assigned as Ga62:-1e
Species of peak 4 assigned as Ti46:O16:-1e
Species of peak 5 assigned as Sc46:O16:-1e
Species of peak 6 assigned as Ca43:F19:-1e
Comment of peak 2 was changed to:  Non-isobaric
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False None None None None None None None None None None None None None
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False None None None None None None None None None None None None None
2 61.934369 ? Non-isobaric nan nan nan nan False None None None None None None None None None None None None None
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False None None None None None None None None None None None None None
4 61.946994 Ti46:O16:-1e - 62 1 61.946992 0.000000 False None None None None None None None None None None None None None
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False None None None None None None None None None None None None None
6 61.956611 Ca43:F19:-1e - 62 1 61.956621 0.000000 False None None None None None None None None None None None None None
7 61.958997 ? - nan nan nan nan False None None None None None None None None None None None None None

Activate hiding of mass values for blind analysis (optional)

By adding peak indeces to the spectrum’s blinded_peaks list, the obtained masses and positions of selected peaks-of-interest can be hidden from the user. This blindfolding can avoid user bias and is automatically lifted once the results are exported.

[5]:
### Optionally turn on blinding of specific peaks of interest to enable blind analysis
spec.set_blinded_peaks([0,3]) # activate blinding for peaks 0 & 3
#spec.set_blinded_peaks([],overwrite=True) # run this to deactivate blinding for all peaks
Blinding is activated for peaks: 0, 3

Select the optimal fit model and perform the peak-shape calibration

Next we need to find both a fit model and a set of model parameters that capture the shape of our peaks as well as possible. In emgfit both of this is achieved with the determine_peak_shape() method. Once the peak-shape calibration has been performed all subsequent fits will be performed with this fixed peak-shape, by only varying the peak centroids, amplitudes and optionally the uniform-baseline parameter bkg_c.

By default determine_peak_shape() performs an automatic model selection in which the shape-calibrant peak is first fitted with a pure Gaussian and then with Hyper-EMG functions with an increasing number of expontential tails on the left and right. The algorithm selects the fit model which yields the smallest \(\chi^2_\text{red}\) without having any of the tail weight parameters \(\eta\) compatible with zero within their uncertainty. Alternatively, the auto-model selection can be turned off with the argument vary_tail_order=False and the fit model can be selected manually with the fit_model argument.

Once the best fit model has been selected the method proceeds with the determination of the peak-shape parameters and shows a detailed report with the fit results.

Some recommendations:

  • It is recommended to do the peak-shape calibration with a chi-squared fit (default) since this yields more robust results and more trusworthy parameter uncertainty estimates. Check the method docs for info on performing the shape calibration with binned maximum likelihood estimation.

  • Ideally the peak-shape calibration is performed on a well-separated peak with high statistics. In this example, the Ca43:F19:-1e peak was selected as peak-shape calibrant. Since the default fit range includes a smaller peak on the right, the range was manually reduced to 0.45u with the x_fit_range argument. If unavoidable, the peak-shape determination can also be performed on partially overlapping peaks since emgfit ensures identical shape parameters for all peaks being fitted.

[6]:
## Peak-shape calibration with default settings, including automatic model selection:
#spec.determine_peak_shape(species_shape_calib='Ca43:F19:-1e')

## Peak-shape calibration with user-defined fit range:
spec.determine_peak_shape(species_shape_calib='Ca43:F19:-1e', x_fit_range=0.0045)

## Peak-shape calibration with user-defined fit model:
#spec.determine_peak_shape(species_shape_calib='Ca43:F19:-1e', fit_model='emg12',
#                           vary_tail_order=False)

##### Determine optimal tail order #####--------------------------------------------------------------------------


### Fitting data with Gaussian ###

_images/emgfit_tutorial_13_1.png
_images/emgfit_tutorial_13_2.png

Gaussian-fit yields reduced chi-square of: 28.42 +- 0.13


### Fitting data with emg01 ###

_images/emgfit_tutorial_13_4.png
_images/emgfit_tutorial_13_5.png

emg01-fit yields reduced chi-square of: 10.57 +- 0.13


### Fitting data with emg10 ###

_images/emgfit_tutorial_13_7.png
_images/emgfit_tutorial_13_8.png

emg10-fit yields reduced chi-square of: 28.72 +- 0.13


### Fitting data with emg11 ###

_images/emgfit_tutorial_13_10.png
_images/emgfit_tutorial_13_11.png

emg11-fit yields reduced chi-square of: 2.6 +- 0.13


### Fitting data with emg12 ###

_images/emgfit_tutorial_13_13.png
_images/emgfit_tutorial_13_14.png

emg12-fit yields reduced chi-square of: 1.22 +- 0.13


### Fitting data with emg21 ###

_images/emgfit_tutorial_13_16.png
_images/emgfit_tutorial_13_17.png
WARNING: p6_tau_m1  = 3.0e-06 +- 8.3e-04 is compatible with zero within uncertainty.
WARNING: p6_eta_m2  = 0.100 +- 0.505 is compatible with zero within uncertainty.
             This tail order is likely overfitting the data and will be excluded from selection.

emg21-fit yields reduced chi-square of: 1.47 +- 0.13


### Fitting data with emg22 ###

_images/emgfit_tutorial_13_19.png
_images/emgfit_tutorial_13_20.png

emg22-fit yields reduced chi-square of: 0.96 +- 0.13


### Fitting data with emg23 ###

_images/emgfit_tutorial_13_22.png
_images/emgfit_tutorial_13_23.png
WARNING: p6_eta_m1  = 0.891 +- 2.155 is compatible with zero within uncertainty.
WARNING: p6_eta_m2  = 0.109 +- 2.155 is compatible with zero within uncertainty.
WARNING: p6_eta_p1  = 0.443 +- 6.534 is compatible with zero within uncertainty.
WARNING: p6_tau_p1  = 3.2e-05 +- 3.9e-05 is compatible with zero within uncertainty.
WARNING: p6_eta_p2  = 0.458 +- 4.860 is compatible with zero within uncertainty.
WARNING: p6_eta_p3  = 0.099 +- 1.674 is compatible with zero within uncertainty.
             This tail order is likely overfitting the data and will be excluded from selection.

emg23-fit yields reduced chi-square of: 0.98 +- 0.13


### Fitting data with emg32 ###

_images/emgfit_tutorial_13_25.png
_images/emgfit_tutorial_13_26.png
WARNING: p6_tau_m2  = 1.5e-04 +- 1.7e-04 is compatible with zero within uncertainty.
WARNING: p6_eta_m3  = 0.007 +- 0.107 is compatible with zero within uncertainty.
WARNING: p6_tau_m3  = 3.2e-04 +- 1.6e-03 is compatible with zero within uncertainty.
             This tail order is likely overfitting the data and will be excluded from selection.

emg32-fit yields reduced chi-square of: 0.98 +- 0.13


### Fitting data with emg33 ###

_images/emgfit_tutorial_13_28.png
_images/emgfit_tutorial_13_29.png
WARNING: p6_eta_m1  = 0.872 +- 3.734 is compatible with zero within uncertainty.
WARNING: p6_eta_m2  = 0.120 +- 3.475 is compatible with zero within uncertainty.
WARNING: p6_tau_m2  = 1.6e-04 +- 2.2e-04 is compatible with zero within uncertainty.
WARNING: p6_eta_m3  = 0.009 +- 0.305 is compatible with zero within uncertainty.
WARNING: p6_tau_m3  = 3.5e-04 +- 2.4e-03 is compatible with zero within uncertainty.
WARNING: p6_eta_p1  = 0.444 +- 10.384 is compatible with zero within uncertainty.
WARNING: p6_tau_p1  = 3.1e-05 +- 4.0e-05 is compatible with zero within uncertainty.
WARNING: p6_eta_p2  = 0.459 +- 7.830 is compatible with zero within uncertainty.
WARNING: p6_eta_p3  = 0.096 +- 2.554 is compatible with zero within uncertainty.
             This tail order is likely overfitting the data and will be excluded from selection.

emg33-fit yields reduced chi-square of: 0.99 +- 0.14


##### RESULT OF AUTOMATIC MODEL SELECTION: #####

    Best fit model determined to be:  emg22
    Corresponding chi²-reduced:  0.96 +- 0.13


##### Peak-shape determination #####------------------------------------------------------------------------------
_images/emgfit_tutorial_13_31.png
_images/emgfit_tutorial_13_32.png

Fit Result

Model: (EMGModel(constant, prefix='bkg_') + EMGModel(emg22, prefix='p6_'))

Determine A_stat_emg for subsequent stat. error estimations (optional)

By default, the statistical uncertainties of Hyper-EMG fits are estimated using the equation:

\(\sigma_{stat} = A_{stat,emg} \cdot \frac{\mathrm{FWHM}}{\sqrt{N_{counts}}}\)

where \(\mathrm{FWHM}\) and \(N_{counts}\) refer to the full width at half maximum and the number of counts in the respective peak. This step can be skipped when the statistical uncertainties are estimated using the get_errors_from_resampling() method (see “Perform parametric bootstrap to get refined statistical uncertainties” section below).

By default a of value \(A_{stat,emg} = 0.52\) will be used for Hyper-EMG models (for Gaussians \(A_{stat,G}=0.425\)).

However, \(A_{stat,emg}\) depends on the peak-shape and can deviate from the default value. Therefore, the determine_A_stat_emg() method can be used to estimate \(A_{stat,emg}\) for the specific peak shape in the spectrum. This is done by fitting many simulated spectra created via bootstrap re-sampling from a reference peak in the spectrum. The reference peak should be well-separated and have decent statistics (e.g. the peak-shape calibrant). For details on how \(A_{stat,emg}\) is estimated see the linked docs of determine_A_stat_emg().

This method will typically run for ~10 minutes if N_spectra=1000 (default) is used. For demonstration purposes here the number of bootstrapped spectra generated for each data point (N_spectra argument) was reduced to 10 to get a quicker run time. This is also the reason for the large scatter of the data points below.

In practice it is convenient to skip this method for the first processing of the spectrum since this will only affect the statistical uncertainties but no other fit properties. Once reasonable fits have been achieved for all peaks of interest in the cells below, the exact uncertainties can be obtained by returning to this cell to execute determine_A_stat_emg() with a decent value for N_spectra and then re-runnning the cells below (then with the update value for determine_A_stat_emg()). The latter is conveniently done by using the Run All Below option in the Cell panel of the Jupyter Notebook.

[7]:
# Determine A_stat_emg and save the resulting plot
# In actual practice N_spectra >= 1000 should be used
spec.determine_A_stat_emg(species='Ca43:F19:-1e', x_range=0.004,
                          plot_filename='outputs/'+filename+'_MLE', N_spectra=10)
Creating synthetic spectra via bootstrap re-sampling and fitting them for A_stat determination.
Depending on the choice of `N_spectra` this can take a few minutes. Interrupt kernel if this takes too long.
Done!

[[Model]]
    Model(powerlaw)
[[Fit Statistics]]
    # fitting method   = leastsq
    # function evals   = 5
    # data points      = 8
    # variables        = 1
    chi-square         = 9.6476e-09
    reduced chi-square = 1.3782e-09
    Akaike info crit   = -162.287958
    Bayesian info crit = -162.208517
    R-squared          = -9.02537707
[[Variables]]
    amplitude:  1.3258e-04 +/- 1.3114e-05 (9.89%) (init = 1)
    exponent:  -0.5 (fixed)
_images/emgfit_tutorial_15_3.png
A_stat of a Gaussian model: 0.425
Default A_stat_emg for Hyper-EMG models: 0.52
A_stat_emg for this spectrum's emg22 fit model: 0.554 +- 0.055

Fit all peaks, perform mass re-calibration & obtain final mass values

The following code fits all peaks in the spectrum, performs the mass (re-)calibration, determines the peak-shape uncertainties and updates the peak properties list with the results including the final mass values and their uncertainties.

The simultaneous mass recalibration is optional and only invoked when the species_mass_calib (or the index_mass_calib) argument are specified. If this feature is not used, the fit_peaks() method requires a pre-existing mass calibration (see Alternative 1 section below). In contrast to determine_peak_shape(), by default fit_peaks() performs a binned maximum likelihood fit (‘MLE’). For chi-square fitting with fit_peaks() see Alternative 2 section below. Fits with fit_peaks() can be restricted to a user-defined mass range or to groups of neighbouring peaks selected by index (see the commented-out lines below).

[8]:
# Maximum likelihood fit of all peaks in the spectrum
spec.fit_peaks(species_mass_calib='Ti46:O16:-1e')

# Maximum likelihood fit of peaks in a user-defined mass range
#spec.fit_peaks(species_mass_calib='Ti46:O16:-1e', x_fit_cen=61.9455, x_fit_range=0.01)

# Maximum likelihood fit of peaks specified by index
#spec.fit_peaks(species_mass_calib='Ti46:O16:-1e', peak_indeces=[3,4,5])
##### Fit peaks of interest #####
_images/emgfit_tutorial_17_1.png
_images/emgfit_tutorial_17_2.png

##### Mass recalibration #####

Relative literature error of mass calibrant:     1.6e-09
Relative statistical error of mass calibrant:    1.2e-08

Recalibration factor:    0.999999703 = 1 -3.0e-07
Relative recalibration error:    1.2e-08


##### Peak-shape uncertainty evaluation #####

Determining centroid shifts of mass calibrant.

All centroid shifts below are corrected for the corresponding shifts of the mass calibrant peak.

Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  0 by  -0.23  /   0.20 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  1 by   0.44  /  -0.44 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  2 by  -0.31  /   0.26 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  3 by  -0.04  /  -0.03 μu/z  & its area by   0.7 /   0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  5 by  -0.04  /  -0.01 μu/z  & its area by   0.5 /  -0.6 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  6 by  -0.02  /  -0.03 μu/z  & its area by   0.2 /  -0.5 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  7 by  -0.12  /  -0.02 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  0 by  -0.54  /   0.47 μu/z  & its area by  -0.2 /   0.2 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  1 by  -1.12  /   1.01 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  2 by  -0.42  /   0.29 μu/z  & its area by  -0.8 /   0.5 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  3 by  -0.22  /   0.08 μu/z  & its area by  -0.5 /   1.1 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  5 by  -0.27  /   0.22 μu/z  & its area by  -0.1 /   0.0 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  6 by  -0.14  /   0.04 μu/z  & its area by  -3.2 /   1.9 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  7 by  -1.99  /   1.86 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  0 by   0.68  /  -0.58 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  1 by   0.32  /  -0.29 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  2 by   0.34  /  -0.42 μu/z  & its area by  -2.0 /   1.2 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  3 by   0.26  /  -0.29 μu/z  & its area by  -0.7 /   1.0 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  5 by   0.24  /  -0.22 μu/z  & its area by  -3.9 /   3.0 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  6 by   0.00  /  -0.10 μu/z  & its area by  -4.5 /   2.7 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  7 by   0.56  /  -0.51 μu/z  & its area by  -0.5 /   0.4 counts.

Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  0 by  -0.81  /   0.77 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  1 by   0.38  /  -0.29 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  2 by  -0.51  /   0.41 μu/z  & its area by   0.3 /  -0.3 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  3 by  -0.14  /   0.02 μu/z  & its area by   0.7 /  -0.1 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  5 by  -0.15  /   0.08 μu/z  & its area by   0.9 /  -0.8 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  6 by  -0.03  /  -0.01 μu/z  & its area by  -0.5 /  -0.7 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  7 by  -0.93  /   0.85 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  0 by  -0.20  /   0.35 μu/z  & its area by   0.3 /  -0.4 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  1 by  -0.10  /   0.05 μu/z  & its area by   0.4 /  -0.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  2 by  -0.09  /   0.17 μu/z  & its area by   1.9 /  -2.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  3 by  -0.06  /   0.16 μu/z  & its area by   0.4 /  -0.8 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  5 by  -0.06  /   0.09 μu/z  & its area by   3.2 /  -3.8 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  6 by   0.00  /  -0.02 μu/z  & its area by   4.6 /  -6.6 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  7 by   0.15  /  -0.11 μu/z  & its area by   0.4 /  -0.5 counts.

Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  0 by  -0.10  /   0.04 μu/z  & its area by  -0.3 /   0.2 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  1 by  -0.33  /   0.15 μu/z  & its area by  -0.5 /   0.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  2 by  -0.23  /   0.18 μu/z  & its area by  -1.8 /   1.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  3 by  -0.17  /   0.13 μu/z  & its area by  -1.7 /   1.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  5 by  -0.29  /   0.16 μu/z  & its area by   0.3 /  -0.8 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  6 by   0.00  /  -0.03 μu/z  & its area by  -6.7 /   5.9 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  7 by  -4.05  /   3.08 μu/z  & its area by   1.0 /  -1.0 counts.

Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  0 by   0.09  /  -0.11 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  1 by   0.35  /  -0.59 μu/z  & its area by   0.3 /  -0.4 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  2 by   0.30  /  -0.42 μu/z  & its area by   0.3 /  -0.5 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  3 by   0.13  /  -0.22 μu/z  & its area by   0.9 /  -0.2 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  5 by   0.02  /  -0.03 μu/z  & its area by   1.0 /  -1.2 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  6 by  -0.01  /   0.00 μu/z  & its area by   0.8 /  -1.0 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  7 by   0.46  /  -0.51 μu/z  & its area by   0.3 /  -0.4 counts.

Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  0 by  -0.03  /  -0.00 μu/z  & its area by   0.4 /  -0.5 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  1 by   0.33  /  -0.43 μu/z  & its area by   0.5 /  -0.7 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  2 by   0.04  /  -0.08 μu/z  & its area by   2.7 /  -3.0 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  3 by   0.05  /  -0.10 μu/z  & its area by   2.3 /  -2.6 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  5 by   0.43  /  -0.51 μu/z  & its area by  -7.8 /   5.2 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  6 by  -0.02  /   0.02 μu/z  & its area by  11.3 /  -8.7 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  7 by   4.33  /  -9.19 μu/z  & its area by  -4.4 /   3.5 counts.

Relative peak-shape error of peak  0:  2.1e-08
Relative peak-shape error of peak  1:  2.5e-08
Relative peak-shape error of peak  2:  1.6e-08
Relative peak-shape error of peak  3:  8.4e-09
Relative peak-shape error of peak  5:  1.2e-08
Relative peak-shape error of peak  6:  3.0e-09
Relative peak-shape error of peak  7:  1.7e-07
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False emg22 MLE 1.210000 38.360000 6.816128 blinded 3.46e-07 1.17e-08 2.05e-08 3.46e-07 blinded 19.986000 blinded
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False emg22 MLE 1.210000 27.670000 6.085935 61.932050 4.07e-07 1.17e-08 2.46e-08 4.08e-07 -62784.026000 23.530000 3.518000
2 61.934369 ? Non-isobaric nan nan nan nan False emg22 MLE 1.210000 3870.630000 62.728006 61.934369 3.44e-08 1.17e-08 1.60e-08 3.97e-08 nan 2.290000 nan
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False emg22 MLE 1.210000 953.360000 31.394839 blinded 6.93e-08 1.17e-08 8.38e-09 7.08e-08 blinded 4.085000 blinded
4 61.946994 Ti46:O16:-1e mass calibrant 62 1 61.946992 0.000000 False emg22 MLE 1.210000 33986.590000 185.140000 61.946992 1.16e-08 1.17e-08 nan nan -48865.272000 nan 0.000000
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False emg22 MLE 1.210000 1572.110000 41.266207 61.949539 5.40e-08 1.17e-08 1.15e-08 5.64e-08 -46492.857000 3.256000 5.788000
6 61.956611 Ca43:F19:-1e shape calibrant 62 1 61.956621 0.000000 False emg22 MLE 1.210000 25957.150000 162.039372 61.956621 1.33e-08 1.17e-08 2.99e-09 1.80e-08 -39896.121000 1.036000 0.197000
7 61.958997 ? - nan nan nan nan False emg22 MLE 1.210000 27.790000 7.937959 61.959016 4.06e-07 1.17e-08 1.66e-07 4.39e-07 nan 25.327000 nan
The values for chi-squared as well as the parameter uncertainties and correlations reported by lmfit below should be taken with caution when your MLE fit includes bins with low statistics. For details see Notes section in the spectrum.peakfit() method documentation.

Fit Result

Model: ((((((((EMGModel(constant, prefix='bkg_') + EMGModel(emg22, prefix='p0_')) + EMGModel(emg22, prefix='p1_')) + EMGModel(emg22, prefix='p2_')) + EMGModel(emg22, prefix='p3_')) + EMGModel(emg22, prefix='p4_')) + EMGModel(emg22, prefix='p5_')) + EMGModel(emg22, prefix='p6_')) + EMGModel(emg22, prefix='p7_'))

Plot the fit curve zoomed to a region of interest (optional)

For more detailed inspection of the fit, a zoom to peaks or regions of interest can be shown with the plot_fit_zoom() method.

[9]:
spec.plot_fit_zoom(peak_indeces=[3,4]) # zoom to region around peaks 3 and 4
_images/emgfit_tutorial_19_0.png
_images/emgfit_tutorial_19_1.png

Perform parametric bootstrap to get refined statistical uncertainties (optional)

The A_stat_emg determination with determine_A_stat_emg() relies on fits of bootstrapped subspectra of a single reference peak. The obtained A_stat_emg factor is then used to estimate the statistical uncertainties of all peaks.

As an alternative, the statistical uncertainty can be estimated for each peak individually using the get_errors_from_resampling() method. In this method synthetic spectra are created for all peaks of interest by resampling events from the best-fit curve (“parametric bootstrap”). As opposed to the non-parametric bootstrap of determine_A_stat_emg(), this technique is also applicable to low statistics peaks (assuming that the fit model describes the data well). The fits of the peaks of interest are re-performed using a large number of synthetic spectra (by default: N_spectra=1000) and the statistical mass and area uncertainties are estimated using the standard deviations of the obtained fit results. Finally, the original statistical mass and area uncertainties in the peak properties table are overwritten with the new values.

[10]:
# NOTE: For quicker run time in this demo, the number of synthetic spectra to fit
# (`N_spectra`) was manually reduced. For reliable results, run this method with
# at least the default value of `N_spectra=1000`.

spec.get_errors_from_resampling(N_spectra=20) # arguments adapted for demonstration
#spec.get_errors_from_resampling() # typical execution with default arguments
Fitting 20 simulated spectra to determine statistical mass and peak area errors.
/opt/hostedtoolcache/Python/3.8.18/x64/lib/python3.8/site-packages/emgfit/sample.py:384: UserWarning: 1 simulated events fell outside the specified sampling range and were discarded. Peak areas and area ratios might deviate from expectation.
  warnings.warn(msg, UserWarning)
Updated the statistical and peak area uncertainties of peak(s) 0, 1, 2, 3, 4, 5, 6, 7.

Re-calculated mass recalibration error from updated statistical uncertainty of mass calibrant.
Updated total mass errors of peaks 0, 1, 2, 3, 5, 6, 7.

Updated peak properties table:
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False emg22 MLE 1.210000 38.360000 6.970207 blinded 3.68e-07 9.67e-09 2.05e-08 3.69e-07 blinded 21.283000 blinded
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False emg22 MLE 1.210000 27.670000 5.536694 61.932050 5.91e-07 9.67e-09 2.46e-08 5.92e-07 -62784.026000 34.125000 3.518000
2 61.934369 ? Non-isobaric nan nan nan nan False emg22 MLE 1.210000 3870.630000 61.573439 61.934369 2.74e-08 9.67e-09 1.60e-08 3.31e-08 nan 1.911000 nan
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False emg22 MLE 1.210000 953.360000 39.030752 blinded 5.61e-08 9.67e-09 8.38e-09 5.75e-08 blinded 3.319000 blinded
4 61.946994 Ti46:O16:-1e mass calibrant 62 1 61.946992 0.000000 False emg22 MLE 1.210000 33986.590000 114.513158 61.946992 9.55e-09 9.67e-09 nan nan -48865.272000 nan 0.000000
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False emg22 MLE 1.210000 1572.110000 36.286418 61.949539 5.90e-08 9.67e-09 1.15e-08 6.08e-08 -46492.857000 3.511000 5.788000
6 61.956611 Ca43:F19:-1e shape calibrant 62 1 61.956621 0.000000 False emg22 MLE 1.210000 25957.150000 109.852711 61.956621 1.20e-08 9.67e-09 2.99e-09 1.57e-08 -39896.121000 0.904000 0.197000
7 61.958997 ? - nan nan nan nan False emg22 MLE 1.210000 27.790000 6.695576 61.959016 5.34e-07 9.67e-09 1.66e-07 5.60e-07 nan 32.295000 nan
         stat. errors from resampling

Get refined peak-shape uncertainties using MCMC parameter samples (optional)

The default peak-shape error estimation in fit_peaks() relies on some simplifying assumptions:

  1. The posterior distributions of the shape parameters follow normal distributions.

  2. The shape parameters are uncorrelated.

In many cases, at least one of those assumptions is violated. Therefore, a refined way of estimating the peak-shape uncertainties has been added to emgfit: get_MC_peakshape_errors(). This method uses Markov-Chain Monte Carlo (MCMC) sampling to estimate the posterior distributions of the shape parameters. The sampling results are compiled in a corner plot/”covariance map” which includes both 1D-histograms of the parameter posteriors and 2D-histograms of the parameter correlations. By randomly drawing shape parameters sets from the obtained MCMC samples one obtains a representation of all peak shapes supported by the data. The calibrant and the peaks of interest are then re-fitted with all drawn shape parameter sets. Refined peak-shape uncertainties are obtained from the RMS deviation of the resulting mass values and peak areas from the best-fit values obtained with fit_peaks(). Usually, accounting for parameter correlations results in significantly smaller peak-shape errors.

The MCMC sampling can also already be performed during the peak-shape calibration using the map_par_covar option of determine_peak_shape(). The corner plot of the parameter covariances can be used to assess whether get_MC_peakshape_errors() should be run. For more details on the MC peak-shape uncertainty estimation see docs of get_MC_peakshape_errors().

[11]:
# NOTE: For quicker run time in this demo, the length of the sampling chain, the thinning
# interval and the number of shape parameter sets to perform fits with were manually
# reduced with the `steps`, `thin` and `N_samples` arguments, respectively. This triggers
# a warning about the insufficient MCMC chain length. For reasonable results those
# parameters should be increased to at least their default values. For these specific
# data, decent results are obtained using the following:
# `steps = 16000, thin = 280, N_samples = 1000`

spec.get_MC_peakshape_errors(steps=1000, thin=20, N_samples=50) # arguments adapted for demonstration
#spec.get_MC_peakshape_errors() # typical execution with default arguments

##### Evaluating posterior PDFs using MCMC sampling #####

Number of varied parameters:                 ndim = 11
Number of MCMC walkers:                  nwalkers = 220
Number of MCMC steps:                       steps = 1000
Number of initial steps to discard:          burn = 500
Length of thinning interval:                 thin = 20
Final number of parameter samples:  N_samples_tot = 5500
Number of CPU cores used:                 n_cores = 4
The chain is shorter than 50 times the integrated autocorrelation time for 11 parameter(s). Use this estimate with caution and run a longer chain!
N/50 = 20;
tau: [72.26995054 64.92317622 70.35294599 69.94821616 70.15518428 74.95706831
 69.53185493 74.55876255 73.11989272 72.48170904 76.70771468]
_images/emgfit_tutorial_23_3.png
_images/emgfit_tutorial_23_4.png
The chain is shorter than 50 times the integrated autocorrelation time for 11 parameter(s). Use this estimate with caution and run a longer chain!
N/50 = 20;
tau: [72.26995054 64.92317622 70.35294599 69.94821616 70.15518428 74.95706831
 69.53185493 74.55876255 73.11989272 72.48170904 76.70771468]
/opt/hostedtoolcache/Python/3.8.18/x64/lib/python3.8/site-packages/emgfit/spectrum.py:2305: UserWarning: Thinning interval `thin` is less than the integrated autocorrelation time for at least one parameter. Consider increasing `thin` MCMC keyword argument to ensure independent parameter samples.
  warnings.warn("Thinning interval `thin` is less than the "
Autocorrelation time for the parameters:
----------------------------------------
     bkg_c: 72.27 steps
    p6_amp: 64.92 steps
     p6_mu: 70.35 steps
  p6_sigma: 69.95 steps
  p6_theta: 70.16 steps
 p6_eta_m1: 74.96 steps
 p6_tau_m1: 69.53 steps
 p6_tau_m2: 74.56 steps
 p6_eta_p1: 73.12 steps
 p6_tau_p1: 72.48 steps
 p6_tau_p2: 76.71 steps

Total number of MCMC parameter sets after discarding burn-in and thinning: 5500


Covariance map for peak 6 with 0.16, 0.50 & 0.84 quantiles (dashed lines) and best-fit values (blue lines):
_images/emgfit_tutorial_23_8.png

##### MC Peak-shape uncertainty evaluation for peaks 0,1,2,3,4,5,6,7 #####

Determining MC recalibration factors from shifted centroids of mass calibrant.

All mass uncertainties below take into account the corresponding mass shifts of the calibrant peak.

Fitting peaks with 50 different MCMC-shape-parameter sets to determine refined peak-shape errors.
### Results ###

         Relative peak-shape (mass) uncertainty    Peak-shape uncertainty of peak areas
         from +-1σ variation / from MC samples     from +-1σ variation / from MC samples
Peak  0:           2.05e-08  /  1.01e-08                          0.9  /    0.3 counts
Peak  1:           2.46e-08  /  9.54e-09                          1.2  /    0.4 counts
Peak  2:           1.60e-08  /  4.80e-09                          4.9  /    2.2 counts
Peak  3:           8.38e-09  /  1.86e-09                          3.8  /    1.2 counts
Peak  4:           0.00e+00  /  0.00e+00                          0.0  /   11.8 counts
Peak  5:           1.15e-08  /  4.56e-09                          9.7  /    4.9 counts
Peak  6:           2.99e-09  /  7.67e-10                         15.8  /    9.1 counts
Peak  7:           1.66e-07  /  6.15e-08                          4.6  /    2.6 counts

Updated area error, peak-shape error and (total) mass error of peak(s) 0,1,2,3,4,5,6,7.


Peak properties table after MC peak-shape error evaluation:
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False emg22 MLE 1.210000 38.360000 6.923146 blinded 3.68e-07 9.67e-09 1.01e-08 3.69e-07 blinded 21.258000 blinded
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False emg22 MLE 1.210000 27.670000 5.415125 61.932050 5.91e-07 9.67e-09 9.54e-09 5.91e-07 -62784.026000 34.100000 3.518000
2 61.934369 ? Non-isobaric nan nan nan nan False emg22 MLE 1.210000 3870.630000 61.422191 61.934369 2.74e-08 9.67e-09 4.80e-09 2.94e-08 nan 1.697000 nan
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False emg22 MLE 1.210000 953.360000 38.868467 blinded 5.61e-08 9.67e-09 1.86e-09 5.69e-08 blinded 3.285000 blinded
4 61.946994 Ti46:O16:-1e mass calibrant 62 1 61.946992 0.000000 False emg22 MLE 1.210000 33986.590000 115.120988 61.946992 9.55e-09 9.67e-09 nan nan -48865.272000 nan 0.000000
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False emg22 MLE 1.210000 1572.110000 35.302795 61.949539 5.90e-08 9.67e-09 4.56e-09 5.99e-08 -46492.857000 3.457000 5.788000
6 61.956611 Ca43:F19:-1e shape calibrant 62 1 61.956621 0.000000 False emg22 MLE 1.210000 25957.150000 109.093094 61.956621 1.20e-08 9.67e-09 7.67e-10 1.54e-08 -39896.121000 0.889000 0.197000
7 61.958997 ? - nan nan nan nan False emg22 MLE 1.210000 27.790000 5.501455 61.959016 5.34e-07 9.67e-09 6.15e-08 5.38e-07 nan 31.038000 nan
         stat. errors from resampling     Monte Carlo peakshape errors

Export fit results

Once all peaks have been fitted, the results can be exported to an XLSX Excel file and a separate file with the peak-shape calibration parameters by using the save_results() method. The XLSX file contains three sheets:

  1. General properties of the spectrum object, e.g. input filename, used versions of emgfit and all relevant dependencies

  2. The peak properties table with all fit results along with linear and logarithmic plots of all peak fits

  3. The mass shifts obtained via \(\pm1\sigma\) variation of the shape parameters in the default peak-shape uncertainty evaluation

[12]:
spec.save_results('outputs/'+filename+' fitting MLE')
Fit results saved to file: outputs/2019-09-13_004-_006 SUMMED High stats 62Ga fitting MLE.xlsx

Peak-shape calibration saved to file: outputs/2019-09-13_004-_006 SUMMED High stats 62Ga fitting MLE_peakshape_calib.txt

That’s it! In principle we’re be done with the fitting at this point. Next we would probably take a look at the output file and proceed with some post-processing in Excel (e.g. combining mass values from different spectra etc.).

However, since emgfit gives the user a large amount of freedom, there’s are a number of things that could have been done differently depending on your preferences. So here is some possible…

Alternative procedures:

The above steps represent a full spectrum analysis. However, emgfit gives the user the freedom to take many different routes in processing the spectrum. Some of the possible alternatives are presented in the following:

Alternative 1: Performing the mass recalibration separately before the ion-of-interest fits

All steps up to the final peak fit are identical. For breviety here we simply create an exact clone of the above spectrum object:

[13]:
import copy
spec2 = copy.deepcopy(spec) # create a clone of the spectrum object

First obtain the recalibration factor from a fit of the mass calibrant

[14]:
spec2.fit_calibrant(species_mass_calib='Ti46:O16:-1e', show_fit_report=False)
##### Calibrant fit #####
_images/emgfit_tutorial_31_1.png
_images/emgfit_tutorial_31_2.png

##### Mass recalibration #####

Relative literature error of mass calibrant:     1.6e-09
Relative statistical error of mass calibrant:    1.2e-08

Recalibration factor:    0.999999703 = 1 -3.0e-07
Relative recalibration error:    1.2e-08

Fit all peaks and use the mass recalibration factor from above to calculate the final mass values

[15]:
spec2.fit_peaks(show_fit_report=False)
##### Fit peaks of interest #####
_images/emgfit_tutorial_33_1.png
_images/emgfit_tutorial_33_2.png

##### Peak-shape uncertainty evaluation #####

Determining centroid shifts of mass calibrant.

All centroid shifts below are corrected for the corresponding shifts of the mass calibrant peak.

Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  0 by  -0.23  /   0.21 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  1 by   0.45  /  -0.43 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  2 by  -0.30  /   0.26 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  3 by  -0.04  /  -0.02 μu/z  & its area by   0.7 /   0.1 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  5 by  -0.04  /  -0.00 μu/z  & its area by   0.5 /  -0.6 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  6 by  -0.01  /  -0.02 μu/z  & its area by   0.2 /  -0.5 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  7 by  -0.11  /  -0.01 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  0 by  -0.53  /   0.48 μu/z  & its area by  -0.2 /   0.2 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  1 by  -1.11  /   1.01 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  2 by  -0.42  /   0.30 μu/z  & its area by  -0.8 /   0.5 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  3 by  -0.22  /   0.08 μu/z  & its area by  -0.5 /   1.1 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  5 by  -0.26  /   0.23 μu/z  & its area by  -0.1 /   0.0 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  6 by  -0.14  /   0.04 μu/z  & its area by  -3.2 /   1.9 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  7 by  -1.98  /   1.87 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  0 by   0.69  /  -0.57 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  1 by   0.33  /  -0.28 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  2 by   0.34  /  -0.41 μu/z  & its area by  -2.0 /   1.2 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  3 by   0.27  /  -0.28 μu/z  & its area by  -0.7 /   1.0 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  5 by   0.25  /  -0.21 μu/z  & its area by  -3.9 /   3.0 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  6 by   0.01  /  -0.09 μu/z  & its area by  -4.5 /   2.7 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  7 by   0.57  /  -0.51 μu/z  & its area by  -0.5 /   0.4 counts.

Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  0 by  -0.81  /   0.78 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  1 by   0.39  /  -0.28 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  2 by  -0.50  /   0.41 μu/z  & its area by   0.3 /  -0.3 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  3 by  -0.13  /   0.03 μu/z  & its area by   0.7 /  -0.1 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  5 by  -0.14  /   0.08 μu/z  & its area by   0.9 /  -0.8 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  6 by  -0.03  /  -0.01 μu/z  & its area by  -0.5 /  -0.7 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  7 by  -0.92  /   0.85 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  0 by  -0.19  /   0.35 μu/z  & its area by   0.3 /  -0.4 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  1 by  -0.09  /   0.05 μu/z  & its area by   0.4 /  -0.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  2 by  -0.09  /   0.18 μu/z  & its area by   1.9 /  -2.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  3 by  -0.05  /   0.17 μu/z  & its area by   0.4 /  -0.8 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  5 by  -0.05  /   0.09 μu/z  & its area by   3.2 /  -3.8 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  6 by   0.01  /  -0.02 μu/z  & its area by   4.6 /  -6.6 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  7 by   0.15  /  -0.11 μu/z  & its area by   0.4 /  -0.5 counts.

Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  0 by  -0.10  /   0.05 μu/z  & its area by  -0.3 /   0.2 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  1 by  -0.32  /   0.15 μu/z  & its area by  -0.5 /   0.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  2 by  -0.23  /   0.19 μu/z  & its area by  -1.8 /   1.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  3 by  -0.17  /   0.14 μu/z  & its area by  -1.7 /   1.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  5 by  -0.28  /   0.16 μu/z  & its area by   0.3 /  -0.8 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  6 by   0.01  /  -0.03 μu/z  & its area by  -6.7 /   5.9 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  7 by  -4.04  /   3.08 μu/z  & its area by   1.0 /  -1.0 counts.

Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  0 by   0.10  /  -0.10 μu/z  & its area by   0.1 /  -0.1 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  1 by   0.36  /  -0.58 μu/z  & its area by   0.3 /  -0.4 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  2 by   0.31  /  -0.42 μu/z  & its area by   0.3 /  -0.5 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  3 by   0.14  /  -0.22 μu/z  & its area by   0.9 /  -0.2 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  5 by   0.03  /  -0.03 μu/z  & its area by   1.0 /  -1.2 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  6 by   0.00  /   0.01 μu/z  & its area by   0.8 /  -1.0 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  7 by   0.47  /  -0.50 μu/z  & its area by   0.3 /  -0.4 counts.

Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  0 by  -0.03  /   0.00 μu/z  & its area by   0.4 /  -0.5 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  1 by   0.34  /  -0.42 μu/z  & its area by   0.5 /  -0.7 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  2 by   0.05  /  -0.08 μu/z  & its area by   2.7 /  -3.0 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  3 by   0.06  /  -0.09 μu/z  & its area by   2.3 /  -2.6 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  5 by   0.43  /  -0.51 μu/z  & its area by  -7.8 /   5.2 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  6 by  -0.02  /   0.02 μu/z  & its area by  11.3 /  -8.7 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  7 by   4.33  /  -9.18 μu/z  & its area by  -4.4 /   3.5 counts.

Relative peak-shape error of peak  0:  2.0e-08
Relative peak-shape error of peak  1:  2.5e-08
Relative peak-shape error of peak  2:  1.6e-08
Relative peak-shape error of peak  3:  8.2e-09
Relative peak-shape error of peak  5:  1.1e-08
Relative peak-shape error of peak  6:  2.8e-09
Relative peak-shape error of peak  7:  1.7e-07
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False emg22 MLE 1.210000 38.360000 6.816128 blinded 3.46e-07 1.17e-08 2.05e-08 3.46e-07 blinded 19.985000 blinded
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False emg22 MLE 1.210000 27.670000 6.085935 61.932050 4.07e-07 1.17e-08 2.46e-08 4.08e-07 -62784.032000 23.529000 3.511000
2 61.934369 ? Non-isobaric nan nan nan nan False emg22 MLE 1.210000 3870.630000 62.728006 61.934369 3.44e-08 1.17e-08 1.57e-08 3.96e-08 nan 2.285000 nan
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False emg22 MLE 1.210000 953.360000 31.394839 blinded 6.93e-08 1.17e-08 8.17e-09 7.08e-08 blinded 4.084000 blinded
4 61.946994 Ti46:O16:-1e mass calibrant 62 1 61.946992 0.000000 False emg22 MLE 1.750000 33933.850000 186.180000 61.946992 1.16e-08 1.17e-08 nan nan -48865.272000 nan 0.000000
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False emg22 MLE 1.210000 1572.110000 41.266207 61.949539 5.40e-08 1.17e-08 1.14e-08 5.64e-08 -46492.864000 3.254000 5.782000
6 61.956611 Ca43:F19:-1e shape calibrant 62 1 61.956621 0.000000 False emg22 MLE 1.210000 25957.150000 162.039372 61.956621 1.33e-08 1.17e-08 2.79e-09 1.79e-08 -39896.127000 1.035000 0.191000
7 61.958997 ? - nan nan nan nan False emg22 MLE 1.210000 27.790000 7.937959 61.959016 4.06e-07 1.17e-08 1.66e-07 4.39e-07 nan 25.324000 nan

Alternative 2: Chi-square instead of MLE fitting

All steps up to the final peak fit are identical. For breviety here we simply create an exact clone of the above spectrum object and re-use the above peak-shape calibration (obtained with chi-square fitting):

[16]:
import copy
spec_chi_sq = copy.deepcopy(spec) # create a clone of the spectrum object

# Use Pearson's chi-squared statistic for A_stat_emg determination
spec_chi_sq.determine_A_stat_emg(species='Ca43:F19:-1e', x_range=0.004,
                                 cost_func='chi-square', N_spectra=10,
                                 plot_filename='outputs/'+filename+'_chi-square')
Creating synthetic spectra via bootstrap re-sampling and fitting them for A_stat determination.
Depending on the choice of `N_spectra` this can take a few minutes. Interrupt kernel if this takes too long.
Done!

[[Model]]
    Model(powerlaw)
[[Fit Statistics]]
    # fitting method   = leastsq
    # function evals   = 5
    # data points      = 8
    # variables        = 1
    chi-square         = 9.9655e-09
    reduced chi-square = 1.4236e-09
    Akaike info crit   = -162.028635
    Bayesian info crit = -161.949193
    R-squared          = -9.47838719
[[Variables]]
    amplitude:  1.4311e-04 +/- 1.3956e-05 (9.75%) (init = 1)
    exponent:  -0.5 (fixed)
_images/emgfit_tutorial_35_3.png
A_stat of a Gaussian model: 0.425
Default A_stat_emg for Hyper-EMG models: 0.52
A_stat_emg for this spectrum's emg22 fit model: 0.598 +- 0.059
[17]:
# Fit with Pearson's chi-squared statistic as cost function
spec_chi_sq.fit_peaks(species_mass_calib='Ti46:O16:-1e', cost_func='chi-square',
                      show_fit_report=False)
##### Fit peaks of interest #####
_images/emgfit_tutorial_36_1.png
_images/emgfit_tutorial_36_2.png

##### Mass recalibration #####

Relative literature error of mass calibrant:     1.6e-09
Relative statistical error of mass calibrant:    1.3e-08

Recalibration factor:    0.999999702 = 1 -3.0e-07
Relative recalibration error:    1.3e-08


##### Peak-shape uncertainty evaluation #####

Determining centroid shifts of mass calibrant.

All centroid shifts below are corrected for the corresponding shifts of the mass calibrant peak.

Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  0 by   0.51  /  -0.48 μu/z  & its area by   0.3 /  -0.3 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  1 by   0.67  /  -0.57 μu/z  & its area by   0.2 /  -0.2 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  2 by  -0.32  /   0.33 μu/z  & its area by  18.1 / -16.6 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  3 by  -0.11  /   0.10 μu/z  & its area by   1.8 /  -1.4 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  5 by  -0.05  /   0.05 μu/z  & its area by   2.1 /  -1.4 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  6 by   0.01  /  -0.01 μu/z  & its area by   7.3 /   6.7 counts.
Re-fitting with sigma  =  8.34e-05 +/- 2.96e-06 shifts peak  7 by  -1.88  /   1.73 μu/z  & its area by   0.3 /  -0.3 counts.

Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  0 by  -0.66  /   0.56 μu/z  & its area by  -0.4 /   0.3 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  1 by  -1.30  /   1.12 μu/z  & its area by  -0.5 /   0.4 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  2 by  -0.63  /   0.57 μu/z  & its area by -17.6 /  18.0 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  3 by  -0.26  /   0.23 μu/z  & its area by  -2.4 /   2.5 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  5 by  -0.16  /   0.16 μu/z  & its area by   0.1 /   0.5 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  6 by  -0.02  /   0.02 μu/z  & its area by   4.6 /   9.3 counts.
Re-fitting with theta  =  7.25e-01 +/- 2.35e-02 shifts peak  7 by  -2.33  /   2.13 μu/z  & its area by   0.2 /  -0.2 counts.

Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  0 by   0.48  /  -0.43 μu/z  & its area by  -0.6 /   0.4 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  1 by   0.06  /  -0.09 μu/z  & its area by  -0.6 /   0.4 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  2 by   0.51  /  -0.42 μu/z  & its area by  -8.5 /   9.4 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  3 by   0.47  /  -0.34 μu/z  & its area by  -2.6 /   2.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  5 by   0.24  /  -0.15 μu/z  & its area by  -4.2 /   3.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  6 by   0.03  /  -0.02 μu/z  & its area by   5.6 /  11.3 counts.
Re-fitting with eta_m1 =  9.20e-01 +/- 2.31e-02 shifts peak  7 by   0.85  /  -0.68 μu/z  & its area by  -0.7 /   0.5 counts.

Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  0 by  -0.62  /   0.61 μu/z  & its area by   0.4 /  -0.4 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  1 by   0.25  /  -0.21 μu/z  & its area by   0.3 /  -0.3 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  2 by  -0.61  /   0.60 μu/z  & its area by  18.6 / -16.6 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  3 by  -0.28  /   0.26 μu/z  & its area by   2.5 /  -1.9 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  5 by  -0.25  /   0.24 μu/z  & its area by   2.3 /  -1.3 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  6 by  -0.04  /   0.03 μu/z  & its area by   9.9 /   8.5 counts.
Re-fitting with tau_m1 =  4.49e-05 +/- 5.88e-06 shifts peak  7 by  -2.28  /   2.37 μu/z  & its area by   0.4 /  -0.4 counts.

Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  0 by   0.07  /  -0.09 μu/z  & its area by   0.4 /  -0.7 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  1 by   0.18  /  -0.40 μu/z  & its area by   0.5 /  -0.8 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  2 by  -0.11  /   0.21 μu/z  & its area by   4.9 /  -5.3 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  3 by  -0.17  /   0.27 μu/z  & its area by   1.8 /  -2.3 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  5 by  -0.01  /   0.08 μu/z  & its area by   3.4 /  -4.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  6 by   0.01  /  -0.01 μu/z  & its area by  11.8 /  -2.5 counts.
Re-fitting with tau_m2 =  1.77e-04 +/- 2.48e-05 shifts peak  7 by   0.05  /   0.10 μu/z  & its area by   0.5 /  -0.8 counts.

Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  0 by  -0.39  /   0.16 μu/z  & its area by  -0.6 /   0.4 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  1 by  -0.71  /   0.37 μu/z  & its area by  -0.8 /   0.6 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  2 by  -0.31  /   0.25 μu/z  & its area by  -9.7 /   9.5 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  3 by  -0.21  /   0.15 μu/z  & its area by  -3.1 /   2.7 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  5 by  -0.22  /   0.17 μu/z  & its area by   0.5 /  -0.8 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  6 by   0.02  /  -0.01 μu/z  & its area by   0.0 /  12.9 counts.
Re-fitting with eta_p1 =  7.90e-01 +/- 4.49e-02 shifts peak  7 by  -7.46  /   5.13 μu/z  & its area by   1.1 /  -1.1 counts.

Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  0 by  -0.29  /   0.25 μu/z  & its area by   0.2 /  -0.3 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  1 by   0.32  /  -0.65 μu/z  & its area by   0.4 /  -0.4 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  2 by   0.54  /  -0.68 μu/z  & its area by  13.7 / -14.2 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  3 by   0.23  /  -0.29 μu/z  & its area by   2.2 /  -2.0 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  5 by  -0.04  /   0.03 μu/z  & its area by   0.9 /  -0.0 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  6 by  -0.00  /  -0.00 μu/z  & its area by   9.5 /   8.9 counts.
Re-fitting with tau_p1 =  1.24e-04 +/- 1.33e-05 shifts peak  7 by   0.30  /  -0.41 μu/z  & its area by   0.3 /  -0.4 counts.

Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  0 by   0.42  /  -0.60 μu/z  & its area by   0.5 /  -0.8 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  1 by   0.65  /  -0.94 μu/z  & its area by   0.7 /  -0.9 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  2 by   0.02  /  -0.08 μu/z  & its area by   5.3 /  -6.3 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  3 by   0.05  /  -0.11 μu/z  & its area by   2.6 /  -3.2 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  5 by   0.47  /  -0.57 μu/z  & its area by  -7.8 /   5.4 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  6 by  -0.02  /   0.02 μu/z  & its area by  17.4 /  -9.3 counts.
Re-fitting with tau_p2 =  4.08e-04 +/- 5.48e-05 shifts peak  7 by   8.42  / -22.14 μu/z  & its area by  -4.3 /   4.3 counts.

Relative peak-shape error of peak  0:  2.2e-08
Relative peak-shape error of peak  1:  3.3e-08
Relative peak-shape error of peak  2:  2.1e-08
Relative peak-shape error of peak  3:  1.2e-08
Relative peak-shape error of peak  5:  1.2e-08
Relative peak-shape error of peak  6:  9.9e-10
Relative peak-shape error of peak  7:  3.8e-07
  x_pos species comment A z m_AME m_AME_error extrapolated fit_model cost_func red_chi area area_error m_ion rel_stat_error rel_recal_error rel_peakshape_error rel_mass_error atomic_ME_keV mass_error_keV m_dev_keV
0 61.927800 Ni62:-1e - 62 1 61.927796 0.000000 False emg22 chi-square 0.990000 40.920000 7.644423 blinded 3.61e-07 1.26e-08 2.24e-08 3.62e-07 blinded 20.892000 blinded
1 61.932021 Cu62:-1e? - 62 1 61.932046 0.000001 False emg22 chi-square 0.990000 28.220000 6.880458 61.932066 4.35e-07 1.26e-08 3.30e-08 4.36e-07 -62769.065000 25.177000 18.479000
2 61.934369 ? Non-isobaric nan nan nan nan False emg22 chi-square 0.990000 3990.650000 74.395022 61.934370 3.66e-08 1.26e-08 2.13e-08 4.42e-08 nan 2.548000 nan
3 61.943618 Ga62:-1e - 62 1 61.943641 0.000001 False emg22 chi-square 0.990000 957.830000 32.298791 blinded 7.47e-08 1.26e-08 1.24e-08 7.67e-08 blinded 4.427000 blinded
4 61.946994 Ti46:O16:-1e mass calibrant 62 1 61.946992 0.000000 False emg22 chi-square 0.990000 34014.370000 184.700000 61.946992 1.25e-08 1.26e-08 nan nan -48865.272000 nan 0.000000
5 61.949527 Sc46:O16:-1e - 62 1 61.949533 0.000001 False emg22 chi-square 0.990000 1573.740000 41.749844 61.949539 5.82e-08 1.26e-08 1.17e-08 6.07e-08 -46492.884000 3.505000 5.762000
6 61.956611 Ca43:F19:-1e shape calibrant 62 1 61.956621 0.000000 False emg22 chi-square 0.990000 25965.460000 164.684606 61.956621 1.43e-08 1.26e-08 9.89e-10 1.91e-08 -39896.168000 1.104000 0.151000
7 61.958997 ? - nan nan nan nan False emg22 chi-square 0.990000 27.460000 8.370657 61.958998 4.41e-07 1.26e-08 3.82e-07 5.84e-07 nan 33.685000 nan