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
// Copyright 2014-2016 Johannes Köster.
// Licensed under the MIT license (http://opensource.org/licenses/MIT)
// This file may not be copied, modified, or distributed
// except according to those terms.
//! # Rust-bio, a bioinformatics library for Rust.
//! This library provides implementations of many algorithms and data structures
//! that are useful for bioinformatics.
//! All provided implementations are rigorously tested via continuous
//! integration.
//! For installation instructions and a general overview, visit
//! https://rust-bio.github.io.
//!
//! Currently, rust-bio provides
//!
//! * most major pattern matching algorithms,
//! * a convenient alphabet implementation,
//! * pairwise alignment,
//! * suffix arrays,
//! * BWT and FM-Index,
//! * FMD-Index for finding supermaximal exact matches,
//! * a q-gram index,
//! * an orf research algorithm,
//! * a rank/select data structure,
//! * FASTQ and FASTA and BED readers and writers,
//! * helper functions for combinatorics and dealing with log probabilities,
//! * an implementation of Hidden Markov Model and related algorithms.
//!
//! # Example
//!
//! ```rust
//! use bio::alphabets;
//! use bio::data_structures::suffix_array::suffix_array;
//! use bio::data_structures::bwt::{bwt, less, Occ};
//! use bio::data_structures::fmindex::{FMIndex, FMIndexable};
//!
//! let text = b"ACGGATGCTGGATCGGATCGCGCTAGCTA$";
//! let pattern = b"ACCG";
//!
//! // Create an FM-Index for a given text.
//! let alphabet = alphabets::dna::iupac_alphabet();
//! let sa = suffix_array(text);
//! let bwt = bwt(text, &sa);
//! let less = less(&bwt, &alphabet);
//! let occ = Occ::new(&bwt, 3, &alphabet);
//! let fmindex = FMIndex::new(&bwt, &less, &occ);
//!
//! let interval = fmindex.backward_search(pattern.iter());
//! let positions = interval.occ(&sa);
//! ```
//!
//! # Multithreaded Example
//!
//! ```rust
//! use bio::alphabets;
//! use bio::data_structures::suffix_array::suffix_array;
//! use bio::data_structures::bwt::{bwt, less, Occ};
//! use bio::data_structures::fmindex::{FMIndex, FMIndexable};
//! use std::sync::Arc;
//! use std::thread;
//!
//! let text = b"ACGGATGCTGGATCGGATCGCGCTAGCTA$";
//! let patterns = vec![b"ACCG", b"TGCT"];
//!
//! // Create an FM-Index for a given text.
//! let alphabet = alphabets::dna::iupac_alphabet();
//! let sa = suffix_array(text);
//! let bwt = Arc::new(bwt(text, &sa));
//! let less = Arc::new(less(bwt.as_ref(), &alphabet));
//! let occ = Arc::new(Occ::new(bwt.as_ref(), 3, &alphabet));
//! let fmindex = Arc::new(FMIndex::new(bwt, less, occ));
//!
//! // Spawn threads to perform backward searches for each interval
//! let interval_calculators = patterns.into_iter().map(|pattern| {
//! let fmindex = fmindex.clone();
//! thread::spawn(move ||
//! fmindex.backward_search(pattern.iter())
//! )
//! }).collect::<Vec<_>>();
//!
//! // Loop through the results, extracting the positions array for each pattern
//! for interval_calculator in interval_calculators {
//! let positions = interval_calculator.join().unwrap().occ(&sa);
//! }
//! ```
//!
//! Documentation and further examples for each module can be found in the module descriptions below.
extern crate bio_types;
#[macro_use]
extern crate approx;
extern crate bit_set;
extern crate bytecount;
extern crate csv;
#[macro_use]
extern crate custom_derive;
extern crate itertools;
extern crate itertools_num;
#[macro_use]
extern crate lazy_static;
extern crate multimap;
extern crate ndarray;
#[macro_use]
extern crate newtype_derive;
extern crate num_integer;
extern crate num_traits;
extern crate ordered_float;
#[macro_use]
extern crate quick_error;
extern crate regex;
extern crate serde;
#[macro_use]
extern crate serde_derive;
extern crate bv;
extern crate fnv;
extern crate statrs;
extern crate vec_map;
pub mod alignment;
pub mod alphabets;
pub mod data_structures;
pub mod io;
pub mod pattern_matching;
pub mod scores;
pub mod seq_analysis;
pub mod stats;
pub mod utils;