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
// Copyright 2015 Vadim Nazarov.
// Licensed under the MIT license (http://opensource.org/licenses/MIT)
// This file may not be copied, modified, or distributed
// except according to those terms.

//! Various subroutines for computing a distance between sequences.
//! Complexity: O(n) for strings of length n for the Hamming distance;
//! O(n * m) for strings of length n and m for the Levenshtein (or edit) distance.
//!
//! # Example
//!
//! ```
//! use bio::alignment::distance::*;
//!
//! let x = b"GTCTGCATGCG";
//! let y = b"TTTAGCTAGCG";
//! // GTCTGCATGCG
//! //  |  ||  |||
//! // TTTAGCTAGCG
//! match hamming(x, y) {
//!     Ok(s)  => println!("Score is: {}", s),  // Score is 5
//!     Err(e) => println!("Error: {}", e),     // No error in this example
//! }
//!
//! let x = b"ACCGTGGAT";
//! let y = b"AAAAACCGTTGAT";
//! // ----ACCGTGGAT
//! //     ||||| |||
//! // AAAAACCGTTGAT
//! let l_score = levenshtein(x, y);  // Score is 5
//! assert_eq!(l_score, 5);
//! ```


use itertools::Zip;
use std::cmp::min;
use std::result::Result;

use utils::TextSlice;


/// Compute the Hamming distance between two strings with `hamming`. If returns the `Result<u32, &str>` type
/// with the first element corresponding to the distance between two strings (a number of mismatches) and the second one to the error message
/// when two strings are not of equal sizes.
///
/// # Example
///
/// ```
/// use bio::alignment::distance::*;
///
/// let x = b"GTCTGCATGCG";
/// let y = b"TTTAGCTAGCG";
/// // GTCTGCATGCG
/// //  |  ||  |||
/// // TTTAGCTAGCG
/// match hamming(x, y) {
///     Ok(s)  => println!("Score is: {}", s),  // Score is 5
///     Err(e) => println!("Error: {}", e),     // No error in this example
/// }
/// assert!(hamming(x, y).is_ok());
/// assert_eq!(hamming(x, y).unwrap(), 5);
/// ```
///
/// In case of the error:
///
/// ```
/// use bio::alignment::distance::*;
///
/// let x = b"GACTATATCGA";
/// let y = b"TTTAGCTC";
/// match hamming(x, y) {
///     Ok(s)  => println!("Score is: {}", s),  // No score because strings are of unequal sizes.
///     Err(e) => println!("Error: {}", e),     // Will print "Error: hamming distance: strings are of unequal sizes!"
/// }
/// assert!(hamming(x, y).is_err());
/// ```
pub fn hamming(alpha: TextSlice, beta: TextSlice) -> Result<u32, &'static str> {
    if alpha.len() == beta.len() {
        let mut score: u32 = 0;
        for (a, b) in Zip::new((alpha, beta)) {
            if a != b {
                score += 1;
            }
        }
        Ok(score)
    } else {
        Err("hamming distance: strings are of unequal sizes!")
    }
}


/// Compute the Levenshtein (or Edit) distance between two strings with `levenshtein`. It returns a distance between two strings,
/// i.e. minimal number of mismatches, insertions and deletions between two strings.
///
/// # Example
///
/// ```
/// use bio::alignment::distance::*;
///
/// let x = b"ACCGTGGAT";
/// let y = b"AAAAACCGTTGAT";
/// // ----ACCGTGGAT
/// //     ||||| |||
/// // AAAAACCGTTGAT
/// let l_score = levenshtein(x, y);  // Score is 5
/// assert_eq!(l_score, 5);
/// ```
#[allow(unused_assignments)]
pub fn levenshtein(alpha: TextSlice, beta: TextSlice) -> u32 {
    let mut columns = [vec!(0u32; alpha.len() + 1), vec!(0u32; alpha.len() + 1)];
    let mut i_prev = 0;
    let mut i_cur = 1;

    for i in 0..columns[0].len() {
        columns[0][i] = i as u32;
    }

    for (j, item) in beta.iter().enumerate() {
        i_cur %= 2;
        i_prev = 1 - i_cur;

        columns[i_cur][0] = 1 + j as u32;
        for i in 1..columns[0].len() {
            columns[i_cur][i] = min(columns[i_prev][i - 1] +
                                    if alpha[i - 1] == *item {
                                        0
                                    } else {
                                        1
                                    },
                                    min(columns[i_cur][i - 1] + 1, columns[i_prev][i] + 1));
        }

        i_cur += 1;
    }

    columns[i_cur - 1][columns[0].len() - 1]
}


#[cfg(test)]
mod tests {
    use super::*;

    #[test]
    fn test_hamming_dist_good() {
        let x = b"GTCTGCATGCG";
        let y = b"TTTAGCTAGCG";
        // GTCTGCATGCG
        //  |  ||  |||
        // TTTAGCTAGCG
        assert!(hamming(x, y).is_ok());
        assert_eq!(hamming(x, y).unwrap(), 5);
    }

    #[test]
    fn test_hamming_dist_bad() {
        let x = b"GACTATATCGA";
        let y = b"TTTAGCTC";
        assert!(hamming(x, y).is_err());
    }

    #[test]
    fn test_levenshtein_dist() {
        let x = b"ACCGTGGAT";
        let y = b"AAAAACCGTTGAT";
        // ----ACCGTGGAT
        //     ||||| |||
        // AAAAACCGTTGAT
        assert_eq!(levenshtein(x, y), 5);
        assert_eq!(levenshtein(x, y), levenshtein(y, x));
        assert_eq!(levenshtein(b"AAA", b"TTTT"), 4);
        assert_eq!(levenshtein(b"TTTT", b"AAA"), 4);
    }
}