Struct nucleic_acid::FMIndex
[−]
[src]
pub struct FMIndex { /* fields omitted */ }
Ferragina-Manzini index (or Full-text index in Minute space) for finding occurrences of substrings in O(1) time.
use nucleic_acid::FMIndex; let text = String::from("GCGTGCCCAGGGCACTGCCGCTGCAGGCGTAGGCATCGCATCACACGCGT"); let index = FMIndex::new(text.as_bytes()); // count the occurrences assert_eq!(0, index.count("CCCCC")); assert_eq!(3, index.count("TG")); // ... or get their positions assert_eq!(index.search("GCGT"), vec![46, 26, 0]);
The current implementation of FM-index is a memory killer, since it stores positions of all bytes in the given data. For the human genome (~3 GB), it consumed ~27 GB of RAM to build the index (in ~4 mins).
That said, it still returns the match results in a few microseconds.
Methods
impl FMIndex
[src]
fn new(data: &[u8]) -> FMIndex
Generate an FM-index for the input data.
fn bwt(&self) -> &[u8]
Get the reference to the inner BWT data.
Note that the length of BWT is one more than the length of the actual text, since it has a null byte to indicate empty string.
fn new_from_bwt(bwt_data: Vec<u8>) -> FMIndex
Generate the FM-index from the BWT data.
It's not a good idea to generate FM-index from scratch all the time, especially for large inputs. This would be very useful when your data is large and remains constant for a while.
FM-index internally uses BWT, and BWT is generated from the suffix array, which takes a lot of time.
If your input doesn't change, then it's better to get the BWT data (using bwt
method), write it
to a file and generate the index from that in the future.
fn nearest(&self, idx: usize, ch: u8) -> usize
Get the nearest position of a character in the internal BWT data.
The count
and search
methods rely on this method for finding occurrences.
For example, we can do soemthing like this,
use nucleic_acid::FMIndex; let fm = FMIndex::new(b"Hello, Hello, Hello" as &[u8]); // initially, the range should be the length of the BWT let mut top = 0; let mut bottom = fm.bwt().len(); let query = b"llo"; // feed the characters in the reverse for ch in query.iter().rev() { top = fm.nearest(top, *ch); bottom = fm.nearest(bottom, *ch); if top >= bottom { return } } // If we get a valid range, then everything in that range is a valid match. // This way, we can get both the count and positions... assert_eq!(3, bottom - top); assert_eq!(vec![17, 10, 3], (top..bottom).map(|i| fm[i]).collect::<Vec<_>>());
This is backward searching. As you feed in the characters along with a position, nearest
will
give you a new position in the index. Once the range becomes invalid (which happens when the
substring doesn't exist), we can bail out. On the contrary, if the range remains valid after
you've fed in all the characters, then every value within in that range is an occurrence.
So, this is useful when you want to cache the repeating ranges. With this, you can build your own count/search functions with caching. It's also useful for making custom approximate matching functions by backtracking whenever there's an invalid range.
fn count(&self, query: &str) -> usize
Count the occurrences of the substring in the original data.
fn search(&self, query: &str) -> Vec<usize>
Get the positions of occurrences of substring in the original data.
Trait Implementations
impl Clone for FMIndex
[src]
fn clone(&self) -> FMIndex
Returns a copy of the value. Read more
fn clone_from(&mut self, source: &Self)
1.0.0
Performs copy-assignment from source
. Read more