Skip to main content

Ebwt

Struct Ebwt 

Source
pub struct Ebwt { /* private fields */ }
Expand description

Extended Burrows-Wheeler Transform (FM-index)

This is the main data structure used by RustAlign for ultrafast sequence alignment. It is memory-mapped from .rai index files.

Implementations§

Source§

impl Ebwt

Source

pub fn load<P: AsRef<Path>>(path_base: P) -> AlignResult<Self>

Load an existing index from disk

§Arguments
  • path_base - Path to the index (without .1.rai, .2.rai, etc. extensions)
§Errors

Returns an error if the index files don’t exist or are corrupted.

Source

pub fn load_large<P: AsRef<Path>>(path_base: P) -> AlignResult<Self>

Load a “large” index variant

Source

pub fn params(&self) -> &EbwtParams

Get the index parameters

Source

pub fn path(&self) -> &Path

Get the path to this index

Source

pub fn is_large(&self) -> bool

Check if this is a large index

Source

pub fn refnames(&self) -> &[String]

Get the chromosome names

Source

pub fn n_pat(&self) -> u32

Get the number of chromosomes

Source

pub fn n_frag(&self) -> u32

Get the number of fragments

Source

pub fn rstarts(&self) -> &[u32]

Get the rstarts array

Source

pub fn plen(&self) -> &[u32]

Get chromosome lengths

Source

pub fn z_off(&self) -> u32

Get the zOff value (position of $ in BWT)

Source

pub fn get_offset(&self, row: u64) -> AlignResult<u64>

Convert a BWT row to a genome offset using suffix array samples and LF-mapping.

This implements RustAlign’s getOffset() algorithm:

  1. If row is zOff, return 0
  2. If row is at a sampled position, return the sampled offset
  3. Otherwise, walk back using LF-mapping until hitting zOff or a sampled row
Source

pub fn joined_to_text_off( &self, qlen: u64, off: u64, fw: bool, ) -> AlignResult<(u32, u64, u64)>

Convert a genome offset to chromosome index and position.

This implements RustAlign’s joinedToTextOff() algorithm using binary search on the rstarts array.

Returns (chromosome_index, position_in_chromosome, chromosome_length)

Source

pub fn resolve_position( &self, row: u64, qlen: u64, fw: bool, ) -> AlignResult<(String, u64, u64)>

Resolve a BWT row to chromosome name and position.

This combines getOffset() and joinedToTextOff() to convert a BWT row directly to a chromosome name and position.

Returns (chromosome_name, position, chromosome_length)

Trait Implementations§

Source§

impl FmIndex for Ebwt

Source§

fn exact_match(&self, pattern: &[Nuc]) -> AlignResult<Vec<usize>>

Exact match search using backward search on FM-index

This implements the classic FM-index backward search algorithm:

  1. Initialize range using ftab for the last ftabChars nucleotides
  2. Extend left through remaining characters using LF-mapping
  3. Return all positions in the final range
Source§

fn len(&self) -> usize

Get the length of the indexed text
Source§

fn is_empty(&self) -> bool

Check if the index is empty
Source§

fn bwt(&self, pos: usize) -> AlignResult<NucPair>

Get forward BWT character at position
Source§

fn rank(&self, c: Nuc, _pos: usize) -> AlignResult<usize>

Rank operation - count occurrences of character up to position
Source§

fn select(&self, c: Nuc, _k: usize) -> AlignResult<usize>

Select operation - find position of k-th occurrence of character

Auto Trait Implementations§

§

impl Freeze for Ebwt

§

impl RefUnwindSafe for Ebwt

§

impl Send for Ebwt

§

impl Sync for Ebwt

§

impl Unpin for Ebwt

§

impl UnsafeUnpin for Ebwt

§

impl UnwindSafe for Ebwt

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.