This convenience function removes non-standard, mitochondrial, and/or sex chromosomes from any GRanges object.
tidyChromosomes(
  gr,
  keep.X = TRUE,
  keep.Y = TRUE,
  keep.M = FALSE,
  keep.nonstandard = FALSE,
  genome = NULL
)Any GRanges object, or any another object with associated
seqinfo (or a Seqinfo object itself). The object should
typically have a standard genome associated with it, e.g. genome(gr)
<- "hg38". gr can also be a list of such GRanges objects.
Logicals indicating which non-autosomes should be kept. By default, sex chromosomes are kept, but mitochondrial and non-standard chromosomes are removed.
An optional string that, if supplied, will be used to set the
genome of gr.
A GRanges object in which both ranges and seqinfo associated
  with trimmed chromosomes have been removed.
Standard chromosomes are defined using the
  standardChromosomes function
  from the GenomeInfoDb package.
# make a GRanges
chrom <- c("chr2", "chr3", "chrX", "chrY", "chrM", "junk")
gr <- GRanges(seqnames = chrom,
              ranges = IRanges(start = 2*(1:6), end = 3*(1:6)),
              strand = "+",
              seqinfo = Seqinfo(chrom))
genome(gr) <- "hg38"
gr
#> GRanges object with 6 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2       2-3      +
#>   [2]     chr3       4-6      +
#>   [3]     chrX       6-9      +
#>   [4]     chrY      8-12      +
#>   [5]     chrM     10-15      +
#>   [6]     junk     12-18      +
#>   -------
#>   seqinfo: 6 sequences from hg38 genome; no seqlengths
tidyChromosomes(gr)
#> GRanges object with 4 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2       2-3      +
#>   [2]     chr3       4-6      +
#>   [3]     chrX       6-9      +
#>   [4]     chrY      8-12      +
#>   -------
#>   seqinfo: 4 sequences from hg38 genome; no seqlengths
tidyChromosomes(gr, keep.M = TRUE)
#> GRanges object with 5 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2       2-3      +
#>   [2]     chr3       4-6      +
#>   [3]     chrX       6-9      +
#>   [4]     chrY      8-12      +
#>   [5]     chrM     10-15      +
#>   -------
#>   seqinfo: 5 sequences from hg38 genome; no seqlengths
tidyChromosomes(gr, keep.M = TRUE, keep.Y = FALSE)
#> GRanges object with 4 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2       2-3      +
#>   [2]     chr3       4-6      +
#>   [3]     chrX       6-9      +
#>   [4]     chrM     10-15      +
#>   -------
#>   seqinfo: 4 sequences from hg38 genome; no seqlengths
tidyChromosomes(gr, keep.nonstandard = TRUE)
#> GRanges object with 5 ranges and 0 metadata columns:
#>       seqnames    ranges strand
#>          <Rle> <IRanges>  <Rle>
#>   [1]     chr2       2-3      +
#>   [2]     chr3       4-6      +
#>   [3]     chrX       6-9      +
#>   [4]     chrY      8-12      +
#>   [5]     junk     12-18      +
#>   -------
#>   seqinfo: 5 sequences from hg38 genome; no seqlengths