rampy API

Subpackages

Submodules

rampy.baseline module

rampy.baseline.als_baseline(y: ndarray, lam: float, p: float, niter: int) ndarray

Asymmetric Least Squares baseline fitting.

rampy.baseline.arPLS_baseline(y: ndarray, lam: float, ratio: float) ndarray

Asymmetrically Reweighted Penalized Least Squares baseline fitting.

rampy.baseline.baseline(x_input: ndarray, y_input: ndarray, roi=array([[0, 100]]), method: str = 'poly', **kwargs) tuple

Subtracts a baseline from an x-y spectrum using various methods.

This function performs baseline subtraction on spectroscopic data by fitting a model to the background signal. It supports multiple correction methods, including polynomial fitting, splines, and advanced penalized least squares techniques.

Parameters:
  • x_input (ndarray) – The x-axis values (e.g., wavelength or wavenumber).

  • y_input (ndarray) – The y-axis values (e.g., intensity or absorbance).

  • roi (ndarray or list of lists, optional) –

    Regions of interest (ROI) for baseline fitting. Must be an n x 2 array or list of lists where each row specifies the start and end of a region. Default is np.array([[0, 100]]). .. rubric:: Example

    • Array: np.array([[100., 200.], [500., 600.]])

    • List: [[100., 200.], [500., 600.]]

    Note: Some methods (e.g., “als”, “arPLS”) do not use the roi parameter but require it to be specified.

  • method (str, optional) –

    The method used for baseline fitting. Default is “poly”. Available options:

    • ”poly”: Polynomial fitting. Requires polynomial_order.

    • ”unispline”: Spline fitting using Scipy’s UnivariateSpline. Requires s.

    • ”gcvspline”: Spline fitting using GCVSpline. Requires s and optionally ese_y.

    • ”gaussian”: Gaussian function fitting. Requires p0_gaussian.

    • ”exp”: Exponential background fitting. Requires p0_exp.

    • ”log”: Logarithmic background fitting. Requires p0_log.

    • ”rubberband”: Rubberband correction using convex hull interpolation.

    • ”als”: Asymmetric Least Squares correction. Requires lam, p, and niter.

    • ”arPLS”: Asymmetrically Reweighted Penalized Least Squares smoothing. Requires lam and ratio.

    • ”drPLS”: Doubly Reweighted Penalized Least Squares smoothing. Requires niter, lam, eta, and ratio.

    • ”whittaker”: Whittaker smoothing with weights applied to regions of interest. Requires lam.

    • ”GP”: Gaussian process method using a rational quadratic kernel.

  • **kwargs

    Additional parameters specific to the chosen method. - polynomial_order (int): Degree of the polynomial for the “poly” method. Default is 1. - s (float): Spline smoothing coefficient for “unispline” and “gcvspline”. Default is 2.0. - ese_y (ndarray): Errors associated with y_input for the “gcvspline” method.

    Defaults to sqrt(abs(y_input)) if not provided.

    • lam (float): Smoothness parameter for methods like “als”, “arPLS”, and others. Default is 1e5.

    • p (float): Weighting parameter for ALS method. Recommended values are between 0.001 and 0.1. Default is 0.01.

    • ratio (float): Convergence ratio parameter for arPLS/drPLS methods. Default is 0.01 for arPLS and 0.001 for drPLS.

    • niter (int): Number of iterations for ALS/drPLS methods. Defaults are 10 for ALS and 100 for drPLS.

    • eta (float): Roughness parameter for drPLS method, between 0 and 1. Default is 0.5.

    • p0_gaussian (list): Initial parameters [a, b, c] for Gaussian fitting: (y = a cdot exp(-log(2) cdot ((x-b)/c)^2)). Default is [1., 1., 1.].

    • p0_exp (list): Initial parameters [a, b, c] for exponential fitting: (y = a cdot exp(b cdot (x-c))). Default is [1., 1., 0.].

    • p0_log (list): Initial parameters [a, b, c, d] for logarithmic fitting: (y = a cdot log(-b cdot (x-c)) - d cdot x^2). Default is [1., 1., 1., 1.].

Returns:

corrected_signal (ndarray): The signal after baseline subtraction. baseline_fitted (ndarray): The fitted baseline.

Return type:

tuple

Raises:

ValueError – If the specified method is not recognized or invalid parameters are provided.

Notes

The input data is standardized before fitting to improve numerical stability during optimization. The fitted baseline is transformed back to the original scale before subtraction.

Examples

Example with polynomial baseline correction:

>>> import numpy as np
>>> x = np.linspace(50, 500, nb_points)
>>> noise = 2.0 * np.random.normal(size=nb_points)
>>> background = 5 * np.sin(x / 50) + 0.1 * x
>>> peak = rampy.gaussian(x, 100.0, 250.0, 7.0)
>>> y = peak + background + noise
>>> roi = np.array([[0, 200], [300, 500]])
>>> corrected_signal, baseline = baseline(x, y, method="poly", polynomial_order=5)

Example with GCVSpline algorithm:

>>> corrected_signal, baseline = baseline(x, y, method="gcvspline", s=2.0)

Example with Whittaker smoothing:

>>> corrected_signal, baseline = baseline(x, y, roi=roi, method="whittaker")

Example with Gaussian process:

>>> corrected_signal, baseline = baseline(x, y, roi=roi, method="GP")

Example with rubberband correction:

>>> corrected_signal, baseline = baseline(x, y, method="rubberband")

Example with ALS algorithm:

>>> corrected_signal, baseline = baseline(x, y, method="als", lam=1e5, p=0.01)
rampy.baseline.drPLS_baseline(y: ndarray, niter: int, lam: float, eta: float, ratio: float) ndarray

Doubly Reweighted Penalized Least Squares baseline fitting.

rampy.baseline.extract_signal(x: ndarray, y: ndarray, roi) ndarray

Extracts the signal from specified regions of interest (ROI) in the x-y data.

This function selects and extracts portions of the input x-y data based on the specified regions of interest (ROI) provided in roi. Each region is defined by a lower and upper bound along the x-axis.

Parameters:
  • x (ndarray) – The x-axis values (e.g., time, wavelength, or other independent variables).

  • y (ndarray) – The y-axis values corresponding to the x-axis values (e.g., signal intensity).

  • roi (ndarray or list of lists) –

    Regions of interest (ROI) where the signal should be extracted. Must be an n x 2 array or a list of lists, where n is the number of regions to extract. Each sublist or row should contain two elements:

    • The lower bound of the region (inclusive).

    • The upper bound of the region (inclusive).

    Example

    • Array: np.array([[10, 20], [50, 70]])

    • List: [[10, 20], [50, 70]]

Returns:

A 2-column array containing the extracted x-y signals from the specified regions.

The first column contains the x values, and the second column contains the corresponding y values.

Return type:

ndarray

Raises:

ValueError – If roi is not a valid n x 2 array or list of lists, or if any region in roi falls outside the range of x.

Notes

  • Overlapping regions in roi are not merged; they are processed as separate regions.

  • If no valid regions are found within roi, an empty array is returned.

Examples

Extracting signal from two regions in an x-y dataset:

>>> import numpy as np
>>> x = np.linspace(0, 100, 101)
>>> y = np.sin(x / 10) + np.random.normal(0, 0.1, size=x.size)
>>> roi = [[10, 20], [50, 70]]
>>> extracted_signal = extract_signal(x, y, roi)
>>> print(extracted_signal)
rampy.baseline.fit_baseline(x: ndarray, y: ndarray, roi: ndarray, yafit: ndarray, method: str, **kwargs) ndarray

Fit the baseline using the specified method.

rampy.baseline.get_portion_interest(x: ndarray, y: ndarray, roi) ndarray
rampy.baseline.rubberband_baseline(x: ndarray, y: ndarray) ndarray

Perform rubberband baseline correction.

Parameters:
  • x (ndarray) – The x-axis values.

  • y (ndarray) – The y values corresponding to x.

Returns:

baseline_fitted – The fitted baseline.

Return type:

ndarray

Raises:

ValueError – If the convex hull does not have enough points for interpolation.

rampy.baseline.standardize_data(x: ndarray, y: ndarray)

Standardize the data using sklearn’s StandardScaler.

rampy.baseline.validate_input(x: ndarray, y: ndarray, roi)

Validate the input arrays.

rampy.filters module

rampy.filters.smooth(x: ndarray, y: ndarray, method: str = 'whittaker', **kwargs) ndarray

Smooths the provided signal using various smoothing algorithms.

This function applies different smoothing techniques to the input signal, including Whittaker smoothing, Savitzky-Golay filtering, spline-based methods, and window-based filters. Each method is designed to reduce noise while preserving signal features.

Parameters:
  • x (np.ndarray) – The x-axis values (e.g., time or wavelength). For window-based methods, these values should be equally spaced.

  • y (np.ndarray) –

    The y-axis values to be smoothed. method (str): The smoothing method to apply. Default is “whittaker”.

    Available options:
    • ”whittaker”: Whittaker smoother (Eilers 2003).

    • ”savgol”: Savitzky-Golay filter.

    • ”GCVSmoothedNSpline”: Spline with generalized cross-validation.

    • ”MSESmoothedNSpline”: Spline with mean squared error criterion (requires the gcvspline library).

    • ”DOFSmoothedNSpline”: Spline with degrees of freedom criterion (requires the gcvspline library).

    • ”flat”: Moving average window.

    • ”hanning”, “hamming”, “bartlett”, “blackman”: Various window filters.

  • **kwargs

    Additional parameters specific to the chosen method. - window_length (int): Length of the filter window for window-based methods

    and Savitzky-Golay filter. Must be a positive odd integer. Default is 5.

    • polyorder (int): Polynomial order for Savitzky-Golay filter. Must be less than window_length. Default is 2.

    • Lambda (float): Smoothing parameter for the Whittaker filter. Higher values produce smoother fits. Default is (10^5).

    • d (int): Difference order in the Whittaker filter. Default is 2.

    • ese_y (float or np.ndarray): Errors associated with y values for spline methods. Default is 1.0.

