irfpy.ica.baseline

Created on Wed Dec 30 16:05:47 2015

@author: M Wieser

LICENCE This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with this program. If not, see <http://www.gnu.org/licenses/>

irfpy.ica.baseline.als(y, lam=1000000.0, p=0.1, itermax=10)[source]

Implements an Asymmetric Least Squares Smoothing baseline correction algorithm (P. Eilers, H. Boelens 2005)

Baseline Correction with Asymmetric Least Squares Smoothing based on https://github.com/vicngtor/BaySpecPlots

Baseline Correction with Asymmetric Least Squares Smoothing Paul H. C. Eilers and Hans F.M. Boelens October 21, 2005

Description from the original documentation:

Most baseline problems in instrumental methods are characterized by a smooth baseline and a superimposed signal that carries the analytical information: a series of peaks that are either all positive or all negative. We combine a smoother with asymmetric weighting of deviations from the (smooth) trend get an effective baseline estimator. It is easy to use, fast and keeps the analytical peak signal intact. No prior information about peak shapes or baseline (polynomial) is needed by the method. The performance is illustrated by simulation and applications to real data.

Inputs:
y:

input data (i.e. chromatogram of spectrum)

lam:

parameter that can be adjusted by user. The larger lambda is, the smoother the resulting background, z

p:

wheighting deviations. 0.5 = symmetric, <0.5: negative deviations are stronger suppressed

itermax:

number of iterations to perform

Output:

the fitted background vector

irfpy.ica.baseline.arpls(y, lam=10000.0, ratio=0.05, itermax=100)[source]

Baseline correction using asymmetrically reweighted penalized least squares smoothing Sung-June Baek, Aaron Park, Young-Jin Ahna and Jaebum Choo, Analyst, 2015, 140, 250 (2015)

Abstract

Baseline correction methods based on penalized least squares are successfully applied to various spectral analyses. The methods change the weights iteratively by estimating a baseline. If a signal is below a previously fitted baseline, large weight is given. On the other hand, no weight or small weight is given when a signal is above a fitted baseline as it could be assumed to be a part of the peak. As noise is distributed above the baseline as well as below the baseline, however, it is desirable to give the same or similar weights in either case. For the purpose, we propose a new weighting scheme based on the generalized logistic function. The proposed method estimates the noise level iteratively and adjusts the weights correspondingly. According to the experimental results with simulated spectra and measured Raman spectra, the proposed method outperforms the existing methods for baseline correction and peak height estimation.

Inputs:
y:

input data (i.e. chromatogram of spectrum)

lam:

parameter that can be adjusted by user. The larger lambda is, the smoother the resulting background, z

ratio:

wheighting deviations: 0 < ratio < 1, smaller values allow less negative values

itermax:

number of iterations to perform

Output:

the fitted background vector

irfpy.ica.baseline.WhittakerSmooth(x, w, lam, differences=1)[source]

Penalized least squares algorithm for background fitting

input
x:

input data (i.e. chromatogram of spectrum)

w:

binary masks (value of the mask is zero if a point belongs to peaks and one otherwise)

lam:

parameter that can be adjusted by user. The larger lambda is, the smoother the resulting background

differences:

integer indicating the order of the difference of penalties

output:

the fitted background vector

irfpy.ica.baseline.airpls(x, lam=100, porder=1, itermax=100)[source]

airpls.py Copyright 2014 Renato Lombardo - renato.lombardo@unipa.it Baseline correction using adaptive iteratively reweighted penalized least squares

This program is a translation in python of the R source code of airPLS version 2.0 by Yizeng Liang and Zhang Zhimin - https://code.google.com/p/airpls Reference: Z.-M. Zhang, S. Chen, and Y.-Z. Liang, Baseline correction using adaptive iteratively reweighted penalized least squares. Analyst 135 (5), 1138-1146 (2010).

Description from the original documentation:

Baseline drift always blurs or even swamps signals and deteriorates analytical results, particularly in multivariate analysis. It is necessary to correct baseline drift to perform further data analysis. Simple or modified polynomial fitting has been found to be effective in some extent. However, this method requires user intervention and prone to variability especially in low signal-to-noise ratio environments. The proposed adaptive iteratively reweighted Penalized Least Squares (airPLS) algorithm doesn’t require any user intervention and prior information, such as detected peaks. It iteratively changes weights of sum squares errors (SSE) between the fitted baseline and original signals, and the weights of SSE are obtained adaptively using between previously fitted baseline and original signals. This baseline estimator is general, fast and flexible in fitting baseline.

Adaptive iteratively reweighted penalized least squares for baseline fitting

input
x:

input data (i.e. chromatogram of spectrum)

lam:

parameter that can be adjusted by user. The larger lambda is, the smoother the resulting background, z

porder:

integer indicating the order of the difference of penalties

output:

the fitted background vector