[][src]Module rust_htslib::bcf

Module for working with VCF and BCF files.

Performance Remarks

Note that BCF corresponds to the in-memory representation of BCF/VCF records in Htslib itself. Thus, it comes without a runtime penalty for parsing, in contrast to reading VCF files.

Example

  • Obtaining 0-based locus index of the VCF record.
  • Obtaining alleles of the VCF record.
  • calculate alt-allele dosage in a mutli-sample VCF / BCF
use crate::rust_htslib::bcf::{Reader, Read};
use std::convert::TryFrom;

let path = &"test/test_string.vcf";
let mut bcf = Reader::from_path(path).expect("Error opening file.");
// iterate through each row of the vcf body.
for (i, record_result) in bcf.records().enumerate() {
    let mut record = record_result.expect("Fail to read record");
    let mut s = String::new();
     for allele in record.alleles() {
         for c in allele {
             s.push(char::from(*c))
         }
         s.push(' ')
     }
    // 0-based position and the list of alleles
    println!("Locus: {}, Alleles: {}", record.pos(), s);
    // number of sample in the vcf
    let sample_count = usize::try_from(record.sample_count()).unwrap();

    // Counting ref, alt and missing alleles for each sample
    let mut n_ref = vec![0; sample_count];
    let mut n_alt = vec![0; sample_count];
    let mut n_missing = vec![0; sample_count];
    let gts = record.genotypes().expect("Error reading genotypes");
    for sample_index in 0..sample_count {
        // for each sample
        for gta in gts.get(sample_index).iter() {
            // for each allele
            match gta.index() {
                Some(0) => n_ref[sample_index] += 1,  // reference allele
                Some(_) => n_alt[sample_index] += 1,  // alt allele
                None => n_missing[sample_index] += 1, // missing allele
            }
        }
    }
}

Re-exports

pub use crate::bcf::header::Header;
pub use crate::bcf::header::HeaderRecord;
pub use crate::bcf::record::Record;

Modules

buffer
header

Module for working with VCF or BCF headers.

record
synced

This module contains the SyncedReader class and related code.

Structs

IndexedReader

An indexed VCF/BCF reader.

Reader

A VCF/BCF reader.

Records
Writer

A VCF/BCF writer.

Enums

Format

Traits

Read

A trait for a BCF reader with a read method.

Functions

set_threads

Safety