Returns:

The smoothed signal sampled on x.

Return type:

np.ndarray

Raises:
  • ValueError – If the input vector is smaller than the window size or if an invalid smoothing method is specified.

  • ImportError – If gcvspline is not installed when using MSESmoothedNSpline or DOFSmoothedNSpline.

Notes

  • The “whittaker” method implements the perfect smoother described by Eilers (2003).

  • The “savgol” method uses scipy.signal.savgol_filter().

  • Spline methods require scipy >= 1.10.0 or the gcvspline package.

  • Window-based methods are implemented following the SciPy Cookbook.

Examples

Smooth a noisy signal using Savitzky-Golay filtering:

>>> import numpy as np
>>> x = np.linspace(0, 100, 101)
>>> y = np.sin(x / 10) + np.random.normal(0, 0.1, size=x.size)
>>> smoothed_signal = smooth(x, y, method="savgol", window_length=7, polyorder=3)

Apply Whittaker smoothing:

>>> smoothed_signal = smooth(x, y, method="whittaker", Lambda=1e6)

Use a moving average window:

>>> smoothed_signal = smooth(x, y, method="flat", window_length=5)
rampy.filters.spectrafilter(spectre: ndarray, filtertype: str, fq: int, numtaps: int, columns: ndarray) ndarray

Applies a Butterworth filter to specific columns of spectral data.

This function filters specific frequencies in the provided spectral data using a Butterworth filter. It supports low-pass, high-pass, band-pass, and band-stop filtering. The input spectra must be provided as an array where the first column represents the x-axis (e.g., wavelength or time), and subsequent columns represent the y-axis values of multiple spectra.

Parameters:
  • spectre (np.ndarray) – A 2D array of spectral data. The first column contains x-axis values (e.g., time or wavelength), and subsequent columns contain y-axis values for multiple spectra.

  • filtertype (str) – The type of filter to apply. Choose from: - ‘low’: Low-pass filter. - ‘high’: High-pass filter. - ‘bandpass’: Band-pass filter. - ‘bandstop’: Band-stop filter.

  • fq (float or np.ndarray) – The cutoff frequency/frequencies for the filter. For ‘low’ or ‘high’ filters, this is a single float value. For ‘bandpass’ or ‘bandstop’, this must be a sequence of two values [low_cutoff, high_cutoff].

  • numtaps (int) – The order of the Butterworth filter. Higher values result in sharper transitions but may introduce more ringing artifacts.

  • columns (np.ndarray) – An array specifying which columns in the input spectre to apply the filter to. Each value should correspond to a valid column index.

Returns:

A 2D array with the same shape as spectre, where the specified columns have been filtered.

Return type:

np.ndarray

Raises:

ValueError – If an invalid filtertype is provided or if fq is improperly specified for ‘bandpass’ or ‘bandstop’ filters.

Notes

  • The x-axis values (first column) are copied directly to the output without modification.

  • The cutoff frequencies (fq) are normalized based on the Nyquist frequency, which is calculated from the sampling rate inferred from the x-axis spacing.

  • The Butterworth filter is applied using zero-phase filtering (scipy.signal.filtfilt) to avoid phase distortion.

Examples

Apply a low-pass filter to a single spectrum:

>>> import numpy as np
>>> import rampy
>>> x = np.linspace(0, 10, 1000)  # Time axis
>>> y = np.sin(2 * np.pi * 5 * x) + np.sin(2 * np.pi * 50 * x)  # Signal with two frequencies
>>> spectre = np.column_stack((x, y))
>>> filtered_spectre = rampy.spectrafilter(spectre, filtertype='low', fq=10, numtaps=4, columns=np.array([1]))

Apply a band-pass filter to multiple spectra:

>>> spectre = np.column_stack((x, y, y * 0.5))  # Two spectra
>>> filtered_spectre = spectrafilter(spectre, filtertype='bandpass', fq=[5, 15], numtaps=4, columns=np.array([1, 2]))

References

rampy.filters.whittaker(y: ndarray, weights: ndarray = None, **kwargs) ndarray

Smooths a signal using the Whittaker smoother.

This function applies Whittaker smoothing to reduce noise in a signal while preserving its features. It uses penalized least squares optimization with a specified smoothing coefficient.

Parameters:
  • y (np.ndarray) – An array containing the values to smooth. The data should be equally spaced.

  • weights (np.ndarray, optional) – Weights for the signal. Regions with weight 0 will not

  • points. (be smoothed. Default is uniform weights across all) –

  • **kwargs

    Additional parameters for Whittaker smoothing. - Lambda (float): The smoothing coefficient; higher values produce smoother results.

    Default is (10^5).

Returns:

An array containing the smoothed values.

Return type:

np.ndarray

Raises:

ValueError – If y or weights are invalid.

Notes

  • This implementation follows the description provided by Eilers (2003).

References

Eilers, P.H.C., 2003. A Perfect Smoother. Anal. Chem. 75, 3631–3636.

Examples

Smooth a noisy signal using default weights and Lambda:

>>> import numpy as np
>>> y = np.sin(np.linspace(0, 10, 100)) + np.random.normal(0, 0.1, size=100)
>>> smoothed_signal = whittaker(y)

Apply custom weights and Lambda:

>>> weights = np.ones(100)
>>> weights[40:60] = 0  # Do not smooth this region
>>> smoothed_signal = whittaker(y, weights=weights, Lambda=1e6)

rampy.plotting module

rampy.plotting.plot_spectrum(x, y, fig=None, baselines=None, baseline_labels=None, smoothed_signals=None, smoothed_labels=None, label='Spectrum', color=None, roi=None, xaxis_title=None, yaxis_title=None, title=None)

Plots a spectrum with optional baselines, smoothed signals, and regions of interest (bir) using Plotly.

This function adds a single spectrum to a Plotly figure. It supports overlaying baselines, smoothed signals, and highlighting regions of interest (bir). If no figure is provided, a new Plotly figure is created.

Parameters:
  • x (array-like) – The x-axis values for the spectrum.

  • y (array-like) – The y-axis values for the spectrum.

  • fig (plotly.graph_objects.Figure, optional) – The Plotly figure to which the spectrum will be added. If None, a new figure will be created. Default is None.

  • baselines (list of array-like, optional) – List of baseline arrays to overlay on the spectrum. Default is None.

  • baseline_labels (list of str, optional) – Labels for each baseline. Default is None.

  • smoothed_signals (list of array-like, optional) – List of smoothed signal arrays to overlay on the spectrum. Default is None.

  • smoothed_labels (list of str, optional) – Labels for each smoothed signal. Default is None.

  • label (str, optional) – Label for the spectrum. Default is “Spectrum”.

  • color (str or tuple, optional) – Base color for the spectrum. If None, a default color will be chosen. Default is None.

  • bir (list or ndarray, optional) – Regions of interest to highlight on the plot. Should be a list or array of shape (n_regions, 2), where each row specifies the start and end of a region (e.g., [[x_start1, x_end1], [x_start2, x_end2]]). Default is None.

  • xaxis_title (str, optional) – Title for the x-axis. If None, no specific title is set. Default is None.

  • yaxis_title (str, optional) – Title for the y-axis. If None, no specific title is set. Default is None.

  • title (str, optional) – Title for the entire plot. If None, no title is set. Default is None.

Returns:

The updated Plotly figure object containing the plotted spectrum.

Return type:

plotly.graph_objects.Figure

Examples

Plot a simple spectrum:

>>> import numpy as np
>>> import plotly.graph_objects as go
>>> x = np.linspace(0, 10, 500)
>>> y = np.sin(x) + np.random.normal(0, 0.1, len(x))
>>> fig = plot_spectrum(
        x=x,
        y=y,
        label="Sample Spectrum",
        xaxis_title="Wavenumber (cm⁻¹)",
        yaxis_title="Intensity (a.u.)",
        title="Sample Spectrum Analysis"
    )
>>> fig.show()

Plot a spectrum with baseline and smoothed signal:

>>> baseline = 0.5 * np.ones_like(x)
>>> smoothed = np.convolve(y - baseline, np.ones(10)/10, mode='same')
>>> fig = plot_spectrum(
>>>        x=x,
>>>        y=y,
>>>        baselines=[baseline],
>>>        baseline_labels=["Constant Baseline"],
>>>        smoothed_signals=[smoothed],
>>>        smoothed_labels=["Moving Average"],
>>>        xaxis_title="Wavenumber (cm⁻¹)",
>>>        yaxis_title="Intensity (a.u.)"
>>>    )
>>> fig.show()

Notes

  • All input arrays are automatically converted to vectors using np.asarray(...).flatten() to ensure compatibility with Plotly.

  • Baselines are displayed as dashed lines with distinct colors.

  • Smoothed signals are displayed as dotted lines with distinct colors.

  • Regions of interest (BIR) are highlighted with semi-transparent yellow rectangles.

  • The function returns the figure object, allowing for further customization if needed.

