UMI4C | R Documentation |
The UMI4C constructor is the function makeUMI4C
. By using
the arguments listed below, performs the necessary steps to analyze UMI-4C
data and summarize it in an object of class UMI4C.
makeUMI4C(
colData,
viewpoint_name = "Unknown",
grouping = "condition",
normalized = TRUE,
ref_umi4c = NULL,
bait_exclusion = 3000,
bait_expansion = 1e+06,
scales = 5:150,
min_win_factor = 0.02,
sd = 2
)
colData |
Data.frame containing the information for constructing the UMI4C experiment object. Needs to contain the following columns:
|
viewpoint_name |
Character indicating the name for the used viewpoint. |
grouping |
Name of the column in colData used to merge the samples or replicates. Set to NULL for skipping grouping. Default: "condition". |
normalized |
Logical indicating whether UMI-4C profiles should be
normalized to the |
ref_umi4c |
Name of the sample or group to use as reference for
normalization. By default is NULL, which means it will use the sample with
less UMIs in the analyzed region. It should be a named vector, where the name
corresponds to the grouping column from |
bait_exclusion |
Region around the bait (in bp) to be excluded from the analysis. Default: 3000bp. |
bait_expansion |
Number of bp upstream and downstream of the bait to use for the analysis (region centered in bait). Default: 1Mb. |
scales |
Numeric vector containing the scales for calculating the domainogram. |
min_win_factor |
Proportion of UMIs that need to be found in a specific window for adaptative trend calculation |
sd |
Stantard deviation for adaptative trend. |
It returns an object of the class UMI4C.
colData
Data.frame containing the information for constructing the UMI4C experiment object. Needs to contain the following columns:
sampleID
: Unique identifier for the sample.
condition
: Condition for performing differential analysis.
Can be control and treatment, two different cell types, etc.
replicate
: Number or ID for identifying different replicates.
file
: Path to the files outputed by contactsUMI4C
.
rowRanges
GRanges object with the coordinates for the restriction fragment ends, their IDs and other additional annotation columns.
metadata
List containing the following elements:
bait
: GRanges object representing the position
of the bait used for the analysis.
scales
: Numeric vector containing the scales used for
calculating the domainogram.
min_win_factor
: Factor for calculating the minimum molecules
requiered in a window for not merging it with the next one when
calculating the adaptative smoothing trend.
grouping
: Columns in colData
used for the different
sample groupings, accessible through groupsUMI4C
.
normalized
: Logical indicating whether samples/groups are
normalized or not.
region
: GRanges with the coordinates of
the genomic window used for analyzing UMI4C data.
ref_umi4c
: Name of the sample or group used as reference for
normalization.
assays
Matrix where each row represents a restriction fragment site
and columns represent each sample or group defined in grouping
.
After running the makeUMI4C
function, it will contain the
following data:
umis
: Raw number of UMIs detected by contactsUMI4C
.
norm_mat
: Normalization factors for each sample/group and fragment end.
trend
: Adaptative smoothing trend of UMIs.
geo_coords
: Geometric coordinates obtained when performing
the adaptative smoothing.
scale
: Scale selected for the adaptative smoothing.
sd
: Stantard deviation for the adaptative smoothing trend.
dgram
List containing the domainograms for each sample. A domainogram is matrix where columns are different scales selected for merging UMI counts and rows are the restriction fragments.
groupsUMI4C
List of UMI4C
objects with the specific groupings.
results
List containing the results for the differential analysis ran
using fisherUMI4C
.
The UMI4C
class extends the SummarizedExperiment class.
UMI4C-methods
# Load sample processed file paths
files <- list.files(system.file("extdata", "CIITA", "count",
package = "UMI4Cats"
),
pattern = "*_counts.tsv",
full.names = TRUE
)
# Create colData including all relevant information
colData <- data.frame(
sampleID = gsub("_counts.tsv.gz", "", basename(files)),
file = files,
stringsAsFactors = FALSE
)
library(tidyr)
colData <- colData %>%
separate(sampleID,
into = c("condition", "replicate", "viewpoint"),
remove = FALSE
)
# Load UMI-4C data and generate UMI4C object
umi <- makeUMI4C(
colData = colData,
viewpoint_name = "CIITA",
grouping = "condition"
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.