pub struct Record { /* private fields */ }Expand description
Represents a record (a line or a site) in BCF file
Implementations§
Source§impl Record
impl Record
Sourcepub fn read<R>(&mut self, reader: &mut R) -> Result<()>where
R: Read + ReadBytesExt,
pub fn read<R>(&mut self, reader: &mut R) -> Result<()>where
R: Read + ReadBytesExt,
read a record (copy bytes from the reader to the record’s interval buffers), and separate fields
Sourcepub fn chrom(&self) -> i32
pub fn chrom(&self) -> i32
get chromosome offset Example:
use bcf_reader::*;
use std::io::Write;
// read data generated by bcftools
// bcftools query -f '%CHROM\n' test.bcf | bgzip -c > test_chrom.gz
let mut chrom_str = String::new();
smart_reader("testdata/test_chrom.gz")
.unwrap()
.read_to_string(&mut chrom_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let s = read_header(&mut f).unwrap();
let header = Header::from_string(&s).unwrap();
let mut record = Record::default();
let mut chrom_str2 = Vec::<u8>::new();
while let Ok(_) = record.read(&mut f) {
write!(
chrom_str2,
"{}\n",
header.get_chrname(record.chrom() as usize)
)
.unwrap();
}
let chrom_str2 = String::from_utf8(chrom_str2).unwrap();
// compare bcftools results and bcf-reader results
assert_eq!(chrom_str, chrom_str2);pub fn n_allele(&self) -> u16
Sourcepub fn fmt_gt(&self, header: &Header) -> NumericValueIter<'_> ⓘ
pub fn fmt_gt(&self, header: &Header) -> NumericValueIter<'_> ⓘ
Returns an iterator over the genotype values in the record’s FORMAT field. If no FORMAT/GT field available, the returned iterator will have items. Example:
use bcf_reader::*;
use std::io::Write;
// read data generated by bcftools
// bcftools query -f '[\t%GT]\n' test.bcf | bgzip -c > test_gt.gz
let mut gt_str = String::new();
smart_reader("testdata/test_gt.gz")
.unwrap()
.read_to_string(&mut gt_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let s = read_header(&mut f).unwrap();
let header = Header::from_string(&s).unwrap();
let mut record = Record::default();
let mut gt_str2 = Vec::<u8>::new();
while let Ok(_) = record.read(&mut f) {
for (i, bn) in record.fmt_gt(&header).enumerate() {
let (noploidy, dot, phased, allele) = bn.unwrap().gt_val();
assert_eq!(noploidy, false); // missing ploidy
let mut sep = '\t';
if i % 2 == 1 {
if phased {
sep = '|';
} else {
sep = '/';
}
}
if dot {
write!(gt_str2, "{sep}.").unwrap();
} else {
write!(gt_str2, "{sep}{allele}").unwrap();
}
}
write!(gt_str2, "\n").unwrap();
}
let gt_str2 = String::from_utf8(gt_str2).unwrap();
// compare bcftools results and bcf-reader results
for (a, b) in gt_str
.split(|c| (c == '\n') || (c == '\t'))
.zip(gt_str2.split(|c| (c == '\n') || (c == '\t')))
{
assert_eq!(a, b);
}Sourcepub fn fmt_field(&self, fmt_key: usize) -> NumericValueIter<'_> ⓘ
pub fn fmt_field(&self, fmt_key: usize) -> NumericValueIter<'_> ⓘ
Returns an iterator over all values for a field in the record’s FORMATs (indiv).
Example:
use bcf_reader::*;
use std::io::Write;
// read data generated by bcftools
// bcftools query -f '[\t%AD]\n' test.bcf | bgzip -c > test_ad.gz
let mut ad_str = String::new();
smart_reader("testdata/test_ad.gz")
.unwrap()
.read_to_string(&mut ad_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let s = read_header(&mut f).unwrap();
let header = Header::from_string(&s).unwrap();
let mut record = Record::default();
let mut ad_str2 = Vec::<u8>::new();
let ad_filed_key = header.get_idx_from_dictionary_str("FORMAT", "AD").unwrap();
while let Ok(_) = record.read(&mut f) {
for (i, val) in record.fmt_field(ad_filed_key).enumerate() {
let val = val.unwrap();
if i % record.n_allele() as usize == 0 {
if ad_str2.last().map(|c| *c == b',') == Some(true) {
ad_str2.pop(); // trim last allele separator
}
ad_str2.push(b'\t'); // sample separator
}
match val.int_val() {
None => {}
Some(ad) => {
write!(ad_str2, "{ad},").unwrap(); // allele separator
}
}
}
// site separator
*ad_str2.last_mut().unwrap() = b'\n'; // sample separator
}
let ad_str2 = String::from_utf8(ad_str2).unwrap();
// compare bcftools results and bcf-reader results
for (a, b) in ad_str
.split(|c| (c == '\n') || (c == '\t'))
.zip(ad_str2.split(|c| (c == '\n') || (c == '\t')))
{
assert_eq!(a, b);
}Sourcepub fn pos(&self) -> i32
pub fn pos(&self) -> i32
get 0-based position (bp) value Example:
use bcf_reader::*;
// read data generated by bcftools
// bcftools query -f '%POS\n' test.bcf | bgzip -c > test_pos.gz
let mut pos_str = String::new();
smart_reader("testdata/test_pos.gz")
.unwrap()
.read_to_string(&mut pos_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let _s = read_header(&mut f).unwrap();
let mut record = Record::default();
let mut pos_str2 = Vec::<u8>::new();
use std::io::Write;
while let Ok(_) = record.read(&mut f) {
write!(pos_str2, "{}\n", record.pos() + 1).unwrap();
}
let pos_str2 = String::from_utf8(pos_str2).unwrap();
// compare bcftools results and bcf-reader results
assert_eq!(pos_str, pos_str2);Sourcepub fn id(&self) -> Result<&str>
pub fn id(&self) -> Result<&str>
get variant ID as &str Example:
use bcf_reader::*;
// read data generated by bcftools
// bcftools query -f '%ID\n' test.bcf | bgzip -c > test_id.gz
let mut id_str = String::new();
smart_reader("testdata/test_id.gz")
.unwrap()
.read_to_string(&mut id_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let _s = read_header(&mut f).unwrap();
let mut record = Record::default();
let mut id_str2 = Vec::<u8>::new();
use std::io::Write;
while let Ok(_) = record.read(&mut f) {
let record_id = record.id().unwrap();
let id = if record_id.is_empty() { "." } else { record_id };
write!(id_str2, "{}\n", id).unwrap();
}
let id_str2 = String::from_utf8(id_str2).unwrap();
// compare bcftools results and bcf-reader results
assert_eq!(id_str, id_str2);Sourcepub fn alleles(&self) -> &[Range<usize>]
pub fn alleles(&self) -> &[Range<usize>]
Returns the ranges of bytes in buf_shared for all alleles in the record. Example:
use bcf_reader::*;
// read data generated by bcftools
// bcftools query -f '%REF,%ALT\n' test.bcf | bgzip -c > test_allele.gz
let mut allele_str = String::new();
smart_reader("testdata/test_allele.gz")
.unwrap()
.read_to_string(&mut allele_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test.bcf").unwrap();
let _s = read_header(&mut f).unwrap();
let mut record = Record::default();
let mut allele_str2 = Vec::<u8>::new();
while let Ok(_) = record.read(&mut f) {
for rng in record.alleles().iter() {
let slice = &record.buf_shared()[rng.clone()];
allele_str2.extend(slice);
allele_str2.push(b',');
}
*allele_str2.last_mut().unwrap() = b'\n';
}
let allele_str2 = String::from_utf8(allele_str2).unwrap();
// compare bcftools results and bcf-reader results
assert_eq!(allele_str, allele_str2);Sourcepub fn info_field_numeric(&self, info_key: usize) -> NumericValueIter<'_> ⓘ
pub fn info_field_numeric(&self, info_key: usize) -> NumericValueIter<'_> ⓘ
Return an iterator of numeric values for an INFO/xxx field. If the key is not found, the returned iterator will have a zero length.
Example:
use bcf_reader::*;
use std::io::Write;
// read data generated by bcftools
// bcftools query -f '%INFO/AF\n' testdata/test2.bcf | bgzip -c > testdata/test2_info_af.gz
let mut info_af_str = String::new();
smart_reader("testdata/test2_info_af.gz")
.unwrap()
.read_to_string(&mut info_af_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test2.bcf").unwrap();
let s = read_header(&mut f).unwrap();
let header = Header::from_string(&s).unwrap();
let mut record = Record::default();
let mut info_af_str2 = Vec::<u8>::new();
let info_af_key = header.get_idx_from_dictionary_str("INFO", "AF").unwrap();
while let Ok(_) = record.read(&mut f) {
record.info_field_numeric(info_af_key).for_each(|nv| {
let af = nv.unwrap().float_val().unwrap();
write!(info_af_str2, "{af},").unwrap();
});
*info_af_str2.last_mut().unwrap() = b'\n'; // line separators
}
let filter_str2 = String::from_utf8(info_af_str2).unwrap();
assert_eq!(info_af_str, filter_str2);Sourcepub fn info_field_str(&self, info_key: usize) -> Result<Option<&str>>
pub fn info_field_str(&self, info_key: usize) -> Result<Option<&str>>
Return str value for an INFO/xxx field. If the key is not found or data type is not string, then return None.
Sourcepub fn filters(&self) -> NumericValueIter<'_> ⓘ
pub fn filters(&self) -> NumericValueIter<'_> ⓘ
iterate an integer for each filter key. If the length of the iterator is 0, it means no filter label is set. Example:
use bcf_reader::*;
use std::io::Write;
// read data generated by bcftools
// bcftools query -f '%FILTER\n' testdata/test2.bcf | bgzip -c > testdata/test2_filters.gz
let mut filter_str = String::new();
smart_reader("testdata/test2_filters.gz")
.unwrap()
.read_to_string(&mut filter_str)
.unwrap();
// read data via bcf-reader
let mut f = smart_reader("testdata/test2.bcf").unwrap();
let s = read_header(&mut f).unwrap();
let header = Header::from_string(&s).unwrap();
let mut record = Record::default();
let mut filter_str2 = Vec::<u8>::new();
let d = header.dict_strings();
while let Ok(_) = record.read(&mut f) {
record.filters().for_each(|nv| {
let nv = nv.unwrap();
let filter_key = nv.int_val().unwrap() as usize;
let dict_string_map = &d[&filter_key];
let filter_name = &dict_string_map["ID"];
write!(filter_str2, "{filter_name};").unwrap();
});
*filter_str2.last_mut().unwrap() = b'\n'; // line separators
}
let filter_str2 = String::from_utf8(filter_str2).unwrap();
// compare bcftools results and bcf-reader results
assert_eq!(filter_str, filter_str2);Sourcepub fn buf_indiv(&self) -> &[u8] ⓘ
pub fn buf_indiv(&self) -> &[u8] ⓘ
Returns the buffer containing indv (sample-level) information
Returns the buffer containing the shared (site-level) information
Trait Implementations§
Auto Trait Implementations§
impl Freeze for Record
impl RefUnwindSafe for Record
impl Send for Record
impl Sync for Record
impl Unpin for Record
impl UnwindSafe for Record
Blanket Implementations§
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left is true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self into a Left variant of Either<Self, Self>
if into_left(&self) returns true.
Converts self into a Right variant of Either<Self, Self>
otherwise. Read more