rampy.functions module

rampy.functions.constant(x, a)

returns a constant value

Parameters:

x (1D array) –

Returns:

y – array filled with a values

Return type:

1D array

rampy.functions.difffull(x1: float, x2: float, t: float, C0: float, C1: float, D: float) float

Equation for the diffusion into a full slab, see Crank 1975

Here we assume the profil to have 2 surfaces of contact on each side

Parameters:
  • C0 (float) – the concentration in the core

  • C1 (float) – the concentration at the border

  • D (float) – the diffusion coefficient in log10 unit, m^2.s^-1

  • x2 (x1 and) – the profil lengths from beginning and end respectively, in meters

  • t (float) – time in seconds

rampy.functions.diffshort(x: ndarray, t: float, C0: float, C1: float, D: float) ndarray

1D equation for the diffusion into a semi-infinite slab, see Crank 1975

Parameters:
  • x (1D array) – the profil length in meters

  • t (float) – time in seconds

  • C0 (float) – the concentration in the core

  • C1 (float) – the concentration at the border

  • D (float) – the diffusion coefficient in log10 unit, m^2.s^-1

Returns:

Cx – concentration at x

Return type:

1D array

rampy.functions.funexp(x, a, b, c)

exponential baseline function

a*exp(b*(x-c))

rampy.functions.funlog(x, a, b, c, d)

log baseline function

a * ln(-b *(x-c)) - d*x**2

rampy.functions.gauss_lsq(params, x)

predicts a sum of gaussian peaks with parameters params

Parameters:
  • params (1D array) – an array of the parameters of the peaks. The number of peaks is assumed to be equal to len(params)/3. In this array, list intensities first, then all peak positions, then all peak half width at half maximum.

  • x (1D array) – x axis

Returns:

y – y values at position x

Return type:

1D array

rampy.functions.gauss_lsq_lfix(params, x)

predicts a sum of gaussian peaks with parameters params

Assumes that all peaks share the same HWHM.

Parameters:
  • params (1D array) – an array of the parameters of the peaks. The number of peaks is assumed to be equal to len(params)/3. In this array, list intensities first, then all peak positions, then the last element is the peaks’ half width at half maximum.

  • x (1D array) – x axis

Returns:

y – y values at position x

Return type:

1D array

rampy.functions.linear(x, a, b)

returns a + b*x

rampy.functions.linear0(x, a)

returns a*x

rampy.functions.multigaussian(x, params)

old attempt to have a multigaussian function, do not use. Will be removed soon.

rampy.functions.poly2(x, a, b, c)

returns a + b*x + c*x*x

rampy.maps module

class rampy.maps.maps(file, spectrometer_type='horiba', map_type='2D', despiking=False, neigh=4, threshold=3)

Bases: object

treat maps of Raman spectra

Parameters:
  • file (str) – filename, including path

  • spectrometer_type (str) – type of spectrometer, choose between “horiba” or “renishaw”, default: ‘horiba’

  • map_type (str) – type of map, choose between “2D” or “1D”, default: ‘2D’

area(y, region_to_investigate)

get the area under the curve in the region to investigate.

The area is calculated by trapezoidal integration, using np.trapz() Do not forget to smooth the signal if necessary prior to using this.

Parameters:
  • y (object intensities) – the intensities to consider. For instance, pass self.normalised for performing the calculation on normalised spectra.

  • region_to_investigate (1x2 array) – the x values of regions where the area will be measured

Returns:

self.A – maximum to make a nice plot

Return type:

ndarray

area_ratio(y, region_to_investigate)

get the area ratio between two regions of interest.

The areas are calculated by trapezoidal integration, using np.trapz() Do not forget to smooth the signals if necessary prior to using this.

Parameters:
  • y (object intensities) – the intensities to consider. For instance, pass self.normalised for performing the calculation on normalised spectra.

  • region_to_investigate (1x2 array) – the x values of regions where the areas will be measured. The two lines record the two regions of interest.

Returns:

self.A_ratio – Area ratio = area region 1 / area region 2

Return type:

ndarray

background(bir, method='poly', **kwargs)

correct a background from the initial signal I on a map using rampy.baseline

Parameters:
  • bir (ndarray) – arrays of the backgroudn interpolation regions.

  • method (string) – see rampy.baseline documentation for methods available. Default is polynomial

  • there. (All kwargs argument for rampy.baseline() will be forwarded and can be used) –

Return type:

Background and corrected spectra area available at self.background and self.I_corrected

centroid(y, region_to_investigate)

calculate the centroid in a given region of interest

Parameters:
  • y (object intensities) – the intensities to normalise. For instance, pass self.normalised for performing the calculation on normalised spectra.

  • region_to_investigate (1x2 array) – the x values of regions where the centroid will be measured

Returns:

self.centroid_position – centroid position for the map

Return type:

ndarray

intensity(y, region_to_investigate)

get the maximum intensity in the region to investigate.

The intensity maximum is estimated from a simple np.max() search. Do not forget to smooth the signal if necessary prior to using this.

Parameters:
  • y (object intensities) – the intensities to consider. For instance, pass self.normalised for performing the calculation on normalised spectra.

  • region_to_investigate (1x2 array) – the x values of regions where the intensity will be measured

Returns:

self.I_max – Intensity maximum

Return type:

ndarray

intensity_ratio(y, region_to_investigate)

get the intensity ratio between two regions of interest.

The intensity maxima are estimated from a simple np.max() search. Do not forget to smooth the signals if necessary prior to using this. :param y: the intensities to consider. For instance, pass self.normalised for performing the calculation on normalised spectra. :type y: object intensities :param region_to_investigate: the x values of regions where the intensity ratios will be measured. The two lines record the two regions of interest. :type region_to_investigate: 2x2 array

Returns:

self.I_ratio – Intensity ratio

Return type:

ndarray

normalise(y, method='intensity')

normalise the spectra to their max intensity, their area or min-max normalisation

This uses the internals of rampy.normalise.

Parameters:
  • y (object intensities) – the intensities to normalise. For instance, if you want to normalised the background corrected I, pass self.I_corrected.

  • method (string) – method used, choose between area, intensity, minmax

Return type:

The normalised spectra are available at self.I_normalised

smooth(y, method='whittaker', **kwargs)

uses the smooth function of rampy to smooth the signals :param y: the intensities to normalise. For instance, if you want to normalised the background corrected I, pass self.I_corrected. :type y: object intensities :param method: Method for smoothing the signal;

choose between savgol (Savitzky-Golay), GCVSmoothedNSpline, MSESmoothedNSpline, DOFSmoothedNSpline, whittaker, flat, hanning, hamming, bartlett, blackman.

Parameters:
  • window_length (int, optional) – The length of the filter window (i.e. the number of coefficients). window_length must be a positive odd integer.

  • polyorder (int, optional) – The order of the polynomial used to fit the samples. polyorder must be less than window_length.

  • Lambda (float, optional) – smoothing parameter of the Whittaker filter described in Eilers (2003). The higher the smoother the fit.

  • d (int, optional) – d parameter in Whittaker filter, see Eilers (2003).

  • ese_y (ndarray, optional) – errors associated with y (for the gcvspline algorithms)

Returns:

self.y_smoothed – the smoothed signal for the map

Return type:

ndarray

rampy.maps.peak(X, Y, lambdas, intensities, function, Xrange, amp, Xmean, sigma, y0, A)

to fit peaks in a map. Work in progress.

rampy.maps.read_horiba(file, map_type='2D')

read Horiba csv maps (1D, 2D)

Parameters:
  • file (str) – filename, including path

  • map_type (str) – 1D map (line) or 2D map, default: 2D

Returns:

  • X (m by n array) – X positions

  • Y (m by n array) – Y position

  • lambdas (n array) – Raman shift

  • intensities (m by n array) – Intensities

rampy.maps.read_renishaw(file)

read Renishaw csv maps

Parameters:

file (str) – filename, including path

Returns:

  • X (m by n array) – X positions

  • Y (m by n array) – Y position

  • lambdas (m by n array) – Raman shift

  • intensities (m by n array) – Intensities

rampy.mixing module

rampy.mixing.mixing_sp(y_fit: ndarray, ref1: ndarray, ref2: ndarray) ndarray

Mixes two reference spectra to match given experimental signals.

This function calculates the fractions of the first reference spectrum (ref1) in a linear combination of ref1 and ref2 that best matches the provided signals (y_fit). The calculation minimizes the sum of the least absolute values of the objective function: ( ext{obj} = sum left| y_{ ext{fit}} - (F_1 cdot ext{ref1} + (1 - F_1) cdot ext{ref2})

ight| ).

Args:
y_fit (np.ndarray): Array containing the experimental signals with shape (m, n),

where m is the number of data points and n is the number of experiments.

ref1 (np.ndarray): Array containing the first reference signal with shape (m,). ref2 (np.ndarray): Array containing the second reference signal with shape (m,).

Returns:

np.ndarray: Array of shape (n,) containing the fractions of ref1 in the mix. Values range between 0 and 1.

Notes:
  • The calculation is performed using cvxpy for optimization.

  • Ensure that y_fit, ref1, and ref2 have compatible dimensions.

Example:

>>> import numpy as np
>>> y_fit = np.array([[0.5, 0.6], [0.4, 0.5], [0.3, 0.4]])
>>> ref1 = np.array([0.5, 0.4, 0.3])
>>> ref2 = np.array([0.2, 0.3, 0.4])
>>> fractions = mixing_sp(y_fit, ref1, ref2)
>>> print(fractions)

