PROCESSING Module
This module contains functions for the processing of NMR spectra, either in time domain or in frequency domain, and the transition between the two domains.
- klassez.processing.abc(ppm, data, n=5, lims=None, alpha=2.75, qfil=False, qfilp={'s': 10, 'u': 4.7})[source]
Automatic computation of a baseline for a spectrum using a thresholding-based method for the detection of the baseline-only region, followed by a weighted linear least squares optimization with a polynomion of degree n-1. The weights are computed on the absolute value of the first derivative of the spectrum. Set
qfil=Trueif there is a very intense solvent peak that would hamper the computation of the threshold.- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as
thr = alpha * np.std(np.gradient(data))- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- baseline: 1darray
Computed baseline
- klassez.processing.abc_v2(ppm, data, SFO1, n=5, lims=None, alpha=5, winsize=2, qfil=False, qfilp={'s': 10, 'u': 4.7})[source]
Automatic computation of a baseline for a spectrum using a thresholding-based method for the detection of the baseline-only region, followed by a weighted linear least squares optimization with a polynomion of degree n-1. Employs the same method for the detection of signal-free regions of
klassez.processing.apk(). Setqfil=Trueif there is a very intense solvent peak that would hamper the computation of the threshold.- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- SFO1float
Nucleus Larmor frequency /MHz
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as thr = alpha * np.std(np.gradient(data))
- winsizefloat
Minimum allowed window containing signals /Hz
- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- baseline1darray
Computed baseline
- klassez.processing.abs(ppm, data, n=5, lims=None, alpha=2.75, qfil=False, qfilp={'s': 10, 'u': 4.7})[source]
Computes the baseline correction on data using
klassez.processing.abc(), and gives back the subtracted spectrum. The imaginary part of the spectrum is reconstructed usingklassez.processing.hilbert().- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as
thr = alpha * np.std(np.gradient(data))- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- S1darray
Baseline-subtracted spectrum
- klassez.processing.abs2(ppm_f2, data, n=5, lims=None, alpha=2.75, qfil=False, qfilp={'s': 10, 'u': 4.7}, FnMODE='States-TPPI')[source]
Baseline correction for 2D datasets. Computes the baseline correction on
datausingklassez.processing.abc()for each row, and gives back the subtracted spectrum. The imaginary part of the spectrum is reconstructed using eitherklassez.processing.hilbert()orklassez.processing.hilbert2()depending onFnMODE.- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as
thr = alpha * np.std(np.gradient(data))- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- S2darray
Baseline-subtracted spectrum, either complex or hypercomplex
- klassez.processing.abs2_v2(ppm_f2, data, SFO1, n=5, lims=None, alpha=5, winsize=2, qfil=False, qfilp={'s': 10, 'u': 4.7}, FnMODE='States-TPPI')[source]
Baseline correction for 2D datasets, alternative version. Computes the baseline correction on
datausingklassez.processing.abc_v2()for each row, and gives back the subtracted spectrum. The imaginary part of the spectrum is reconstructed using eitherklassez.processing.hilbert()orklassez.processing.hilbert2()depending onFnMODE.- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as
thr = alpha * np.std(np.gradient(data))- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- S2darray
Baseline-subtracted spectrum, either complex or hypercomplex
- klassez.processing.abs_v2(ppm, data, SFO1, n=5, lims=None, alpha=5, winsize=2, qfil=False, qfilp={'s': 10, 'u': 4.7})[source]
Computes the baseline correction on data using
klassez.processing.abc_v2(), and gives back the subtracted spectrum. The imaginary part of the spectrum is reconstructed usingklassez.processing.hilbert().- Parameters:
- ppm1darray
PPM scale of the spectrum
- data1darray
The spectrum to baseline-correct
- nint
Number of coefficients of the polynomial baseline
- limstuple or None
Limits for the region on which to compute the baseline, in ppm
- alphafloat
The threshold will be set as
thr = alpha * np.std(np.gradient(data))- winsizefloat
Minimum allowed window containing signals /Hz
- qfilbool
Choose whether to apply a filter on the solvent region (True) or not (False)
- qfilpdict
Parameters to be used to compute the filter if
qfil=True. Keys:'u'= center of the filter in ppm's'= width of the filter in Hz
- Returns:
- S1darray
Baseline-subtracted spectrum
- klassez.processing.acme(data, m=1, a=5e-05)[source]
Automated phase Correction based on Minimization of Entropy. This algorithm allows for automatic phase correction by minimizing the entropy of the m-th derivative of the spectrum, as explained in detail by L. Chen et. al..
Defined the entropy of h as:
\[S = - \sum_j h[j] \ln( h[j] )\]and
\[h = \frac { | R[j]^{(m)} | }{ \sum_j | R[j]^{(m)} | }\]where
\[R = \Re\{ d \, e^{-i \phi} \}\]and \(R^{(m)}\) is the m-th derivative of R, the objective function to minimize is:
\[S + P(R)\]where P(R) is a penalty function for negative values of the spectrum.
The phase correction is applied using
klassez.processing.ps(). The valuesp0andp1are fitted using Nelder-Mead algorithm.- Parameters:
- data1darray
Spectrum to be phased, complex
- mint
Order of the derivative to be computed
- afloat
Weighting factor for the penalty function
- Returns:
- p0ffloat
Fitted zero-order phase correction, in degrees
- p1ffloat
Fitted first-order phase correction, in degrees
- klassez.processing.align(ppm_scale, data, lims, u_off=0.5, ref_idx=0)[source]
Performs the calibration of a pseudo-2D experiment by circular-shifting the spectra of an appropriate amount. The target function aims to minimize the superimposition between a reference spectrum and the others using a brute-force method.
- Parameters:
- ppm_scale1darray
ppm scale of the spectrum to calibrate
- data2darray
Complex-valued spectrum
- limstuple
(ppm sx, ppm dx) of the calibration region
- u_offfloat
Maximum offset for the circular shift, in ppm
- ref_idxint
Index of the spectrum to be used as reference
- Returns:
- data_roll2darray
Calibrated data
- u_cal_ppmlist
Correction for the ppm scale of each experiment
- klassez.processing.align_old(ppm_scale, data, lims, u_off=0.5, ref_idx=0)[source]
Performs the calibration of a pseudo-2D experiment by circular-shifting the spectra of an appropriate amount. The target function aims to minimize the superimposition between a reference spectrum and the others using a brute-force method.
Error
Old function!! Legacy
- Parameters:
- ppm_scale1darray
ppm scale of the spectrum to calibrate
- data2darray
Complex-valued spectrum
- limstuple
(ppm sx, ppm dx) of the calibration region
- u_offfloat
Maximum offset for the circular shift, in ppm
- ref_idxint
Index of the spectrum to be used as reference
- Returns:
- data_roll2darray
Calibrated data
- u_callist
Number of point of which the spectra have been circular-shifted
- u_cal_ppmlist
Correction for the ppm scale of each experiment
- klassez.processing.apk(ppm, data, SFO1, alpha=3, winsize=50, ap1=True, seethrough=False)[source]
Performs automatic phase correction.
The algorithm starts with the computation of a mask to separate signal from baseline-only regions. This is done via an iterative thresholding, i.e. a point is “signal” if the first derivative of the spectrum in that point is higher than its standard deviation by
alphatimes:\[d'[k] > \alpha \, std(d') \implies k \in \text{signal region}\]The selection is further refined by repeating the same procedure on the original data. Then, the regions separated by less than winsize are joined together, and the presence of actual peaks in the region is checked with a pick-picker.
At this point, each region is phased independently with only phase 0. The phase angle is tested in a brute-force manner. The cost function minimizes the area below the straight line that connects the borders of the window. The first-order phase correction is calculated with a weighted linear regression, where the weights are the integrals of the magnitude of each region.
Tip
The choice of
alphaandwinsizecan be important for the outcome. Higheralphavalues make the detection of peaks more stringent -> for spectra with high SNR and sharp peaks a suitable value is 4-6. Decreasingwinsizemakes the algorithm to estimate more regions.- Parameters:
- ppm1darray
ppm scale of the spectrum
- data1darray
Spectrum
- SFO1float
Nucleus’ Larmor frequency /MHz
- alphafloat
Factor that multiplies the std of the spectrum to set the threshold
- winsizefloat
Minimum size of the window that can contain peaks /Hz
- ap1bool
True to adjust both zero and first order, False for only phase zero
- seethroughbool
If True, draws a series of diagnostic figures to see what the algorithm is doing
- Returns:
- datap1darray
Phased data
- valuestuple
Found values
(p0, p1)
- klassez.processing.apodf(size, wf)[source]
Generates a function to be used as apodization function on the basis of the
wfdictionary. The behavior is controlled by themodekey. Eachmodereads the attributes of the corresponding positions, e.g.:mode: em, qsin, qsin lb: 5, 0, 0 ssb: 0, 2, 3
will compute the function:
processing.em(lb=5) * processing.qsin(ssb=2) * processing.qsin(ssb=3)
- Parameters:
- sizetuple
Dimension of data to be windowed
- wfdict
Dictionary of window functions modes and parameters
- Returns:
- apod_funcnp.ndarray
Custom apodization function of dimension
size
- klassez.processing.blp(data, pred=1, order=8)[source]
Applies backward linear prediction by calling
klassez.processing.lp()withmode='b'.- Parameters:
- data1darray
FID to be linear-predicted
- predint
Number of points to predict
- orderint
Number of coefficients to use for the prediction
- Returns:
- lpdata1darray
FID with linear prediction applied.
See also
- klassez.processing.cadzow(data, n, nc, print_head=True)[source]
Performs Cadzow denoising on data, which is a 1D array of
Npoints. The algorithm works as follows:Transform data in a Hankel matrix
Hof dimensions(N-n, n)Make SVD on \(H = U S V\)
Keep only the first
ncsingular values, and put all the rest to 0 (S -> S’)Rebuild \(H' = U S' V\)
Average the antidiagonals to rebuild the Hankel-type structure, then make 1D array
Set
print_head=Trueto display the fancy heading.- Parameters:
- data1darray
Input data
- nint
Number of columns of the Hankel matrix.
- ncint
Number of singular values to keep.
- print_headbool
Set it to True to display the fancy heading.
- Returns:
- datap1darray
Denoised data
- klassez.processing.cadzow_2D(data, n, nc, i=True, f=0.005, itermax=100, print_time=True)[source]
Performs the Cadzow denoising method on a 2D spectrum, one transient at the time. This function calls either Cadzow or iterCadzow, depending on the parameter
i: True foriterCadzow, False for normalcadzow.- Parameters:
- data2darray
Input data
- nint
Number of columns of the Hankel matrix.
- ncint
Number of singular values to keep.
- ibool
Calls
processing.cadzowifi=False, orprocessing.iterCadzowifi=True.- itermaxint
Maximum number of iterations allowed.
- ffloat
Factor for the arrest criterion.
- print_timebool
Set it to True to display the time spent.
- Returns:
- datap2darray
Denoised data
- klassez.processing.calc_nc(data, s_n)[source]
Calculates the optimal number of components to be used for either MCR or SVD, given the standard deviation of the noise. The threshold value is calculated as stated in Theorem 1 of this article .
- Parameters:
- data2darray
Input data
- s_nfloat
Noise standard deviation
- Returns:
- n_cint
Number of components
- klassez.processing.calibration(ppmscale, S, ref=None)[source]
Scroll the ppm scale of spectrum to make calibration. The interface offers two guidelines: the red one, labelled ‘reference signal’ remains fixed, whereas the green one (‘calibration value’) moves with the ppm scale. The ideal calibration procedure consists in placing the red line on the signal you want to use as reference, and the green line on the ppm value that the reference signal must assume in the calibrated spectrum. Then, scroll with the mouse until the two lines are superimposed. You can use another spectrum as reference.
- Parameters:
- ppmscale1darray
The ppm scale to be calibrated
- S1darray
The spectrum to calibrate
- reflist of 1darray or Spectrum_1D object
Reference spectrum to be used for calibration. If list,
[ppm scale, spectrum]
- Returns:
- offsetfloat
Difference between original scale and new scale. This must be summed up to the original ppm scale to calibrate the spectrum.
- klassez.processing.convdta(data, grpdly=0, scaling=1)[source]
Removes the digital filtering to obtain a spectrum similar to the command CONVDTA performed by TopSpin. However, they will differ a little bit because of the digitization. These differences are not invisible to human’s eye.
Note
Since TopSpin version 4.0 the algorithm changed somehow, and the result is different.
- Parameters:
- datandarray
FID with digital filter
- grpdlyint
Number of points that the digital filter consists of. Key $GRPDLY in acqus file
- scalingfloat
Scaling factor of the resulting FID. Needed to match TopSpin’s intensities.
- Returns:
- data_inndarray
FID without the digital filter. It will have
grpdlypoints less thandata.
- klassez.processing.convolve(in1, in2)[source]
Perform the convolution of the two array by multiplying their inverse Fourier transform. The two arrays must have the same dimension.
- Parameters:
- in1ndarray
First array
- in2ndarray
Second array
- Returns:
- cnvndarray
Convolved array
- klassez.processing.eae(data)[source]
Shuffles data if the spectrum is acquired with
FnMODE = 'Echo-Antiecho'. NOTE: introduces -90° phase shift in F1, to be corrected after the processingpdata = np.zeros_like(data) pdata[::2] = (data[::2].real - data[1::2].real) + 1j*(data[::2].imag - data[1::2].imag) pdata[1::2] = -(data[::2].imag + data[1::2].imag) + 1j*(data[::2].real + data[1::2].real)
- Parameters:
- data2darray
FID in echo-antiecho format
- Returns:
- pdata2darray
FID in States-TPPI format
- klassez.processing.em(data, lb, sw)[source]
Exponential apodization
\[a(x) = \exp\biggl[-\pi \frac{LB}{2\,SW} x \biggr]\]- Parameters:
- data: ndarray
Input data
- lb: float
Lorentzian broadening. It should be positive.
- sw: float
Spectral width /Hz
- Returns:
- datap: ndarray
Apodized data
See also
- klassez.processing.extend_taq(old_taq, newsize=None)[source]
Extend the acquisition timescale to a longer size, using the same dwell time
- Parameters:
- old_taq1darray
Old timescale
- newsizeint
New size of acqusition timescale, in points
- Returns:
- new_taq1darray
Extended timescale
- klassez.processing.fp(data, wf=None, zf=None, fcor=0.5, tdeff=0)[source]
Performs the full processing of a 1D NMR FID (
data). Also works for pseudo-2D: only applies the processing on the last dimension.- Parameters:
- datandarray
Input data
- wfdict
{'mode': function to be used, 'parameters': different from each function}- zfint
final size of spectrum
- fcorfloat
weighting factor for the FID first point
- tdeffint
number of points of the FID to be used for the processing.
- Returns:
- datapndarray
Processed data
- klassez.processing.ft(data0, alt=False, fcor=0.5)[source]
Fourier transform in NMR sense. This means it returns the reversed spectrum.
- Parameters:
- data0ndarray
Array to Fourier-transform
- altbool
negates the sign of the odd points, then take their complex conjugate. Required for States-TPPI processing.
- fcorfloat
weighting factor for FID 1st point. Default value (0.5) prevents baseline offset
- Returns:
- dataftndarray
Transformed data
- klassez.processing.gm(data, lb, gb, gc, sw)[source]
Gaussian apodization. The parameter
lbcontrols the sharpening factor of a rising exponential, and behaves exactly as inprocessing.em. In contrast,gbcontrols the gaussian decay factor.gcmoves the central point of the gaussian filter.\[a(x) = \exp\biggl\{ -\pi \frac{LB}{SW} x - \biggl[ \pi \frac{GB}{SW} \biggl(GC (N-1) -x \biggr) \biggr]^2 \biggr\}\]Apply this function VERY CAREFULLY. Choose the right values through the interactive processing.
- Parameters:
- datandarray
Input data
- lbfloat
Lorentzian sharpening /Hz. It should be negative.
- gbfloat
Gaussian broadening. It should be positive.
- gcfloat
Gaussian center, relatively to the FID length: \(0 \leq g_c \leq 1\)
- swfloat
Spectral width /Hz
- Returns:
- pdatandarray
Processed data
See also
- klassez.processing.gmb(data, lb, gb, sw)[source]
Bruker-style Gaussian apodization.
\[a(x) = \exp\biggl[ -N \frac{\pi}{SW} \biggl( LB\,x - \frac{LB}{2\,GB}x^2 \biggr) \biggr]\]Apply this function VERY CAREFULLY. Choose the right values through the interactive processing.
- Parameters:
- datandarray
Input data
- lbfloat
Lorentzian sharpening /Hz. It should be negative.
- gbfloat
Gaussian broadening. It should be positive.
- swfloat
Spectral width /Hz
- Returns:
- pdatandarray
Processed data
See also
- klassez.processing.hilbert(f, axis=-1)[source]
Computes the Hilbert transform of real vector
fin order to retrieve its imaginary part. The algorithm computes the convolution by means of FT, as follows: #. make IFT off:a#. computeh = [1j for x in range(N) if x < N/2 else -1j]#. Computeb = h * a#. Buildd = a + 1j*b#. make FT ofd:F#. replace the real part ofFwithfImportant
Make sure that the original spectrum was zero-filled to at least twice the original size of the FID.
- Parameters:
- fndarray
Array of which you want to compute the imaginary part
- axisint
Axis along which to compute the Hilbert transform. Default is -1 (last axis).
- Returns:
- f_cplxndarray
Complex version of
f
See also
- klassez.processing.hilbert2(data)[source]
Retrieve the imaginary parts of a hypercomplex dataset, when you only have the rr part.
Given
Ht = klassez.processing.hilbert:rr = Ht(rr).real ir = Ht(rr).imag ri = - Ht(rr.T).T.imag ii = Ht(ri).imag
- Parameters:
- data2darray
rr part
- Returns:
- rr2darray
Real part in f2, real part in f1
- ir2darray
Imaginary part in f2, real part in f1
- ri2darray
Real part in f2, imaginary part in f1
- ii2darray
Imaginary part in f2, imaginary part in f1
- klassez.processing.ift(data0, alt=False, fcor=0.5)[source]
Inverse Fourier transform in NMR sense. This means that the input dataset is reversed before to do iFT.
- Parameters:
- data0ndarray
Array to Fourier-transform
- altbool
negates the sign of the odd points, then take their complex conjugate. Required for States-TPPI processing.
- fcorfloat
weighting factor for FID 1st point. Default value (0.5) prevents baseline offset
- Returns:
- dataftndarray
Transformed data
- klassez.processing.integral(fx, x=None, dx=None, lims=None, use_bas=False)[source]
Calculates the primitive of
fx. Iffxis a multidimensional array, the integrals are computed along the last dimension.- Parameters:
- fxndarray
Function (array) to integrate
- x1darray or None
Independent variable. Determines the integration step. If None, it is the point scale
- dxfloat or None
Integration step. If
None, computes it from the resolution ofx- limstuple or None
Integration range in the
xscale. If None, the whole function is integrated.- use_basbool
Subtracts the straight line that connects the limit window before the integration (
True) or not (False)
- Returns:
- Fxndarray
Integrated function.
See also
- klassez.processing.integrate(fx, x=None, dx=None, lims=None, use_bas=False)[source]
Calculates the definite integral of
fxasI = F[-1] - F[0](basically it applies the fundamental theorem of calculus). Iffxis a multidimensional array, the integrals are computed along the last dimension.- Parameters:
- fxndarray
Function (array) to integrate
- x1darray or None
Independent variable. Can determine the integration step. If
None, it is the point scale.- dxfloat or None
Integration step. If
None, computes it from the resolution ofx- limstuple or None
Integration range according to
x. IfNone, the whole function is integrated.- use_basbool
Subtracts the straight line that connects the limit window before the integration (
True) or not (False)
- Returns:
- integfloat
Integrated function.
See also
- klassez.processing.inv_convolve(in1, in2)[source]
Perform the inverse-convolution of the two array by dividing their inverse Fourier transform. The two arrays must have the same dimension.
Important
This operation involves a division!!! Might give unexpected and unpleasant results!
- Parameters:
- in1ndarray
First array
- in2ndarray
Second array
- Returns:
- cnvndarray
Deconvolved array
- klassez.processing.inv_fp(data, wf=None, size=None, fcor=0.5)[source]
Performs the full inverse processing of a 1D NMR spectrum (
data).- Parameters:
- data1darray
Spectrum
- wfdict
{'mode': function to be used, 'parameters': different from each function}- sizeint
initial size of the FID
- fcorfloat
weighting factor for the FID first point
- Returns:
- pdata1darray
FID
- klassez.processing.inv_xfb(data, wf=[None, None], size=(None, None), fcor=[0.5, 0.5], FnMODE='States-TPPI')[source]
Reverts the full processing of a 2D NMR FID (data).
- Parameters:
- data2darray
Input data, complex or hypercomplex
- wflist of dict
list of two entries [F1, F2]. Each entry is a dictionary of window functions
- sizelist of int
Initial size of FID
- fcorlist of float
first fid point weighting factor [F1, F2]
- FnMODEstr
Acquisition mode in F1
- Returns:
- data2darray
Processed data
- klassez.processing.iterCadzow(data, n, nc, itermax=100, f=0.005, print_head=True, print_time=True)[source]
Performs Cadzow denoising on data, which is a 1D array of
Npoints, in an iterative manner. The algorithm works as follows:Transform data in a Hankel matrix H of dimensions
(N-n, n)Make SVD on \(H = U S V\)
Keep only the first
ncsingular values, and put all the rest to 0 (S -> S’)Rebuild \(H' = U S' V\)
Average the antidiagonals to rebuild the Hankel-type structure, then make 1D array
Check arrest criterion: if it is not reached, go to 1, else exit.
The arrest criterion is, at the k-th step:
\[\frac{S_{(k-1)}[nc-1] }{ S_{(k-1)}[0] } - \frac{S_{(k)}[nc-1] }{ S_{(k)}[0] } < f \frac{S_{(0)}[nc-1] }{ S_{(0)}[0] }\]- Parameters:
- data1darray
Data to be processed
- nint
Number of columns of the Hankel matrix
- ncint
Number of singular values to preserve
- itermaxint
max number of iterations allowed
- ffloat
factor that appears in the arrest criterion
- print_timebool
set it to True to show the time it took
- print_headbool
set it to True to display the fancy heading.
- Returns:
- datap1darray
Denoised data
See also
- klassez.processing.lp(data, pred=1, order=8, mode='b')[source]
Apply linear prediction on the dataset. This method solves the linear system
\[D a = d\]where a is the array of linear prediction coefficients.
- Parameters:
- data1darray
FID to be linear-predicted
- predint
Number of points to predict
- orderint
Number of coefficients to use for the prediction
- modestr
'f'for forward linear prediction,'b'for backward linear prediction
- Returns:
- newdata1darray
FID with linear prediction applied.
- klassez.processing.lrd(data, nc)[source]
Denoising method based on Low-Rank Decomposition. The algorithm performs a singular value decomposition on data, then keeps only the first
ncsingular values while setting all the others to 0. Finally, rebuilds the data matrix using the modified singular values.- Parameters:
- data2darray
Data to be denoised
- ncint
Number of components, i.e. number of singular values to keep
- Returns:
- data_out2darray
Denoised data
- klassez.processing.make_scale(size, dw, rev=True)[source]
Computes the frequency scale of the NMR spectrum, given the number of points and the employed dwell time (the REAL one, not the TopSpin one!).
rev=Trueis required for the correct frequency arrangement in the NMR sense.- Parameters:
- size: int
Number of points of the frequency scale
- dwfloat
Time spacing in the time dimension
- rev: bool
Reverses the scale
- Returns:
- fqscale: 1darray
The computed frequency scale.
- klassez.processing.mask_sgn_basl(ppm, data, SFO1, alpha=3, winsize=50)[source]
Given an NMR spectrum, this function estimates the signal and baseline regions, and return a list of slices to cut the spectrum accordingly.
- Parameters:
- ppm1darray
ppm scale of the spectrum
- data1darray
Spectrum
- SFO1float
Nucleus’ Larmor frequency /MHz
- alphafloat
Factor that multiplies the std of the spectrum to set the threshold
- winsizefloat
Minimum size of the window that can contain peaks /Hz
- Returns:
- peak_sliceslist of slices
Slices that trim the data in the signal-only regions
- basl_sliceslist of slices
Slices that trim the data in the baseline-only regions
- klassez.processing.mcr(input_data, nc, f=10, tol=0.001, itermax=10000.0, P='H', oncols=True)[source]
This is an implementation of Multivariate Curve Resolution for the denoising of 2D NMR data. Let us consider a matrix D, of dimensions
(m, n), where the starting data are stored. The final purpose of MCR is to decompose the D matrix as follows:\[D = CS + E\]where C and S are matrices of dimension
(m, nc)and(nc, n), respectively, and E contains the part of the data that are not reproduced by the factorization. Being D the FID of a NMR spectrum, C will contain time evolutions of the indirect dimension, and S will contain transients in the direct dimension.The total MCR workflow can be separated in two parts: a first algorithm that produces an initial guess for the three matrices C, S and E (
simplisma), and an optimization step that aims at the removal of the unwanted features of the data by iteratively filling the E matrix (mcr_als). This function returns the denoised datasets, CS, and the single C and S matrices.- Parameters:
- input_data2darray or 3darray
a 3D array containing the set of 2D NMR datasets to be coprocessed stacked along the first dimension. A single 2D array can be passed, if the denoising of a single dataset is desired.
- ncint
number of purest components to be looked for;
- ffloat
percentage of allowed noise;
- tolfloat
tolerance for the arrest criterion;
- itermaxint
maximum number of allowed iterations
- Pstr or 2darray
'H'for horizontal stacking,'V'for vertical stacking, or custom matrix as explained in the description ofmcr_stack- oncolsbool
True to estimate
Swithprocessing.simplisma, False to estimateC.
- Returns:
- CS_f2darray or 3darray
Final denoised data matrix
- C_f2darray or 3darray
Final C matrix
- S_f2darray or 3darray
Final S matrix
- klassez.processing.mcr_als(D, C, S, itermax=10000, tol=1e-05)[source]
Performs alternating least squares to get the final
CandSmatrices. Being the fundamental MCR equation:\[D = CS + E\]At the k-th step of the iterative cycle:
\(C_{(k)} = D S^+_{(k-1)}\)
\(S_{(k)} = C^+_{(k)} D\)
\(E_{(k)} = D - C_{(k)} S_{(k)}\)
Defined
rCandrSas the Frobenius norm of the difference ofCandSmatrices between two subsequent steps:\[rC = || C_{(k)} - C_{(k-1)} || \qquad rS = || S_{(k)} - S_{(k-1)} ||\]The convergence is reached when both
CandSchange less thantoltimes the first iteration:\[\frac{ rC_{(k)} }{rC_{1}} \leq tol \quad \text{and} \quad \frac{ rS_{(k)} }{rS_{1}} \leq tol\]- Parameters:
- D2darray
Input data, of dimensions
(m, n)- C2darray
Estimation of the
Cmatrix, of dimensions(m, nc).- S2darray
Estimation of the
Smatrix, of dimensions(nc, n).- itermaxint
Maximum number of iterations
- tolfloat
Threshold for the arrest criterion.
- Returns:
- C2darray
Optimized C matrix, of dimensions
(m, nc).- S2darray
Optimized S matrix, of dimensions
(nc, n).
- klassez.processing.mcr_stack(input_data, P='H')[source]
Performs matrix augmentation by assembling input_data according to the positioning matrix
P.Phas two default modes:'H'= horizontal stacking;'V'= vertical stacking. Otherwise, a customPmatrix can be given as follows. The entries of thePmatrix are the indices of the data in input_data. The shape of the matrix determines the final arrangement. See example. If each dataset ininput_datahas dimensions(m, n)andPhas dimensions(u,v), then the returned data matrix will have dimensions(m*u, n*v).- Parameters:
- input_data3darray or list of 2darray
Contains the spectra to be stacked together. The index that runs on the datasets must be the first one.
- Pstr or 2darray
'H'for horizontal stacking,'V'for vertical stacking, or custom matrix as explained in the description
- Returns:
- data2darray
Augmented data matrix.
Examples
If
input_data = [a, b, c, d, e, f], and one wants to obtain[[a, b], [d,c], [f, e]], the correspondantPmatrix is:P = [ [0, 1], [3, 2], [5, 4] ]
See also
- klassez.processing.mcr_unpack(C, S, nds, P='H')[source]
Reverts matrix augmentation of
klassez.processing.mcr_stack(). The denoised spectra can be calculated by matrix multiplication:for k in range(nds): D[k] = C_f[k] @ S_f[k]
- Parameters:
- C2darray
MCR C matrix
- S2darray
MCR S matrix
- ndsint
number of experiments
- Pstr or 2darray
'H'for horizontal stacking,'V'for vertical stacking, or custom matrix as explained in the description ofmcr_stack
- Returns:
- C_flist of 2darray
Disassembled MCR C matrix
- S_flist of 2darray
Disassembled MCR C matrix
See also
- klassez.processing.pknl(data, grpdly=0, onfid=False)[source]
Compensate for the Bruker group delay at the beginning of FID through a first-order phase correction of
p1 = - 360 * GRPDLYdegrees. If applied on the FID (onfid=True), it is equivalent to a left circular shift ofGRPDLYpoints. However, in order to accomodate for also non-integerGRPDLY, it is computed by doing the Fourier transform on the fly.- Parameters:
- datandarray
Input data. Be sure it is complex!
- grpdlyint
Number of points that make the group delay.
- onfidbool
If it is True, performs FT before to apply the phase correction, and IFT after.
- Returns:
- datapndarray
Corrected data
- klassez.processing.print(*args, c=None, s=None, sep=' ', end='\n', file=None, flush=False)[source]
Contains a series of processing functions for different purposes
- klassez.processing.ps(data, ppmscale=None, p0=None, p1=None, pivot=None, interactive=False, reference=None)[source]
Applies phase correction on the last dimension of data. The pivot is set at the center of the spectrum by default. Missing parameters will be inserted interactively.
- Parameters:
- datandarray
Input data
- ppmscale1darray or None
PPM scale of the spectrum. Required for pivot and interactive phase correction
- p0float
Zero-order phase correction angle /degrees
- p1float
First-order phase correction angle /degrees
- pivotfloat or None.
First-order phase correction pivot /ppm. If None, it is the center of the spectrum.
- interactivebool
If True, all the parameters will be ignored and the interactive phase correction panel will be opened.
- referencelist of 1darray or Spectrum_1D object
Reference spectrum to be used for phasing. Can be also given as
[ppm, spectrum]
- Returns:
- datapndarray
Phased data
- final_valuestuple
Employed values of the phase correction.
(p0, p1, pivot)
See also
- klassez.processing.qfil(ppm, data, u, s, SFO1)[source]
Suppress signals in the spectrum using a gaussian filter.
- Parameters:
- ppm1darray
ppm scale of the spectrum
- datandarray
Data to be processed. The filter is applied on the last dimension
- ufloat
Position of the filter /ppm
- sfloat
Width of the filter (FWHM) /Hz
- SFO1float
Spectrometer larmor frequency
- Returns:
- pdatandarray
Filtered data
- klassez.processing.qpol(fid)[source]
Fits the FID with a 4-th degree polynomion, then subtracts it from the original FID. The real and imaginary channels are treated separately.
- Parameters:
- fidndarray
Self-explanatory.
- Returns:
- fid_corrndarray
Processed FID
- klassez.processing.qsin(data, ssb)[source]
Sine-squared apodization.
\[a(x) = \sin^2 \biggl[\frac{\pi}{SSB} + \pi \biggl(1 - \frac{1}{SSB}\biggr) x \biggr]\]- Parameters:
- data: ndarray
FID to be processed
- ssb: int
Sine bell shift.
- Returns:
- datap: ndarray
Apodized data
See also
- klassez.processing.quad(fid)[source]
Subtracts from the FID the arithmetic mean of its last quarter. The real and imaginary channels are treated separately.
- Parameters:
- fidndarray
Self-explanatory.
- Returns:
- fidndarray
Processed FID.
- klassez.processing.repack_2D(rr, ir, ri, ii)[source]
Renconstruct hypercomplex 2D NMR data given the 4 components
- Parameters:
- rr2darray
Real F2, Real F1
- ir2darray
Imaginary F2, Real F1
- ri2darray
Real F2, Imaginary F1
- ii2darray
Imaginary F2, Imaginary F1
- Returns:
- data2darray
Hypecomplex matrix
- klassez.processing.rndc(data)[source]
Robust Noise Derivative Calculation, reference . Employed coefficients:
(42, 48, 27, 8, 1)/512Used to compute the first derivative of a function.- Parameters:
- data1darray
Input data
- Returns:
- dy1darray
First derivative of data. First and last 5 points are set to zero.
- klassez.processing.roll_dirac(y, x, off=0, onfid=False)[source]
Perform a circular shift on
yby convolution with a Dirac delta, centered atoffin thextimescale.- Parameters:
- yndarray
Data to shift. The shift will be applied on the last dimension.
- x1darray
Scale on which to consider
off.- offfloat
Shift value on
x- onfidbool
Set it to True if
yis a FID. The shift then is done by multiplying, instead of convolving.
- Returns:
- y_rollndarray
Shifted data
- klassez.processing.rpbc(data, split_imag=False, n=5, basl_method='huber', basl_thresh=0.2, basl_itermax=2000, **phase_kws)[source]
Reversed Phase and Baseline Correction. Allows for the automatic phase correction and baseline subtraction of NMR spectra. It is called “reversed” because the baseline is actually computed and subtracted before to perform the phase correction.
The baseline is computed using a low-order polynomion, built on a scale that goes from -1 to 1, whose coefficients are obtained minimizing a non-quadratic cost function. It is recommended to use either
"tq"(truncated quadratic, much faster) or"huber"(Huber function, slower but sometimes more accurate). The user is requested to choose between separating the real and imaginary channel in this step. The order of the polynomion and the threshold value are the key parameters for obtaining a good baseline. The used function isklassez.processing.polyn_basl()The phase correction is computed on the baseline-subtracted complex data as described in the SINC algorithm. The default parameters are generally fine, but in case of data with poor SNR (approximately SNR < 10) better results can be obtained by increasing the value of the
e1parameter. The employed function isklassez.fit.SINC_phase()Note
Not excellent results. Computation might be slow
- Parameters:
- data1darray
Data to be processed, complex-valued
- split_imagbool
If True, computes the baseline on the real and imaginary part separately; else, the set of polynomion coefficients are forced to be the same for both
- nint
Number of coefficients of the polynomion, i.e. it will be of degree n-1
- basl_methodstr
Cost function to be minimized for the baseline computation. Look for
fit.CostFunc,"method"attribute- basl_threshfloat
Relative threshold value for the non-quadratic behaviour of the cost function. Look for
fit.CostFunc,"s"attribute- basl_itermaxint
Maximun number of iterations allowed during the baseline fitting procedure
- phase_kwskeyworded arguments
Optional arguments for the phase correction. Look for
fit.SINC_phasekeyworded arguments for details.
- Returns:
- y1darray
Processed data
- p0float
Zero-order phase correction angle, in degrees
- p1float
First-order phase correction angle, in degrees
- c1darray
Set of coefficients to be used for the baseline computation, starting from the 0-order coefficient
- klassez.processing.simplisma(D, nc, f=10, oncols=True)[source]
Finds the first
ncpurest components of matrixDusing the simplisma algorithm, proposed by Windig and Guilment . Ifoncols=True, this function estimatesSwith simplisma, then calculates \(C = D S^+\). Ifoncols=False, this function estimatesCwith simplisma, then calculates \(S = C^+ D\).fdefines the percentage of allowed noise.- Parameters:
- D2darray
Input data, of dimensions
(m, n)- ncint
Number of components to be found. This determines the final size of the
CandSmatrices.- ffloat
Percentage of allowed noise.
- oncolsbool
If True, simplisma estimates the
Smatrix, otherwise estimatesC.
- Returns:
- C2darray
Estimation of the
Cmatrix, of dimensions(m, nc).- S2darray
Estimation of the
Smatrix, of dimensions(nc, n).
- klassez.processing.sin(data, ssb)[source]
Sine apodization.
\[a(x) = \sin \biggl[\frac{\pi}{SSB} + \pi \biggl(1 - \frac{1}{SSB}\biggr) x \biggr]\]- Parameters:
- data: ndarray
FID to be processed
- ssb: int
Sine bell shift.
- Returns:
- datap: ndarray
Apodized data
See also
- klassez.processing.sl_bas(x, y, lims=None)[source]
Computes the straight line that connects
ybetweenlimson thexscale. Ifyis a 2darray, a 2darray of lines is returned, one correspondant to each row ofy.- Parameters:
- y1darray or 2darray
Array for which to compute the connecting lines
- x1darray
Scale for the referencing of
lims- limstuple
(left, right)for selecting the endpoints
- Returns:
- bas_overlay1darray or 2darray
Computed connecting lines, with the same shape of
y. Of note, this is squeezed before returning!
- klassez.processing.sl_bas_onidx(y, x_idx)[source]
Computes the straight line that connects
y[min(x_idx)]withy[max(x_idx)]. Ifyis a 2darray, a 2darray of lines is returned, one correspondant to each row ofy. The lines will bemax(x_idx) - min(x_idx)points long.- Parameters:
- y1darray or 2darray
Array for which to compute the connecting lines
- x_idxtuple
Points indices of the endpoints of the connecting lines
- Returns:
- bas1darray or 2darray
Computed connecting lines. Of note, this is squeezed before returning!
See also
- klassez.processing.smooth_g(d, m)[source]
Apply a smoothing with a gaussian filter by convolution. The width of the filter is
1/m.- Parameters:
- d1darray
Data to be smoothed
- mfloat
Inverse width of the filter /pt
- Returns:
- yc1darray
Smoothed data
- klassez.processing.split_echo_train(datao, n, n_echoes, i_p=0)[source]
Separate a CPMG echo-train FID into echoes so to be processed separately. The first decay, i.e. the native FID, is extracted, and corresponds to echo number 0. Then, for each echo, the left side (reversed) is summed up to its right part.
- Parameters:
- dataondarray
FID with an echo train on its last dimension
- nint
number of points that separate one echo from the next
- n_echoesint
number of echoes to extract. If it is 0, extracts only the first decay
- i_pint
Number of offset points
- Returns:
- data_p(n+1)darray
Separated echoes
See also
- klassez.processing.splitcomb(data, taq, J=53.8)[source]
Applies the processing required for the IPAP virtual decoupling scheme in the direct dimension. The data structure must be with the IP in the first half of the direct dimension, and with AP in the second half. The default J is 53.8 Hz, correspondant to the CO-Ca coupling.
- Parameters:
- data2darray
FID of the spectrum to process
- taq1darray
Acquisition timescale
- Jfloat
Scalar coupling constant of the coupling to suppress, in Hz
- Returns:
- datap2darray
Decoupled data. The direct dimension is halved with respect to the original FID
- klassez.processing.stack_fids(*fids, filename=None)[source]
Stacks together FIDs in order to create a pseudo-2D experiment. This function can handle either arrays or Spectrum_1D objects.
- Parameters:
- fidssequence of 1darrays or Spectrum_1D objects
Input data.
- filenamestr
Location for a .npy file to be saved. If None, no file is created.
- Returns:
- p2d2darray
Stacked FIDs.
- klassez.processing.sum_echo_train(datao, n, n_echoes, i_p=0)[source]
Sum up a CPMG echo-train FID into echoes so to be enchance the SNR. This function calls
processing.split_echo_trainwith the same parameters.- Parameters:
- dataondarray
FID with an echo train on its last dimension
- nint
number of points that separate one echo from the next
- n_echoesint
number of echoes to sum
- i_pint
Number of offset points
- Returns:
- data_pndarray
Summed echoes
See also
- klassez.processing.tabula_rasa(data, lvl=0.05, cmap=<matplotlib.colors.LinearSegmentedColormap object>)[source]
This function is to be used in SIFT algorithm. Allows interactive selection using a Lasso widget of the region of the spectrum which contain signal. Returns a masking matrix, of the same shape as data, whose entries are 1 inside the selection and 0 outside.
- klassez.processing.td_eff(data, tdeff)[source]
Uses only the first
tdeffpoints of data. Equivalent to box-like apodization.tdeffmust be a list as long as the dimensions:tdeff = [F1, F2, ..., Fn]
- Parameters:
- datandarray
Data to be trimmed
- tdefflist of int
Number of points to be used in each dimension
- Returns:
- datainndarray
Trimmed data
- klassez.processing.tp_hyper(data)[source]
Computes the hypercomplex transpose of data. Needed for the processing of data acquired in a phase-sensitive manner in the indirect dimension.
- Parameters:
- data2darray
Hypercomplex data to be transposed
- Returns:
- datap2darray
Transposed data
- klassez.processing.unpack_2D(data)[source]
Separates hypercomplex data into 4 distinct components
- Parameters:
- data2darray
Hypercomplex matrix
- Returns:
- rr2darray
Real F2, Real F1
- ir2darray
Imaginary F2, Real F1
- ri2darray
Real F2, Imaginary F1
- ii2darray
Imaginary F2, Imaginary F1
- klassez.processing.whittaker_smoother(data, n=2, s_f=1, w=None)[source]
Adapted from P.H.C. Eilers, Anal. Chem 2003, 75, 3631-3636. Implementation of the smoothing algorithm proposed by Whittaker in 1923.
- Parameters:
- data1darray
Data to be smoothed
- nint
Order of the difference to be computed
- s_ffloat
Smoothing factor
- w1darray or None
Array of weights. If None, no weighting is applied.
- Returns:
- z1darray
Smoothed data
- klassez.processing.xfb(data, wf=[None, None], zf=[None, None], fcor=[0.5, 0.5], tdeff=[0, 0], u=True, FnMODE='States-TPPI')[source]
Performs the full processing of a 2D NMR FID (
data). The returned values depend onu: if it is True, returns a sequence of 2darrays depending on FnMODE, otherwise just the complex/hypercomplex data after FT in both dimensions.- Parameters:
- data2darray
Input data
- wfsequence of dict
(F1, F2);
{'mode': function to be used, 'parameters': different from each function}- zfsequence of int
final size of spectrum, (F1, F2)
- fcorsequence of float
weighting factor for the FID first point, (F1, F2)
- tdeffsequence of int
number of points of the FID to be used for the processing, (F1, F2)
- ubool
choose if to unpack the hypercomplex spectrum into separate arrays or not
- FnMODEstr
Acquisition mode in F1
- Returns:
- datap2darray or tuple of 2darray
Processed data or tuple of 2darray
- klassez.processing.zf(data, size)[source]
Zero-filling of
dataup tosizein its last dimension.- Parameters:
- datandarray
Array to be zero-filled
- sizeint
Number of points of the last dimension after zero-filling
- Returns:
- datazfndarray
Zero-filled data
See also