1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793
// Copyright 2014-2018 Johannes Köster, Christopher Schröder, Henning Timm.
// Licensed under the MIT license (http://opensource.org/licenses/MIT)
// This file may not be copied, modified, or distributed
// except according to those terms.
//! Structs and trait to read and write files in FASTA format.
//!
//! # Example
//!
//! ## Read
//!
//! In this example, we parse a fasta file from stdin and compute some statistics
//!
//! ```
//! use bio::io::fasta;
//! use std::io;
//!
//! let mut reader = fasta::Reader::new(io::stdin());
//!
//! let mut nb_reads = 0;
//! let mut nb_bases = 0;
//!
//! for result in reader.records() {
//! let record = result.expect("Error during fasta record parsing");
//! println!("{}", record.id());
//!
//! nb_reads += 1;
//! nb_bases += record.seq().len();
//! }
//!
//! println!("Number of reads: {}", nb_reads);
//! println!("Number of bases: {}", nb_bases);
//! ```
//!
//! We can also use a `while` loop to iterate over records.
//! This is slightly faster than the `for` loop.
//! ```
//! use bio::io::fasta;
//! use std::io;
//! let mut records = fasta::Reader::new(io::stdin()).records();
//!
//! let mut nb_reads = 0;
//! let mut nb_bases = 0;
//!
//! while let Some(Ok(record)) = records.next() {
//! nb_reads += 1;
//! nb_bases += record.seq().len();
//! }
//!
//! println!("Number of reads: {}", nb_reads);
//! println!("Number of bases: {}", nb_bases);
//! ```
//!
//! ## Write
//!
//! In this example we generate 10 random sequences with length 100 and write them to stdout.
//!
//! ```
//! use std::io;
//! use bio::io::fasta;
//!
//! let mut seed = 42;
//!
//! let nucleotides = [b'A', b'C', b'G', b'T'];
//!
//! let mut writer = fasta::Writer::new(io::stdout());
//!
//! for _ in 0..10 {
//! let seq = (0..100).map(|_| {
//! seed = ((seed ^ seed << 13) ^ seed >> 7) ^ seed << 17; // don't use this random generator
//! nucleotides[seed % 4]
//! }).collect::<Vec<u8>>();
//!
//! writer.write("random", None, seq.as_slice()).expect("Error writing record.");
//! }
//! ```
//!
//! ## Read and Write
//!
//! In this example we filter reads from stdin on sequence length and write them to stdout
//!
//! ```
//! use bio::io::fasta;
//! use bio::io::fasta::FastaRead;
//! use std::io;
//!
//! let mut reader = fasta::Reader::new(io::stdin());
//! let mut writer = fasta::Writer::new(io::stdout());
//! let mut record = fasta::Record::new();
//!
//! while let Ok(()) = reader.read(&mut record) {
//! if record.is_empty() {
//! break;
//! }
//!
//! if record.seq().len() > 100 {
//! writer
//! .write_record(&record)
//! .ok()
//! .expect("Error writing record.");
//! }
//! }
//! ```
//!
//! ## Index
//!
//! Random access to FASTA files is facilitated by [`Index`] and [`IndexedReader`]. The FASTA files
//! must already be indexed with [`samtools faidx`](https://www.htslib.org/doc/faidx.html).
//!
//! In this example, we read in the first 10 bases of the sequence named "chr1".
//!
//! ```rust
//! use bio::io::fasta::IndexedReader;
//! // create dummy files
//! const FASTA_FILE: &[u8] = b">chr1\nGTAGGCTGAAAA\nCCCC";
//! const FAI_FILE: &[u8] = b"chr1\t16\t6\t12\t13";
//!
//! let seq_name = "chr1";
//! let start: u64 = 0; // start is 0-based, inclusive
//! let stop: u64 = 10; // stop is 0-based, exclusive
//! // load the index
//! let mut faidx = IndexedReader::new(std::io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
//! // move the pointer in the index to the desired sequence and interval
//! faidx
//! .fetch(seq_name, start, stop)
//! .expect("Couldn't fetch interval");
//! // read the subsequence defined by the interval into a vector
//! let mut seq = Vec::new();
//! faidx.read(&mut seq).expect("Couldn't read the interval");
//! assert_eq!(seq, b"GTAGGCTGAA");
//! ```
use std::cmp::min;
use std::collections;
use std::convert::AsRef;
use std::fs;
use std::io;
use std::io::prelude::*;
use std::path::Path;
use crate::utils::{Text, TextSlice};
use anyhow::Context;
use std::fmt;
/// Maximum size of temporary buffer used for reading indexed FASTA files.
const MAX_FASTA_BUFFER_SIZE: usize = 512;
/// Trait for FASTA readers.
pub trait FastaRead {
fn read(&mut self, record: &mut Record) -> io::Result<()>;
}
/// A FASTA reader.
#[derive(Default, Clone, Eq, PartialEq, Ord, PartialOrd, Hash, Debug, Serialize, Deserialize)]
pub struct Reader<B> {
reader: B,
line: String,
}
impl Reader<io::BufReader<fs::File>> {
/// Read FASTA from given file path.
pub fn from_file<P: AsRef<Path> + std::fmt::Debug>(path: P) -> anyhow::Result<Self> {
fs::File::open(&path)
.map(Reader::new)
.with_context(|| format!("Failed to read fasta from {:#?}", path))
}
/// Read FASTA from give file path and a capacity
pub fn from_file_with_capacity<P: AsRef<Path> + std::fmt::Debug>(
capacity: usize,
path: P,
) -> anyhow::Result<Self> {
fs::File::open(&path)
.map(|file| Reader::with_capacity(capacity, file))
.with_context(|| format!("Failed to read fasta from {:#?}", path))
}
}
impl<R> Reader<io::BufReader<R>>
where
R: io::Read,
{
/// Create a new Fasta reader given an instance of `io::Read`.
///
/// # Example
/// ```rust
/// # use std::io;
/// # use bio::io::fasta::Reader;
/// # fn main() {
/// # const fasta_file: &'static [u8] = b">id desc
/// # AAAA
/// # ";
/// let reader = Reader::new(fasta_file);
/// # }
/// ```
pub fn new(reader: R) -> Self {
Reader {
reader: io::BufReader::new(reader),
line: String::new(),
}
}
/// Create a new Fasta reader given a capacity and an instance of `io::Read`.
///
/// # Example
/// ```rust
/// # use std::io;
/// # use bio::io::fasta::Reader;
/// # fn main() {
/// # const fasta_file: &'static [u8] = b">id desc
/// # AAAA
/// # ";
/// let reader = Reader::with_capacity(16384, fasta_file);
/// # }
/// ```
pub fn with_capacity(capacity: usize, reader: R) -> Self {
Reader {
reader: io::BufReader::with_capacity(capacity, reader),
line: String::new(),
}
}
}
impl<B> Reader<B>
where
B: io::BufRead,
{
/// Create a new Fasta reader with an object that implements `io::BufRead`.
///
/// # Example
/// ```rust
/// # use std::io;
/// # use bio::io::fasta::Reader;
/// # fn main() {
/// # const fasta_file: &'static [u8] = b">id desc
/// # AAAA
/// # ";
/// let buffer = io::BufReader::with_capacity(16384, fasta_file);
/// let reader = Reader::from_bufread(buffer);
/// # }
/// ```
pub fn from_bufread(bufreader: B) -> Self {
Reader {
reader: bufreader,
line: String::new(),
}
}
/// Return an iterator over the records of this Fasta file.
///
/// # Example
/// ```rust
/// # use std::io;
/// # use bio::io::fasta::Reader;
/// # use bio::io::fasta::Record;
/// # fn main() {
/// # const fasta_file: &'static [u8] = b">id desc
/// # AAAA
/// # ";
/// # let reader = Reader::new(fasta_file);
/// for record in reader.records() {
/// let record = record.unwrap();
/// assert_eq!(record.id(), "id");
/// assert_eq!(record.desc().unwrap(), "desc");
/// assert_eq!(record.seq().to_vec(), b"AAAA");
/// }
/// # }
/// ```
pub fn records(self) -> Records<B> {
Records {
reader: self,
error_has_occured: false,
}
}
}
impl<B> FastaRead for Reader<B>
where
B: io::BufRead,
{
/// Read the next FASTA record into the given `Record`.
/// An empty record indicates that no more records can be read.
///
/// Use this method when you want to read records as fast as
/// possible because it allows the reuse of a `Record` allocation.
///
/// The [records](Reader::records) iterator provides a more ergonomic
/// approach to accessing FASTA records.
///
/// # Errors
///
/// This function will return an error if the record is incomplete,
/// syntax is violated or any form of I/O error is encountered.
///
/// # Example
///
/// ```rust
/// use bio::io::fasta::Record;
/// use bio::io::fasta::{FastaRead, Reader};
///
/// const fasta_file: &'static [u8] = b">id desc
/// AAAA
/// ";
/// let mut reader = Reader::new(fasta_file);
/// let mut record = Record::new();
///
/// // Check for errors parsing the record
/// reader
/// .read(&mut record)
/// .expect("fasta reader: got an io::Error or could not read_line()");
///
/// assert_eq!(record.id(), "id");
/// assert_eq!(record.desc().unwrap(), "desc");
/// assert_eq!(record.seq().to_vec(), b"AAAA");
/// ```
fn read(&mut self, record: &mut Record) -> io::Result<()> {
record.clear();
if self.line.is_empty() {
self.reader.read_line(&mut self.line)?;
if self.line.is_empty() {
return Ok(());
}
}
if !self.line.starts_with('>') {
return Err(io::Error::new(
io::ErrorKind::Other,
"Expected > at record start.",
));
}
let mut header_fields = self.line[1..].trim_end().splitn(2, char::is_whitespace);
record.id = header_fields.next().map(|s| s.to_owned()).unwrap();
record.desc = header_fields.next().map(|s| s.to_owned());
loop {
self.line.clear();
self.reader.read_line(&mut self.line)?;
if self.line.is_empty() || self.line.starts_with('>') {
break;
}
record.seq.push_str(self.line.trim_end());
}
Ok(())
}
}
/// A FASTA index as created by SAMtools (.fai).
#[derive(Default, Clone, Eq, PartialEq, Debug, Serialize, Deserialize)]
pub struct Index {
inner: Vec<IndexRecord>,
name_to_rid: collections::HashMap<String, usize>,
}
impl Index {
/// Open a FASTA index from a given `io::Read` instance.
pub fn new<R: io::Read>(fai: R) -> csv::Result<Self> {
let mut inner = vec![];
let mut name_to_rid = collections::HashMap::new();
let mut fai_reader = csv::ReaderBuilder::new()
.delimiter(b'\t')
.has_headers(false)
.from_reader(fai);
for (rid, row) in fai_reader.deserialize().enumerate() {
let record: IndexRecord = row?;
name_to_rid.insert(record.name.clone(), rid);
inner.push(record);
}
Ok(Index { inner, name_to_rid })
}
/// Open a FASTA index from a given file path.
pub fn from_file<P: AsRef<Path> + std::fmt::Debug>(path: &P) -> anyhow::Result<Self> {
fs::File::open(&path)
.map_err(csv::Error::from)
.and_then(Self::new)
.with_context(|| format!("Failed to read fasta index from {:#?}", path))
}
/// Open a FASTA index given the corresponding FASTA file path.
/// That is, for ref.fasta we expect ref.fasta.fai.
pub fn with_fasta_file<P: AsRef<Path>>(fasta_path: &P) -> anyhow::Result<Self> {
let mut fai_path = fasta_path.as_ref().as_os_str().to_owned();
fai_path.push(".fai");
Self::from_file(&fai_path)
}
/// Return a vector of sequences described in the index.
pub fn sequences(&self) -> Vec<Sequence> {
// sort kv pairs by rid to preserve order
self.inner
.iter()
.map(|record| Sequence {
name: record.name.clone(),
len: record.len,
})
.collect()
}
}
/// A FASTA reader with an index as created by SAMtools (.fai).
#[derive(Debug)]
pub struct IndexedReader<R: io::Read + io::Seek> {
reader: io::BufReader<R>,
pub index: Index,
fetched_idx: Option<IndexRecord>,
start: Option<u64>,
stop: Option<u64>,
}
impl IndexedReader<fs::File> {
/// Read from a given file path. This assumes the index ref.fasta.fai to be
/// present for FASTA ref.fasta.
pub fn from_file<P: AsRef<Path> + std::fmt::Debug>(path: &P) -> anyhow::Result<Self> {
let index = Index::with_fasta_file(path)?;
fs::File::open(&path)
.map(|f| Self::with_index(f, index))
.map_err(csv::Error::from)
.with_context(|| format!("Failed to read fasta from {:#?}", path))
}
}
impl<R: io::Read + io::Seek> IndexedReader<R> {
/// Read from a FASTA and its index, both given as `io::Read`. FASTA has to
/// be `io::Seek` in addition.
pub fn new<I: io::Read>(fasta: R, fai: I) -> csv::Result<Self> {
let index = Index::new(fai)?;
Ok(IndexedReader {
reader: io::BufReader::new(fasta),
index,
fetched_idx: None,
start: None,
stop: None,
})
}
/// Read from a FASTA and its index, the first given as `io::Read`, the
/// second given as index object.
pub fn with_index(fasta: R, index: Index) -> Self {
IndexedReader {
reader: io::BufReader::new(fasta),
index,
fetched_idx: None,
start: None,
stop: None,
}
}
/// Fetch an interval from the sequence with the given name for reading.
///
/// `start` and `stop` are 0-based and `stop` is exclusive - i.e. `[start, stop)`
///
/// # Example
///
/// ```rust
/// use bio::io::fasta::IndexedReader;
/// // create dummy files
/// const FASTA_FILE: &[u8] = b">chr1\nGTAGGCTGAAAA\nCCCC";
/// const FAI_FILE: &[u8] = b"chr1\t16\t6\t12\t13";
///
/// let seq_name = "chr1";
/// let start: u64 = 0; // start is 0-based, inclusive
/// let stop: u64 = 10; // stop is 0-based, exclusive
/// // load the index
/// let mut faidx = IndexedReader::new(std::io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
/// // move the pointer in the index to the desired sequence and interval
/// faidx
/// .fetch(seq_name, start, stop)
/// .expect("Couldn't fetch interval");
/// // read the subsequence defined by the interval into a vector
/// let mut seq = Vec::new();
/// faidx.read(&mut seq).expect("Couldn't read the interval");
/// assert_eq!(seq, b"GTAGGCTGAA");
/// ```
///
/// # Errors
/// If the `seq_name` does not exist within the index.
pub fn fetch(&mut self, seq_name: &str, start: u64, stop: u64) -> io::Result<()> {
let idx = self.idx(seq_name)?;
self.start = Some(start);
self.stop = Some(stop);
self.fetched_idx = Some(idx);
Ok(())
}
/// Fetch an interval from the sequence with the given record index for reading.
///
/// `start` and `stop` are 0-based and `stop` is exclusive - i.e. `[start, stop)`
///
/// # Example
///
/// ```rust
/// use bio::io::fasta::IndexedReader;
/// // create dummy files
/// const FASTA_FILE: &[u8] = b">chr1\nGTAGGCTGAAAA\nCCCC";
/// const FAI_FILE: &[u8] = b"chr1\t16\t6\t12\t13";
///
/// let rid: usize = 0;
/// let start: u64 = 0; // start is 0-based, inclusive
/// let stop: u64 = 10; // stop is 0-based, exclusive
/// // load the index
/// let mut faidx = IndexedReader::new(std::io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
/// // move the pointer in the index to the desired sequence and interval
/// faidx
/// .fetch_by_rid(rid, start, stop)
/// .expect("Couldn't fetch interval");
/// // read the subsequence defined by the interval into a vector
/// let mut seq = Vec::new();
/// faidx.read(&mut seq).expect("Couldn't read the interval");
/// assert_eq!(seq, b"GTAGGCTGAA");
/// ```
///
/// # Errors
/// If `rid` does not exist within the index.
pub fn fetch_by_rid(&mut self, rid: usize, start: u64, stop: u64) -> io::Result<()> {
let idx = self.idx_by_rid(rid)?;
self.start = Some(start);
self.stop = Some(stop);
self.fetched_idx = Some(idx);
Ok(())
}
/// Fetch the whole sequence with the given name for reading.
pub fn fetch_all(&mut self, seq_name: &str) -> io::Result<()> {
let idx = self.idx(seq_name)?;
self.start = Some(0);
self.stop = Some(idx.len);
self.fetched_idx = Some(idx);
Ok(())
}
/// Fetch the whole sequence with the given record index for reading.
pub fn fetch_all_by_rid(&mut self, rid: usize) -> io::Result<()> {
let idx = self.idx_by_rid(rid)?;
self.start = Some(0);
self.stop = Some(idx.len);
self.fetched_idx = Some(idx);
Ok(())
}
/// Read the fetched sequence into the given vector.
pub fn read(&mut self, seq: &mut Text) -> io::Result<()> {
let idx = self.fetched_idx.clone();
match (idx, self.start, self.stop) {
(Some(idx), Some(start), Some(stop)) => self.read_into_buffer(idx, start, stop, seq),
_ => Err(io::Error::new(
io::ErrorKind::Other,
"No sequence fetched for reading.",
)),
}
}
/// Return an iterator yielding the fetched sequence.
pub fn read_iter(&mut self) -> io::Result<IndexedReaderIterator<'_, R>> {
let idx = self.fetched_idx.clone();
match (idx, self.start, self.stop) {
(Some(idx), Some(start), Some(stop)) => self.read_into_iter(idx, start, stop),
_ => Err(io::Error::new(
io::ErrorKind::Other,
"No sequence fetched for reading.",
)),
}
}
fn read_into_buffer(
&mut self,
idx: IndexRecord,
start: u64,
stop: u64,
seq: &mut Text,
) -> io::Result<()> {
if stop > idx.len {
return Err(io::Error::new(
io::ErrorKind::Other,
"FASTA read interval was out of bounds",
));
} else if start > stop {
return Err(io::Error::new(
io::ErrorKind::Other,
"Invalid query interval",
));
}
let mut bases_left = stop - start;
let mut line_offset = self.seek_to(&idx, start)?;
seq.clear();
while bases_left > 0 {
bases_left -= self.read_line(&idx, &mut line_offset, bases_left, seq)?;
}
Ok(())
}
fn read_into_iter(
&mut self,
idx: IndexRecord,
start: u64,
stop: u64,
) -> io::Result<IndexedReaderIterator<'_, R>> {
if stop > idx.len {
return Err(io::Error::new(
io::ErrorKind::Other,
"FASTA read interval was out of bounds",
));
} else if start > stop {
return Err(io::Error::new(
io::ErrorKind::Other,
"Invalid query interval",
));
}
let bases_left = stop - start;
let line_offset = self.seek_to(&idx, start)?;
let capacity = min(
MAX_FASTA_BUFFER_SIZE,
min(bases_left, idx.line_bases) as usize,
);
Ok(IndexedReaderIterator {
reader: self,
record: idx,
bases_left,
line_offset,
buf: Vec::with_capacity(capacity),
buf_idx: 0,
})
}
/// Return the IndexRecord for the given sequence name or io::Result::Err
fn idx(&self, seqname: &str) -> io::Result<IndexRecord> {
match self.index.name_to_rid.get(seqname) {
Some(rid) => self.idx_by_rid(*rid),
None => Err(io::Error::new(
io::ErrorKind::Other,
format!("Unknown sequence name: {}.", seqname),
)),
}
}
/// Return the IndexRecord for the given record index or io::Result::Err
fn idx_by_rid(&self, rid: usize) -> io::Result<IndexRecord> {
match self.index.inner.get(rid) {
Some(record) => Ok(record.clone()),
None => Err(io::Error::new(
io::ErrorKind::Other,
"Invalid record index in fasta file.",
)),
}
}
/// Seek to the given position in the specified FASTA record. The position
/// of the cursor on the line that the seek ended on is returned.
fn seek_to(&mut self, idx: &IndexRecord, start: u64) -> io::Result<u64> {
assert!(start <= idx.len);
let line_offset = start % idx.line_bases;
let line_start = start / idx.line_bases * idx.line_bytes;
let offset = idx.offset + line_start + line_offset;
self.reader.seek(io::SeekFrom::Start(offset))?;
Ok(line_offset)
}
/// Tries to read up to `bases_left` bases from the current line into `buf`,
/// returning the actual number of bases read. Depending on the amount of
/// whitespace per line, the current `line_offset`, and the amount of bytes
/// returned from `BufReader::fill_buf`, this function may return Ok(0)
/// multiple times in a row.
fn read_line(
&mut self,
idx: &IndexRecord,
line_offset: &mut u64,
bases_left: u64,
buf: &mut Vec<u8>,
) -> io::Result<u64> {
let (bytes_to_read, bytes_to_keep) = {
let src = self.reader.fill_buf()?;
if src.is_empty() {
return Err(io::Error::new(
io::ErrorKind::UnexpectedEof,
"FASTA file is truncated.",
));
}
let bases_on_line = idx.line_bases - min(idx.line_bases, *line_offset);
let bases_in_buffer = min(src.len() as u64, bases_on_line);
let (bytes_to_read, bytes_to_keep) = if bases_in_buffer <= bases_left {
let bytes_to_read = min(src.len() as u64, idx.line_bytes - *line_offset);
(bytes_to_read, bases_in_buffer)
} else {
(bases_left, bases_left)
};
buf.extend_from_slice(&src[..bytes_to_keep as usize]);
(bytes_to_read, bytes_to_keep)
};
self.reader.consume(bytes_to_read as usize);
assert!(bytes_to_read > 0);
*line_offset += bytes_to_read;
if *line_offset >= idx.line_bytes {
*line_offset = 0;
}
Ok(bytes_to_keep)
}
}
/// Record of a FASTA index.
#[derive(Clone, Eq, PartialEq, Debug, Serialize, Deserialize)]
struct IndexRecord {
name: String,
len: u64,
offset: u64,
line_bases: u64,
line_bytes: u64,
}
/// A sequence record returned by the FASTA index.
#[derive(Default, Clone, Eq, PartialEq, Ord, PartialOrd, Hash, Debug, Serialize, Deserialize)]
pub struct Sequence {
pub name: String,
pub len: u64,
}
#[derive(Debug)]
pub struct IndexedReaderIterator<'a, R: io::Read + io::Seek> {
reader: &'a mut IndexedReader<R>,
record: IndexRecord,
bases_left: u64,
line_offset: u64,
buf: Vec<u8>,
buf_idx: usize,
}
impl<'a, R: io::Read + io::Seek + 'a> IndexedReaderIterator<'a, R> {
fn fill_buffer(&mut self) -> io::Result<()> {
assert!(self.bases_left > 0);
self.buf.clear();
let bases_to_read = min(self.buf.capacity() as u64, self.bases_left);
// May loop one or more times; see IndexedReader::read_line.
while self.buf.is_empty() {
self.bases_left -= self.reader.read_line(
&self.record,
&mut self.line_offset,
bases_to_read,
&mut self.buf,
)?;
}
self.buf_idx = 0;
Ok(())
}
}
impl<'a, R: io::Read + io::Seek + 'a> Iterator for IndexedReaderIterator<'a, R> {
type Item = io::Result<u8>;
fn next(&mut self) -> Option<Self::Item> {
if self.buf_idx < self.buf.len() {
let item = Some(Ok(self.buf[self.buf_idx]));
self.buf_idx += 1;
item
} else if self.bases_left > 0 {
if let Err(e) = self.fill_buffer() {
self.bases_left = 0;
self.buf_idx = self.buf.len();
return Some(Err(e));
}
self.buf_idx = 1;
Some(Ok(self.buf[0]))
} else {
None
}
}
fn size_hint(&self) -> (usize, Option<usize>) {
let hint = self.bases_left as usize + (self.buf.len() - self.buf_idx);
(hint, Some(hint))
}
}
/// A Fasta writer.
#[derive(Debug)]
pub struct Writer<W: io::Write> {
writer: io::BufWriter<W>,
}
impl Writer<fs::File> {
/// Write to the given file path.
#[allow(clippy::wrong_self_convention)]
pub fn to_file<P: AsRef<Path>>(path: P) -> io::Result<Self> {
fs::File::create(path).map(Writer::new)
}
/// Write to the given file path and a buffer capacity
pub fn to_file_with_capacity<P: AsRef<Path>>(capacity: usize, path: P) -> io::Result<Self> {
fs::File::create(path).map(|file| Writer::with_capacity(capacity, file))
}
}
impl<W: io::Write> Writer<W> {
/// Create a new Fasta writer.
pub fn new(writer: W) -> Self {
Writer {
writer: io::BufWriter::new(writer),
}
}
/// Create a new Fasta writer with a capacity of write buffer
pub fn with_capacity(capacity: usize, writer: W) -> Self {
Writer {
writer: io::BufWriter::with_capacity(capacity, writer),
}
}
/// Create a new Fasta writer with a given BufWriter
pub fn from_bufwriter(bufwriter: io::BufWriter<W>) -> Self {
Writer { writer: bufwriter }
}
/// Directly write a [`fasta::Record`](struct.Record.html).
///
/// # Errors
/// If there is an issue writing to the `Writer`.
///
/// # Examples
/// ```rust
/// use bio::io::fasta::{Record, Writer};
/// use std::fs;
/// use std::io;
/// use std::path::Path;
///
/// let path = Path::new("test.fa");
/// let file = fs::File::create(path).unwrap();
/// {
/// let handle = io::BufWriter::new(file);
/// let mut writer = Writer::new(handle);
/// let record = Record::with_attrs("id", Some("desc"), b"ACGT");
///
/// let write_result = writer.write_record(&record);
/// assert!(write_result.is_ok());
/// }
///
/// let actual = fs::read_to_string(path).unwrap();
/// let expected = ">id desc\nACGT\n";
///
/// assert!(fs::remove_file(path).is_ok());
/// assert_eq!(actual, expected)
/// ```
pub fn write_record(&mut self, record: &Record) -> io::Result<()> {
self.write(record.id(), record.desc(), record.seq())
}
/// Write a Fasta record with given id, optional description and sequence.
pub fn write(&mut self, id: &str, desc: Option<&str>, seq: TextSlice<'_>) -> io::Result<()> {
self.writer.write_all(b">")?;
self.writer.write_all(id.as_bytes())?;
if let Some(desc) = desc {
self.writer.write_all(b" ")?;
self.writer.write_all(desc.as_bytes())?;
}
self.writer.write_all(b"\n")?;
self.writer.write_all(seq)?;
self.writer.write_all(b"\n")?;
Ok(())
}
/// Flush the writer, ensuring that everything is written.
pub fn flush(&mut self) -> io::Result<()> {
self.writer.flush()
}
}
/// A FASTA record.
#[derive(Default, Clone, Eq, PartialEq, Ord, PartialOrd, Hash, Debug, Serialize, Deserialize)]
pub struct Record {
id: String,
desc: Option<String>,
seq: String,
}
impl Record {
/// Create a new instance.
pub fn new() -> Self {
Record {
id: String::new(),
desc: None,
seq: String::new(),
}
}
/// Create a new `Record` from given attributes.
///
/// # Examples
/// ```rust
/// use bio::io::fasta::Record;
///
/// let read_id = "read1";
/// let description = Some("sampleid=foobar");
/// let sequence = b"ACGT";
/// let record = Record::with_attrs(read_id, description, sequence);
///
/// assert_eq!(">read1 sampleid=foobar\nACGT\n", record.to_string())
/// ```
pub fn with_attrs(id: &str, desc: Option<&str>, seq: TextSlice<'_>) -> Self {
let desc = desc.map(|desc| desc.to_owned());
Record {
id: id.to_owned(),
desc,
seq: String::from_utf8(seq.to_vec()).unwrap(),
}
}
/// Check if record is empty.
pub fn is_empty(&self) -> bool {
self.id.is_empty() && self.desc.is_none() && self.seq.is_empty()
}
/// Check validity of Fasta record.
pub fn check(&self) -> Result<(), &str> {
if self.id().is_empty() {
return Err("Expecting id for Fasta record.");
}
if !self.seq.is_ascii() {
return Err("Non-ascii character found in sequence.");
}
Ok(())
}
/// Return the id of the record.
pub fn id(&self) -> &str {
self.id.as_ref()
}
/// Return descriptions if present.
pub fn desc(&self) -> Option<&str> {
match self.desc.as_ref() {
Some(desc) => Some(desc),
None => None,
}
}
/// Return the sequence of the record.
pub fn seq(&self) -> TextSlice<'_> {
self.seq.as_bytes()
}
/// Clear the record.
fn clear(&mut self) {
self.id.clear();
self.desc = None;
self.seq.clear();
}
}
impl fmt::Display for Record {
/// Allows for using `Record` in a given formatter `f`. In general this is for
/// creating a `String` representation of a `Record` and, optionally, writing it to
/// a file.
///
/// # Errors
/// Returns [`std::fmt::Error`](https://doc.rust-lang.org/std/fmt/struct.Error.html)
/// if there is an issue formatting to the stream.
///
/// # Examples
///
/// Read in a Fasta `Record` and create a `String` representation of it.
///
/// ```rust
/// use bio::io::fasta::Reader;
/// use std::fmt::Write;
/// // create a "fake" fasta file
/// let fasta: &'static [u8] = b">id comment1 comment2\nACGT\n";
/// let mut records = Reader::new(fasta).records().map(|r| r.unwrap());
/// let record = records.next().unwrap();
///
/// let mut actual = String::new();
/// // populate `actual` with a string representation of our record
/// write!(actual, "{}", record).unwrap();
///
/// let expected = std::str::from_utf8(fasta).unwrap();
///
/// assert_eq!(actual, expected)
/// ```
fn fmt(&self, f: &mut fmt::Formatter<'_>) -> Result<(), fmt::Error> {
let header = match self.desc() {
Some(d) => format!("{} {}", self.id().to_owned(), d),
None => self.id().to_owned(),
};
write!(
f,
">{}\n{}\n",
header,
std::str::from_utf8(self.seq()).unwrap(),
)
}
}
/// An iterator over the records of a Fasta file.
#[derive(Default, Clone, Eq, PartialEq, Ord, PartialOrd, Hash, Debug, Serialize, Deserialize)]
pub struct Records<B>
where
B: io::BufRead,
{
reader: Reader<B>,
error_has_occured: bool,
}
impl<B> Iterator for Records<B>
where
B: io::BufRead,
{
type Item = io::Result<Record>;
fn next(&mut self) -> Option<io::Result<Record>> {
if self.error_has_occured {
None
} else {
let mut record = Record::new();
match self.reader.read(&mut record) {
Ok(()) if record.is_empty() => None,
Ok(()) => Some(Ok(record)),
Err(err) => {
self.error_has_occured = true;
Some(Err(err))
}
}
}
}
}
#[cfg(test)]
mod tests {
use super::*;
use std::fmt::Write as FmtWrite;
use std::io;
const FASTA_FILE: &[u8] = b">id desc
ACCGTAGGCTGA
CCGTAGGCTGAA
CGTAGGCTGAAA
GTAGGCTGAAAA
CCCC
>id2
ATTGTTGTTTTA
ATTGTTGTTTTA
ATTGTTGTTTTA
GGGG
";
const FAI_FILE: &[u8] = b"id\t52\t9\t12\t13
id2\t40\t71\t12\t13
";
const TRUNCATED_FASTA: &[u8] = b">id desc\nACCGTAGGCTGA";
const FASTA_FILE_CRLF: &[u8] = b">id desc\r
ACCGTAGGCTGA\r
CCGTAGGCTGAA\r
CGTAGGCTGAAA\r
GTAGGCTGAAAA\r
CCCC\r
>id2\r
ATTGTTGTTTTA\r
ATTGTTGTTTTA\r
ATTGTTGTTTTA\r
GGGG\r
";
const FAI_FILE_CRLF: &[u8] = b"id\t52\t10\t12\t14\r
id2\t40\t78\t12\t14\r
";
const FASTA_FILE_NO_TRAILING_LF: &[u8] = b">id desc
GTAGGCTGAAAA
CCCC";
const FAI_FILE_NO_TRAILING_LF: &[u8] = b"id\t16\t9\t12\t13";
const WRITE_FASTA_FILE: &[u8] = b">id desc
ACCGTAGGCTGA
>id2
ATTGTTGTTTTA
";
struct ReaderMock {
seek_fails: bool,
read_fails: bool,
}
impl Read for ReaderMock {
fn read(&mut self, buf: &mut [u8]) -> io::Result<usize> {
if self.read_fails {
Err(io::Error::new(io::ErrorKind::Other, "Read set to fail"))
} else {
Ok(buf.len())
}
}
}
impl Seek for ReaderMock {
fn seek(&mut self, pos: io::SeekFrom) -> io::Result<u64> {
if let io::SeekFrom::Start(pos) = pos {
if self.seek_fails {
Err(io::Error::new(io::ErrorKind::Other, "Seek set to fail"))
} else {
Ok(pos)
}
} else {
unimplemented!();
}
}
}
#[test]
fn test_reader() {
let reader = Reader::new(FASTA_FILE);
let ids = ["id", "id2"];
let descs = [Some("desc"), None];
let seqs: [&[u8]; 2] = [
b"ACCGTAGGCTGACCGTAGGCTGAACGTAGGCTGAAAGTAGGCTGAAAACCCC",
b"ATTGTTGTTTTAATTGTTGTTTTAATTGTTGTTTTAGGGG",
];
for (i, r) in reader.records().enumerate() {
let record = r.expect("Error reading record");
assert_eq!(record.check(), Ok(()));
assert_eq!(record.id(), ids[i]);
assert_eq!(record.desc(), descs[i]);
assert_eq!(record.seq(), seqs[i]);
}
let reader = Reader::with_capacity(100, FASTA_FILE);
for (i, r) in reader.records().enumerate() {
let record = r.expect("Error reading record");
assert_eq!(record.check(), Ok(()));
assert_eq!(record.id(), ids[i]);
assert_eq!(record.desc(), descs[i]);
assert_eq!(record.seq(), seqs[i]);
}
let reader = Reader::from_bufread(io::BufReader::new(FASTA_FILE));
for (i, r) in reader.records().enumerate() {
let record = r.expect("Error reading record");
assert_eq!(record.check(), Ok(()));
assert_eq!(record.id(), ids[i]);
assert_eq!(record.desc(), descs[i]);
assert_eq!(record.seq(), seqs[i]);
}
}
#[test]
fn test_faread_trait() {
let path = "genome.fa.gz";
let mut fa_reader: Box<dyn FastaRead> = match path.ends_with(".gz") {
true => Box::new(Reader::new(io::BufReader::new(FASTA_FILE))),
false => Box::new(Reader::new(FASTA_FILE)),
};
// The read method can be called, since it is implemented by
// FQRead. Right now, the records method would not work.
let mut record = Record::new();
fa_reader.read(&mut record).unwrap();
// Check if the returned result is correct.
assert_eq!(record.check(), Ok(()));
assert_eq!(record.id(), "id");
assert_eq!(record.desc(), Some("desc"));
assert_eq!(
record.seq().to_vec(),
b"ACCGTAGGCTGACCGTAGGCTGAACGTAGGCTGAAAGTAGGCTGAAAACCCC".to_vec()
);
}
#[test]
fn test_reader_wrong_header() {
let mut reader = Reader::new(&b"!test\nACGTA\n"[..]);
let mut record = Record::new();
assert!(
reader.read(&mut record).is_err(),
"read() should return Err if FASTA header is malformed"
);
}
#[test]
fn test_reader_no_id() {
let mut reader = Reader::new(&b">\nACGTA\n"[..]);
let mut record = Record::new();
reader.read(&mut record).unwrap();
assert!(
record.check().is_err(),
"check() should return Err if FASTA header is empty"
);
}
#[test]
fn test_reader_non_ascii_sequence() {
let mut reader = Reader::new(&b">id\nACGTA\xE2\x98\xB9AT\n"[..]);
let mut record = Record::new();
reader.read(&mut record).unwrap();
assert!(
record.check().is_err(),
"check() should return Err if FASTA sequence is not ASCII"
);
}
#[test]
fn test_reader_read_fails() {
let mut reader = Reader::new(ReaderMock {
seek_fails: false,
read_fails: true,
});
let mut record = Record::new();
assert!(
reader.read(&mut record).is_err(),
"read() should return Err if Read::read fails"
);
}
#[test]
fn test_reader_read_fails_iter() {
let reader = Reader::new(ReaderMock {
seek_fails: false,
read_fails: true,
});
let mut records = reader.records();
assert!(
records.next().unwrap().is_err(),
"next() should return Err if Read::read fails"
);
assert!(
records.next().is_none(),
"next() should return None after error has occurred"
);
}
#[test]
fn test_reader_from_file_path_doesnt_exist_returns_err() {
let path = Path::new("/I/dont/exist.fasta");
let error = Reader::from_file(path)
.unwrap_err()
.downcast::<String>()
.unwrap();
assert_eq!(&error, "Failed to read fasta from \"/I/dont/exist.fasta\"")
}
#[test]
fn test_record_with_attrs_without_description() {
let record = Record::with_attrs("id_str", None, b"ATGCGGG");
assert_eq!(record.id(), "id_str");
assert_eq!(record.desc(), None);
assert_eq!(record.seq(), b"ATGCGGG");
}
#[test]
fn test_record_with_attrs_with_description() {
let record = Record::with_attrs("id_str", Some("desc"), b"ATGCGGG");
assert_eq!(record.id(), "id_str");
assert_eq!(record.desc(), Some("desc"));
assert_eq!(record.seq(), b"ATGCGGG");
}
#[test]
fn test_index_sequences() {
let reader = IndexedReader::new(io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
let sequences = reader.index.sequences();
assert_eq!(sequences.len(), 2);
assert_eq!(
sequences[0],
Sequence {
name: "id".into(),
len: 52,
}
);
assert_eq!(
sequences[1],
Sequence {
name: "id2".into(),
len: 40,
}
);
}
#[test]
fn test_indexed_reader() {
_test_indexed_reader(&FASTA_FILE, &FAI_FILE, _read_buffer);
_test_indexed_reader_truncated(_read_buffer);
_test_indexed_reader_extreme_whitespace(_read_buffer);
}
#[test]
fn test_indexed_reader_crlf() {
_test_indexed_reader(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_buffer);
}
#[test]
fn test_indexed_reader_iter() {
_test_indexed_reader(&FASTA_FILE, &FAI_FILE, _read_iter);
_test_indexed_reader_truncated(_read_iter);
_test_indexed_reader_extreme_whitespace(_read_iter);
}
#[test]
fn test_indexed_reader_iter_crlf() {
_test_indexed_reader(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_iter);
}
fn _test_indexed_reader<'a, F>(fasta: &'a [u8], fai: &'a [u8], read: F)
where
F: Fn(&mut IndexedReader<io::Cursor<&'a [u8]>>, &str, u64, u64) -> io::Result<Vec<u8>>,
{
let mut reader = IndexedReader::new(io::Cursor::new(fasta), fai).unwrap();
// Test reading various substrings of the sequence
assert_eq!(read(&mut reader, "id", 1, 5).unwrap(), b"CCGT");
assert_eq!(
read(&mut reader, "id", 1, 31).unwrap(),
b"CCGTAGGCTGACCGTAGGCTGAACGTAGGC"
);
assert_eq!(read(&mut reader, "id", 13, 23).unwrap(), b"CGTAGGCTGA");
assert_eq!(
read(&mut reader, "id", 36, 52).unwrap(),
b"GTAGGCTGAAAACCCC"
);
assert_eq!(
read(&mut reader, "id2", 12, 40).unwrap(),
b"ATTGTTGTTTTAATTGTTGTTTTAGGGG"
);
assert_eq!(read(&mut reader, "id2", 12, 12).unwrap(), b"");
assert_eq!(read(&mut reader, "id2", 12, 13).unwrap(), b"A");
// Minimal sequence spanning new-line
assert_eq!(read(&mut reader, "id", 11, 13).unwrap(), b"AC");
assert!(read(&mut reader, "id2", 12, 11).is_err());
assert!(read(&mut reader, "id2", 12, 1000).is_err());
assert!(read(&mut reader, "id3", 0, 1).is_err());
}
fn _test_indexed_reader_truncated<'a, F>(read: F)
where
F: Fn(&mut IndexedReader<io::Cursor<&'a [u8]>>, &str, u64, u64) -> io::Result<Vec<u8>>,
{
let mut reader = IndexedReader::new(io::Cursor::new(TRUNCATED_FASTA), FAI_FILE).unwrap();
assert_eq!(read(&mut reader, "id", 0, 12).unwrap(), b"ACCGTAGGCTGA");
assert!(read(&mut reader, "id", 0, 13).is_err()); // read past EOF
assert!(read(&mut reader, "id", 36, 52).is_err()); // seek and read past EOF
assert!(read(&mut reader, "id2", 12, 40).is_err()); // seek and read past EOF
}
fn _test_indexed_reader_extreme_whitespace<F>(read: F)
where
F: Fn(&mut IndexedReader<io::Cursor<Vec<u8>>>, &str, u64, u64) -> io::Result<Vec<u8>>,
{
// Test to exercise the case where we cannot consume all whitespace at once. More than
// DEFAULT_BUF_SIZE (a non-public constant set to 8 * 1024) whitespace is used to ensure
// that it can't all fit in the BufReader at once.
let mut seq = Vec::new();
seq.push(b'A');
seq.resize(10000, b' ');
seq.push(b'B');
let fasta = io::Cursor::new(seq);
let fai = io::Cursor::new(Vec::from(&b"id\t2\t0\t1\t10000"[..]));
let mut reader = IndexedReader::new(fasta, fai).unwrap();
assert_eq!(read(&mut reader, "id", 0, 2).unwrap(), b"AB");
}
fn _read_buffer<T>(
reader: &mut IndexedReader<T>,
seqname: &str,
start: u64,
stop: u64,
) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch(seqname, start, stop)?;
reader.read(&mut seq)?;
Ok(seq)
}
fn _read_iter<T>(
reader: &mut IndexedReader<T>,
seqname: &str,
start: u64,
stop: u64,
) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch(seqname, start, stop)?;
for nuc in reader.read_iter()? {
seq.push(nuc?);
}
Ok(seq)
}
#[test]
fn test_indexed_reader_all() {
_test_indexed_reader_all(&FASTA_FILE, &FAI_FILE, _read_buffer_all);
}
#[test]
fn test_indexed_reader_crlf_all() {
_test_indexed_reader_all(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_buffer_all);
}
#[test]
fn test_indexed_reader_iter_all() {
_test_indexed_reader_all(&FASTA_FILE, &FAI_FILE, _read_iter_all);
}
#[test]
fn test_indexed_reader_iter_crlf_all() {
_test_indexed_reader_all(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_iter_all);
}
fn _test_indexed_reader_all<'a, F>(fasta: &'a [u8], fai: &'a [u8], read: F)
where
F: Fn(&mut IndexedReader<io::Cursor<&'a [u8]>>, &str) -> io::Result<Vec<u8>>,
{
let mut reader = IndexedReader::new(io::Cursor::new(fasta), fai).unwrap();
assert_eq!(
read(&mut reader, "id").unwrap(),
&b"ACCGTAGGCTGACCGTAGGCTGAACGTAGGCTGAAAGTAGGCTGAAAACCCC"[..]
);
assert_eq!(
read(&mut reader, "id2").unwrap(),
&b"ATTGTTGTTTTAATTGTTGTTTTAATTGTTGTTTTAGGGG"[..]
);
}
fn _read_buffer_all<T>(reader: &mut IndexedReader<T>, seqname: &str) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch_all(seqname)?;
reader.read(&mut seq)?;
Ok(seq)
}
fn _read_iter_all<T>(reader: &mut IndexedReader<T>, seqname: &str) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch_all(seqname)?;
for nuc in reader.read_iter()? {
seq.push(nuc?);
}
Ok(seq)
}
#[test]
fn test_indexed_reader_by_rid_all() {
_test_indexed_reader_by_rid_all(&FASTA_FILE, &FAI_FILE, _read_buffer_by_rid_all);
}
#[test]
fn test_indexed_reader_crlf_by_rid_all() {
_test_indexed_reader_by_rid_all(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_buffer_by_rid_all);
}
#[test]
fn test_indexed_reader_iter_by_rid_all() {
_test_indexed_reader_by_rid_all(&FASTA_FILE, &FAI_FILE, _read_iter_by_rid_all);
}
#[test]
fn test_indexed_reader_iter_crlf_by_rid_all() {
_test_indexed_reader_by_rid_all(&FASTA_FILE_CRLF, &FAI_FILE_CRLF, _read_iter_by_rid_all);
}
fn _test_indexed_reader_by_rid_all<'a, F>(fasta: &'a [u8], fai: &'a [u8], read: F)
where
F: Fn(&mut IndexedReader<io::Cursor<&'a [u8]>>, usize) -> io::Result<Vec<u8>>,
{
let mut reader = IndexedReader::new(io::Cursor::new(fasta), fai).unwrap();
assert_eq!(
read(&mut reader, 0).unwrap(),
&b"ACCGTAGGCTGACCGTAGGCTGAACGTAGGCTGAAAGTAGGCTGAAAACCCC"[..]
);
assert_eq!(
read(&mut reader, 1).unwrap(),
&b"ATTGTTGTTTTAATTGTTGTTTTAATTGTTGTTTTAGGGG"[..]
);
}
fn _read_buffer_by_rid_all<T>(
reader: &mut IndexedReader<T>,
seq_index: usize,
) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch_all_by_rid(seq_index)?;
reader.read(&mut seq)?;
Ok(seq)
}
fn _read_iter_by_rid_all<T>(
reader: &mut IndexedReader<T>,
seq_index: usize,
) -> io::Result<Vec<u8>>
where
T: Seek + Read,
{
let mut seq = vec![];
reader.fetch_all_by_rid(seq_index)?;
for nuc in reader.read_iter()? {
seq.push(nuc?);
}
Ok(seq)
}
#[test]
fn test_indexed_reader_iter_size_hint() {
let mut reader = IndexedReader::new(io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
reader.fetch("id", 2, 4).unwrap();
let mut iterator = reader.read_iter().unwrap();
assert_eq!(iterator.size_hint(), (2, Some(2)));
assert_eq!(iterator.next().unwrap().unwrap(), b'C');
assert_eq!(iterator.size_hint(), (1, Some(1)));
assert_eq!(iterator.next().unwrap().unwrap(), b'G');
assert_eq!(iterator.size_hint(), (0, Some(0)));
assert!(iterator.next().is_none());
assert_eq!(iterator.size_hint(), (0, Some(0)));
}
#[test]
fn test_indexed_reader_reused_buffer() {
let mut reader = IndexedReader::new(io::Cursor::new(FASTA_FILE), FAI_FILE).unwrap();
let mut seq = Vec::new();
reader.fetch("id", 1, 5).unwrap();
reader.read(&mut seq).unwrap();
assert_eq!(seq, b"CCGT");
reader.fetch("id", 13, 23).unwrap();
reader.read(&mut seq).unwrap();
assert_eq!(seq, b"CGTAGGCTGA");
}
#[test]
fn test_indexed_reader_no_trailing_lf() {
let mut reader = IndexedReader::new(
io::Cursor::new(FASTA_FILE_NO_TRAILING_LF),
FAI_FILE_NO_TRAILING_LF,
)
.unwrap();
let mut seq = Vec::new();
reader.fetch("id", 0, 16).unwrap();
reader.read(&mut seq).unwrap();
assert_eq!(seq, b"GTAGGCTGAAAACCCC");
}
#[test]
fn test_indexed_reader_bad_reader() {
let bad_reader = ReaderMock {
seek_fails: false,
read_fails: false,
};
let mut reader = IndexedReader::new(bad_reader, FAI_FILE).unwrap();
let mut seq = Vec::new();
reader.fetch("id", 0, 10).unwrap();
assert!(reader.read(&mut seq).is_ok())
}
#[test]
fn test_indexed_reader_read_seek_fails() {
let bad_reader = ReaderMock {
seek_fails: true,
read_fails: false,
};
let mut reader = IndexedReader::new(bad_reader, FAI_FILE).unwrap();
let mut seq = Vec::new();
reader.fetch("id", 0, 10).unwrap();
assert!(reader.read(&mut seq).is_err());
}
#[test]
fn test_indexed_reader_read_read_fails() {
let bad_reader = ReaderMock {
seek_fails: false,
read_fails: true,
};
let mut reader = IndexedReader::new(bad_reader, FAI_FILE).unwrap();
let mut seq = Vec::new();
reader.fetch("id", 0, 10).unwrap();
assert!(reader.read(&mut seq).is_err());
}
#[test]
fn test_indexed_reader_iter_seek_fails() {
let bad_reader = ReaderMock {
seek_fails: true,
read_fails: false,
};
let mut reader = IndexedReader::new(bad_reader, FAI_FILE).unwrap();
reader.fetch("id", 0, 10).unwrap();
assert!(reader.read_iter().is_err());
}
#[test]
fn test_indexed_reader_iter_read_fails() {
let bad_reader = ReaderMock {
seek_fails: false,
read_fails: true,
};
let mut reader = IndexedReader::new(bad_reader, FAI_FILE).unwrap();
reader.fetch("id", 0, 10).unwrap();
let mut iterator = reader.read_iter().unwrap();
assert!(iterator.next().unwrap().is_err());
assert!(
iterator.next().is_none(),
"next() should return none after error has occurred"
);
}
#[test]
fn test_indexed_reader_no_fetch_read_fails() {
let reader = ReaderMock {
seek_fails: false,
read_fails: false,
};
let mut reader = IndexedReader::new(reader, FAI_FILE).unwrap();
let mut seq = vec![];
assert!(reader.read(&mut seq).is_err());
}
#[test]
fn test_indexed_reader_no_fetch_read_iter_fails() {
let reader = ReaderMock {
seek_fails: false,
read_fails: false,
};
let mut reader = IndexedReader::new(reader, FAI_FILE).unwrap();
assert!(reader.read_iter().is_err());
}
#[test]
fn test_writer() {
let mut writer = Writer::new(Vec::new());
writer.write("id", Some("desc"), b"ACCGTAGGCTGA").unwrap();
writer.write("id2", None, b"ATTGTTGTTTTA").unwrap();
writer.flush().unwrap();
assert_eq!(writer.writer.get_ref(), &WRITE_FASTA_FILE);
let mut writer = Writer::with_capacity(100, Vec::new());
writer.write("id", Some("desc"), b"ACCGTAGGCTGA").unwrap();
writer.write("id2", None, b"ATTGTTGTTTTA").unwrap();
writer.flush().unwrap();
assert_eq!(writer.writer.get_ref(), &WRITE_FASTA_FILE);
let mut writer = Writer::from_bufwriter(std::io::BufWriter::with_capacity(100, Vec::new()));
writer.write("id", Some("desc"), b"ACCGTAGGCTGA").unwrap();
writer.write("id2", None, b"ATTGTTGTTTTA").unwrap();
writer.flush().unwrap();
assert_eq!(writer.writer.get_ref(), &WRITE_FASTA_FILE);
}
#[test]
fn test_display_record_no_desc_id_without_space_after() {
let fasta: &'static [u8] = b">id\nACGT\n";
let mut records = Reader::new(fasta).records().map(|r| r.unwrap());
let record = records.next().unwrap();
let mut actual = String::new();
write!(actual, "{}", record).unwrap();
let expected = std::str::from_utf8(fasta).unwrap();
assert_eq!(actual, expected)
}
#[test]
fn test_display_record_with_desc_id_has_space_between_id_and_desc() {
let fasta: &'static [u8] = b">id comment1 comment2\nACGT\n";
let mut records = Reader::new(fasta).records().map(|r| r.unwrap());
let record = records.next().unwrap();
let mut actual = String::new();
write!(actual, "{}", record).unwrap();
let expected = std::str::from_utf8(fasta).unwrap();
assert_eq!(actual, expected)
}
#[test]
fn test_index_record_idx_by_rid_invalid_index_returns_error() {
let reader = ReaderMock {
seek_fails: false,
read_fails: false,
};
let index_reader = IndexedReader::new(reader, FAI_FILE).unwrap();
let actual = index_reader.idx_by_rid(99999).unwrap_err();
let expected = io::Error::new(io::ErrorKind::Other, "Invalid record index in fasta file.");
assert_eq!(actual.kind(), expected.kind());
assert_eq!(actual.to_string(), expected.to_string())
}
#[test]
fn test_index_record_fetch_by_rid_second_index_returns_second_record() {
let reader = ReaderMock {
seek_fails: false,
read_fails: false,
};
let mut index_reader = IndexedReader::new(reader, FAI_FILE).unwrap();
let actual = index_reader.fetch_by_rid(1, 1, 3);
assert!(actual.is_ok());
assert_eq!(
index_reader.fetched_idx,
Some(IndexRecord {
name: String::from("id2"),
len: 40,
offset: 71,
line_bases: 12,
line_bytes: 13
})
)
}
#[test]
fn test_writer_to_file_dir_doesnt_exist_returns_err() {
let path = Path::new("/I/dont/exist.fa");
let actual = Writer::to_file(path).unwrap_err();
let expected = io::Error::new(io::ErrorKind::NotFound, "foo");
assert_eq!(actual.kind(), expected.kind());
}
#[test]
fn test_writer_to_file_dir_exists_returns_ok() {
let file = tempfile::NamedTempFile::new().expect("Could not create temp file");
let path = file.path();
assert!(Writer::to_file(path).is_ok());
assert!(Writer::to_file_with_capacity(100, path).is_ok());
}
#[test]
fn test_write_record() {
let path = Path::new("test.fa");
let file = fs::File::create(path).unwrap();
{
let handle = io::BufWriter::new(file);
let mut writer = Writer { writer: handle };
let record = Record::with_attrs("id", Some("desc"), b"ACGT");
let write_result = writer.write_record(&record);
assert!(write_result.is_ok());
}
let actual = fs::read_to_string(path).unwrap();
let expected = ">id desc\nACGT\n";
assert!(fs::remove_file(path).is_ok());
assert_eq!(actual, expected)
}
}