rampy.ml_classification module

class rampy.ml_classification.mlclassificator(x, y, **kwargs)

Bases: object

Perform automatic classification of spectral data using scikit-learn machine learning algorithms.

This class supports various classification algorithms and allows customization of hyperparameters. It also handles scaling and splitting of training and testing datasets.

x

Training spectra organized in rows (1 row = one spectrum).

Type:

np.ndarray

y

Target labels for training data.

Type:

np.ndarray

X_test

Testing spectra organized in rows.

Type:

np.ndarray

y_test

Target labels for testing data.

Type:

np.ndarray

algorithm

Machine learning algorithm to use. Options: “Nearest Neighbors”, “Linear SVM”, “RBF SVM”, “Gaussian Process”, “Decision Tree”, “Random Forest”, “Neural Net”, “AdaBoost”, “Naive Bayes”, “QDA”.

Type:

str

scaling

Whether to scale the data during fitting and prediction.

Type:

bool

scaler

Type of scaler to use (“MinMaxScaler” or “StandardScaler”).

Type:

str

test_size

Fraction of the dataset to use as a testing dataset if X_test and y_test are not provided.

Type:

float

rand_state

Random seed for reproducibility. Default is 42.

Type:

int

params_

Hyperparameters for the selected algorithm.

Type:

dict

model

Scikit-learn model instance.

X_scaler

Scikit-learn scaler instance for X values.

fit(params_: dict = None)

Scale data and train or re-train the model with the specified algorithm.

This method initializes and trains the model if it hasn’t been trained yet. If a model already exists (from a previous fit), it reuses the existing model and optionally updates its hyperparameters.

Parameters:

params (dict, optional) – Hyperparameters for the selected algorithm. If provided, these parameters will override any previously set parameters.

Raises:

ValueError – If an invalid algorithm is specified or if scaling is inconsistent.

predict(X)

Predict target values using the trained model.

Parameters:

X (np.ndarray) – Samples to predict with shape (n_samples, n_features).

Returns:

Predicted target values with shape (n_samples,).

Return type:

np.ndarray

Notes

  • If scaling is enabled, input samples will be scaled before prediction.

Raises:

ValueError – If the model has not been fitted yet.

refit()

Re-train a model previously trained with fit()

scale_data()

Scale training and testing data.

rampy.ml_exploration module

class rampy.ml_exploration.mlexplorer(x, **kwargs)

Bases: object

use machine learning algorithms from scikit learn to explore spectroscopic datasets

Performs automatic scaling and train/test split before NMF or PCA fit.

x

Spectra; n_features = n_frequencies.

Type:

{array-like, sparse matrix}, shape = (n_samples, n_features)

X_test

spectra organised in rows (1 row = one spectrum) that you want to use as a testing dataset. THose spectra should not be present in the x (training) dataset. The spectra should share a common X axis.

Type:

{array-like, sparse matrix}, shape = (n_samples, n_features)

algorithm

“PCA”, “NMF”, default = “PCA”

Type:

String,

scaling

True or False. If True, data will be scaled prior to fitting (see below),

Type:

Bool

scaler

the type of scaling performed. Choose between MinMaxScaler or StandardScaler, see http://scikit-learn.org/stable/modules/preprocessing.html for details. Default = “MinMaxScaler”.

Type:

String

test_size

the fraction of the dataset to use as a testing dataset; only used if X_test and y_test are not provided.

Type:

float

rand_state

the random seed that is used for reproductibility of the results. Default = 42.

Type:

Float64

model

A Scikit Learn object model, see scikit learn library documentation.

Type:

Scikit learn model

Remarks
-------
For details on hyperparameters of each algorithms, please directly consult the documentation of SciKit Learn at
http
Type:

//scikit-learn.org/stable/

Results for machine learning algorithms can vary from run to run. A way to solve that is to fix the random_state.

Example

Given an array X of n samples by m frequencies, and Y an array of n x 1 concentrations

>>> explo = rampy.mlexplorer(X) # X is an array of signals built by mixing two partial components
>>> explo.algorithm = 'NMF' # using Non-Negative Matrix factorization
>>> explo.nb_compo = 2 # number of components to use
>>> explo.test_size = 0.3 # size of test set
>>> explo.scaler = "MinMax" # scaler
>>> explo.fit() # fitting!
>>> W = explo.model.transform(explo.X_train_sc) # getting the mixture array
>>> H = explo.X_scaler.inverse_transform(explo.model.components_) # components in the original space
>>> plt.plot(X,H.T) # plot the two components
fit()

Train the model with the indicated algorithm.

Do not forget to tune the hyperparameters.

predict(X)

Predict using the model.

Parameters:

X ({array-like, sparse matrix}, shape = (n_samples, n_features)) – Samples.

Returns:

  • C (array, shape = (n_samples,)) – Returns predicted values.

  • Remark

  • ——

  • if self.scaling == “yes”, scaling will be performed on the input X.

refit()

Train the model with the indicated algorithm.

Do not forget to tune the hyperparameters.

rampy.ml_regressor module

rampy.ml_regressor.chemical_splitting(Pandas_DataFrame, target, split_fraction=0.3, rand_state=42)

split datasets depending on their chemistry

Parameters:
  • Pandas_DataFrame (Pandas DataFrame) – The input DataFrame with in the first row the names of the different data compositions

  • label (string) – The target in the DataFrame according to which we will split the dataset

  • split_fraction (float, between 0 and 1) – This is the amount of splitting you want, in reference to the second output dataset (see OUTPUTS).

  • rand_state (float64) – the random seed that is used for reproductibility of the results. Default = 42.

Returns:

  • frame1 (Pandas DataFrame) – A DataSet with (1-split_fraction) datas from the input dataset with unique chemical composition / names

  • frame2 (Pandas DataFrame) – A DataSet with split_fraction datas from the input dataset with unique chemical composition / names

  • frame1_idx (ndarray) – Contains the indexes of the data picked in Pandas_DataFrame to construct frame1

  • frame2_idx (ndarray) – Contains the indexes of the data picked in Pandas_DataFrame to construct frame2

Notes

This function avoids the same chemical dataset to be found in different training/testing/validating datasets that are used in ML.

Indeed, it is worthless to put data from the same original dataset / with the same chemical composition in the training / testing / validating datasets. This creates a initial bias in the splitting process…

Another way of doing that would be to write:

>>> grouped = Pandas_DataFrame.groupby(by='label')
>>> k = [i for i in grouped.groups.keys()]
>>> k_train, k_valid = model_selection.train_test_split(np.array(k),test_size=0.40,random_state=100)
>>> train = Pandas_DataFrame.loc[Pandas_DataFrame['label'].isin(k_train)]
>>> valid = Pandas_DataFrame.loc[Pandas_DataFrame['label'].isin(k_valid)]

(results will vary slightly as variable k is sorted but not variable names in the function below)

class rampy.ml_regressor.mlregressor(x, y, **kwargs)

Bases: object

use machine learning algorithms from scikit learn to perform regression between spectra and an observed variable.

x

Spectra; n_features = n_frequencies.

Type:

{array-like, sparse matrix}, shape = (n_samples, n_features)

y

Returns predicted values.

Type:

array, shape = (n_samples,)

X_test

spectra organised in rows (1 row = one spectrum) that you want to use as a testing dataset. THose spectra should not be present in the x (training) dataset. The spectra should share a common X axis.

Type:

{array-like, sparse matrix}, shape = (n_samples, n_features)

y_test

the target that you want to use as a testing dataset. Those targets should not be present in the y (training) dataset.

Type:

array, shape = (n_samples,)

algorithm

“KernelRidge”, “SVM”, “LinearRegression”, “Lasso”, “ElasticNet”, “NeuralNet”, “BaggingNeuralNet”, default = “SVM”

Type:

String,

scaling

True or False. If True, data will be scaled during fitting and prediction with the requested scaler (see below),

Type:

Bool

scaler

the type of scaling performed. Choose between MinMaxScaler or StandardScaler, see http://scikit-learn.org/stable/modules/preprocessing.html for details. Default = “MinMaxScaler”.

Type:

String

test_size

the fraction of the dataset to use as a testing dataset; only used if X_test and y_test are not provided.

Type:

float

rand_state

the random seed that is used for reproductibility of the results. Default = 42.

Type:

Float64

param_kr

contain the values of the hyperparameters that should be provided to KernelRidge and GridSearch for the Kernel Ridge regression algorithm.

Type:

Dictionary

param_svm

containg the values of the hyperparameters that should be provided to SVM and GridSearch for the Support Vector regression algorithm.

Type:

Dictionary

param_neurons

contains the parameters for the Neural Network (MLPregressor model in sklearn). Default= dict(hidden_layer_sizes=(3,),solver = ‘lbfgs’,activation=’relu’,early_stopping=True)

Type:

Dictionary

param_bagging

contains the parameters for the BaggingRegressor sklearn function that uses a MLPregressor base method. Default= dict(n_estimators=100, max_samples=1.0, max_features=1.0, bootstrap=True,

bootstrap_features=False, oob_score=False, warm_start=False, n_jobs=1, random_state=rand_state, verbose=0)

Type:

Dictionary

prediction_train

the predicted target values for the training y dataset.

Type:

Array{Float64}

prediction_test

the predicted target values for the testing y_test dataset.

Type:

Array{Float64}

