Peak_detection: Validate nmr_dataset_peak_table objects

Description Usage Arguments Value Parallelization See Also Examples

Description

The function detects peaks on an nmr_dataset_1D object, using speaq::detectSpecPeaks. detectSpecPeaks divides the whole spectra into smaller segments and uses MassSpecWavelet::peakDetectionCWT for peak detection.

This function is based on speaq::dohCluster.

The function allows the integration of a given ppm vector with a specific width.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
validate_nmr_dataset_peak_table(nmr_dataset_peak_table)

nmr_detect_peaks(
  nmr_dataset,
  nDivRange_ppm = 0.1,
  scales = seq(1, 16, 2),
  baselineThresh = NULL,
  SNR.Th = 3
)

nmr_detect_peaks_plot(nmr_dataset, peak_data, NMRExperiment, ...)

nmr_detect_peaks_tune_snr(
  ds,
  NMRExperiment = NULL,
  SNR_thresholds = seq(from = 2, to = 6, by = 0.1)
)

nmr_align(
  nmr_dataset,
  peak_data,
  NMRExp_ref = NULL,
  maxShift_ppm = 0.0015,
  acceptLostPeak = FALSE
)

nmr_integrate_peak_positions(
  samples,
  peak_pos_ppm,
  peak_width_ppm = 0.006,
  ...
)

get_integration_with_metadata(integration_object)

Arguments

nmr_dataset_peak_table

An nmr_dataset_peak_table object

nmr_dataset

An nmr_dataset_1D

nDivRange_ppm

Segment size, in ppms, to divide the spectra and search for peaks.

scales

The parameter of peakDetectionCWT function of MassSpecWavelet package, look it up in the original function.

baselineThresh

It will remove all peaks under an intensity set by baselineThresh. If you set it to 'NULL', nmr_detect_peaks will automatically compute an aproximate value considering baseline between 9.5 and 10.0 ppm (automatically calculation using baselineThresh = NULL will not work if spectra were not interpolated up to 10.0 ppm)

SNR.Th

The parameter of peakDetectionCWT function of MassSpecWavelet package, look it up in the original function. If you set -1, the function will itself re-compute this value.

peak_data

The detected peak data given by nmr_detect_peaks.

NMRExperiment

A string with the single NMRExperiment used explore the SNR thresholds. If not given, use the first one.

...

Arguments passed on to nmr_integrate_regions

regions

A named list. Each element of the list is a region, given as a named numeric vector of length two with the range to integrate. The name of the region will be the name of the column

ds

An nmr_dataset_1D dataset

SNR_thresholds

A numeric vector with the SNR thresholds to explore

NMRExp_ref

NMRExperiment of the reference to use for alignment

maxShift_ppm

The maximum shift allowed, in ppm

acceptLostPeak

This is an option for users, TRUE is the default value. If the users believe that all the peaks in the peak list are true positive, change it to FALSE.

samples

A nmr_dataset object

peak_pos_ppm

The peak positions, in ppm

peak_width_ppm

The peak widths (or a single peak width for all peaks)

integration_object

A nmr_dataset object

Value

The nmr_dataset_peak_table unchanged

This function is useful for its side-effects: Stopping in case of error

A data frame with the NMRExperiment, the sample index, the position in ppm and index and the peak intensity

Plot peak detection results

A list with the following elements:

An nmr_dataset_1D, with the spectra aligned

Integrate peak positions

Get integrals with metadata from integrate peak positions

integration dataframe

Parallelization

This function accepts parallellization with future strategies. You can use plan(multiprocess) or plan(sequential) before calling this function to determine if it should be parallellized or not.

See Also

nmr_align for peak alignment with the detected peak table

