#include "config.h"
#include <stdlib.h>
#if defined(_MSC_VER)
#include <intrin.h>
#else
#include <emmintrin.h>
#endif
#include "parasail.h"
#include "parasail/memory.h"
#include "parasail/internal_sse.h"
static inline __m128i _mm_blendv_epi8_rpl(__m128i a, __m128i b, __m128i mask) {
a = _mm_andnot_si128(mask, a);
a = _mm_or_si128(a, _mm_and_si128(mask, b));
return a;
}
#ifdef PARASAIL_TABLE
static inline void arr_store_si128(
int *array,
__m128i vWH,
int32_t i,
int32_t s1Len,
int32_t j,
int32_t s2Len)
{
if (0 <= i+0 && i+0 < s1Len && 0 <= j-0 && j-0 < s2Len) {
array[1LL*(i+0)*s2Len + (j-0)] = (int16_t)_mm_extract_epi16(vWH, 7);
}
if (0 <= i+1 && i+1 < s1Len && 0 <= j-1 && j-1 < s2Len) {
array[1LL*(i+1)*s2Len + (j-1)] = (int16_t)_mm_extract_epi16(vWH, 6);
}
if (0 <= i+2 && i+2 < s1Len && 0 <= j-2 && j-2 < s2Len) {
array[1LL*(i+2)*s2Len + (j-2)] = (int16_t)_mm_extract_epi16(vWH, 5);
}
if (0 <= i+3 && i+3 < s1Len && 0 <= j-3 && j-3 < s2Len) {
array[1LL*(i+3)*s2Len + (j-3)] = (int16_t)_mm_extract_epi16(vWH, 4);
}
if (0 <= i+4 && i+4 < s1Len && 0 <= j-4 && j-4 < s2Len) {
array[1LL*(i+4)*s2Len + (j-4)] = (int16_t)_mm_extract_epi16(vWH, 3);
}
if (0 <= i+5 && i+5 < s1Len && 0 <= j-5 && j-5 < s2Len) {
array[1LL*(i+5)*s2Len + (j-5)] = (int16_t)_mm_extract_epi16(vWH, 2);
}
if (0 <= i+6 && i+6 < s1Len && 0 <= j-6 && j-6 < s2Len) {
array[1LL*(i+6)*s2Len + (j-6)] = (int16_t)_mm_extract_epi16(vWH, 1);
}
if (0 <= i+7 && i+7 < s1Len && 0 <= j-7 && j-7 < s2Len) {
array[1LL*(i+7)*s2Len + (j-7)] = (int16_t)_mm_extract_epi16(vWH, 0);
}
}
#endif
#ifdef PARASAIL_ROWCOL
static inline void arr_store_rowcol(
int *row,
int *col,
__m128i vWH,
int32_t i,
int32_t s1Len,
int32_t j,
int32_t s2Len)
{
if (i+0 == s1Len-1 && 0 <= j-0 && j-0 < s2Len) {
row[j-0] = (int16_t)_mm_extract_epi16(vWH, 7);
}
if (j-0 == s2Len-1 && 0 <= i+0 && i+0 < s1Len) {
col[(i+0)] = (int16_t)_mm_extract_epi16(vWH, 7);
}
if (i+1 == s1Len-1 && 0 <= j-1 && j-1 < s2Len) {
row[j-1] = (int16_t)_mm_extract_epi16(vWH, 6);
}
if (j-1 == s2Len-1 && 0 <= i+1 && i+1 < s1Len) {
col[(i+1)] = (int16_t)_mm_extract_epi16(vWH, 6);
}
if (i+2 == s1Len-1 && 0 <= j-2 && j-2 < s2Len) {
row[j-2] = (int16_t)_mm_extract_epi16(vWH, 5);
}
if (j-2 == s2Len-1 && 0 <= i+2 && i+2 < s1Len) {
col[(i+2)] = (int16_t)_mm_extract_epi16(vWH, 5);
}
if (i+3 == s1Len-1 && 0 <= j-3 && j-3 < s2Len) {
row[j-3] = (int16_t)_mm_extract_epi16(vWH, 4);
}
if (j-3 == s2Len-1 && 0 <= i+3 && i+3 < s1Len) {
col[(i+3)] = (int16_t)_mm_extract_epi16(vWH, 4);
}
if (i+4 == s1Len-1 && 0 <= j-4 && j-4 < s2Len) {
row[j-4] = (int16_t)_mm_extract_epi16(vWH, 3);
}
if (j-4 == s2Len-1 && 0 <= i+4 && i+4 < s1Len) {
col[(i+4)] = (int16_t)_mm_extract_epi16(vWH, 3);
}
if (i+5 == s1Len-1 && 0 <= j-5 && j-5 < s2Len) {
row[j-5] = (int16_t)_mm_extract_epi16(vWH, 2);
}
if (j-5 == s2Len-1 && 0 <= i+5 && i+5 < s1Len) {
col[(i+5)] = (int16_t)_mm_extract_epi16(vWH, 2);
}
if (i+6 == s1Len-1 && 0 <= j-6 && j-6 < s2Len) {
row[j-6] = (int16_t)_mm_extract_epi16(vWH, 1);
}
if (j-6 == s2Len-1 && 0 <= i+6 && i+6 < s1Len) {
col[(i+6)] = (int16_t)_mm_extract_epi16(vWH, 1);
}
if (i+7 == s1Len-1 && 0 <= j-7 && j-7 < s2Len) {
row[j-7] = (int16_t)_mm_extract_epi16(vWH, 0);
}
if (j-7 == s2Len-1 && 0 <= i+7 && i+7 < s1Len) {
col[(i+7)] = (int16_t)_mm_extract_epi16(vWH, 0);
}
}
#endif
#ifdef PARASAIL_TABLE
#define FNAME parasail_nw_stats_table_diag_sse2_128_16
#else
#ifdef PARASAIL_ROWCOL
#define FNAME parasail_nw_stats_rowcol_diag_sse2_128_16
#else
#define FNAME parasail_nw_stats_diag_sse2_128_16
#endif
#endif
parasail_result_t* FNAME(
const char * const restrict _s1, const int _s1Len,
const char * const restrict _s2, const int s2Len,
const int open, const int gap, const parasail_matrix_t *matrix)
{
int32_t N = 0;
int32_t PAD = 0;
int32_t PAD2 = 0;
int32_t s1Len = 0;
int32_t s1Len_PAD = 0;
int32_t s2Len_PAD = 0;
int16_t * restrict s1 = NULL;
int16_t * restrict s2B = NULL;
int16_t * restrict _H_pr = NULL;
int16_t * restrict _HM_pr = NULL;
int16_t * restrict _HS_pr = NULL;
int16_t * restrict _HL_pr = NULL;
int16_t * restrict _F_pr = NULL;
int16_t * restrict _FM_pr = NULL;
int16_t * restrict _FS_pr = NULL;
int16_t * restrict _FL_pr = NULL;
int16_t * restrict s2 = NULL;
int16_t * restrict H_pr = NULL;
int16_t * restrict HM_pr = NULL;
int16_t * restrict HS_pr = NULL;
int16_t * restrict HL_pr = NULL;
int16_t * restrict F_pr = NULL;
int16_t * restrict FM_pr = NULL;
int16_t * restrict FS_pr = NULL;
int16_t * restrict FL_pr = NULL;
parasail_result_t *result = NULL;
int32_t i = 0;
int32_t j = 0;
int32_t end_query = 0;
int32_t end_ref = 0;
int16_t NEG_LIMIT = 0;
int16_t POS_LIMIT = 0;
int16_t score = 0;
int16_t matches = 0;
int16_t similar = 0;
int16_t length = 0;
__m128i vNegLimit;
__m128i vPosLimit;
__m128i vSaturationCheckMin;
__m128i vSaturationCheckMax;
__m128i vNegInf;
__m128i vOpen;
__m128i vGap;
__m128i vZero;
__m128i vNegInf0;
__m128i vOne;
__m128i vN;
__m128i vGapN;
__m128i vNegOne;
__m128i vI;
__m128i vJreset;
__m128i vMaxH;
__m128i vMaxM;
__m128i vMaxS;
__m128i vMaxL;
__m128i vILimit;
__m128i vILimit1;
__m128i vJLimit;
__m128i vJLimit1;
__m128i vIBoundary;
PARASAIL_CHECK_NULL(_s2);
PARASAIL_CHECK_GT0(s2Len);
PARASAIL_CHECK_GE0(open);
PARASAIL_CHECK_GE0(gap);
PARASAIL_CHECK_NULL(matrix);
if (matrix->type == PARASAIL_MATRIX_TYPE_PSSM) {
PARASAIL_CHECK_NULL_PSSM_STATS(_s1);
}
else {
PARASAIL_CHECK_NULL(_s1);
PARASAIL_CHECK_GT0(_s1Len);
}
N = 8;
PAD = N-1;
PAD2 = PAD*2;
s1Len = matrix->type == PARASAIL_MATRIX_TYPE_SQUARE ? _s1Len : matrix->length;
s1Len_PAD = s1Len+PAD;
s2Len_PAD = s2Len+PAD;
i = 0;
j = 0;
end_query = s1Len-1;
end_ref = s2Len-1;
NEG_LIMIT = (-open < matrix->min ? INT16_MIN + open : INT16_MIN - matrix->min) + 1;
POS_LIMIT = INT16_MAX - matrix->max - 1;
score = NEG_LIMIT;
matches = NEG_LIMIT;
similar = NEG_LIMIT;
length = NEG_LIMIT;
vNegLimit = _mm_set1_epi16(NEG_LIMIT);
vPosLimit = _mm_set1_epi16(POS_LIMIT);
vSaturationCheckMin = vPosLimit;
vSaturationCheckMax = vNegLimit;
vNegInf = _mm_set1_epi16(NEG_LIMIT);
vOpen = _mm_set1_epi16(open);
vGap = _mm_set1_epi16(gap);
vZero = _mm_set1_epi16(0);
vNegInf0 = _mm_insert_epi16(vZero, NEG_LIMIT, 7);
vOne = _mm_set1_epi16(1);
vN = _mm_set1_epi16(N);
vGapN = _mm_set1_epi16(gap*N);
vNegOne = _mm_set1_epi16(-1);
vI = _mm_set_epi16(0,1,2,3,4,5,6,7);
vJreset = _mm_set_epi16(0,-1,-2,-3,-4,-5,-6,-7);
vMaxH = vNegInf;
vMaxM = vNegInf;
vMaxS = vNegInf;
vMaxL = vNegInf;
vILimit = _mm_set1_epi16(s1Len);
vILimit1 = _mm_subs_epi16(vILimit, vOne);
vJLimit = _mm_set1_epi16(s2Len);
vJLimit1 = _mm_subs_epi16(vJLimit, vOne);
vIBoundary = _mm_set_epi16(
-open-0*gap,
-open-1*gap,
-open-2*gap,
-open-3*gap,
-open-4*gap,
-open-5*gap,
-open-6*gap,
-open-7*gap);
#ifdef PARASAIL_TABLE
result = parasail_result_new_table3(s1Len, s2Len);
#else
#ifdef PARASAIL_ROWCOL
result = parasail_result_new_rowcol3(s1Len, s2Len);
#else
result = parasail_result_new_stats();
#endif
#endif
if (!result) return NULL;
result->flag |= PARASAIL_FLAG_NW | PARASAIL_FLAG_DIAG
| PARASAIL_FLAG_STATS
| PARASAIL_FLAG_BITS_16 | PARASAIL_FLAG_LANES_8;
#ifdef PARASAIL_TABLE
result->flag |= PARASAIL_FLAG_TABLE;
#endif
#ifdef PARASAIL_ROWCOL
result->flag |= PARASAIL_FLAG_ROWCOL;
#endif
s1 = parasail_memalign_int16_t(16, s1Len+PAD);
s2B = parasail_memalign_int16_t(16, s2Len+PAD2);
_H_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_HM_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_HS_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_HL_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_F_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_FM_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_FS_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
_FL_pr = parasail_memalign_int16_t(16, s2Len+PAD2);
s2 = s2B+PAD;
H_pr = _H_pr+PAD;
HM_pr = _HM_pr+PAD;
HS_pr = _HS_pr+PAD;
HL_pr = _HL_pr+PAD;
F_pr = _F_pr+PAD;
FM_pr = _FM_pr+PAD;
FS_pr = _FS_pr+PAD;
FL_pr = _FL_pr+PAD;
if (!s1) return NULL;
if (!s2B) return NULL;
if (!_H_pr) return NULL;
if (!_HM_pr) return NULL;
if (!_HS_pr) return NULL;
if (!_HL_pr) return NULL;
if (!_F_pr) return NULL;
if (!_FM_pr) return NULL;
if (!_FS_pr) return NULL;
if (!_FL_pr) return NULL;
for (i=0; i<s1Len; ++i) {
s1[i] = matrix->mapper[(unsigned char)_s1[i]];
}
for (i=s1Len; i<s1Len_PAD; ++i) {
s1[i] = 0;
}
for (j=0; j<s2Len; ++j) {
s2[j] = matrix->mapper[(unsigned char)_s2[j]];
}
for (j=-PAD; j<0; ++j) {
s2[j] = 0;
}
for (j=s2Len; j<s2Len_PAD; ++j) {
s2[j] = 0;
}
for (j=0; j<s2Len; ++j) {
H_pr[j] = -open - j*gap;
HM_pr[j] = 0;
HS_pr[j] = 0;
HL_pr[j] = 0;
F_pr[j] = NEG_LIMIT;
FM_pr[j] = 0;
FS_pr[j] = 0;
FL_pr[j] = 0;
}
for (j=-PAD; j<0; ++j) {
H_pr[j] = 0;
HM_pr[j] = 0;
HS_pr[j] = 0;
HL_pr[j] = 0;
F_pr[j] = 0;
FM_pr[j] = 0;
FS_pr[j] = 0;
FL_pr[j] = 0;
}
for (j=s2Len; j<s2Len+PAD; ++j) {
H_pr[j] = 0;
HM_pr[j] = 0;
HS_pr[j] = 0;
HL_pr[j] = 0;
F_pr[j] = 0;
FM_pr[j] = 0;
FS_pr[j] = 0;
FL_pr[j] = 0;
}
H_pr[-1] = 0;
for (i=0; i<s1Len; i+=N) {
__m128i case1 = vZero;
__m128i case2 = vZero;
__m128i vNH = vZero;
__m128i vNM = vZero;
__m128i vNS = vZero;
__m128i vNL = vZero;
__m128i vWH = vZero;
__m128i vWM = vZero;
__m128i vWS = vZero;
__m128i vWL = vZero;
__m128i vE = vNegInf0;
__m128i vE_opn = vNegInf;
__m128i vE_ext = vNegInf;
__m128i vEM = vZero;
__m128i vES = vZero;
__m128i vEL = vZero;
__m128i vF = vNegInf0;
__m128i vF_opn = vNegInf;
__m128i vF_ext = vNegInf;
__m128i vFM = vZero;
__m128i vFS = vZero;
__m128i vFL = vZero;
__m128i vJ = vJreset;
__m128i vs1 = _mm_set_epi16(
s1[i+0],
s1[i+1],
s1[i+2],
s1[i+3],
s1[i+4],
s1[i+5],
s1[i+6],
s1[i+7]);
__m128i vs2 = vNegInf;
const int * const restrict matrow0 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+0] : ((i+0 >= s1Len) ? s1Len-1 : i+0))];
const int * const restrict matrow1 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+1] : ((i+1 >= s1Len) ? s1Len-1 : i+1))];
const int * const restrict matrow2 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+2] : ((i+2 >= s1Len) ? s1Len-1 : i+2))];
const int * const restrict matrow3 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+3] : ((i+3 >= s1Len) ? s1Len-1 : i+3))];
const int * const restrict matrow4 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+4] : ((i+4 >= s1Len) ? s1Len-1 : i+4))];
const int * const restrict matrow5 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+5] : ((i+5 >= s1Len) ? s1Len-1 : i+5))];
const int * const restrict matrow6 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+6] : ((i+6 >= s1Len) ? s1Len-1 : i+6))];
const int * const restrict matrow7 = &matrix->matrix[matrix->size * ((matrix->type == PARASAIL_MATRIX_TYPE_SQUARE) ? s1[i+7] : ((i+7 >= s1Len) ? s1Len-1 : i+7))];
vNH = _mm_insert_epi16(vNH, H_pr[-1], 7);
vWH = _mm_insert_epi16(vWH, -open - i*gap, 7);
H_pr[-1] = -open - (i+N)*gap;
for (j=0; j<s2Len+PAD; ++j) {
__m128i vMat;
__m128i vNWH = vNH;
__m128i vNWM = vNM;
__m128i vNWS = vNS;
__m128i vNWL = vNL;
vNH = _mm_srli_si128(vWH, 2);
vNH = _mm_insert_epi16(vNH, H_pr[j], 7);
vNM = _mm_srli_si128(vWM, 2);
vNM = _mm_insert_epi16(vNM, HM_pr[j], 7);
vNS = _mm_srli_si128(vWS, 2);
vNS = _mm_insert_epi16(vNS, HS_pr[j], 7);
vNL = _mm_srli_si128(vWL, 2);
vNL = _mm_insert_epi16(vNL, HL_pr[j], 7);
vF = _mm_srli_si128(vF, 2);
vF = _mm_insert_epi16(vF, F_pr[j], 7);
vFM = _mm_srli_si128(vFM, 2);
vFM = _mm_insert_epi16(vFM, FM_pr[j], 7);
vFS = _mm_srli_si128(vFS, 2);
vFS = _mm_insert_epi16(vFS, FS_pr[j], 7);
vFL = _mm_srli_si128(vFL, 2);
vFL = _mm_insert_epi16(vFL, FL_pr[j], 7);
vF_opn = _mm_subs_epi16(vNH, vOpen);
vF_ext = _mm_subs_epi16(vF, vGap);
vF = _mm_max_epi16(vF_opn, vF_ext);
case1 = _mm_cmpgt_epi16(vF_opn, vF_ext);
vFM = _mm_blendv_epi8_rpl(vFM, vNM, case1);
vFS = _mm_blendv_epi8_rpl(vFS, vNS, case1);
vFL = _mm_blendv_epi8_rpl(vFL, vNL, case1);
vFL = _mm_adds_epi16(vFL, vOne);
vE_opn = _mm_subs_epi16(vWH, vOpen);
vE_ext = _mm_subs_epi16(vE, vGap);
vE = _mm_max_epi16(vE_opn, vE_ext);
case1 = _mm_cmpgt_epi16(vE_opn, vE_ext);
vEM = _mm_blendv_epi8_rpl(vEM, vWM, case1);
vES = _mm_blendv_epi8_rpl(vES, vWS, case1);
vEL = _mm_blendv_epi8_rpl(vEL, vWL, case1);
vEL = _mm_adds_epi16(vEL, vOne);
vs2 = _mm_srli_si128(vs2, 2);
vs2 = _mm_insert_epi16(vs2, s2[j], 7);
vMat = _mm_set_epi16(
matrow0[s2[j-0]],
matrow1[s2[j-1]],
matrow2[s2[j-2]],
matrow3[s2[j-3]],
matrow4[s2[j-4]],
matrow5[s2[j-5]],
matrow6[s2[j-6]],
matrow7[s2[j-7]]
);
vNWH = _mm_adds_epi16(vNWH, vMat);
vWH = _mm_max_epi16(vNWH, vE);
vWH = _mm_max_epi16(vWH, vF);
case1 = _mm_cmpeq_epi16(vWH, vNWH);
case2 = _mm_cmpeq_epi16(vWH, vF);
vWM = _mm_blendv_epi8_rpl(
_mm_blendv_epi8_rpl(vEM, vFM, case2),
_mm_adds_epi16(vNWM,
_mm_and_si128(
_mm_cmpeq_epi16(vs1,vs2),
vOne)),
case1);
vWS = _mm_blendv_epi8_rpl(
_mm_blendv_epi8_rpl(vES, vFS, case2),
_mm_adds_epi16(vNWS,
_mm_and_si128(
_mm_cmpgt_epi16(vMat,vZero),
vOne)),
case1);
vWL = _mm_blendv_epi8_rpl(
_mm_blendv_epi8_rpl(vEL, vFL, case2),
_mm_adds_epi16(vNWL, vOne), case1);
{
__m128i cond = _mm_cmpeq_epi16(vJ,vNegOne);
vWH = _mm_blendv_epi8_rpl(vWH, vIBoundary, cond);
vWM = _mm_andnot_si128(cond, vWM);
vWS = _mm_andnot_si128(cond, vWS);
vWL = _mm_andnot_si128(cond, vWL);
vE = _mm_blendv_epi8_rpl(vE, vNegInf, cond);
vEM = _mm_andnot_si128(cond, vEM);
vES = _mm_andnot_si128(cond, vES);
vEL = _mm_andnot_si128(cond, vEL);
}
if (j > PAD) {
vSaturationCheckMin = _mm_min_epi16(vSaturationCheckMin, vWH);
vSaturationCheckMax = _mm_max_epi16(vSaturationCheckMax, vWH);
vSaturationCheckMax = _mm_max_epi16(vSaturationCheckMax, vWM);
vSaturationCheckMax = _mm_max_epi16(vSaturationCheckMax, vWS);
vSaturationCheckMax = _mm_max_epi16(vSaturationCheckMax, vWL);
vSaturationCheckMax = _mm_max_epi16(vSaturationCheckMax, vWL);
}
#ifdef PARASAIL_TABLE
arr_store_si128(result->stats->tables->score_table, vWH, i, s1Len, j, s2Len);
arr_store_si128(result->stats->tables->matches_table, vWM, i, s1Len, j, s2Len);
arr_store_si128(result->stats->tables->similar_table, vWS, i, s1Len, j, s2Len);
arr_store_si128(result->stats->tables->length_table, vWL, i, s1Len, j, s2Len);
#endif
#ifdef PARASAIL_ROWCOL
arr_store_rowcol(result->stats->rowcols->score_row, result->stats->rowcols->score_col, vWH, i, s1Len, j, s2Len);
arr_store_rowcol(result->stats->rowcols->matches_row, result->stats->rowcols->matches_col, vWM, i, s1Len, j, s2Len);
arr_store_rowcol(result->stats->rowcols->similar_row, result->stats->rowcols->similar_col, vWS, i, s1Len, j, s2Len);
arr_store_rowcol(result->stats->rowcols->length_row, result->stats->rowcols->length_col, vWL, i, s1Len, j, s2Len);
#endif
H_pr[j-7] = (int16_t)_mm_extract_epi16(vWH,0);
HM_pr[j-7] = (int16_t)_mm_extract_epi16(vWM,0);
HS_pr[j-7] = (int16_t)_mm_extract_epi16(vWS,0);
HL_pr[j-7] = (int16_t)_mm_extract_epi16(vWL,0);
F_pr[j-7] = (int16_t)_mm_extract_epi16(vF,0);
FM_pr[j-7] = (int16_t)_mm_extract_epi16(vFM,0);
FS_pr[j-7] = (int16_t)_mm_extract_epi16(vFS,0);
FL_pr[j-7] = (int16_t)_mm_extract_epi16(vFL,0);
{
__m128i cond_valid_I = _mm_cmpeq_epi16(vI, vILimit1);
__m128i cond_valid_J = _mm_cmpeq_epi16(vJ, vJLimit1);
__m128i cond_all = _mm_and_si128(cond_valid_I, cond_valid_J);
vMaxH = _mm_blendv_epi8_rpl(vMaxH, vWH, cond_all);
vMaxM = _mm_blendv_epi8_rpl(vMaxM, vWM, cond_all);
vMaxS = _mm_blendv_epi8_rpl(vMaxS, vWS, cond_all);
vMaxL = _mm_blendv_epi8_rpl(vMaxL, vWL, cond_all);
}
vJ = _mm_adds_epi16(vJ, vOne);
}
vI = _mm_adds_epi16(vI, vN);
vIBoundary = _mm_subs_epi16(vIBoundary, vGapN);
}
for (i=0; i<N; ++i) {
int16_t value;
value = (int16_t) _mm_extract_epi16(vMaxH, 7);
if (value > score) {
score = value;
matches = (int16_t) _mm_extract_epi16(vMaxM, 7);
similar = (int16_t) _mm_extract_epi16(vMaxS, 7);
length= (int16_t) _mm_extract_epi16(vMaxL, 7);
}
vMaxH = _mm_slli_si128(vMaxH, 2);
vMaxM = _mm_slli_si128(vMaxM, 2);
vMaxS = _mm_slli_si128(vMaxS, 2);
vMaxL = _mm_slli_si128(vMaxL, 2);
}
if (_mm_movemask_epi8(_mm_or_si128(
_mm_cmplt_epi16(vSaturationCheckMin, vNegLimit),
_mm_cmpgt_epi16(vSaturationCheckMax, vPosLimit)))) {
result->flag |= PARASAIL_FLAG_SATURATED;
score = 0;
matches = 0;
similar = 0;
length = 0;
end_query = 0;
end_ref = 0;
}
result->score = score;
result->end_query = end_query;
result->end_ref = end_ref;
result->stats->matches = matches;
result->stats->similar = similar;
result->stats->length = length;
parasail_free(_FL_pr);
parasail_free(_FS_pr);
parasail_free(_FM_pr);
parasail_free(_F_pr);
parasail_free(_HL_pr);
parasail_free(_HS_pr);
parasail_free(_HM_pr);
parasail_free(_H_pr);
parasail_free(s2B);
parasail_free(s1);
return result;
}