model

A Scikit Learn object model, see scikit learn library documentation.

Type:

Scikit learn model

X_scaler

A Scikit Learn scaler object for the x values.

Y_scaler

A Scikit Learn scaler object for the y values.

Example

Given an array X of n samples by m frequencies, and Y an array of n x 1 concentrations

>>> model = rampy.mlregressor(X,y)
>>> model.algorithm("SVM")
>>> model.user_kernel = 'poly'
>>> model.fit()
>>> y_new = model.predict(X_new)

Remarks

For details on hyperparameters of each algorithms, please directly consult the documentation of SciKit Learn at:

http://scikit-learn.org/stable/

For Support Vector and Kernel Ridge regressions, mlregressor performs a cross_validation search with using 5 KFold cross validators.

If the results are poor with Support Vector and Kernel Ridge regressions, you will have to tune the param_grid_kr or param_grid_svm dictionnary that records the hyperparameter space to investigate during the cross validation.

Results for machine learning algorithms can vary from run to run. A way to solve that is to fix the random_state. For neural nets, results from multiple neural nets (bagging technique) may also generalise better, such that it may be better to use the BaggingNeuralNet function.

fit()

Scale data and train the model with the indicated algorithm.

Do not forget to tune the hyperparameters.

Parameters:

algorithm (String,) – algorithm to use. Choose between “KernelRidge”, “SVM”, “LinearRegression”, “Lasso”, “ElasticNet”, “NeuralNet”, “BaggingNeuralNet”, default = “SVM”

predict(X)

Predict using the model.

Parameters:

X ({array-like, sparse matrix}, shape = (n_samples, n_features)) – Samples.

Returns:

  • C (array, shape = (n_samples,)) – Returns predicted values.

  • Remark

  • ——

  • if self.scaling == “yes”, scaling will be performed on the input X.

refit()

Re-train a model previously trained with fit()

rampy.peak_area module

rampy.peak_area.gaussianarea(amp, HWHM, **options)

returns the area of a Gaussian peak

Parameters:
  • amp (float or ndarray) – amplitude of the peak

  • HWHM (float or ndarray) – half-width at half-maximum

  • eseAmplitude (float or ndarray, optional) – standard deviation on amp; Default = None

  • eseHWHM (float or ndarray, optional) – standard deviation on HWHM; Default = None

Returns:

peak area esearea (float or ndarray): error on peak area; will be None if no errors on amp and HWHM were provided.

Return type:

area (float or ndarray)

rampy.peak_area.peakarea(shape: str, amp: float, HWHM: float, pos: float = None, a3: float = None, L_ratio: float = None, ese_amp: float = None, ese_HWHM: float = None) tuple[float, float | None]

Computes the area of a peak given its shape and parameters.

For Gaussian peaks, the area is calculated analytically. For other shapes (Lorentzian, pseudo-Voigt, Pearson VII), the area is calculated using numerical integration (trapezoidal rule).

Parameters:
  • shape (str) – The shape of the peak. Must be one of ‘gaussian’, ‘lorentzian’, ‘pseudovoigt’, or ‘pearson7’.

  • amp (float, list of floats, or np.ndarray) – Amplitude of the peak.

  • HWHM (float, list of floats, or np.ndarray) – Half-width at half-maximum of the peak.

  • pos (float, list of floats, np.ndarray, optional) – Peak position (required for non-Gaussian shapes). Default is None.

  • a3 (float, list of floats, np.ndarray, optional) – Shape parameter for Pearson VII peaks. Required if shape=’pearson7’.

  • L_ratio (float, list of floats, np.ndarray, optional) – Lorentzian component ratio for pseudo-Voigt peaks. Must be between 0 and 1. Required if shape=’pseudovoigt’. Default is None.

  • ese_amp (float, list of floats, np.ndarray, optional) – Standard deviation of the amplitude. Used to calculate error on area. Default is None.

  • ese_HWHM (float, list of floats, np.ndarray, optional) – Standard deviation of HWHM. Used to calculate error on area. Default is None.

Returns:

A tuple containing:
  • area (float): The computed area of the peak.

  • ese_area (float or None): The error on the computed area if ese_amp and ese_HWHM are provided; otherwise, None.

Return type:

tuple[float, float | None]

Raises:
  • ValueError – If required parameters are missing or invalid, or if the list of floats or the np.ndarray for the peak parameters have different sizes.

  • NotImplementedError – If an unsupported shape is specified.

Notes

  • For Gaussian peaks, the formula used is: ( ext{area} = sqrt{pi / ln(2)} cdot ext{amp} cdot ext{HWHM} ).

  • For other shapes, numerical integration is performed over a range of ( pm 10 cdot ext{HWHM} ).

Example

Compute the area of a Gaussian peak:

>>> amp = 10
>>> HWHM = 2
>>> area, ese_area = peakarea("gaussian", amp=amp, HWHM=HWHM)
>>> print(area)

Compute the area of a Lorentzian peak with numerical integration:

>>> amp = 10
>>> HWHM = 2
>>> pos = 5
>>> area, _ = peakarea("lorentzian", amp=amp, HWHM=HWHM, pos=pos)
>>> print(area)

Compute the area of several Lorentzian peaks with numerical integration (beware that in this case, the lists should have equal sizes):

>>> amp = [10., 5., 3.]
>>> HWHM = 2 # they share a common width
>>> pos = [2., 5., 5.5]
>>> area, _ = peakarea("lorentzian", amp=amp, HWHM=HWHM, pos=pos)
>>> print(area)

rampy.peak_shapes module

rampy.peak_shapes.gaussian(x: ndarray, amp, freq, HWHM) ndarray

Computes a Gaussian peak.

Parameters:
  • x (np.ndarray) – Positions at which the signal should be sampled.

  • amp (float or np.ndarray) – Amplitude of the Gaussian peak.

  • freq (float or np.ndarray) – Frequency/position of the Gaussian component.

  • HWHM (float or np.ndarray) – Half-width at half-maximum.

Returns:

The computed Gaussian signal.

Return type:

np.ndarray

Notes

Formula used: ( ext{amp} cdot exp(-log(2) cdot ((x - ext{freq}) / ext{HWHM})^2) ).

Examples

You can create a single peak like:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.gaussian(x, 10., 3., 1.0)

You can also create an array with several peaks in columns using arrays for the peak parameters:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.gaussian(x.reshape(-1, 1), np.array([[10, 10.]]), np.array([[3, 7.]]), np.array([[1., 1.]]))

In this case, y will be an array of shape (len(x), 2) with one peak per columns. Your peaks can even share parameters, here the HWHM:

>>> y = rp.gaussian(x.reshape(-1, 1), np.array([[10, 10.]]), np.array([[3, 7.]]), 2.0)

The composite signal is simply given by

>>> y_sum = np.sum(y, axis=1)
rampy.peak_shapes.lorentzian(x: ndarray, amp, freq, HWHM) ndarray

Computes a Lorentzian peak.

Parameters:
  • x (np.ndarray) – Positions at which the signal should be sampled.

  • amp (float or np.ndarray) – Amplitude of the Lorentzian peak.

  • freq (float or np.ndarray) – Frequency/position of the Lorentzian component.

  • HWHM (float or np.ndarray) – Half-width at half-maximum.

Returns:

The computed Lorentzian signal.

Return type:

np.ndarray

Notes

Formula used: ( ext{amp} / (1 + ((x - ext{freq}) / ext{HWHM})^2) ).

Examples

You can create a single peak like:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.lorentzian(x, 10., 3., 1.0)

You can also create an array with several peaks in columns using arrays for the peak parameters:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.lorentzian(x.reshape(-1, 1), np.array([[10, 10.]]), np.array([[3, 7.]]), np.array([[1., 1.]]))

In this case, y will be an array of shape (len(x), 2) with one peak per columns. Your peaks can even share parameters, here the HWHM:

>>> y = rp.lorentzian(x.reshape(-1, 1), np.array([[10, 10.]]), np.array([[3, 7.]]), 2.0)

The composite signal is simply given by

>>> y_sum = np.sum(y, axis=1)
rampy.peak_shapes.pearson7(x, a0, a1, a2, a3)

Computes a Pearson VII peak.

Parameters:
  • x (np.ndarray) – Positions at which the signal should be sampled.

  • a0 (float or np.ndarray) – Amplitude parameter.

  • a1 (float or np.ndarray) – Frequency/position parameter.

  • a2 (float or np.ndarray) – Width parameter.

  • a3 (float or np.ndarray) – Shape parameter.

Returns:

The computed Pearson VII signal.

Return type:

np.ndarray

Notes

Formula used: ( a_0 / ((1 + ((x - a_1) / a_2)^2 cdot (2^{(1/a_3)} - 1))^{a_3}) ).

Examples

You can create a single peak like:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.pearson7(x, 10., 3., 1.0, 0.5)

You can also create an array with several peaks in columns using arrays for the peak parameters:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.pearson7(x.reshape(-1, 1), np.array([[10, 10]]), np.array([[3, 7]]), np.array([[1., 1.]]), np.array([[0.5, 0.5]]))

In this case, y will be an array of shape (len(x), 2) with one peak per columns. Your peaks can even share parameters, here the a3:

>>> y = rp.pearson7(x.reshape(-1, 1), np.array([[10, 10]]), np.array([[3, 7]]), np.array([[1., 1.]]), 0.5)

The composite signal is simply given by