Other nmr_dataset_peak_table functions: [.nmr_dataset_peak_table(), format.nmr_dataset_peak_table(), is.nmr_dataset_peak_table(), load_and_save_functions, new_nmr_dataset_peak_table(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), print.nmr_dataset_peak_table()

Other class helper functions: format.nmr_dataset_1D(), format.nmr_dataset_peak_table(), format.nmr_dataset(), is.nmr_dataset_1D(), is.nmr_dataset_peak_table(), new_nmr_dataset_1D(), new_nmr_dataset_peak_table(), new_nmr_dataset(), print.nmr_dataset_1D(), print.nmr_dataset_peak_table(), print.nmr_dataset(), validate_nmr_dataset_family(), validate_nmr_dataset()

Other peak detection functions: Pipelines, nmr_baseline_threshold(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions(), regions_from_peak_table()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Other peak detection functions: Pipelines, nmr_baseline_threshold(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions(), regions_from_peak_table()

Other peak detection functions: Pipelines, nmr_baseline_threshold(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions(), regions_from_peak_table()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Other peak detection functions: Pipelines, nmr_baseline_threshold(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions(), regions_from_peak_table()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Other alignment functions: Pipelines, nmr_align_find_ref()

Other peak alignment functions: nmr_align_find_ref()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Other peak integration functions: Pipelines, computes_peak_width_ppm(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Other peak integration functions: Pipelines, computes_peak_width_ppm(), nmr_identify_regions_blood(), nmr_identify_regions_cell(), nmr_identify_regions_urine(), nmr_integrate_regions()

Other nmr_dataset_1D functions: [.nmr_dataset_1D(), computes_peak_width_ppm(), file_lister(), files_to_rDolphin(), format.nmr_dataset_1D(), is.nmr_dataset_1D(), load_and_save_functions, new_nmr_dataset_1D(), nmr_align_find_ref(), nmr_baseline_removal(), nmr_baseline_threshold(), nmr_exclude_region(), nmr_integrate_regions(), nmr_interpolate_1D(), nmr_meta_add(), nmr_meta_export(), nmr_meta_get_column(), nmr_meta_get(), nmr_normalize(), nmr_pca_build_model(), nmr_pca_outliers_filter(), nmr_pca_outliers_plot(), nmr_pca_outliers_robust(), nmr_pca_outliers(), nmr_ppm_resolution(), plot.nmr_dataset_1D(), plot_webgl(), print.nmr_dataset_1D(), rdCV_PLS_RF_ML(), rdCV_PLS_RF(), save_files_to_rDolphin(), to_ChemoSpec(), validate_nmr_dataset()

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
dir_to_demo_dataset <- system.file("dataset-demo", package = "AlpsNMR")
nmr_dataset <- nmr_read_samples_dir(dir_to_demo_dataset)
dataset_1D <- nmr_interpolate_1D(nmr_dataset, axis = c(min = -0.5, max = 10, by = 2.3E-4))

sample_10 <- filter(dataset_1D, NMRExperiment == "10")

# 1.Peak detection in the dataset.
peak_data <- nmr_detect_peaks(dataset_1D,
                              nDivRange_ppm = 0.1, # Size of detection segments
                              scales = seq(1, 16, 2),
                              baselineThresh = 0, # Minimum peak intensity
                              SNR.Th = 4) # Signal to noise ratio

#nmr_detect_peaks_plot(sample_10, peak_data, "NMRExp_ref")

peaks_detected <- nmr_detect_peaks_tune_snr(sample_10, 
                                            SNR_thresholds = seq(from = 2, 
                                            to = 3, by = 0.5))


# 2.Find the reference spectrum to align with.
NMRExp_ref <- nmr_align_find_ref(dataset_1D, peak_data)

# 3.Spectra alignment using the ref spectrum and a maximum alignment shift
nmr_dataset <- nmr_align(dataset_1D, # the dataset
                         peak_data, # detected peaks
                         NMRExp_ref = NMRExp_ref, # ref spectrum
                         maxShift_ppm = 0.0015, # max alignment shift
                         acceptLostPeak = FALSE) # lost peaks
                         
# 4.PEAK INTEGRATION (please, consider previous normalization step).
# First we take the peak table from the reference spectrum
peak_data_ref <- filter(peak_data, NMRExperiment == NMRExp_ref)

# Then we integrate spectra considering the peaks from the ref spectrum
nmr_peak_table <- nmr_integrate_peak_positions(
                      samples = nmr_dataset,
                      peak_pos_ppm = peak_data_ref$ppm,
                      peak_width_ppm = NULL)

validate_nmr_dataset_peak_table(nmr_peak_table)

#If you wanted the final peak table before machine learning you can run
nmr_peak_table_completed <- get_integration_with_metadata(nmr_peak_table)

AlpsNMR documentation built on April 1, 2021, 6:02 p.m.