ramantools package
ramantools.ramantools module
- ramantools.ramantools.bgsubtract(x_data, y_data, polyorder=1, toplot=False, fitmask=None, hmin=50, hmax=10000, wmin=4, wmax=60, prom=10, exclusion_factor=6, peak_pos=None)[source]
Takes Raman shift and Raman intensity data and automatically finds peaks in the spectrum, using
scipy.find_peaks
. These peaks are then used to define the areas of the background signal. In the areas with the peaks removed, the background is fitted by a polynomial of order given by the optional argument:polyorder
. The fit is performed byscipy.optimize.curve_fit
. The function returns the Raman intensity counts with the background removed, the background polinomial values themselves and the coefficients of the background fit results, as used bynumpy.polyval
.In cases, where the automatic peak find is not functioning as expected, one can pass the values in
x_data
, at which peaks appear. In this case, thewmin
option determines the width of all peaks.If a
fitmask
is supplied for fitting, the fitmask is not calculated and only a polynomial fit is performed. This can decrease the runtime.- Parameters:
x_data (
numpy
array) – Raman shift valuesy_data (
numpy
array) – Raman intensity valuespolyorder (int, optional) – order of polynomial used to fit the background, defaults to 1
toplot (bool, optional) – if True a plot of: the fit, the background used and positions of the peaks is shown, defaults to False
fitmask (
numpy
array) – Fitmask to be used for polynomial fitting.hmin (float, optional) – minimum height of the peaks passed to
scipy.signal.find_peaks
, defaults to 50hmax (float, optional) – maximum height of the peaks passed to
scipy.signal.find_peaks
, defaults to 10000wmin (float, optional) – minimum width of the peaks, passed to
scipy.signal.find_peaks
, defaults to 4wmax (float, optional) – maximum width of the peaks passed to
scipy.signal.find_peaks
, defaults to 60prom (float, optional) – prominence of the peaks, passed to
scipy.signal.find_peaks
, defaults to 10exclusion_factor (float, optional) – this parameter multiplies the width of the peaks found by
scipy.signal.find_peaks
, or specified bywmin
if the peak positions are passed by hand, usingpeak_pos
, defaults to 6peak_pos (list of floats, optional) – list of the peak positions in
x_data
values used for exclusion, defaults to None
- Returns:
y_data_nobg
,bg_values
,coeff
,params_used_at_run
,mask
,covar
- Return type:
tuple: (
numpy
array,numpy
array,numpy
array, dictionary,numpy
array,numpy
array)
y_data_nobg
: Raman counts, with the background subtracted,bg_values
: the polynomial values of the fit, at the Raman shift positions,coeff
: coefficients of the polynomial fit, as used by:numpy.polyval
,params_used_at_run
: parameters used at runtimemask
: the calculated fitmaskcovar
: covariance of the fit parameters
Note
Using the option:
peak_pos
, awmin*exclusion_factor/2
region (measured in datapoints) on both sides of the peaks is excluded from the background fit. If automatic peak finding is used, the exclusion area is calculated in a similar way, but the width of the individual peaks are used, as determined byscipy.signal.find_peaks
.
- ramantools.ramantools.gaussian(x, x0=1580, ampl=10, width=15, offset=0)[source]
Gaussian function. Width and amplitude parameters have the same meaning as for
lorentz()
.- Parameters:
x (float,
numpy
array, etc.) – values for the x coordinatex0 (float) – shift along the x corrdinate
ampl (float) – amplitude of the peak
width (float) – FWHM of the peak
offset (float) – offset along the function value
- Returns:
values of a Gaussian function
- Return type:
float,
numpy
array, etc.
- ramantools.ramantools.lorentz(x, x0=1580, ampl=1, width=14, offset=0)[source]
Single Lorentz function
- Parameters:
x (float,
numpy
array, etc.) – values for the x coordinatex0 (float) – shift along the x corrdinate
ampl (float) – amplitude of the peak
width (float) – FWHM of the peak
offset (float) – offset along the function value
- Returns:
values of a single Lorentz function
- Return type:
float,
numpy
array, etc.
Note
The area of the peak can be given by:
area = np.pi * amplitude * width / 2
- ramantools.ramantools.lorentz2(x, x01=2700, ampl1=1, width1=15, x02=2730, ampl2=1, width2=15, offset=0)[source]
Double Lorentz function
- Returns:
values of a double Lorentz function
- Return type:
float,
numpy
array, etc.- Parameters:
x (float,
numpy
array, etc.) – values for the x coordinatex0 (float) – shift along the x corrdinate
area (float) – area of the peak
width (float) – width (FWHM) of the peak
offset (float) – offset along the function value
- ramantools.ramantools.peakfit(xrobj, func=<function lorentz>, fitresult=None, stval=None, bounds=None, toplot=False, width=None, height=None, **kwargs)[source]
Fitting a function to peaks in the data contained in
xrobj
, which can be a single spectrum, a map or a selected spectrum from a map.- Parameters:
xrobj (
xarray
DataArray) –xarray
DataArray, of a single spectrum or a map.func (function, optional) – function to be used for fitting, defaults to lorentz
fitresult (
xarray
Dataset, optional) – anxarray
Dataset of a previous fit calculation, with matching dimensions. If this is passed topeakfit()
, the fit calculation in skipped and the passed Dataset is used.stval (dictionary of
func
parameters, optional) – starting values for the fit parameters offunc
. You are free to specify only some of the values, the rest will be filled by defaults. Defaults are given in the starting values for keyword arguments infunc
.bounds (dictionary of
func
parameters, with tuples containing lower and upper bounds, optional) – bounds for the fit parameters, used byxarray.curvefit
. Simlar dictionary, likestval
, but the values area a list, with lower and upper components. Defaults to Nonetoplot (boolean, optional) – plot the fit result, defaults to
False
width (int or float, optional) – width parameter of an
xarray
map to be used in conjunction withtoplot = True
height (int or float, optional) – height parameter of an
xarray
map to be used in conjunction withtoplot = True
- Returns:
fitted parameters of
func
and covariances in a Dataset- Return type:
xarray
Dataset- Example:
import ramantools as rt map_path = r'data path on you machine' info_path = r'metadata path on your machine' # use raw strings, starting with `r'` to escape special characters, such as backslash # Load a map map = rt.ramanmap(map_path, info_path) # Creating the dictionary for the starting values and the bounds # The default function is `lorentz`, with parameters: p = {'x0': 2724, 'ampl': 313, 'width': 49, 'offset': 0} # passing the starting values contained in `p` and bounds: `b` to the `peakfit()` method. b = {'x0': [2500, 2900], 'ampl': [0, 900], 'width': [20, 100], 'offset': [-10, 50]} mapfit = rt.peakfit(m_nobg.mapxr, stval = p, bounds = b, toplot = True)
Note
Use
toplot
= True to tweak the starting values. Iftoplot
= True, in case of a map, if nowidth
andheight
are specified, the middle of the map is used for plotting.Passing a
bounds
dictionary topeakfit()
seems to increase the fitting time significantly. This might be an issue withxarray.DataArray.curvefit
.By passing a previous fit result, using the optional parameter
fitresult
, we can just plot the fit result at multiple regions of the map.In case of using double Lorentz fitting, the names of the parameters change! See:
lorentz2()
.
See also
It is good practice, to crop the data to the vicinity of the peak you want to fit to.
- ramantools.ramantools.polynomial_fit(order, x_data, y_data)[source]
Polinomial fit to x_data, y_data
- Parameters:
order (int) – order of the polinomial to be fit
x_data (
numpy
array) – x coordinate of the data, this would be Raman shifty_data (
numpy
array) – y coordinate of the data, this would be Raman intensity
- Returns:
coefficients of the polinomial
coeff
, as used bynumpy.polyval
, covariance matrixcovar
, as returned byscipy.optimize.curve_fit
- Return type:
tuple: (
numpy
array,numpy
array)
- class ramantools.ramantools.ramanmap(map_path, info_path)[source]
Bases:
object
Container for Raman maps, imported from a text file. The text file needs to be exported as a “table” from Witec Project or Witec Control. Additional info also needs to be exported, containing the metadata for the measurement. This is the text next to the map data in the Witec software.
- Returns:
object containing the data and metadata
- Return type:
ramanmap
instance- Parameters:
map_path (str) – Path to the text file, containing the Raman map, exported from Witec
info_path (str) – Path to the info file, containing the metadata, exported from Witec
Most important variables of the
ramanmap
instance:- Variables:
mapxr – (type
xarray
DataArray) all data, coordinates and metadatamap – (type
numpy
array) Raman intensity valuesramanshift – (type
numpy
array) Raman shift values for the datapoints stored in mapmask – (type:
numpy
array) A boolean array of the same length as theramanshift
. It’s only available ifsinglespec.remove_bg()
orramanmap.remove_bg()
is called.samplename – (type: str) name of the sample, as shown in the Witec software.
mapname – (type: str) contains the name of the Raman map, as shown in the Witec software.
For a compete list see example below.
- Example:
import ramantools as rt map_path = r'data path on you machine' info_path = r'metadata path on your machine' # use raw strings, starting with `r'` to escape special characters, such as backslash map = rt.ramanmap(map_path, info_path) # list of the variables stored in the `ramanmap` instance print(list(map.__dict__))
- calibrate(peakshift, calibfactor=0, width=None, height=None, **kwargs)[source]
Calibrating a Raman map. It uses
peakfit()
to find the Raman peak closest topeakshift
and returns a newramanmap
instance, with the ramanshift coordinate offset to have the position of the peak atpeakshift
. If the optional argumentcalibfactor
is passed,peakshift
is ignored and the data is shifted by the given value. All possible keyword arguments ofpeakfit()
can be passed.- Parameters:
peakshift (float) – expected position, in 1/cm, of the peak used for calibration
calibfactor (int, optional) – If the calibration factor is known it can be passed directly. In this case
peakshift
is ignored, defaults to 0width (float, optional) – width coordinate of the spectrum, which will be used for calibration. Defaults to the middle of the map.
height (float, optional) – height coordinate of the spectrum, which will be used for calibration. Defaults to the middle of the map.
- Returns:
calibrated
ramanmap
instance- Return type:
- crr(cutoff=2, window=2, **kwargs)[source]
Tool for removing cosmic rays from a spectroscopy maps. The CRR peaks are determined as the standard deviation of the data: std times the cutoff value, in the window sized vicinity of each pixel.
- Parameters:
cutoff (int, optional) – multiplication factor for the data’s standard deviation; defaults to 2.
window (int, optional) – size of the neighborhood to consider; defaults to 2.
- Returns:
instance of the
ramanmap
class with the cosmic-ray peaks removed.- Return type:
Note
If CRR is not satisfactory, keep reducing the cutoff value and compare to the original data.
- history()[source]
Display the notes accumulated in the ‘comments’ attribute of the
ramanmap.mapxr
xarray variable.
- normalize(peakshift, width=None, height=None, mode='const', **kwargs)[source]
Normalize the Raman spectrum to the peak at
peakshift
. Returns a normalizedramanmap
instance. An exception will be raised if the background has not been removed. It usespeakfit()
to find the amplitude of the peak to be normalized. It accepts all keyword arguments accepted bypeakfit()
.- Parameters:
peakshift (float) – rough position of the peak in
ramanmap.mapxr.ramanshift
dimensionmode (str, optional) – Has two modes: ‘const’ and ‘individual’. defaults to ‘const’.
width (float, optional) – width coordinate of the spectrum, which will be used for normalization in ‘const’ mode, defaults to the middle of the map.
height (float, optional) – height coordinate of the spectrum, which will be used for normalization in ‘const’ mode, defaults to the middle of the map.
- Returns:
normalized
ramanmap
instance- Return type:
- Raises:
ValueError – Background needs to be removed for normalization to make sense.
ValueError – mode parameter must be either: ‘const’ or ‘individual’.
Note
Attributes of
ramanmap.mapxr
are updated to reflect the fact that the normalized peak intensities are dimensionless, with a new long_name.In
mode == 'individual'
, each spectrum in the map will be normalized to the local peak amplitude. Inmode == 'const'
, the peak at the position specified bywidth
andheight
is used for normalization. Ifmode == 'individual'
, thewidth
andheight
parameters are ignored.
- peakmask(peakpos, cutoff=0.1, width=None, height=None, **kwargs)[source]
Create a boolean mask for the map, where the mean Raman intensity of the peak at
peakpos
is larger than the peak mean in the selected spectrum by thecutoff
value. The method also returns the croppedxarray
DataArray, with the values that are cropped replaced by NaNs. The method needs a reference spectrum for determining the “typical” mean of the peak amplitude. This is also used to determine the cutoff value for the rest of the map as:'selected spectrum mean value' * cutoff
The optional
width
andheight
parameters can be passed, which selects that spectrum for reference. If these are not passed, the spectrum in the middle of the map is taken as reference.- Parameters:
peakpos (float) – position in 1/cm of the peak we want to create the mask for
cutoff (float, optional) – cutoff value, interpreted as a percentage. Values between 0 and 1. Defaults to 0.1
width (float, optional) – width parameter of the spectrum in the map we want to have as a reference, defaults to None
height (float, optional) – height parameter of the spectrum in the map we want to have as a reference, defaults to None
- Returns:
mapmasked
:ramanmap
instance containing the cropped mappeakmask
:xarray
DataArray containing the mask
- Return type:
tuple: (
ramanmap
,xarray
)
- plotspec(width, height, shift)[source]
Plots a Raman map at a given Raman shift and displays alongside a selected spectrum at a specified width and height. Needs width and height coordinates for the single spectrum and the Raman shift where we want to plot the Raman intensity in the map.
- Returns:
none
- Parameters:
width (float) – ‘width’ coordinate in um, from
ramanmap.mapxr
height (float) – ‘height’ coordinate in um, from
ramanmap.mapxr
shift (float) – ‘ramanshift’ coordinate in um, from
ramanmap.mapxr
- print_metadata()[source]
Prints the metadata of the
ramanmap
instance, imported from the info file.- Returns:
none
- remove_bg(mode='const', fitmask=None, height=None, width=None, **kwargs)[source]
Remove the background of Raman maps. It takes the same optional arguments as
bgsubtract()
. Default fit function is a first order polynomial. This can be changed by thepolyorder
parameter.There are several modes of background fitting, as a function of the optional variables:
mode
andfitmask
:mode
= ‘const’: (this is the default) Subtracts the same background from all spectra in the map. The background is determined by either a fitmask calculated bybgsubtract()
or passed to the method using thefitmask
parameter. Additionally, aheight
andwidth
parameter can be passed in which case, the spectrum at those coordinates is used to determine the background. If a fitmask is supplied, it is used instead of the coordinates.mode
= ‘individual’: An individual background is subtracted from all spectra.If
fitmask
= None,bgsubtract()
does a peak search for each Raman spectrum and the fitmask is determined based on the parameters passed tobgsubtract()
.If a
fitmask
is passed toremove_bg()
,bgsubtract()
only does the polynomial fit.
- Parameters:
mode (str, optional, default: ‘const’) – Values can be: ‘const’, ‘individual’
fitmask (
numpy
array, optional) – fitmask as returned bybgsubtract()
or as contained in the singlespec.mask variable of thesinglespec
class, or inramanmap
height (float, optional) – height coordinate of the
ramanmap.mapxr
xarraywidth (float, optional) – width coordinate of the
ramanmap.mapxr
xarray
- Returns:
map_mod
: newramanmap
instance, containing the data with background removed.coeff
: the coefficients andcovar
: covariance of the polynomial fit, as supplied bypolynomial_fit()
.
- Return type:
tuple: (
ramanmap
,numpy
,numpy
)
Note
If the peaks and background are complicated, it is advised to test the background fit, by selecting the most difficult spectrum from the map and tweaking the fit parameters directly, using
bgsubtract()
.Metadata is copied over to the returned
ramanmap
instance, because there are no unit changes, in removing the background. After running, the ‘comments’ attribute of the new xarray instance is updated with the background fit information.See also
mode = 'individual'
is not implemented yet.- Example:
import ramantools as rt spec_path = r'data path on you machine' info_path = r'metadata path on your machine' # use raw strings, starting with `r'` to escape special characters, such as backslash # Loading a single spectrum from files m = rt.ramanmap(spec_path, info_path) # create a new `ramanmap` instance `new_m`, with the background removed new_m, coeff, covar = m.remove_bg(fitmask = mask, toplot = True) # plot the spectrum in the middle of the map new_m.mapxr.sel(width = self.size_x/2, height = self.size_y/2, method = 'nearest').plot()
- class ramantools.ramantools.singlespec(spec_path, info_path)[source]
Bases:
object
Container for Raman single spectra, imported from a text file. The text file needs to be exported as a “table” from Witec Project or Witec Control Additional info also needs to be exported, containing the metadata for the measurement. This is the text next to the map data in the Witec software. It takes two arguments, the first is the path to the file containing the spectroscopy data, the second is the path to the metadata.
- Parameters:
spec_path (str) – Path to the text file, containing the Raman spectrum, exported from Witec
info_path (str) – Path to the info file, containing the metadata, exported from Witec
Most important variables of the
singlespec
instance:- Variables:
ssxr – (type
xarray
DataArray) all data, coordinates and metadatamask – (type:
numpy
array) A boolean array of the same length as theramanshift
. It’s only available ifsinglespec.remove_bg()
is called.counts – (type
numpy
array) Raman intensity valuesramanshift – (type
numpy
array) Raman shift values for the datapoints stored in mapsamplename – (type: str) name of the sample, as shown in the Witec software.
specname – (type: str) contains the name of the Raman single spectrum, as shown in the Witec software.
For a complete list see example below.
- Returns:
singlespec
instance containing the data and metadata- Return type:
singlespec
instance- Example:
import ramantools as rt spec_path = r'data path on you machine' info_path = r'metadata path on your machine' # use raw strings, starting with `r'` to escape special characters, such as backslash single_spectrum = rt.singlespec(spec_path, info_path) # list of variables stored in the `singlespec` instance print(list(single_spectrum.__dict__))
- calibrate(peakshift, calibfactor=0, **kwargs)[source]
Calibrating a single Raman spectrum. It uses
peakfit()
to find the Raman peak closest topeakshift
and returns a newsinglespec
instance, with the ramanshift coordinate offset to have the position of the peak atpeakshift
. If the optional argumentcalibfactor
is passed,peakshift
is ignored and the data is shifted by the given value. All possible keyword arguments ofpeakfit()
can be passed.- Parameters:
peakshift (float) – expected position, in 1/cm, of the peak used for calibration
calibfactor (int, optional) – If the calibration factor is known it can be passed directly. In this case
peakshift
is ignored, defaults to 0
- Returns:
calibrated
singlespec
instance- Return type:
- crr(cutoff=2, window=2, **kwargs)[source]
Tool for removing cosmic rays from a single spectrum. The CRR peaks are determined as the standard deviation of the data: std times the cutoff value, in the window sized vicinity of each pixel.
- Parameters:
cutoff (int, optional) – multiplication factor for the data’s standard deviation; defaults to 2.
window (int, optional) – size of the neighborhood to consider; defaults to 2.
- Returns:
instance of the
singlespec
class with the cosmic-ray peaks removed.- Return type:
Note
If CRR is not satisfactory, keep reducing the cutoff value and compare to the original data.
- history()[source]
Display the notes accumulated in the ‘comments’ attribute of the
singlespec.ssxr
xarray variable.
- normalize(peakshift, **kwargs)[source]
Normalize the Raman spectrum to the peak at
peakshift
. Returns a normalizedsinglespec
instance. An exception will be raised if the background has not been removed. It usespeakfit()
to find the amplitude of the peak to be normalized. It accepts all keyword arguments accepted bypeakfit()
.- Parameters:
peakshift (float) – rough position of the peak in
singlespec.ssxr.ramanshift
dimension- Returns:
normalized
singlespec
instance- Return type:
- Raises:
ValueError – Background needs to be removed for normalization to make sense.
Note
Attributes of
singlespec.ssxr
are updated to reflect the fact that the normalized peak intensities are dimensionless, with a new long_name.
- print_metadata()[source]
Prints the metadata of the
singlespec
instance, imported from the info file.- Returns:
none
- remove_bg(**kwargs)[source]
Remove the background of Raman spectra. The
xarray
variable of thesinglespec
instance is updated with the new dataset, with the background removed.It takes the same optional arguments as
bgsubtract()
. Default fit function is a first order polynomial. This can be changed by thepolyorder
parameter.- Returns:
singlesp_mod
: newsinglespec
instance, containing the data with background removed.coeff
: the coefficients andcovar
: covariance of the polynomial fit, as supplied bypolynomial_fit()
.
- Return type:
tuple: (
singlespec
,numpy
,numpy
)
Note
Metadata is copied over to the returned
singlespec
instance, because there are no unit changes, in removing the background. After running, the ‘comments’ attribute of the new xarray instance is updated with the background fit information.- Example:
import ramantools as rt spec_path = r'data path on you machine' info_path = r'metadata path on your machine' # use raw strings, starting with `r'` to escape special characters, such as backslash # Loading a single spectrum from files single_spectrum = rt.singlespec(spec_path, info_path) # Using remove_bg to fit and remove the background # In this example, we let remove_bg() find the peaks automatically. In this case, if no options are passed, the fit is returned. new_ss, coeff, covar = single_spectrum.remove_bg() # `new_ss` is the new `singlespec` instance containing the ssxr `xarray` object, with the background removed # In this example, we also want to plot the result and we select the peaks by hand, by using `peak_pos`. new_ss, coeff, covar = single_spectrum.remove_bg(toplot = True, peak_pos = [520, 1583, 2700], wmin = 15) # Fitting a third order polynomial new_ss, coeff, covar = single_spectrum.remove_bg(polyorder = 3) # Replot the `xarray` DataArray, which has the background removed new_ss.ssxr.plot()