1 module parasail; 2 import std.meta:AliasSeq; 3 import std.stdio:FILE; 4 5 extern(C): 6 7 //Struct declarations 8 struct parasail_result_t { 9 int score; 10 int end_query; 11 int end_ref; 12 int flag; 13 void *extra; 14 } 15 struct parasail_matrix_t; 16 struct parasail_cigar_t{ 17 uint *seq; 18 int len; 19 int beg_query; 20 int beg_ref; 21 } 22 //Matrix functions 23 parasail_matrix_t* parasail_matrix_lookup(char * matrix); 24 parasail_matrix_t* parasail_matrix_create(const char *alphabet, const int match, const int mismatch); 25 void parasail_matrix_set_value(parasail_matrix_t *matrix, int row, int col, int value); 26 void parasail_matrix_free(parasail_matrix_t* matrix); 27 28 29 alias ALIGNMENT_TYPES = AliasSeq!("nw","sg","sg_qb","sg_qe","sg_qx","sg_db","sg_de","sg_dx","sg_qb_de","sg_qe_db","sw"); 30 alias OPTIONS = AliasSeq!("_stats","_table","_rowcol","_scan"); 31 alias VEC_TYPES = AliasSeq!("_striped","_scan","_diag"); 32 alias VEC_BITS = AliasSeq!("_8","_16","_32","_64","_sat"); 33 alias VEC_OPTIONS = AliasSeq!("_sse2_128","_sse41_128","_avx2_256","_altivec_128","_neon_128"); 34 35 enum FUNC_SIG="(const char * s1,const int s1Len,"~ 36 "const char *s2,const int s2Len,"~ 37 "const int open,const int gap,"~ 38 "const parasail_matrix_t* matrix);\n"; 39 40 enum FUNC_PRE="parasail_result_t* parasail_"; 41 //Meta-programming to generate all alignment functions 42 //Non-vectorized, reference implementations. 43 //parasail_ {nw,sg,sg_qb,sg_qe,sg_qx,sg_db,sg_de,sg_dx,sg_qb_de,sg_qe_db,sw} [_stats] [{_table,_rowcol}] [_scan] 44 string generateNonVecFuns(){ 45 string funcs; 46 foreach (aln; ALIGNMENT_TYPES) 47 { 48 funcs~=FUNC_PRE~aln~FUNC_SIG; 49 } 50 foreach (aln; ALIGNMENT_TYPES) 51 { 52 funcs~=FUNC_PRE~aln~OPTIONS[0]~FUNC_SIG; 53 funcs~=FUNC_PRE~aln~OPTIONS[1]~FUNC_SIG; 54 funcs~=FUNC_PRE~aln~OPTIONS[2]~FUNC_SIG; 55 funcs~=FUNC_PRE~aln~OPTIONS[3]~FUNC_SIG; 56 57 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[1]~FUNC_SIG; 58 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[2]~FUNC_SIG; 59 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[3]~FUNC_SIG; 60 61 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[1]~OPTIONS[3]~FUNC_SIG; 62 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[2]~OPTIONS[3]~FUNC_SIG; 63 64 funcs~=FUNC_PRE~aln~OPTIONS[1]~OPTIONS[3]~FUNC_SIG; 65 funcs~=FUNC_PRE~aln~OPTIONS[2]~OPTIONS[3]~FUNC_SIG; 66 } 67 return funcs; 68 } 69 mixin(generateNonVecFuns); 70 71 //Non-vectorized, traceback-capable reference implementations. 72 //parasail_ {nw,sg,sg_qb,sg_qe,sg_qx,sg_db,sg_de,sg_dx,sg_qb_de,sg_qe_db,sw} _trace [_scan] 73 string generateNonVecTraceFuns(){ 74 string funcs; 75 foreach (aln; ALIGNMENT_TYPES) 76 { 77 funcs~=FUNC_PRE~aln~"_trace"~FUNC_SIG; 78 } 79 foreach (aln; ALIGNMENT_TYPES) 80 { 81 funcs~=FUNC_PRE~aln~"_trace"~"_scan"~FUNC_SIG; 82 } 83 return funcs; 84 } 85 mixin(generateNonVecTraceFuns); 86 87 //Vectorized. 88 //parasail_ {nw,sg,sg_qb,sg_qe,sg_qx,sg_db,sg_de,sg_dx,sg_qb_de,sg_qe_db,sw} [_stats] [{_table,_rowcol}] {_striped,_scan,_diag} [{_sse2_128,_sse41_128,_avx2_256,_altivec_128,_neon_128}] {_8,_16,_32,_64,_sat} 89 string generateVecFuns(){ 90 string funcs; 91 foreach (aln; ALIGNMENT_TYPES) 92 { 93 foreach (type; VEC_TYPES) 94 { 95 foreach (bit; VEC_BITS) 96 { 97 funcs~=FUNC_PRE~aln~type~bit~FUNC_SIG; 98 } 99 } 100 } 101 foreach (aln; ALIGNMENT_TYPES) 102 { 103 foreach (type; VEC_TYPES) 104 { 105 foreach (bit; VEC_BITS) 106 { 107 foreach (vec; VEC_OPTIONS) 108 { 109 funcs~=FUNC_PRE~aln~type~vec~bit~FUNC_SIG; 110 funcs~=FUNC_PRE~aln~OPTIONS[0]~type~vec~bit~FUNC_SIG; 111 funcs~=FUNC_PRE~aln~OPTIONS[1]~type~vec~bit~FUNC_SIG; 112 funcs~=FUNC_PRE~aln~OPTIONS[2]~type~vec~bit~FUNC_SIG; 113 114 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[1]~type~vec~bit~FUNC_SIG; 115 funcs~=FUNC_PRE~aln~OPTIONS[0]~OPTIONS[2]~type~vec~bit~FUNC_SIG; 116 } 117 } 118 } 119 } 120 return funcs; 121 } 122 mixin(generateVecFuns); 123 124 //Vectorized, traceback-capable. 125 //parasail_ {nw,sg,sg_qb,sg_qe,sg_qx,sg_db,sg_de,sg_dx,sg_qb_de,sg_qe_db,sw} _trace {_striped,_scan,_diag} [{_sse2_128,_sse41_128,_avx2_256,_altivec_128,_neon_128}] {_8,_16,_32,_64,_sat} 126 string generateVecTraceFuns(){ 127 string funcs; 128 foreach (aln; ALIGNMENT_TYPES) 129 { 130 foreach (type; VEC_TYPES) 131 { 132 foreach (bit; VEC_BITS) 133 { 134 funcs~=FUNC_PRE~aln~"_trace"~type~bit~FUNC_SIG; 135 } 136 } 137 } 138 foreach (aln; ALIGNMENT_TYPES) 139 { 140 foreach (type; VEC_TYPES) 141 { 142 foreach (bit; VEC_BITS) 143 { 144 foreach (vec; VEC_OPTIONS) 145 { 146 funcs~=FUNC_PRE~aln~"_trace"~type~vec~bit~FUNC_SIG; 147 } 148 } 149 } 150 } 151 return funcs; 152 } 153 mixin(generateVecTraceFuns); 154 155 mixin(FUNC_PRE~"_nw_banded"~FUNC_SIG); 156 157 void parasail_result_free(parasail_result_t *result); 158 159 //Cigar Functions 160 void parasail_cigar_free(parasail_cigar_t *cigar); 161 char parasail_cigar_decode_op(uint cigar_int); 162 uint parasail_cigar_decode_len(uint cigar_int); 163 char* parasail_cigar_decode(parasail_cigar_t *cigar); 164 parasail_cigar_t* parasail_result_get_cigar( 165 parasail_result_t *result, 166 const char *seqA, int lena, 167 const char *seqB, int lenb, 168 const parasail_matrix_t *matrix); 169 170 int parasail_result_get_score(parasail_result_t *result); 171 int parasail_result_get_end_query(parasail_result_t *result); 172 int parasail_result_get_end_ref(parasail_result_t *result); 173 int parasail_result_get_matches(parasail_result_t *result); 174 int parasail_result_get_similar(parasail_result_t *result); 175 int parasail_result_get_length(parasail_result_t *result); 176 int* parasail_result_get_score_table(parasail_result_t *result); 177 int* parasail_result_get_matches_table(parasail_result_t *result); 178 int* parasail_result_get_similar_table(parasail_result_t *result); 179 int* parasail_result_get_length_table(parasail_result_t *result); 180 int* parasail_result_get_score_row(parasail_result_t *result); 181 int* parasail_result_get_matches_row(parasail_result_t *result); 182 int* parasail_result_get_similar_row(parasail_result_t *result); 183 int* parasail_result_get_length_row(parasail_result_t *result); 184 int* parasail_result_get_score_col(parasail_result_t *result); 185 int* parasail_result_get_matches_col(parasail_result_t *result); 186 int* parasail_result_get_similar_col(parasail_result_t *result); 187 int* parasail_result_get_length_col(parasail_result_t *result); 188 int* parasail_result_get_trace_table(parasail_result_t *result); 189 int* parasail_result_get_trace_ins_table(parasail_result_t *result); 190 int* parasail_result_get_trace_del_table(parasail_result_t *result); 191 192 int parasail_result_is_nw(parasail_result_t *result); 193 int parasail_result_is_sg(parasail_result_t *result); 194 int parasail_result_is_sw(parasail_result_t *result); 195 int parasail_result_is_saturated(parasail_result_t *result); 196 int parasail_result_is_banded(parasail_result_t *result); 197 int parasail_result_is_scan(parasail_result_t *result); 198 int parasail_result_is_striped(parasail_result_t *result); 199 int parasail_result_is_diag(parasail_result_t *result); 200 int parasail_result_is_blocked(parasail_result_t *result); 201 int parasail_result_is_stats(parasail_result_t *result); 202 int parasail_result_is_stats_table(parasail_result_t *result); 203 int parasail_result_is_stats_rowcol(parasail_result_t *result); 204 int parasail_result_is_table(parasail_result_t *result); 205 int parasail_result_is_rowcol(parasail_result_t *result); 206 int parasail_result_is_trace(parasail_result_t *result); 207 208 void parasail_traceback_generic( 209 const char *seqA, int lena, 210 const char *seqB, int lenb, 211 const char *nameA, 212 const char *nameB, 213 const parasail_matrix_t *matrix, 214 parasail_result_t *result, 215 char match, /* character to use for a match */ 216 char pos, /* character to use for a positive-value mismatch */ 217 char neg, /* character to use for a negative-value mismatch */ 218 int width, /* width of traceback to display before wrapping */ 219 int name_width, 220 int use_stats); /* if 0, don't display stats, if non-zero, summary stats displayed */ 221 222 void parasail_traceback_generic_extra( 223 const char *seqA, int lena, 224 const char *seqB, int lenb, 225 const char *nameA, 226 const char *nameB, 227 const parasail_matrix_t *matrix, 228 parasail_result_t *result, 229 char match, /* character to use for a match */ 230 char pos, /* character to use for a positive-value mismatch */ 231 char neg, /* character to use for a negative-value mismatch */ 232 int width, /* width of traceback to display before wrapping */ 233 int name_width, 234 int use_stats, /* if 0, don't display stats, if non-zero, summary stats displayed */ 235 int int_width, /* width used for reference and query indexes */ 236 FILE *stream); /* to print to custom file stream */