>>> y_sum = np.sum(y, axis=1)
rampy.peak_shapes.pseudovoigt(x: ndarray, amp, freq, HWHM, L_ratio) ndarray

Computes a pseudo-Voigt peak.

Parameters:
  • x (np.ndarray) – Positions at which the signal should be sampled. Can be a vector or array.

  • amp (float) – Amplitude of the pseudo-Voigt peak.

  • freq (float or np.ndarray) – Frequency/position of the Gaussian component.

  • HWHM (float or np.ndarray) – Half-width at half-maximum.

  • L_ratio (float) – Ratio of the Lorentzian component. Must be between 0 and 1 inclusive.

Returns:

The computed pseudo-Voigt signal.

Return type:

np.ndarray

Raises:

ValueError – If L_ratio is not between 0 and 1.

Notes

Formula used: ( (1 - L_{ ext{ratio}}) cdot ext{gaussian}(x, ext{amp}, ext{freq}, ext{HWHM}) +

L_{ ext{ratio}} cdot ext{lorentzian}(x, ext{amp}, ext{freq}, ext{HWHM}) ).

Examples

You can create a single peak like:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.pseudovoigt(x, 10., 3., 1.0, 0.5)

You can also create an array with several peaks in columns using arrays for the peak parameters:

>>> x = np.linspace(0, 10, 1000)
>>> y = rp.pseudovoigt(x.reshape(-1, 1), np.array([[10, 10]]), np.array([[3, 7]]), np.array([[1., 1.]]), np.array([[0.5, 0.5]]))

In this case, y will be an array of shape (len(x), 2) with one peak per columns. Your peaks can even share parameters, here the L_ratio:

>>> y = rp.pseudovoigt(x.reshape(-1, 1), np.array([[10, 10]]), np.array([[3, 7]]), np.array([[1., 1.]]), 0.5)

The composite signal is simply given by

>>> y_sum = np.sum(y, axis=1)

rampy.rameau module

rampy.rameau.DG2017_calibrate(dictio)

Fit a calibration by optimizing the K coefficient in the DG2017 method

Parameters:

dictio (dictionary) – dictionary with arrays named “feo”, “rws” and “water”.

Returns:

popt – The optimize a and b parameters of the equation K = a * [FeO wt%] + b.

Return type:

ndarray

rampy.rameau.DG2017_predict(dictio, a=0.096, b=0.663)

Calculate the K coefficient for the DG2017 method.

Parameters:
  • dictio (dict) – a dictionary with ndarrays named “feo” and “rws”

  • b (a and) – factors in the equation: K = a * [FeO wt%] + b; default values from Di Genova et al. (2017)

Returns:

H2O (wt %) – The water content of the glasses calculated as Rws * (a * [FeO wt%] + b)

Return type:

ndarray

rampy.rameau.LL2012_calibrate(dictio)

Fit a calibration line following equations (2) and (3) from Le Losq et al. (2012)

Parameters:

dictio – dictionary with arrays named “feo”, “rws” and “water”.

Returns:

A – The parameter in the equation (3) of Le Losq et al. (2012).

Return type:

float

rampy.rameau.LL2012_predict(dictio, A=0.007609)

Predict the water content using the equation (3) from Le Losq et al. (2012)

Parameters:

dictio (dict) – a dictionary with ndarray named “rws”

Returns:

The glass water contents in wt%

Return type:

H2O

rampy.rameau.fit_spectra(data_liste, method='LL2012', delim='\t', path_in='./raw/', laser=514.532, spline_coeff=0.001, poly_coeff=3)

Calculate the ratios of water and silicate signals from Raman spectra

Parameters:
  • data_liste (Pandas DataFrame) – Contains the list of spectra, see provided file as an example

  • method (string) – The used method. LL2012: Le Losq et al. (2012); DG2017: Di Genova et al. (2017). See references.

  • delim (string) – File delimiter. Use ‘ ‘ for tabulated text or ‘,’ for comma separated text.

  • path_in (string) – Path for the spectra

  • laser (float) – Laser line wavelength in nm

  • spline_coeff (float) – Smoothing coefficient for the spline baseline. An array of size len(data_liste) can be provided. Default = 0.001.

  • poly_coeff (int) – Polynomial coefficient for the polynomial baseline function. Default = 3 (DG2017 method). Set to 2 for Behrens et al. (2006) method.

Returns:

  • x (ndarray) – Common x axis.

  • y_all (ndarray) – All raw spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • y_all_corr (ndarray) – All corrected spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • y_all_base (ndarray) – All baselines for spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • rws (ndarray) – The ratio of the water integrated intensity over that of silicate signals.

  • rw (ndarray) – The integrated intensity of water signal.

  • rs (ndarray) – The integrated intensity of silicate signals.

Raises:

IOError – If method is not set to LL2012 or DG2017.

References

  1. Le Losq, D. R. Neuville, R. Moretti, J. Roux, Determination of water content in silicate glasses using Raman spectrometry: Implications for the study of explosive volcanism. American Mineralogist. 97, 779–790 (2012).

  2. Di Genova et al., Effect of iron and nanolites on Raman spectra of volcanic glasses: A reassessment of existing strategies to estimate the water content. Chemical Geology. 475, 76–86 (2017).

class rampy.rameau.rameau(data_liste, path_spectra='./raw/')

Bases: object

treat Raman spectra of glass to retrieve the glass water content

Parameters:
  • data_liste (Pandas dataframe) – A Pandas dataframe containing the data and various meta information.

  • path_spectra (string) – Path for the spectra. Default = ‘./raw/’

x

a 1D array (Nx1) containing the common x axis (wavelength) of the spectra.

Type:

ndarray

y

a NxM array (with M the number of spectra) containing the raw intensities of the spectra.

Type:

ndarray

y_corr

a NxM array (with M the number of spectra) containing the corrected intensities of the spectra.

Type:

ndarray

y_base

a NxM array (with M the number of spectra) containing the backgrounds of the spectra.

Type:

ndarray

rws

a 1D array (Nx1) containing the ratio between the integrated intensities of the water and silicate signals.

Type:

ndarray

rw

a 1D array (Nx1) containing the integrated intensities of the water signal.

Type:

ndarray

rs

a 1D array (Nx1) containing the integrated intensities of the silicate signals.

water

the known glass water content provided in data_liste (set to 0 if predicting for unknowns)

water_predicted

the predicted glass water content provided in data_liste (set to 0 if predicting for unknowns)

p

calibration coefficient(s) of the LL2012 or DG2017 method

Type:

ndarray

names

filenames indicated in the data_liste input

Type:

pandas dataframe

Notes

Uses either the LL2012 method (Le Losq et al., 2012), the DG2017 method (Di Genova et al., 2017) or the external calibration method (Thomas et al., 2008). See references.

In the LL2012 method, a cubic spline is fitted to the regions of interest provided in self.data_liste (see example). The spline is smoothed by the spline_coeff of the data_reduction method. The water content is calculated following eq. (3) of LL2012, with the A coefficient either provided or calculated by the method self.calibrate().

In the DG2017 method, a third-order polynomial is fitted to the spectra following the instructions of Di Genova et al. (2017). The water content is calculated as wt% H2O = Rws * (a * [FeO wt%] + b) with a and b the coefficients either provided or calculated by the method self.calibrate().

In the external calibration method, a cross-multiplication is computed using a reference spectrum and the proportionality relation between the water content of the glass and the water peak area of its spectrum.

References

LL2102: C. Le Losq, D. R. Neuville, R. Moretti, J. Roux, Determination of water content in silicate glasses using Raman spectrometry: Implications for the study of explosive volcanism. American Mineralogist. 97, 779–790 (2012). DG2017: D. Di Genova et al., Effect of iron and nanolites on Raman spectra of volcanic glasses: A reassessment of existing strategies to estimate the water content. Chemical Geology. 475, 76–86 (2017). External calibration: S.-M. Thomas, R. Thomas, P. Davidson, P. Reichart, M. Koch-Muller, G. Dollinger, Application of Raman Spectroscopy to Quantify Trace Water Concentrations in Glasses and Garnets. American Mineralogist 2008, 93 (10), 1550–1557.

calibrate(data_calib=None, method='LL2012', delim='\t', path_calib='./raw/', laser=488.0, spline_coeff=0.005, poly_coeff=3)

Fit a calibration by optimizing the K coefficient(s) on a chosen dataset of calibration spectra.

Parameters:
  • self (object) – rameau object with treated spectra (see data_reduction method)

  • data_calib (Pandas dataframe) – A Pandas dataframe containing the calibration data (same informations as self.data_liste but for calibration spectra). If data_calib is not given, the calibration is performed on self.data_liste.

  • method (string) – the method used for calibration; choose between “LL2012” or “DG2017”, default = “LL2012”.

  • delim (string) – File delimiter. Use ‘ ‘ for tabulated text or ‘,’ for comma separated text. Default = ‘ ‘.

  • path_calib (string) – Path for the calibration spectra. Default = ‘./raw/’.

  • laser (float) – Laser line wavelength in nm of the calibration spectra. Default = 488.0.

  • spline_coeff (float) – Smoothing coefficient for the spline baseline. An array of size len(data_liste) can be provided. Default = 0.005.

  • poly_coeff (int) – Polynomial coefficient for the polynomial baseline function. Default = 3 (DG2017 method; set to 2 for Behrens et al. (2006) method).

