Title: | Tools for the preprocessing and analysis of the Illumina microarrays on the detector (bead) level |
---|---|
Description: | Package blima includes several algorithms for the preprocessing of Illumina microarray data. It focuses to the bead level analysis and provides novel approach to the quantile normalization of the vectors of unequal lengths. It provides variety of the methods for background correction including background subtraction, RMA like convolution and background outlier removal. It also implements variance stabilizing transformation on the bead level. There are also implemented methods for data summarization. It also provides the methods for performing T-tests on the detector (bead) level and on the probe level for differential expression testing. |
Authors: | Vojtěch Kulvait |
Maintainer: | Vojtěch Kulvait <[email protected]> |
License: | GPL-3 |
Version: | 1.41.0 |
Built: | 2024-11-18 03:16:06 UTC |
Source: | https://github.com/bioc/blima |
Package blima includes several algorithms for the preprocessing of Illumina microarray data. It focuses to the bead level analysis and provides novel approach to the quantile normalization of the vectors of unequal lengths. It provides variety of the methods for background correction including background subtraction, RMA like convolution and background outlier removal. It also implements variance stabilizing transformation on the bead level. There are also implemented methods for data summarization. It provides the methods for performing T-tests on the detector (bead) level and on the probe level for differential expression testing.
The DESCRIPTION file:
Encoding: | UTF-8 |
Package: | blima |
Type: | Package |
Title: | Tools for the preprocessing and analysis of the Illumina microarrays on the detector (bead) level |
Version: | 1.41.0 |
Date: | 2017-09-23 |
Author: | Vojtěch Kulvait |
Maintainer: | Vojtěch Kulvait <[email protected]> |
Description: | Package blima includes several algorithms for the preprocessing of Illumina microarray data. It focuses to the bead level analysis and provides novel approach to the quantile normalization of the vectors of unequal lengths. It provides variety of the methods for background correction including background subtraction, RMA like convolution and background outlier removal. It also implements variance stabilizing transformation on the bead level. There are also implemented methods for data summarization. It also provides the methods for performing T-tests on the detector (bead) level and on the probe level for differential expression testing. |
License: | GPL-3 |
LazyLoad: | yes |
Depends: | R(>= 3.3) |
Imports: | beadarray(>= 2.0.0), Biobase(>= 2.0.0), Rcpp (>= 0.12.8), BiocGenerics, grDevices, stats, graphics |
LinkingTo: | Rcpp |
Suggests: | xtable, blimaTestingData, BiocStyle, illuminaHumanv4.db, lumi, knitr |
URL: | https://bitbucket.org/kulvait/blima |
biocViews: | Microarray, Preprocessing, Normalization, DifferentialExpression, GeneRegulation, GeneExpression |
VignetteBuilder: | knitr |
Config/pak/sysreqs: | libicu-dev libpng-dev libssl-dev |
Repository: | https://bioc.r-universe.dev |
RemoteUrl: | https://github.com/bioc/blima |
RemoteRef: | HEAD |
RemoteSha: | c59eaff0b1196f71ec0ab38c08153ae372d450cd |
Index of help topics:
aggregateAndPreprocess Aggregate data bacgroundCorrect Data background correction. bacgroundCorrectSingleArray Data background correction. backgroundChannelSubtract Background channel subtraction backgroundChannelSubtractSingleArray Background channel subtraction blima-package Package for the preprocessing and analysis of the Illumina microarrays on the detector (bead) level. channelExistsIntegrityWithLogicalVectorList Internal function checkIntegrity Internal function checkIntegrityLogical Internal function checkIntegrityOfListOfBeadLevelDataObjects Internal function checkIntegrityOfSingleBeadLevelDataObject Internal function chipArrayStatistics Statistics of beadLevelData createSummarizedMatrix Summarized value matrix. doAction Internal function doProbeTTests T-test for probe level data. doTTests T-test for bead (detector) level data. filterBg Bg correct vector getNextVector Support probe and beadl level testing. initMeanDistribution initMeanDistribution insertColumn Internal function to support chipArrayStatistics interpolateSortedVector Interpolate sorted vector interpolateSortedVectorRcpp_ interpolateSortedVectorRcpp log2TransformPositive Log2 transform of numbers >1. meanDistribution Produce sorted double vector with mean distribution. nonParametricEstimator INTERNAL FUNCTION Xie background correct. nonPositiveCorrect Correct non positive nonPositiveCorrectSingleArray Correct non positive numberOfDistributionElements Internal performXieCorrection INTERNAL FUNCTION Xie background correct. plotBackgroundImageAfterCorrection Plot background image after correction plotBackgroundImageBeforeCorrection Plot background image before correction quantileNormalize Bead level quantile normalization. readToVector Support doTTests function. selectedChannelTransform Channel transformation selectedChannelTransformSingleArray Channel transformation singleArrayNormalize Bead level quantile normalization. singleChannelExistsIntegrityWithLogicalVector Internal function singleCheckIntegrityLogicalVector Internal function singleNumberOfDistributionElements Internal updateMeanDistribution updateMeanDistribution varianceBeadStabilise Bead level VST. varianceBeadStabiliseSingleArray Bead level VST. vstFromLumi Function from LGPL lumi package 2.16.0 writeBackgroundImages Write Background Images xieBacgroundCorrect Xie background correct. xieBacgroundCorrectSingleArray INTERNAL FUNCTION Xie background correct.
Vojtěch Kulvait Maintainer: Vojtěch Kulvait <[email protected]>
This function is not intended to direct use. It helps perform work of doProbeTTests function. For each probe it prints mean and sd of an quality.
aggregateAndPreprocess(x, quality = "qua", transformation = NULL)
aggregateAndPreprocess(x, quality = "qua", transformation = NULL)
x |
Two column matrix to agregate with columns "ProbeID" and quality. |
quality |
Quality to analyze, default is "qua". |
transformation |
Function of input data trasformation, default is NULL. Any function which for input value returns transformed value may be supplied. T-test then will be evaluated on transformed data, consider use log2TranformPositive. |
Some return value
Vojtěch Kulvait
Background correction procedure selecting beads with background Intensity I_b |mean - I_b | > k*SD(I_bs) for exclusion.
bacgroundCorrect(b, normalizationMod = NULL, channelBackground = "GrnB", k = 3, channelBackgroundFilter = "bgf", channelAndVector = NULL)
bacgroundCorrect(b, normalizationMod = NULL, channelBackground = "GrnB", k = 3, channelBackgroundFilter = "bgf", channelAndVector = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelBackground |
Name of channel to normalize. |
k |
Parameter of method stringency (default is 3). |
channelBackgroundFilter |
Filtered beads will have weight 0 and non filtered weight 1. |
channelAndVector |
Represents vector to bitvise multiple to the channelBackgroundFilter vector. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
Background correction procedure selecting beads with background Intensity I_b |mean - I_b | > k*SD(I_bs) for exclusion, internal.
bacgroundCorrectSingleArray(b, normalizationMod = NULL, channelBackground = "GrnB", k = 3, channelBackgroundFilter = "bgf", channelAndVector = NULL)
bacgroundCorrectSingleArray(b, normalizationMod = NULL, channelBackground = "GrnB", k = 3, channelBackgroundFilter = "bgf", channelAndVector = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelBackground |
Name of channel to normalize. |
k |
Parameter of method stringency (default is 3). |
channelBackgroundFilter |
Filtered beads will have weight 0 and non filtered weight 1. |
channelAndVector |
Represents vector to bitvise multiple to the channelBackgroundFilter vector. |
Vojtěch Kulvait
Function to subtract one channel from another producing new channel. Standard graphic subtraction.
backgroundChannelSubtract(b, normalizationMod = NULL, channelSubtractFrom = "GrnF", channelSubtractWhat = "GrnB", channelResult = "Grn")
backgroundChannelSubtract(b, normalizationMod = NULL, channelSubtractFrom = "GrnF", channelSubtractWhat = "GrnB", channelResult = "Grn")
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for performing on all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelSubtractFrom |
Name of channel to subtract from. |
channelSubtractWhat |
Name of channel to subtract. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
INTERNAL FUNCTION Correction for positive values only
backgroundChannelSubtractSingleArray(b, normalizationMod = NULL, channelSubtractFrom = "GrnF", channelSubtractWhat = "GrnB", channelResult = "Grn")
backgroundChannelSubtractSingleArray(b, normalizationMod = NULL, channelSubtractFrom = "GrnF", channelSubtractWhat = "GrnB", channelResult = "Grn")
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelSubtractFrom |
Name of channel to subtract from. |
channelSubtractWhat |
Name of channel to subtract. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
Vojtěch Kulvait
Test existence of channel slot based on vector list
channelExistsIntegrityWithLogicalVectorList(b, spotsToCheck = NULL, slotToCheck, action = c("returnText", "warn", "error"))
channelExistsIntegrityWithLogicalVectorList(b, spotsToCheck = NULL, slotToCheck, action = c("returnText", "warn", "error"))
b |
List of beadLevelData objects. |
spotsToCheck |
NULL for check all spots from b. Otherwise specifies logical vector of the length equals to the number of arrays in b with TRUE for checking. |
slotToCheck |
Slot name to check |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
Check integrity of the list of beadLevelData objects or single beadLevelData object returns waslist.
checkIntegrity(b, action = c("warn", "error"))
checkIntegrity(b, action = c("warn", "error"))
b |
List of beadLevelData objects or single. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Returns value if the object was list or not before calling this function.
Vojtěch Kulvait
Check integrity of the list of logical objects, internal.
checkIntegrityLogical(xx, b, action = c("returnText", "warn", "error"))
checkIntegrityLogical(xx, b, action = c("returnText", "warn", "error"))
xx |
List of logical objects compatible with a list b. |
b |
List of beadLevelData objects. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
Check integrity of the list of beadLevelData objects, internal.
checkIntegrityOfListOfBeadLevelDataObjects(listb, action = c("returnText", "warn", "error"))
checkIntegrityOfListOfBeadLevelDataObjects(listb, action = c("returnText", "warn", "error"))
listb |
List of beadLevelData objects. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
Check integrity of single beadLevelData object, internal.
checkIntegrityOfSingleBeadLevelDataObject(b, action = c("returnText", "warn", "error"))
checkIntegrityOfSingleBeadLevelDataObject(b, action = c("returnText", "warn", "error"))
b |
beadLevelData object. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
This function returns table with statistics of single beadLevelData object indexed by order of spots. It prints number of beads on each array spot mean foreground intensity and optionally mean background intensity, mean number of beads in probe set and unbiased estimate of standard deviations of these parameters. Optionaly you can also obtain percentage of removed beads within excludedOnSDMultiple multiple of standard deviations from the background value.
chipArrayStatistics(b, includeBeadStatistic = TRUE, channelForeground = "GrnF", channelBackground = "GrnB", includeBackground = TRUE, excludedOnSDMultiple = NA)
chipArrayStatistics(b, includeBeadStatistic = TRUE, channelForeground = "GrnF", channelBackground = "GrnB", includeBackground = TRUE, excludedOnSDMultiple = NA)
b |
Single beadLevelData object. |
includeBeadStatistic |
Include number of beads per probe in output. |
channelForeground |
Name of channel of foreground. |
channelBackground |
Name of channel of background. |
includeBackground |
Whether to output background data. |
excludedOnSDMultiple |
If positive number, print how much percents of the background lies more than excludedOnSDMultiple multipliers of standard deviation estimate away from background mean. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To print basic statistic data about blimatesting[[1]] object. data(blimatesting) array1stats = chipArrayStatistics(blimatesting[[1]], includeBeadStatistic=TRUE, excludedOnSDMultiple=3) array1pheno = pData(blimatesting[[1]]@experimentData$phenoData) array1stats = data.frame(array1pheno$Name, array1stats) colnames(array1stats)[1] <- "Array"; print(array1stats); }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To print basic statistic data about blimatesting[[1]] object. data(blimatesting) array1stats = chipArrayStatistics(blimatesting[[1]], includeBeadStatistic=TRUE, excludedOnSDMultiple=3) array1pheno = pData(blimatesting[[1]]@experimentData$phenoData) array1stats = data.frame(array1pheno$Name, array1stats) colnames(array1stats)[1] <- "Array"; print(array1stats); }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
This function creates summarized matrix of values of certain type.
createSummarizedMatrix(b, spotsToProcess = NULL, quality = "qua", channelInclude = "bgf", annotationTag = NULL)
createSummarizedMatrix(b, spotsToProcess = NULL, quality = "qua", channelInclude = "bgf", annotationTag = NULL)
b |
List of beadLevelData objects (or single object). |
spotsToProcess |
NULL for processing all spots in b. Otherwise specifies logical vector of the length equals to the number of arrays in b. |
quality |
Quality to matrize. |
channelInclude |
This field allows user to set channel with weights which have to be from 0,1. All zero weighted items are excluded from summarization. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to "bgf"). |
annotationTag |
Tag from annotation file which to use in resulting matrix as colname. |
Vojtěch Kulvait
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #Create summarization of nonnormalized data from GrnF column. data(blimatesting) blimatesting = bacgroundCorrect(blimatesting, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). nonnormalized = createSummarizedMatrix(blimatesting, quality="GrnF", channelInclude="bgf", annotationTag="Name") head(nonnormalized) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #Create summarization of nonnormalized data from GrnF column. data(blimatesting) blimatesting = bacgroundCorrect(blimatesting, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). nonnormalized = createSummarizedMatrix(blimatesting, quality="GrnF", channelInclude="bgf", annotationTag="Name") head(nonnormalized) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
Performs action of certain type
doAction(message, action = c("returnText", "warn", "error"))
doAction(message, action = c("returnText", "warn", "error"))
message |
Text message. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
This function does aggregated probe level t-tests on the data provided by the object beadLevelData from package beadarray.
doProbeTTests(b, c1, c2, quality = "qua", channelInclude = "bgf", correction = "BY", transformation = NULL)
doProbeTTests(b, c1, c2, quality = "qua", channelInclude = "bgf", correction = "BY", transformation = NULL)
b |
List of beadLevelData objects (or single object). |
c1 |
List of logical vectors of data to assign to the first group (or single vector). |
c2 |
List of logical vectors of data to assign to the second group (or single vector). |
quality |
Quality to analyze, default is "qua". |
channelInclude |
This field allows user to set channel with weights which have to be 0,1. All zero weighted items are excluded from t-test. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to "bgf"). |
correction |
Multiple testing adjustment method as defined by p.adjust function, default is "BY". |
transformation |
Function of input data trasformation, default is NULL. Any function which for input value returns transformed value may be supplied. T-test then will be evaluated on transformed data, consider use log2TranformPositive. |
Vojtěch Kulvait
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, variance stabilization and quantile normalization then test on probe level, bead level and print top 10 results. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") beadTest = doTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") probeTest = doProbeTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") adrToSymbol <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) adrToSymbol <- adrToSymbol[,c("ArrayAddress", "SymbolReannotated") ] colnames(adrToSymbol) <- c("Array_Address_Id", "Symbol") probeTestID = probeTest[,"ProbeID"] beadTestID = beadTest[,"ProbeID"] probeTestFC = abs(probeTest[,"mean1"]-probeTest[,"mean2"]) beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) probeTestP = probeTest[,"adjustedp"] beadTestP = beadTest[,"adjustedp"] probeTestMeasure = (1-probeTestP)*probeTestFC beadTestMeasure = (1-beadTestP)*beadTestFC probeTest = cbind(probeTestID, probeTestMeasure) beadTest = cbind(beadTestID, beadTestMeasure) colnames(probeTest) <- c("ArrayAddressID", "difexPL") colnames(beadTest) <- c("ArrayAddressID", "difexBL") tocmp <- merge(probeTest, beadTest) tocmp = merge(tocmp, adrToSymbol, by.x="ArrayAddressID", by.y="Array_Address_Id") tocmp = tocmp[, c("ArrayAddressID", "Symbol", "difexPL", "difexBL")] sortPL = sort(-tocmp[,"difexPL"], index.return=TRUE)$ix sortBL = sort(-tocmp[,"difexBL"], index.return=TRUE)$ix beadTop10 = tocmp[sortBL[1:10],] probeTop10 = tocmp[sortPL[1:10],] print(beadTop10) print(probeTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, variance stabilization and quantile normalization then test on probe level, bead level and print top 10 results. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") beadTest = doTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") probeTest = doProbeTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") adrToSymbol <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) adrToSymbol <- adrToSymbol[,c("ArrayAddress", "SymbolReannotated") ] colnames(adrToSymbol) <- c("Array_Address_Id", "Symbol") probeTestID = probeTest[,"ProbeID"] beadTestID = beadTest[,"ProbeID"] probeTestFC = abs(probeTest[,"mean1"]-probeTest[,"mean2"]) beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) probeTestP = probeTest[,"adjustedp"] beadTestP = beadTest[,"adjustedp"] probeTestMeasure = (1-probeTestP)*probeTestFC beadTestMeasure = (1-beadTestP)*beadTestFC probeTest = cbind(probeTestID, probeTestMeasure) beadTest = cbind(beadTestID, beadTestMeasure) colnames(probeTest) <- c("ArrayAddressID", "difexPL") colnames(beadTest) <- c("ArrayAddressID", "difexBL") tocmp <- merge(probeTest, beadTest) tocmp = merge(tocmp, adrToSymbol, by.x="ArrayAddressID", by.y="Array_Address_Id") tocmp = tocmp[, c("ArrayAddressID", "Symbol", "difexPL", "difexBL")] sortPL = sort(-tocmp[,"difexPL"], index.return=TRUE)$ix sortBL = sort(-tocmp[,"difexBL"], index.return=TRUE)$ix beadTop10 = tocmp[sortBL[1:10],] probeTop10 = tocmp[sortPL[1:10],] print(beadTop10) print(probeTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
This function does t-tests on the data provided by the object beadLevelData from package beadarray.
doTTests(b, c1, c2, quality = "qua", channelInclude = "bgf", correction = "BY", transformation = NULL)
doTTests(b, c1, c2, quality = "qua", channelInclude = "bgf", correction = "BY", transformation = NULL)
b |
List of beadLevelData objects (or single object). |
c1 |
List of logical vectors of data to assign to the first group (or single vector). |
c2 |
List of logical vectors of data to assign to the second group (or single vector). |
quality |
Quality to analyze, default is "qua". |
channelInclude |
This field allows user to set channel with weights which have to be 0,1. All zero weighted items are excluded from t-test. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to "bgf"). |
correction |
Multiple testing adjustment method as defined by p.adjust function, default is "BY". |
transformation |
Function of input data trasformation, default is NULL. Any function which for input value returns transformed value may be supplied. T-test then will be evaluated on transformed data, consider use log2TransformPositive. |
Vojtěch Kulvait
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, variance stabilization and quantile normalization then test on probe level, bead level and print top 10 results. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") beadTest = doTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") probeTest = doProbeTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") adrToSymbol <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) adrToSymbol <- adrToSymbol[,c("ArrayAddress", "SymbolReannotated") ] colnames(adrToSymbol) <- c("Array_Address_Id", "Symbol") probeTestID = probeTest[,"ProbeID"] beadTestID = beadTest[,"ProbeID"] probeTestFC = abs(probeTest[,"mean1"]-probeTest[,"mean2"]) beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) probeTestP = probeTest[,"adjustedp"] beadTestP = beadTest[,"adjustedp"] probeTestMeasure = (1-probeTestP)*probeTestFC beadTestMeasure = (1-beadTestP)*beadTestFC probeTest = cbind(probeTestID, probeTestMeasure) beadTest = cbind(beadTestID, beadTestMeasure) colnames(probeTest) <- c("ArrayAddressID", "difexPL") colnames(beadTest) <- c("ArrayAddressID", "difexBL") tocmp <- merge(probeTest, beadTest) tocmp = merge(tocmp, adrToSymbol, by.x="ArrayAddressID", by.y="Array_Address_Id") tocmp = tocmp[, c("ArrayAddressID", "Symbol", "difexPL", "difexBL")] sortPL = sort(-tocmp[,"difexPL"], index.return=TRUE)$ix sortBL = sort(-tocmp[,"difexBL"], index.return=TRUE)$ix beadTop10 = tocmp[sortBL[1:10],] probeTop10 = tocmp[sortPL[1:10],] print(beadTop10) print(probeTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, variance stabilization and quantile normalization then test on probe level, bead level and print top 10 results. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") beadTest = doTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") probeTest = doProbeTTests(blimatesting, groups1Mod, groups2Mod, "qua", "bgf") adrToSymbol <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) adrToSymbol <- adrToSymbol[,c("ArrayAddress", "SymbolReannotated") ] colnames(adrToSymbol) <- c("Array_Address_Id", "Symbol") probeTestID = probeTest[,"ProbeID"] beadTestID = beadTest[,"ProbeID"] probeTestFC = abs(probeTest[,"mean1"]-probeTest[,"mean2"]) beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) probeTestP = probeTest[,"adjustedp"] beadTestP = beadTest[,"adjustedp"] probeTestMeasure = (1-probeTestP)*probeTestFC beadTestMeasure = (1-beadTestP)*beadTestFC probeTest = cbind(probeTestID, probeTestMeasure) beadTest = cbind(beadTestID, beadTestMeasure) colnames(probeTest) <- c("ArrayAddressID", "difexPL") colnames(beadTest) <- c("ArrayAddressID", "difexBL") tocmp <- merge(probeTest, beadTest) tocmp = merge(tocmp, adrToSymbol, by.x="ArrayAddressID", by.y="Array_Address_Id") tocmp = tocmp[, c("ArrayAddressID", "Symbol", "difexPL", "difexBL")] sortPL = sort(-tocmp[,"difexPL"], index.return=TRUE)$ix sortBL = sort(-tocmp[,"difexBL"], index.return=TRUE)$ix beadTop10 = tocmp[sortBL[1:10],] probeTop10 = tocmp[sortPL[1:10],] print(beadTop10) print(probeTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
Background correction procedure selecting beads with background Intensity I_b |mean - I_b | > k*SD(I_bs) for exclusion, internal.
filterBg(x, k = 3)
filterBg(x, k = 3)
x |
Vector to correct |
k |
Parameter of method stringency (default is 3). |
Vojtěch Kulvait
Internal function supporting probe and beadl level testing.
getNextVector(what, from, length)
getNextVector(what, from, length)
what |
Two column sorted matrix with probe values. |
from |
Index to start on |
length |
nrow(what) |
Vojtěch Kulvait
This is internal function not intended to direct use which initializes mean distribution.
initMeanDistribution(srt, prvku)
initMeanDistribution(srt, prvku)
srt |
vector of sorted values |
prvku |
number of items in meanDistribution |
Vojtěch Kulvait
Internal
insertColumn(matrix, column, name)
insertColumn(matrix, column, name)
matrix |
Object to insert column to |
column |
Column to insert |
name |
Name of column to assign. |
Vojtěch Kulvait
Interpolates given sorted vector to the vector of different length. It does not sort input vector thus for unsorted vectors do not guarantee functionality. Internal function.
interpolateSortedVector(vector, newSize)
interpolateSortedVector(vector, newSize)
vector |
Sorted vector to interpolate. |
newSize |
Size of the vector to produce. |
Vojtěch Kulvait
interpolateSortedVectorRcpp_(vector, newSize)
interpolateSortedVectorRcpp_(vector, newSize)
vector |
|
newSize |
Vojtěch Kulvait
Transformation function are popular in beadarray package. Here this is similar concept. This function allow user to perform log transformation before doing t-tests.
log2TransformPositive(x)
log2TransformPositive(x)
x |
Number to transform. |
This function returns logarithm of base 2 for numbers >=1 and zero for numbers <1.
Vojtěch Kulvait
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, quantile normalization and then bead level t-test on log data run. Vst is not performed in this scheme. Top 10 probes is then printed according to certain measure. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = quantileNormalize(blimatesting, normalizationMod=c, channelNormalize="GrnF", channelOutput="qua", channelInclude="bgf") beadTest <- doTTests(blimatesting, groups1Mod, groups2Mod, transformation=log2TransformPositive, quality="qua", channelInclude="bgf") symbol2address <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) symbol2address <- symbol2address[,c("SymbolReannotated", "ArrayAddress") ] colnames(symbol2address) <- c("Symbol", "ArrayAddressID") beadTest = merge(beadTest, symbol2address, by.x="ProbeID", by.y="ArrayAddressID") beadTestID = beadTest[,c("ProbeID", "Symbol")] beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) beadTestP = beadTest[,"adjustedp"] beadTestMeasure = (1-beadTestP)*beadTestFC beadTest = cbind(beadTestID, beadTestMeasure) colnames(beadTest) <- c("ArrayAddressID", "Symbol", "difexBL") sortBL = sort(-beadTest[,"difexBL"], index.return=TRUE)$ix beadTop10 = beadTest[sortBL[1:10],] print(beadTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
if(require("blimaTestingData") && require("illuminaHumanv4.db") && interactive()) { #To perform background correction, quantile normalization and then bead level t-test on log data run. Vst is not performed in this scheme. Top 10 probes is then printed according to certain measure. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() groups1Mod = list() groups2Mod = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) groups1Mod[[i]] = p$Group %in% groups1; groups2Mod[[i]] = p$Group %in% groups2; c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod =c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = quantileNormalize(blimatesting, normalizationMod=c, channelNormalize="GrnF", channelOutput="qua", channelInclude="bgf") beadTest <- doTTests(blimatesting, groups1Mod, groups2Mod, transformation=log2TransformPositive, quality="qua", channelInclude="bgf") symbol2address <- merge(toTable(illuminaHumanv4ARRAYADDRESS), toTable(illuminaHumanv4SYMBOLREANNOTATED)) symbol2address <- symbol2address[,c("SymbolReannotated", "ArrayAddress") ] colnames(symbol2address) <- c("Symbol", "ArrayAddressID") beadTest = merge(beadTest, symbol2address, by.x="ProbeID", by.y="ArrayAddressID") beadTestID = beadTest[,c("ProbeID", "Symbol")] beadTestFC = abs(beadTest[,"mean1"]-beadTest[,"mean2"]) beadTestP = beadTest[,"adjustedp"] beadTestMeasure = (1-beadTestP)*beadTestFC beadTest = cbind(beadTestID, beadTestMeasure) colnames(beadTest) <- c("ArrayAddressID", "Symbol", "difexBL") sortBL = sort(-beadTest[,"difexBL"], index.return=TRUE)$ix beadTop10 = beadTest[sortBL[1:10],] print(beadTop10) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and illuminaHumanv4.db by running BiocManager::install('illuminaHumanv4.db')."); }
This function processes arrays in the object beadLevelData from package beadarray and returns sorted double vector. The vector has length prvku. And the distribution of this vector is a "mean" of all distributions of distributionChannel quantity in arrays. In case that probe numbers are different from prvku it does some averaging.
meanDistribution(b, normalizationMod = NULL, distributionChannel = "Grn", channelInclude = NULL, prvku)
meanDistribution(b, normalizationMod = NULL, distributionChannel = "Grn", channelInclude = NULL, prvku)
b |
Object beadLevelData from package beadarray or list of these objects |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes (defaults to NULL). |
distributionChannel |
Channel to do mean distribution from (defaults to "Grn"). |
channelInclude |
This field allows user to set channel with weights which have to be in 0,1. All zero weighted items are excluded from quantile normalization and the value asigned to such probes is a close to value which would be assigned to them if not being excluded. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to NULL). |
prvku |
Number of items in a resulting double vector. Prvku must not be more than minimal number of indluded items in any distributionChannel. |
Vojtěch Kulvait
INTERNAL This function is not intended for direct use. Background correction according to non parametric estimator in Xie, Yang, Xinlei Wang, and Michael Story. "Statistical Methods of Background Correction for Illumina BeadArray Data." Bioinformatics 25, no. 6 (March 15, 2009): 751-57. doi:10.1093/bioinformatics/btp040. The method is applied on the bead level.
nonParametricEstimator(toCorrectAll, toCorrectNeg)
nonParametricEstimator(toCorrectAll, toCorrectNeg)
toCorrectAll |
|
toCorrectNeg |
Vojtěch Kulvait
Correction for positive values only
nonPositiveCorrect(b, normalizationMod = NULL, channelCorrect = "GrnF", channelBackgroundFilter = "bgf", channelAndVector = NULL)
nonPositiveCorrect(b, normalizationMod = NULL, channelCorrect = "GrnF", channelBackgroundFilter = "bgf", channelAndVector = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelCorrect |
Name of channel to correct. |
channelBackgroundFilter |
Filtered beads will have weight 0 and non filtered weight 1. |
channelAndVector |
Represents vector to bitvise multiple to the channelBackgroundFilter vector. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
INTERNAL FUNCTION Correction for positive values only
nonPositiveCorrectSingleArray(b, normalizationMod = NULL, channelCorrect = "GrnF", channelBackgroundFilter = "bgf", channelAndVector = NULL)
nonPositiveCorrectSingleArray(b, normalizationMod = NULL, channelCorrect = "GrnF", channelBackgroundFilter = "bgf", channelAndVector = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelCorrect |
Name of channel to correct. |
channelBackgroundFilter |
Filtered beads will have weight 0 and non filtered weight 1. |
channelAndVector |
Represents vector to bitvise multiple to the channelBackgroundFilter vector. |
Vojtěch Kulvait
Internal function
numberOfDistributionElements(b, normalizationMod = NULL, channelInclude = NULL)
numberOfDistributionElements(b, normalizationMod = NULL, channelInclude = NULL)
b |
Object beadLevelData from package beadarray or list of these objects |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelInclude |
Vojtěch Kulvait
INTERNAL This function is not intended for direct use. Background correction according to non parametric estimator in Xie, Yang, Xinlei Wang, and Michael Story. "Statistical Methods of Background Correction for Illumina BeadArray Data." Bioinformatics 25, no. 6 (March 15, 2009): 751-57. doi:10.1093/bioinformatics/btp040. ###The method is applied on the bead level.
performXieCorrection(value, alpha, mu, sigma)
performXieCorrection(value, alpha, mu, sigma)
value |
|
alpha |
|
mu |
|
sigma |
Vojtěch Kulvait
This function plots image of background distribution versus to foreground after background subtraction.
plotBackgroundImageAfterCorrection(b, index, channelForeground = "GrnF", channelBackground = "GrnB", SDMultiple = 3, includePearson = FALSE)
plotBackgroundImageAfterCorrection(b, index, channelForeground = "GrnF", channelBackground = "GrnB", SDMultiple = 3, includePearson = FALSE)
b |
Single beadLevelData object. |
index |
Index of spot to generate. |
channelForeground |
Name of channel of foreground. |
channelBackground |
Name of channel of background. |
SDMultiple |
Correct on this level. |
includePearson |
Include Pearson corelation. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #Write background images after correction. This function prints graph for condition D4. Call dev.off() to close. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) index = base::match("D4", p$Name) plotBackgroundImageAfterCorrection(blimatesting[[2]], index) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #Write background images after correction. This function prints graph for condition D4. Call dev.off() to close. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) index = base::match("D4", p$Name) plotBackgroundImageAfterCorrection(blimatesting[[2]], index) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
This function plots image of background distribution versus to foreground before background subtraction.
plotBackgroundImageBeforeCorrection(b, index, channelForeground = "GrnF", channelBackground = "GrnB", includePearson = FALSE)
plotBackgroundImageBeforeCorrection(b, index, channelForeground = "GrnF", channelBackground = "GrnB", includePearson = FALSE)
b |
Single beadLevelData object. |
index |
Index of spot to generate. |
channelForeground |
Name of channel of foreground. |
channelBackground |
Name of channel of background. |
includePearson |
Include Pearson corelation. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #Write background images before correction. This function prints graph for condition D4. Call dev.off() to close. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) index = base::match("D4", p$Name) plotBackgroundImageBeforeCorrection(blimatesting[[2]], index) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #Write background images before correction. This function prints graph for condition D4. Call dev.off() to close. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) index = base::match("D4", p$Name) plotBackgroundImageBeforeCorrection(blimatesting[[2]], index) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
This function does quantile normalization of object beadLevelData from package beadarray.
quantileNormalize(b, normalizationMod = NULL, channelNormalize = "Grn", channelOutput = "qua", channelInclude = NULL, dst)
quantileNormalize(b, normalizationMod = NULL, channelNormalize = "Grn", channelOutput = "qua", channelInclude = NULL, dst)
b |
Object beadLevelData from package beadarray or list of these objects |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelNormalize |
Name of channel to normalize. |
channelOutput |
Name of output normalized channel. |
channelInclude |
This field allows user to set channel with weights which have to be in 0,1. All zero weighted items are excluded from quantile normalization and the value asigned to such probes is a close to value which would be assigned to them if not being excluded. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to NULL). |
dst |
User can specify sorted vector which represents distribution that should be assigned to items. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction, variance stabilization and quantile normalization. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod = processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod = processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction, variance stabilization and quantile normalization. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod = processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod = processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
Internal function supporting doTTests function.
readToVector(what, from, length, quality)
readToVector(what, from, length, quality)
what |
Item to read. |
from |
From index. |
length |
Length of vector. |
quality |
Column. |
Vojtěch Kulvait
Function to transform channel data.
selectedChannelTransform(b, normalizationMod = NULL, channelTransformFrom, channelResult, transformation = NULL)
selectedChannelTransform(b, normalizationMod = NULL, channelTransformFrom, channelResult, transformation = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for performing on all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelTransformFrom |
Name of channel to transform. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
transformation |
Function of input data trasformation, default is NULL. Any function which for input value returns transformed value may be supplied. T-test then will be evaluated on transformed data, consider use log2TranformPositive. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction on blimatesting object for two groups. Background correction is followed by correction for non positive data. Array spots out of selected groups will not be processed. data(blimatesting) #Prepare logical vectors corresponding to conditions A and E. groups1 = "A"; groups2 = "E"; sampleNames = list() c = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) c[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod=c, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod=c, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
Function to transform channel data.
selectedChannelTransformSingleArray(b, normalizationMod = NULL, channelTransformFrom, channelResult, transformation)
selectedChannelTransformSingleArray(b, normalizationMod = NULL, channelTransformFrom, channelResult, transformation)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for performing on all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelTransformFrom |
Name of channel to transform. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
transformation |
Function of input data trasformation, default is NULL. Any function which for input value returns transformed value may be supplied. T-test then will be evaluated on transformed data, consider use log2TranformPositive. |
Vojtěch Kulvait
This function does quantile normalization of object beadLevelData from package beadarray. Internal function not intended to direct use. Please use quantileNormalize.
singleArrayNormalize(b, normalizationMod = NULL, channelNormalize = "Grn", channelOutput = "qua", channelInclude = NULL, dst)
singleArrayNormalize(b, normalizationMod = NULL, channelNormalize = "Grn", channelOutput = "qua", channelInclude = NULL, dst)
b |
Object beadLevelData from package beadarray |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b. |
channelNormalize |
Name of channel to normalize. |
channelOutput |
Name of output normalized channel. |
channelInclude |
This field allows user to set channel with weights which have to be in 0,1. All zero weighted items are excluded from quantile normalization and the value asigned to such probes is a close to value which would be assigned to them if not being excluded. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to NULL). |
dst |
This field must be sorted. It is a distribution of values to assign to ports. By default this distribution is computed using meanDistribution function. |
Vojtěch Kulvait
Test existence of channel slot based on logical list
singleChannelExistsIntegrityWithLogicalVector(b, spotsToCheck = NULL, slotToCheck, action = c("returnText", "warn", "error"))
singleChannelExistsIntegrityWithLogicalVector(b, spotsToCheck = NULL, slotToCheck, action = c("returnText", "warn", "error"))
b |
single beadLevelData object |
spotsToCheck |
NULL for check all spots from b. Otherwise specifies logical vector of the length equals to the number of arrays in b with TRUE for checking. |
slotToCheck |
Slot name to check |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
Check integrity of the logical object, internal.
singleCheckIntegrityLogicalVector(xx, b, action = c("returnText", "warn", "error"))
singleCheckIntegrityLogicalVector(xx, b, action = c("returnText", "warn", "error"))
xx |
Logical object compatible with b. |
b |
Single beadLevelData object. |
action |
What type of action is required in case of invalid object structure. Either return text different from TRUE, warn or error. |
Vojtěch Kulvait
Internal function
singleNumberOfDistributionElements(b, normalizationMod = NULL, channelInclude = NULL)
singleNumberOfDistributionElements(b, normalizationMod = NULL, channelInclude = NULL)
b |
Object beadLevelData from package beadarray |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
channelInclude |
Vojtěch Kulvait
This is internal function not intended to direct use. Updates mean distribution.
updateMeanDistribution(meanDistribution, srt, arraysUsed)
updateMeanDistribution(meanDistribution, srt, arraysUsed)
meanDistribution |
|
srt |
vector of sorted values |
arraysUsed |
number of arrays allready used to create distribution |
Vojtěch Kulvait
This function does variance stabilising step on bead level.
varianceBeadStabilise(b, normalizationMod = NULL, quality = "qua", channelInclude = "bgf", channelOutput = "vst")
varianceBeadStabilise(b, normalizationMod = NULL, quality = "qua", channelInclude = "bgf", channelOutput = "vst")
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equal to the number of arrays in b or list of such vectors if b is a list of beadLevelData classes. |
quality |
Quality to analyze, default is "qua". |
channelInclude |
This field allows user to set channel with weights which have to be in 0,1. All zero weighted items are excluded from t-test. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to "bgf"). |
channelOutput |
Output from VST. |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #To perform background correction, variance stabilization and quantile normalization. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod = processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod = processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #To perform background correction, variance stabilization and quantile normalization. data(blimatesting) #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(c). groups1 = "A"; groups2 = "E"; sampleNames = list() processingMod = list() for(i in 1:length(blimatesting)) { p = pData(blimatesting[[i]]@experimentData$phenoData) processingMod[[i]] = p$Group %in% c(groups1, groups2); sampleNames[[i]] = p$Name } #Background correction and quantile normalization followed by testing including log2TransformPositive transformation. blimatesting = bacgroundCorrect(blimatesting, normalizationMod = processingMod, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, normalizationMod = processingMod, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = varianceBeadStabilise(blimatesting, normalizationMod = processingMod, quality="GrnF", channelInclude="bgf", channelOutput="vst") blimatesting = quantileNormalize(blimatesting, normalizationMod = processingMod, channelNormalize="vst", channelOutput="qua", channelInclude="bgf") }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
This function is not intended to direct use it takes single beadLevelData object and do bead level variance stabilisation.
varianceBeadStabiliseSingleArray(b, normalizationMod = NULL, quality = "qua", channelInclude = "bgf", channelOutput = "vst")
varianceBeadStabiliseSingleArray(b, normalizationMod = NULL, quality = "qua", channelInclude = "bgf", channelOutput = "vst")
b |
Object beadLevelData. |
normalizationMod |
NULL for normalization of all input b. Otherwise specifies logical vector of the length equals to the number of arrays in b. |
quality |
Quality to analyze, default is "qua". |
channelInclude |
This field allows user to set channel with weights which have to be in 0,1. All zero weighted items are excluded from t-test. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to "bgf"). |
channelOutput |
Output from VST. |
Vojtěch Kulvait
This function is derived from copy and paste of lumi::vst function. Since lumi package has extensive imports I decided to hardcode this function to the blima instead of importing lumi package.
vstFromLumi(u, std, nSupport = min(length(u), 500), backgroundStd = NULL, lowCutoff = 1/3)
vstFromLumi(u, std, nSupport = min(length(u), 500), backgroundStd = NULL, lowCutoff = 1/3)
u |
The mean of probe beads |
std |
The standard deviation of the probe beads |
nSupport |
Something for c3 guess. |
backgroundStd |
Estimate the background variance c3. Input should be variance according to article, not SD. |
lowCutoff |
Something for c3 guess. |
authors are Pan Du, Simon Lin, the function was edited by Vojtěch Kulvait
http://www.bioconductor.org/packages/release/bioc/html/lumi.html
This function writes images with background distribution according to foreground before and after background subtraction.
writeBackgroundImages(b, spotsToGenerate = NULL, imageType = c("jpg", "png", "eps"), channelForeground = "GrnF", channelBackground = "GrnB", SDMultiple = 3, includePearson = FALSE, outputDir = getwd(), width = 505, height = 505)
writeBackgroundImages(b, spotsToGenerate = NULL, imageType = c("jpg", "png", "eps"), channelForeground = "GrnF", channelBackground = "GrnB", SDMultiple = 3, includePearson = FALSE, outputDir = getwd(), width = 505, height = 505)
b |
Single beadLevelData object. |
spotsToGenerate |
NULL for generate images for all spots from b. Otherwise specifies logical vector of the length equals to the number of arrays in b with TRUE for images to generate. |
imageType |
Type of images produced, either jpg, png or eps |
channelForeground |
Name of channel of foreground. |
channelBackground |
Name of channel of background. |
SDMultiple |
Correct on this level. |
includePearson |
Include Pearson corelation. |
outputDir |
Directory where to output images. |
width |
Width of image (default 505 fits well for 86mm 150dpi illustration in Bioinformatics journal:) |
height |
Height of image |
Vojtěch Kulvait
if(require("blimaTestingData") && interactive()) { #Write background images before and after correction for background into /tmp directory. This function creates two jpg images for condition D. Output files are /tmp/6898481102_D_CORRECTED.jpg and /tmp/6898481102_D.jpg. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) spotsToGenerate = p$Group %in% "D"; writeBackgroundImages(blimatesting[[2]], imageType="jpg", spotsToGenerate=spotsToGenerate, includePearson=FALSE, outputDir="/tmp", width=505, height=505) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
if(require("blimaTestingData") && interactive()) { #Write background images before and after correction for background into /tmp directory. This function creates two jpg images for condition D. Output files are /tmp/6898481102_D_CORRECTED.jpg and /tmp/6898481102_D.jpg. data(blimatesting) p = pData(blimatesting[[2]]@experimentData$phenoData) spotsToGenerate = p$Group %in% "D"; writeBackgroundImages(blimatesting[[2]], imageType="jpg", spotsToGenerate=spotsToGenerate, includePearson=FALSE, outputDir="/tmp", width=505, height=505) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData')."); }
Background correction according to non parametric estimator in Xie, Yang, Xinlei Wang, and Michael Story. "Statistical Methods of Background Correction for Illumina BeadArray Data." Bioinformatics 25, no. 6 (March 15, 2009): 751-57. doi:10.1093/bioinformatics/btp040.###The method is applied on the bead level.
xieBacgroundCorrect(b, normalizationMod = NULL, negativeArrayAddresses, channelCorrect, channelResult, channelInclude = NULL)
xieBacgroundCorrect(b, normalizationMod = NULL, negativeArrayAddresses, channelCorrect, channelResult, channelInclude = NULL)
b |
List of beadLevelData objects (or single object). |
normalizationMod |
NULL for processing all spots in b. Otherwise specifies logical vector of the length equals to the number of arrays in b. |
negativeArrayAddresses |
Vector of addresses of negative control probes on array |
channelCorrect |
Slot to perform convolution correction. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
channelInclude |
This field allows user to set channel with weights which have to be from 0,1. All zero weighted items are excluded from summarization. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to NULL). |
Vojtěch Kulvait
if(require("blimaTestingData") && exists("annotationHumanHT12V4") && interactive()) { #Create vector of negative array addresses. negAdr = unique(annotationHumanHT12V4$Controls[annotationHumanHT12V4$Controls$Reporter_Group_Name=="negative", "Array_Address_Id"]) #Create summarization of nonnormalized data from GrnF column. data(blimatesting) blimatesting = bacgroundCorrect(blimatesting, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = xieBacgroundCorrect(blimatesting, negativeArrayAddresses=negAdr, channelCorrect="GrnF", channelResult="GrnFXIE", channelInclude="bgf") #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). xiecorrected = createSummarizedMatrix(blimatesting, quality="GrnFXIE", channelInclude="bgf", annotationTag="Name") head(xiecorrected) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and prepare annotationHumanHT12V4 object according to blimaTestingData manual."); }
if(require("blimaTestingData") && exists("annotationHumanHT12V4") && interactive()) { #Create vector of negative array addresses. negAdr = unique(annotationHumanHT12V4$Controls[annotationHumanHT12V4$Controls$Reporter_Group_Name=="negative", "Array_Address_Id"]) #Create summarization of nonnormalized data from GrnF column. data(blimatesting) blimatesting = bacgroundCorrect(blimatesting, channelBackgroundFilter="bgf") blimatesting = nonPositiveCorrect(blimatesting, channelCorrect="GrnF", channelBackgroundFilter="bgf", channelAndVector="bgf") blimatesting = xieBacgroundCorrect(blimatesting, negativeArrayAddresses=negAdr, channelCorrect="GrnF", channelResult="GrnFXIE", channelInclude="bgf") #Prepare logical vectors corresponding to conditions A(groups1Mod), E(groups2Mod) and both(processingMod). xiecorrected = createSummarizedMatrix(blimatesting, quality="GrnFXIE", channelInclude="bgf", annotationTag="Name") head(xiecorrected) }else { print("To run this example, please install blimaTestingData package from bioconductor by running BiocManager::install('blimaTestingData') and prepare annotationHumanHT12V4 object according to blimaTestingData manual."); }
INTERNAL This function is not intended for direct use. Background correction according to non parametric estimator in Xie, Yang, Xinlei Wang, and Michael Story. "Statistical Methods of Background Correction for Illumina BeadArray Data." Bioinformatics 25, no. 6 (March 15, 2009): 751-57. doi:10.1093/bioinformatics/btp040. The method is applied on the bead level.
xieBacgroundCorrectSingleArray(b, normalizationMod = NULL, negativeArrayAddresses, channelCorrect, channelResult, channelInclude = NULL)
xieBacgroundCorrectSingleArray(b, normalizationMod = NULL, negativeArrayAddresses, channelCorrect, channelResult, channelInclude = NULL)
b |
Single beadLevelData object. |
normalizationMod |
NULL for processing all spots in b. Otherwise specifies logical vector of the length equals to the number of arrays in b. |
negativeArrayAddresses |
Vector of addresses of negative control probes on array |
channelCorrect |
Slot to perform convolution correction. |
channelResult |
Result channel, if this channel exists it will be overwritten. |
channelInclude |
This field allows user to set channel with weights which have to be from 0,1. All zero weighted items are excluded from summarization. You can turn this off by setting this NULL. This option may be used together with bacgroundCorrect method or/and with beadarray QC (defaults to NULL). |
Vojtěch Kulvait