MSstatsSummarizeSingleTMP.Rd
Tukey Median Polish summarization for a single protein
MSstatsSummarizeSingleTMP( single_protein, impute, censored_symbol, remove50missing )
single_protein | feature-level data for a single protein |
---|---|
impute | only for summaryMethod = "TMP" and censoredInt = 'NA' or '0'. TRUE (default) imputes 'NA' or '0' (depending on censoredInt option) by Accelated failure model. FALSE uses the values assigned by cutoffCensored |
censored_symbol | Missing values are censored or at random. 'NA' (default) assumes that all 'NA's in 'Intensity' column are censored. '0' uses zero intensities as censored intensity. In this case, NA intensities are missing at random. The output from Skyline should use '0'. Null assumes that all NA intensites are randomly missing. |
remove50missing | only for summaryMethod = "TMP". TRUE removes the runs which have more than 50% missing values. FALSE is default. |
list of two data.tables: one with fitted survival model, the other with protein-level data
raw = DDARawData method = "TMP" cens = "NA" impute = TRUE # currently, MSstats only supports MBimpute = FALSE for linear summarization MSstatsConvert::MSstatsLogsSettings(FALSE) input = MSstatsPrepareForDataProcess(raw, 2, NULL)#> INFO [2021-07-05 20:05:30] ** Features with one or two measurements across runs are removed. #> INFO [2021-07-05 20:05:30] ** Fractionation handled. #> INFO [2021-07-05 20:05:30] ** Updated quantification data to make balanced design. Missing values are marked by NAinput = MSstatsNormalize(input, "EQUALIZEMEDIANS") input = MSstatsMergeFractions(input) input = MSstatsHandleMissing(input, "TMP", TRUE, "NA", 0.999)#> INFO [2021-07-05 20:05:30] ** Log2 intensities under cutoff = 13.456 were considered as censored missing values. #> INFO [2021-07-05 20:05:30] ** Log2 intensities = NA were considered as censored missing values.#> INFO [2021-07-05 20:05:30] ** Use all features that the dataset originally has.input = MSstatsPrepareForSummarization(input, method, impute, cens, FALSE) input_split = split(input, input$PROTEIN) single_protein_summary = MSstatsSummarizeSingleTMP(input_split[[1]], impute, cens, FALSE) head(single_protein_summary[[1]])#> RUN LogIntensities Protein #> 1: 1 21.28437 bovine #> 2: 2 20.85653 bovine #> 3: 3 20.67521 bovine #> 4: 4 21.60443 bovine #> 5: 5 21.82186 bovine #> 6: 6 21.20445 bovine