Returns:

  • popt (ndarray or float) – The optimized parameter(s); if method = “DG2017”, popt=np.array([a,b]), parameters of the equation K = a * [FeO wt%] + b. if method = “LL2017”, popt = A (float), with A parameter in the equation (3) of Le Losq et al. (2012).

  • self.x_calib (ndarray) – Common x axis.

  • self.y_calib (ndarray) – All raw calibration spectra from data_calib in an array of length len(x) and with as many column as spectra.

  • self.y_corr_calib (ndarray) – All corrected calibration spectra from data_calib in an array of length len(x) and with as many column as spectra.

  • self.y_base_calib (ndarray) – All baselines for calibration spectra from data_calib in an array of length len(x) and with as many column as spectra.

  • self.rws_calib (ndarray) – The ratio of the water integrated intensity over that of silicate signals for calibration spectra from data_calib.

  • self.rw_calib (ndarray) – The integrated intensity of water signal for calibration spectra from data_calib.

  • self.rs_calib (ndarray) – The integrated intensity of silicate signals for calibration spectra from data_calib.

data_reduction(method='LL2012', delim='\t', laser=514.532, spline_coeff=0.001, poly_coeff=3)

process Raman spectra of glass to calculate the Rws ratio

Parameters:
  • self (object) – a rameau object that has been initiated.

  • method (string) – The used method. LL2012: Le Losq et al. (2012); DG2017: Di Genova et al. (2017). See references. Default = “LL2012”.

  • delim (string) – File delimiter. Use ‘ ‘ for tabulated text or ‘,’ for comma separated text. Default = ‘ ‘.

  • laser (float) – Laser line wavelength in nm. Default = 514.532.

  • spline_coeff (float) – Smoothing coefficient for the spline baseline. An array of size len(data_liste) can be provided. Default = 0.001.

  • poly_coeff (int) – Polynomial coefficient for the polynomial baseline function. Default = 3 (DG2017 method; set to 2 for Behrens et al. (2006) method).

Returns:

  • self.x (ndarray) – Common x axis.

  • self.y_all (ndarray) – All raw spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • self.y_all_corr (ndarray) – All corrected spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • self.y_all_base (ndarray) – All baselines for spectra from data_liste in an array of length len(x) and with as many column as spectra.

  • self.rws (ndarray) – The ratio of the water integrated intensity over that of silicate signals.

  • self.rw (ndarray) – The integrated intensity of water signal.

  • self.rs (ndarray) – The integrated intensity of silicate signals.

external_calibration(path_ref='./raw/Standards/', roi=array([[2900, 3100], [3700, 3800]]), lb=3200, hb=3750, show_fig=False)

Predict water content using an external calibration and reference spectra. For this method, each spectrum from self.data_liste must have a reference spectrum with a known water concentration (filled in columns “Ref” and “Water Ref” respectively).

Parameters:
  • self (object) – a rameau object that has been initiated.

  • path_ref (string) – Path for reference spectra.

  • roi (array((2,2))) – Region of interest at the beginning and the end of the water peak to fit the baseline.

  • lb (int) – Lower limit for the water peak area calculation.

  • hb (int) – Upper limit for the water peak area calculation.

  • show_fig (boolean) – to show figures of water peak and baseline for each spectrum from self.data_list and its reference spectrum.

Returns:

  • wat_list (ndarray) – Array containing predicted water for spectra from self.data_liste.

  • areas_list (ndarray) – Array containing water peak area for spectra from self.data_liste.

  • areas_lref (ndarray) – Array containing water peak area for reference spectra.

names = []
p = None
predict(method='LL2012')

predict the water content from the Rws

Parameters:
  • self (object) – rameau object with treated spectra (see data_reduction method).

  • method (string) – the method used; choose between “LL2012” or “DG2017”, default = “LL2012”.

Returns:

H2O – The glass water contents in wt%

Return type:

array

rs = []
rw = []
rws = []
water = []
water_predicted = []
x = []
y = []
y_base = []
y_corr = []

rampy.spectranization module

rampy.spectranization.centroid(x: ndarray, y: ndarray, smoothing: bool = False, **kwargs) ndarray

Calculates the centroid of y signal(s).

The centroid is calculated as ( ext{centroid} = sum(y / sum(y) cdot x) ).

Parameters:
  • x (np.ndarray) – A 2D array of shape (m values, n samples) containing the x-axis values.

  • y (np.ndarray) – A 2D array of shape (m values, n samples) containing the y-axis values.

  • smoothing (bool) – Whether to smooth the signals before calculating centroids. If True, smoothing is applied using rampy.smooth with additional arguments passed via kwargs.

Returns:

A 1D array of centroids for each sample in y.

Return type:

np.ndarray

Example

>>> import numpy as np
>>> import rampy as rp
>>> x = np.linspace(0, 10, 100).reshape(-1, 1)
>>> y = rp.gaussian(x, 10., 50., 2.0)
>>> centroids = rp.centroid(x, y)
rampy.spectranization.convert_x_units(x: ndarray, laser_nm: float = 532.0, way: str = 'nm_to_cm-1') ndarray

Converts between nanometers and inverse centimeters for Raman spectroscopy.

Parameters:
  • x (np.ndarray) – Array of x-axis values to convert.

  • laser_nm (float) – Wavelength of the excitation laser in nanometers. Default is 532.0 nm.

  • way (str) – Conversion direction. Options are: - “nm_to_cm-1”: Convert from nanometers to inverse centimeters. - “cm-1_to_nm”: Convert from inverse centimeters to nanometers.

Returns:

Converted x-axis values.

Return type:

np.ndarray

Raises:

ValueError – If an invalid conversion direction is specified.

Example

Convert from nanometers to inverse centimeters:

>>> import rampy as rp
>>> x_nm = np.array([600.0])
>>> x_cm_1 = rp.convert_x_units(x_nm)

Convert from inverse centimeters to nanometers:

>>> x_cm_1 = np.array([1000.0])
>>> x_nm = rp.convert_x_units(x_cm_1, way="cm-1_to_nm")
rampy.spectranization.despiking(x: ndarray, y: ndarray, neigh: int = 4, threshold: int = 3) ndarray

Removes spikes from a 1D signal using a threshold-based approach.

This function identifies spikes in a signal by comparing local residuals to a threshold based on the root mean square error (RMSE). Spikes are replaced with the mean of neighboring points.

Parameters:
  • x (np.ndarray) – A 1D array containing the x-axis values of the signal.

  • y (np.ndarray) – A 1D array containing the y-axis values of the signal to despike.

  • neigh (int) – The number of neighboring points to use for calculating average values during despiking and for smoothing. Default is 4.

  • threshold (int) – The multiplier of RMSE used to identify spikes. Default is 3.

Returns:

A 1D array of the despiked signal.

Return type:

np.ndarray

Example

>>> import numpy as np
>>> import rampy as rp
>>> x = np.linspace(0, 10, 100)
>>> y = rp.gaussian(x, 10., 50., 2.0)
>>> y_despiked = rp.despiking(x, y)
rampy.spectranization.flipsp(sp: ndarray) ndarray

Sorts or flips a spectrum along its row dimension based on x-axis values.

Parameters:

sp (np.ndarray) – A 2D array where the first column contains x-axis values and subsequent columns contain y-values.

Returns:

The input array sorted in ascending order based on the first column.

Return type:

np.ndarray

Notes

  • Uses np.argsort to ensure sorting regardless of initial order.

Example

>>> import numpy as np
>>> import rampy as rp
>>> sp = np.array([[300, 30], [100, 10], [200, 20]])
>>> sorted_sp = rp.flipsp(sp)
rampy.spectranization.invcm_to_nm(x_inv_cm: ndarray, laser_nm: float = 532.0) ndarray

Converts Raman shifts from inverse centimeters (cm⁻¹) to nanometers (nm).

Parameters:
  • x_inv_cm (float or np.ndarray) – Raman shift(s) in inverse centimeters (cm⁻¹).

  • laser_nm (float) – Wavelength of the excitation laser in nanometers. Default is 532.0 nm.

Returns:

Wavelength(s) in nanometers corresponding to the Raman shifts.

Return type:

float or np.ndarray

Example

>>> import rampy as rp
>>> x_inv_cm = 1000.0
>>> wavelength_nm = rp.invcm_to_nm(x_inv_cm)
rampy.spectranization.nm_to_invcm(x: ndarray, laser_nm: float = 532.0) ndarray

Converts wavelengths from nanometers (nm) to Raman shifts in inverse centimeters (cm⁻¹).

Parameters:
  • x (float or np.ndarray) – Wavelength(s) in nanometers.

  • laser_nm (float) – Wavelength of the excitation laser in nanometers. Default is 532.0 nm.

Returns:

Raman shift(s) in inverse centimeters (cm⁻¹).

Return type:

float or np.ndarray

Example

>>> import rampy as rp
>>> wavelength_nm = 600
>>> shift_inv_cm = nm_to_invcm(wavelength_nm)
rampy.spectranization.normalise(y: ndarray, x: ndarray = None, method: str = 'intensity') ndarray

Normalizes the y signal(s) using specified methods.

This function normalizes the input y signal(s) based on the chosen method: by area under the curve, maximum intensity, or min-max scaling.

Parameters:
  • y (np.ndarray) – A 2D array of shape (m values, n samples) containing the y values to normalize.

  • x (np.ndarray, optional) – A 2D array of shape (m values, n samples) containing the x values corresponding to y. Required for area normalization. Default is None.

  • method (str) – The normalization method to use. Options are: - ‘area’: Normalize by the area under the curve. - ‘intensity’: Normalize by the maximum intensity. - ‘minmax’: Normalize using min-max scaling.

