#include "config.h"
#include <stdint.h>
#include <stdlib.h>
#include "parasail.h"
#include "parasail/memory.h"
#include "parasail/internal_neon.h"
#define SG_SUFFIX _scan_neon_128_8
#define SG_SUFFIX_PROF _scan_profile_neon_128_8
#include "sg_helper.h"
#ifdef PARASAIL_TABLE
static inline void arr_store_si128(
int *array,
simde__m128i vH,
int32_t t,
int32_t seglen,
int32_t d,
int32_t dlen)
{
array[1LL*( 0*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 0);
array[1LL*( 1*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 1);
array[1LL*( 2*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 2);
array[1LL*( 3*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 3);
array[1LL*( 4*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 4);
array[1LL*( 5*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 5);
array[1LL*( 6*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 6);
array[1LL*( 7*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 7);
array[1LL*( 8*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 8);
array[1LL*( 9*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 9);
array[1LL*(10*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 10);
array[1LL*(11*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 11);
array[1LL*(12*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 12);
array[1LL*(13*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 13);
array[1LL*(14*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 14);
array[1LL*(15*seglen+t)*dlen + d] = (int8_t)simde_mm_extract_epi8(vH, 15);
}
#endif
#ifdef PARASAIL_ROWCOL
static inline void arr_store_col(
int *col,
simde__m128i vH,
int32_t t,
int32_t seglen)
{
col[ 0*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 0);
col[ 1*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 1);
col[ 2*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 2);
col[ 3*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 3);
col[ 4*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 4);
col[ 5*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 5);
col[ 6*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 6);
col[ 7*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 7);
col[ 8*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 8);
col[ 9*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 9);
col[10*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 10);
col[11*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 11);
col[12*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 12);
col[13*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 13);
col[14*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 14);
col[15*seglen+t] = (int8_t)simde_mm_extract_epi8(vH, 15);
}
#endif
#ifdef PARASAIL_TABLE
#define FNAME parasail_sg_flags_table_scan_neon_128_8
#define PNAME parasail_sg_flags_table_scan_profile_neon_128_8
#else
#ifdef PARASAIL_ROWCOL
#define FNAME parasail_sg_flags_rowcol_scan_neon_128_8
#define PNAME parasail_sg_flags_rowcol_scan_profile_neon_128_8
#else
#define FNAME parasail_sg_flags_scan_neon_128_8
#define PNAME parasail_sg_flags_scan_profile_neon_128_8
#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,
int s1_beg, int s1_end, int s2_beg, int s2_end)
{
parasail_profile_t *profile = NULL;
parasail_result_t *result = NULL;
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_SQUARE) {
PARASAIL_CHECK_NULL(s1);
PARASAIL_CHECK_GT0(s1Len);
}
profile = parasail_profile_create_neon_128_8(s1, s1Len, matrix);
if (!profile) return NULL;
result = PNAME(profile, s2, s2Len, open, gap, s1_beg, s1_end, s2_beg, s2_end);
parasail_profile_free(profile);
return result;
}
parasail_result_t* PNAME(
const parasail_profile_t * const restrict profile,
const char * const restrict s2, const int s2Len,
const int open, const int gap,
int s1_beg, int s1_end, int s2_beg, int s2_end)
{
int32_t i = 0;
int32_t j = 0;
int32_t k = 0;
int32_t end_query = 0;
int32_t end_ref = 0;
int32_t s1Len = 0;
const parasail_matrix_t *matrix = NULL;
int32_t segWidth = 0;
int32_t segLen = 0;
int32_t offset = 0;
int32_t position = 0;
simde__m128i* restrict pvP = NULL;
simde__m128i* restrict pvE = NULL;
int8_t* restrict boundary = NULL;
simde__m128i* restrict pvHt = NULL;
simde__m128i* restrict pvH = NULL;
simde__m128i* restrict pvGapper = NULL;
simde__m128i vGapO;
simde__m128i vGapE;
int8_t NEG_LIMIT = 0;
int8_t POS_LIMIT = 0;
simde__m128i vZero;
int8_t score = 0;
simde__m128i vNegLimit;
simde__m128i vPosLimit;
simde__m128i vSaturationCheckMin;
simde__m128i vSaturationCheckMax;
simde__m128i vMaxH;
simde__m128i vPosMask;
simde__m128i vNegInfFront;
simde__m128i vSegLenXgap;
parasail_result_t *result = NULL;
PARASAIL_CHECK_NULL(profile);
PARASAIL_CHECK_NULL(profile->profile8.score);
PARASAIL_CHECK_NULL(profile->matrix);
PARASAIL_CHECK_GT0(profile->s1Len);
PARASAIL_CHECK_NULL(s2);
PARASAIL_CHECK_GT0(s2Len);
PARASAIL_CHECK_GE0(open);
PARASAIL_CHECK_GE0(gap);
i = 0;
j = 0;
k = 0;
s1Len = profile->s1Len;
end_query = s1Len-1;
end_ref = s2Len-1;
matrix = profile->matrix;
segWidth = 16;
segLen = (s1Len + segWidth - 1) / segWidth;
offset = (s1Len - 1) % segLen;
position = (segWidth - 1) - (s1Len - 1) / segLen;
pvP = (simde__m128i*)profile->profile8.score;
vGapO = simde_mm_set1_epi8(open);
vGapE = simde_mm_set1_epi8(gap);
NEG_LIMIT = (-open < matrix->min ? INT8_MIN + open : INT8_MIN - matrix->min) + 1;
POS_LIMIT = INT8_MAX - matrix->max - 1;
vZero = simde_mm_setzero_si128();
score = NEG_LIMIT;
vNegLimit = simde_mm_set1_epi8(NEG_LIMIT);
vPosLimit = simde_mm_set1_epi8(POS_LIMIT);
vSaturationCheckMin = vPosLimit;
vSaturationCheckMax = vNegLimit;
vMaxH = vNegLimit;
vPosMask = simde_mm_cmpeq_epi8(simde_mm_set1_epi8(position),
simde_mm_set_epi8(0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15));
vNegInfFront = vZero;
vNegInfFront = simde_mm_insert_epi8(vNegInfFront, NEG_LIMIT, 0);
vSegLenXgap = simde_mm_adds_epi8(vNegInfFront,
simde_mm_slli_si128(simde_mm_set1_epi8(-segLen*gap), 1));
#ifdef PARASAIL_TABLE
result = parasail_result_new_table1(segLen*segWidth, s2Len);
#else
#ifdef PARASAIL_ROWCOL
result = parasail_result_new_rowcol1(segLen*segWidth, s2Len);
#else
result = parasail_result_new();
#endif
#endif
if (!result) return NULL;
result->flag |= PARASAIL_FLAG_SG | PARASAIL_FLAG_SCAN
| PARASAIL_FLAG_BITS_8 | PARASAIL_FLAG_LANES_16;
result->flag |= s1_beg ? PARASAIL_FLAG_SG_S1_BEG : 0;
result->flag |= s1_end ? PARASAIL_FLAG_SG_S1_END : 0;
result->flag |= s2_beg ? PARASAIL_FLAG_SG_S2_BEG : 0;
result->flag |= s2_end ? PARASAIL_FLAG_SG_S2_END : 0;
#ifdef PARASAIL_TABLE
result->flag |= PARASAIL_FLAG_TABLE;
#endif
#ifdef PARASAIL_ROWCOL
result->flag |= PARASAIL_FLAG_ROWCOL;
#endif
pvE = parasail_memalign_simde__m128i(16, segLen);
boundary = parasail_memalign_int8_t(16, s2Len+1);
pvHt= parasail_memalign_simde__m128i(16, segLen);
pvH = parasail_memalign_simde__m128i(16, segLen);
pvGapper = parasail_memalign_simde__m128i(16, segLen);
if (!pvE) return NULL;
if (!boundary) return NULL;
if (!pvHt) return NULL;
if (!pvH) return NULL;
if (!pvGapper) return NULL;
{
int32_t index = 0;
for (i=0; i<segLen; ++i) {
int32_t segNum = 0;
simde__m128i h;
simde__m128i e;
for (segNum=0; segNum<segWidth; ++segNum) {
int64_t tmp = s1_beg ? 0 : (-open-gap*(segNum*segLen+i));
h.i8[segNum] = tmp < INT8_MIN ? INT8_MIN : tmp;
tmp = tmp - open;
e.i8[segNum] = tmp < INT8_MIN ? INT8_MIN : tmp;
}
simde_mm_store_si128(&pvH[index], h);
simde_mm_store_si128(&pvE[index], e);
++index;
}
}
{
boundary[0] = 0;
for (i=1; i<=s2Len; ++i) {
int64_t tmp = s2_beg ? 0 : (-open-gap*(i-1));
boundary[i] = tmp < INT8_MIN ? INT8_MIN : tmp;
}
}
{
simde__m128i vGapper = simde_mm_subs_epi8(vZero,vGapO);
for (i=segLen-1; i>=0; --i) {
simde_mm_store_si128(pvGapper+i, vGapper);
vGapper = simde_mm_subs_epi8(vGapper, vGapE);
vSaturationCheckMin = simde_mm_min_epi8(vSaturationCheckMin, vGapper);
}
}
for (j=0; j<s2Len; ++j) {
simde__m128i vE;
simde__m128i vHt;
simde__m128i vF;
simde__m128i vH;
simde__m128i vHp;
simde__m128i *pvW;
simde__m128i vW;
vHp = simde_mm_load_si128(pvH+(segLen-1));
vHp = simde_mm_slli_si128(vHp, 1);
vHp = simde_mm_insert_epi8(vHp, boundary[j], 0);
pvW = pvP + matrix->mapper[(unsigned char)s2[j]]*segLen;
vHt = simde_mm_subs_epi8(vNegLimit, pvGapper[0]);
vF = vNegLimit;
for (i=0; i<segLen; ++i) {
vH = simde_mm_load_si128(pvH+i);
vE = simde_mm_load_si128(pvE+i);
vW = simde_mm_load_si128(pvW+i);
vE = simde_mm_max_epi8(
simde_mm_subs_epi8(vE, vGapE),
simde_mm_subs_epi8(vH, vGapO));
vHp = simde_mm_adds_epi8(vHp, vW);
vF = simde_mm_max_epi8(vF, simde_mm_adds_epi8(vHt, pvGapper[i]));
vHt = simde_mm_max_epi8(vE, vHp);
simde_mm_store_si128(pvE+i, vE);
simde_mm_store_si128(pvHt+i, vHt);
vHp = vH;
}
vHt = simde_mm_slli_si128(vHt, 1);
vHt = simde_mm_insert_epi8(vHt, boundary[j+1], 0);
vF = simde_mm_max_epi8(vF, simde_mm_adds_epi8(vHt, pvGapper[0]));
for (i=0; i<segWidth-2; ++i) {
simde__m128i vFt = simde_mm_slli_si128(vF, 1);
vFt = simde_mm_adds_epi8(vFt, vSegLenXgap);
vF = simde_mm_max_epi8(vF, vFt);
}
vF = simde_mm_slli_si128(vF, 1);
vF = simde_mm_adds_epi8(vF, vNegInfFront);
vH = simde_mm_max_epi8(vHt, vF);
for (i=0; i<segLen; ++i) {
vHt = simde_mm_load_si128(pvHt+i);
vF = simde_mm_max_epi8(
simde_mm_subs_epi8(vF, vGapE),
simde_mm_subs_epi8(vH, vGapO));
vH = simde_mm_max_epi8(vHt, vF);
simde_mm_store_si128(pvH+i, vH);
vSaturationCheckMin = simde_mm_min_epi8(vSaturationCheckMin, vH);
vSaturationCheckMin = simde_mm_min_epi8(vSaturationCheckMin, vF);
vSaturationCheckMax = simde_mm_max_epi8(vSaturationCheckMax, vH);
#ifdef PARASAIL_TABLE
arr_store_si128(result->tables->score_table, vH, i, segLen, j, s2Len);
#endif
}
{
simde__m128i vCompare;
vH = simde_mm_load_si128(pvH + offset);
vCompare = simde_mm_and_si128(vPosMask, simde_mm_cmpgt_epi8(vH, vMaxH));
vMaxH = simde_mm_max_epi8(vH, vMaxH);
if (simde_mm_movemask_epi8(vCompare)) {
end_ref = j;
}
#ifdef PARASAIL_ROWCOL
for (k=0; k<position; ++k) {
vH = simde_mm_slli_si128(vH, 1);
}
result->rowcols->score_row[j] = (int8_t) simde_mm_extract_epi8 (vH, 15);
#endif
}
}
#ifdef PARASAIL_ROWCOL
for (i=0; i<segLen; ++i) {
simde__m128i vH = simde_mm_load_si128(pvH + i);
arr_store_col(result->rowcols->score_col, vH, i, segLen);
}
#endif
if (s2_end)
{
for (k=0; k<position; ++k) {
vMaxH = simde_mm_slli_si128(vMaxH, 1);
}
score = (int8_t) simde_mm_extract_epi8(vMaxH, 15);
end_query = s1Len-1;
}
if (s1_end)
{
int8_t *t = (int8_t*)pvH;
int32_t column_len = segLen * segWidth;
for (i = 0; i<column_len; ++i, ++t) {
int32_t temp = i / segWidth + i % segWidth * segLen;
if (temp >= s1Len) continue;
if (*t > score) {
score = *t;
end_query = temp;
end_ref = s2Len-1;
}
else if (*t == score && end_ref == s2Len-1 && temp < end_query) {
end_query = temp;
}
}
}
if (!s1_end && !s2_end) {
{
simde__m128i vH = simde_mm_load_si128(pvH + offset);
for (k=0; k<position; ++k) {
vH = simde_mm_slli_si128(vH, 1);
}
score = (int8_t) simde_mm_extract_epi8 (vH, 15);
end_ref = s2Len - 1;
end_query = s1Len - 1;
}
}
if (simde_mm_movemask_epi8(simde_mm_or_si128(
simde_mm_cmplt_epi8(vSaturationCheckMin, vNegLimit),
simde_mm_cmpgt_epi8(vSaturationCheckMax, vPosLimit)))) {
result->flag |= PARASAIL_FLAG_SATURATED;
score = 0;
end_query = 0;
end_ref = 0;
}
result->score = score;
result->end_query = end_query;
result->end_ref = end_ref;
parasail_free(pvGapper);
parasail_free(pvH);
parasail_free(pvHt);
parasail_free(boundary);
parasail_free(pvE);
return result;
}
SG_IMPL_ALL
SG_IMPL_PROF_ALL