Returns:

A 2D array of normalized y signals with the same shape as the input y.

Return type:

np.ndarray

Raises:
  • ValueError – If x is not provided when using the ‘area’ normalization method.

  • NotImplementedError – If an invalid normalization method is specified.

Example

>>> import numpy as np
>>> import rampy as rp
>>> x = np.linspace(0, 10, 100)
>>> y = rp.gaussian(x, 10., 50., 2.0)
>>> y_norm = rp.normalise(y, x=x, method="area")
rampy.spectranization.resample(x: ndarray, y: ndarray, x_new: ndarray, **kwargs) ndarray

Resamples a y signal along new x-axis values using interpolation.

Parameters:
  • x (np.ndarray) – Original x-axis values.

  • y (np.ndarray) – Original y-axis values corresponding to x.

  • x_new (np.ndarray) – New x-axis values for resampling.

  • **kwargs

    Additional arguments passed to scipy.interpolate.interp1d.

    • kind (str or int): Type of interpolation (‘linear’, ‘cubic’, etc.). Default is ‘linear’.

    • bounds_error (bool): If True, raises an error when extrapolation is required. Default is False.

    • fill_value (float or str): Value used for out-of-bounds points. Default is NaN or “extrapolate”.

Returns:

Resampled y-values corresponding to x_new.

Return type:

np.ndarray

Example

>>> import numpy as np
>>> import rampy as rp
>>> original_x = np.array([100, 200, 300])
>>> original_y = np.array([10, 20, 30])
>>> new_x = np.linspace(100, 300, 5)
>>> resampled_y = rp.resample(original_x, original_y, new_x)
rampy.spectranization.shiftsp(sp: ndarray, shift: float) ndarray

Shifts the x-axis values of a spectrum by a given amount.

Parameters:
  • sp (np.ndarray) – A 2D array where the first column contains x-axis values (e.g., frequency or wavenumber) and subsequent columns contain y-values.

  • shift (float) – The amount by which to shift the x-axis values. Negative values shift left; positive values shift right.

Returns:

The input array with shifted x-axis values.

Return type:

np.ndarray

Example

>>> import numpy as np
>>> import rampy as rp
>>> sp = np.array([[100, 10], [200, 20], [300, 30]])
>>> shifted_sp = rp.shiftsp(sp, shift=50)
rampy.spectranization.spectraoffset(spectre, oft)

Applies vertical offsets to spectra.

This function offsets the y-values of each spectrum by a specified amount.

Parameters:
  • spectre (np.ndarray) – A 2D array where rows represent x-axis values and columns represent spectra (first column is x-axis).

  • oft (np.ndarray) – A 1D array specifying offset values for each spectrum.

Returns:

A 2D array with the same shape as spectre, where specified columns have been vertically offset.

Return type:

np.ndarray

Raises:

ValueError – If the length of oft does not match the number of spectra.

Example

>>> import numpy as np
>>> import rampy as rp
>>> spectre = np.array([[100, 10], [200, 20], [300, 30]])
>>> offsets = np.array([5])
>>> offset_spectra = rp.spectraoffset(spectre, offsets)
rampy.spectranization.spectrarray(name: ndarray, sh: int, sf: int, x: ndarray) ndarray

Constructs a unified array of spectra with a common x-axis.

This function reads spectral data from multiple files, resamples the y-values to match a common x-axis, and combines them into a single array.

Parameters:
  • name (np.ndarray) – Array of file names containing the spectra.

  • sh (int) – Number of header lines to skip in each file.

  • sf (int) – Number of footer lines to skip in each file.

  • x (np.ndarray) – The common x-axis values to which all spectra will be resampled.

Returns:

A 2D array where the first column contains the common x-axis values and subsequent columns contain the resampled y-values for each spectrum.

Return type:

np.ndarray

Raises:

ValueError – If any file contains invalid or missing data.

Notes

  • The function uses np.genfromtxt to read spectral data and resample for interpolation.

  • NaN values in the input data are automatically removed.

Example

>>> import numpy as np
>>> import rampy as rp
>>> filenames = np.array(["spectrum1.txt", "spectrum2.txt"])
>>> x_common = np.linspace(100, 2000, 500)
>>> spectra_array = rp.spectrarray(filenames, sh=1, sf=1, x=x_common)
rampy.spectranization.spectrataux(spectres: ndarray) ndarray

Calculates the rate of change for each frequency in a set of spectra.

This function fits a second-order polynomial to the y-values at each frequency across multiple spectra and calculates the polynomial coefficients.

Parameters:

spectres (np.ndarray) – A 2D array where the first column contains the common x-axis (frequencies) and subsequent columns contain y-values for multiple spectra.

Returns:

A 2D array where the first column contains the frequencies and subsequent columns contain the polynomial coefficients for each frequency.

Return type:

np.ndarray

Raises:

ValueError – If curve fitting fails or input data is invalid.

Example

>>> import numpy as np
>>> import rampy as rp
>>> spectres = np.array([[100, 10, 12], [200, 20, 22], [300, 30, 32]])
>>> taux = rp.spectrataux(spectres)

rampy.tlcorrection module

rampy.tlcorrection.tlcorrection(x: ndarray, y: ndarray, temperature: float, wavelength: float, **kwargs) tuple

Corrects Raman spectra for temperature and excitation line effects.

This function applies corrections to Raman spectra to account for temperature and laser excitation wavelength effects. It supports multiple correction equations and normalization methods, making it suitable for a variety of materials and experimental conditions.

Parameters:
  • x (np.ndarray) – Raman shifts in cm⁻¹.

  • y (np.ndarray) – Intensity values (e.g., counts).

  • temperature (float) – Temperature in °C.

  • wavelength (float) – Wavelength of the laser that excited the sample, in nm.

  • correction (str, optional) – The correction equation to use. Options are: - ‘long’: Default equation from Galeener and Sen (1978) with a (v_0^3) coefficient correction. - ‘galeener’: Original equation from Galeener and Sen (1978), based on Shuker and Gammon (1970). - ‘hehlen’: Equation from Hehlen et al. (2010), preserving the Boson peak signal. Default is ‘long’.

  • normalisation (str, optional) – Normalization method for the corrected data. Options are: - ‘intensity’: Normalize by maximum intensity. - ‘area’: Normalize by total area under the curve. - ‘no’: No normalization. Default is ‘area’.

  • density (float, optional) – Density of the studied material in kg/m³, used only with the ‘hehlen’ equation. Default is 2210.0 (density of silica).

Returns:

  • x (np.ndarray): Raman shift values after correction.

  • ycorr (np.ndarray): Corrected intensity values.

  • ese_corr (np.ndarray): Propagated errors calculated as (sqrt{y}) on raw intensities.

Return type:

tuple[np.ndarray, np.ndarray, np.ndarray]

Raises:

ValueError – If an invalid correction or normalization method is specified.

Notes

  • The ‘galeener’ equation is a modification of Shuker and Gammon’s formula to account for ((v_0 - v)^4) dependence of Raman intensity.

  • The ‘long’ equation includes a (v_0^3) coefficient to remove cubic meter dimensions, as used in several studies like Mysen et al. (1982).

  • The ‘hehlen’ equation avoids signal suppression below 500 cm⁻¹, preserving features like the Boson peak in glasses.

References

  • Galeener, F.L., & Sen, P.N. (1978). Theory of the first-order vibrational spectra of disordered solids. Physical Review B, 17(4), 1928–1933.

  • Hehlen, B. (2010). Inter-tetrahedra bond angle of permanently densified silicas extracted from their Raman spectra. Journal of Physics: Condensed Matter, 22(2), 025401.

  • Brooker, M.H., Nielsen, O.F., & Praestgaard, E. (1988). Assessment of correction procedures for reduction of Raman spectra. Journal of Raman Spectroscopy, 19(2), 71–78.

  • Mysen, B.O., Finger, L.W., Virgo, D., & Seifert, F.A. (1982). Curve-fitting of Raman spectra of silicate glasses. American Mineralogist, 67(7-8), 686–695.

  • Neuville, D.R., & Mysen, B.O. (1996). Role of aluminium in the silicate network: In situ high-temperature study of glasses and melts on the join SiO₂-NaAlO₂. Geochimica et Cosmochimica Acta, 60(9), 1727–1737.

  • Le Losq, C., Neuville, D.R., Moretti, R., & Roux, J. (2012). Determination of water content in silicate glasses using Raman spectrometry: Implications for the study of explosive volcanism. American Mineralogist, 97(5-6), 779–790.

  • Shuker, R., & Gammon, R.W. (1970). Raman-scattering selection-rule breaking and the density of states in amorphous materials. Physical Review Letters, 25(4), 222.

Examples

Correct a simple spectrum using default parameters:

>>> import numpy as np
>>> x = np.array([100, 200, 300])  # Raman shifts in cm⁻¹
>>> y = np.array([10, 20, 30])     # Intensity values
>>> temperature = 25.0             # Temperature in °C
>>> wavelength = 532.0             # Wavelength in nm
>>> x_corr, y_corr, ese_corr = correct_spectra(x, y, temperature, wavelength)

Use a specific correction equation and normalization method:

>>> x_corr, y_corr, ese_corr = correct_spectra(
        x=x,
        y=y,
        temperature=25,
        wavelength=532,
        correction='hehlen',
        normalisation='intensity',
        density=2500
    )