forked from soedinglab/MMseqs2
-
Notifications
You must be signed in to change notification settings - Fork 0
/
LinsearchIndexReader.cpp
281 lines (252 loc) · 11.6 KB
/
LinsearchIndexReader.cpp
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
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
//
// Created by Martin Steinegger on 2019-01-04.
//
#include "FileUtil.h"
#include "LinsearchIndexReader.h"
#include "PrefilteringIndexReader.h"
#include "Debug.h"
#include "Timer.h"
#include "NucleotideMatrix.h"
#include "SubstitutionMatrix.h"
#include "ReducedMatrix.h"
#include "KmerIndex.h"
#include "kmersearch.h"
#ifndef SIZE_T_MAX
#define SIZE_T_MAX ((size_t) -1)
#endif
template <int TYPE>
size_t LinsearchIndexReader::pickCenterKmer(KmerPosition<short> *hashSeqPair, size_t splitKmerCount) {
size_t writePos = 0;
size_t prevHash = hashSeqPair[0].kmer;
if (TYPE == Parameters::DBTYPE_NUCLEOTIDES) {
prevHash = BIT_SET(prevHash, 63);
}
size_t prevHashStart = 0;
size_t prevSetSize = 0;
for (size_t elementIdx = 0; elementIdx < splitKmerCount + 1; elementIdx++) {
size_t currKmer = hashSeqPair[elementIdx].kmer;
if (TYPE == Parameters::DBTYPE_NUCLEOTIDES) {
currKmer = BIT_SET(currKmer, 63);
}
if (prevHash != currKmer) {
size_t indexToPick = 0;
size_t randIdx = prevHashStart + indexToPick;
size_t kmer = hashSeqPair[randIdx].kmer;
if (TYPE == Parameters::DBTYPE_NUCLEOTIDES) {
kmer = BIT_SET(hashSeqPair[randIdx].kmer, 63);
}
// remove singletones from set
if (kmer != SIZE_T_MAX) {
hashSeqPair[writePos].kmer = hashSeqPair[randIdx].kmer;
hashSeqPair[writePos].pos = hashSeqPair[randIdx].pos;
hashSeqPair[writePos].seqLen = hashSeqPair[randIdx].seqLen;
hashSeqPair[writePos].id = hashSeqPair[randIdx].id;
writePos++;
}
prevHashStart = elementIdx;
}
if (hashSeqPair[elementIdx].kmer == SIZE_T_MAX) {
break;
}
prevSetSize++;
prevHash = hashSeqPair[elementIdx].kmer;
if (TYPE == Parameters::DBTYPE_NUCLEOTIDES) {
prevHash = BIT_SET(prevHash, 63);
}
}
hashSeqPair[writePos].kmer = SIZE_T_MAX;
return writePos;
}
template size_t LinsearchIndexReader::pickCenterKmer<0>(KmerPosition<short> *hashSeqPair, size_t splitKmerCount);
template size_t LinsearchIndexReader::pickCenterKmer<1>(KmerPosition<short> *hashSeqPair, size_t splitKmerCount);
template <int TYPE>
void LinsearchIndexReader::mergeAndWriteIndex(DBWriter & dbw, std::vector<std::string> tmpFiles, int alphSize, int kmerSize) {
KmerIndex kmerIndex(alphSize, kmerSize);
dbw.writeStart(0);
Debug(Debug::INFO) << "Merge splits ... ";
const int fileCnt = tmpFiles.size();
FILE ** files = new FILE*[fileCnt];
KmerPosition<short> **entries = new KmerPosition<short>*[fileCnt];
size_t * entrySizes = new size_t[fileCnt];
size_t * offsetPos = new size_t[fileCnt];
size_t * dataSizes = new size_t[fileCnt];
// init structures
for(size_t file = 0; file < tmpFiles.size(); file++){
files[file] = FileUtil::openFileOrDie(tmpFiles[file].c_str(),"r",true);
size_t dataSize;
entries[file] = (KmerPosition<short>*)FileUtil::mmapFile(files[file], &dataSize);
dataSizes[file] = dataSize;
entrySizes[file] = dataSize/sizeof(KmerPosition<short>);
offsetPos[file] = 0;
}
std::priority_queue<FileKmer, std::vector<FileKmer>, CompareRepSequenceAndIdAndDiag> queue;
// read one entry for each file
for(int file = 0; file < fileCnt; file++ ){
size_t offset = offsetPos[file];
if(offset < entrySizes[file]){
KmerPosition<short> currKmerPosition = entries[file][offset];
size_t currKmer = currKmerPosition.kmer;
bool isReverse = false;
if(TYPE == Parameters::DBTYPE_NUCLEOTIDES){
isReverse = (BIT_CHECK(currKmerPosition.kmer, 63) == false);
currKmer = BIT_CLEAR(currKmer, 63);
}
queue.push(FileKmer(currKmer, currKmerPosition.id, currKmerPosition.pos, currKmerPosition.seqLen, isReverse, file));
}
}
std::string prefResultsOutString;
prefResultsOutString.reserve(100000000);
FileKmer res;
size_t prevKmer = SIZE_T_MAX;
while(queue.empty() == false) {
res = queue.top();
queue.pop();
{
size_t offset = offsetPos[res.file];
if(offset + 1 < entrySizes[res.file]){
size_t currKmer = entries[res.file][offset + 1].kmer;
bool isReverse = false;
if(TYPE == Parameters::DBTYPE_NUCLEOTIDES){
isReverse = (BIT_CHECK(entries[res.file][offset + 1].kmer, 63) == false);
currKmer = BIT_CLEAR(currKmer, 63);
}
queue.push(FileKmer(currKmer, entries[res.file][offset + 1].id,
entries[res.file][offset + 1].pos, entries[res.file][offset + 1].seqLen,
isReverse, res.file));
offsetPos[res.file] = offset + 1;
}
}
if(prevKmer != res.kmer){
if (kmerIndex.needsFlush(res.kmer) == true) {
kmerIndex.flush(dbw);
}
kmerIndex.addElementSorted(res.kmer, res.id, res.pos, res.seqLen, res.reverse);
}
prevKmer = res.kmer;
}
kmerIndex.flush(dbw);
dbw.writeEnd(PrefilteringIndexReader::ENTRIES, 0);
dbw.alignToPageSize();
// clear memory
for(size_t file = 0; file < tmpFiles.size(); file++) {
fclose(files[file]);
FileUtil::munmapData((void*)entries[file], dataSizes[file]);
}
delete [] dataSizes;
delete [] offsetPos;
delete [] entries;
delete [] entrySizes;
delete [] files;
// write index
Debug(Debug::INFO) << "Write ENTRIESOFFSETS (" << PrefilteringIndexReader::ENTRIESOFFSETS << ")\n";
kmerIndex.setupOffsetTable();
dbw.writeData((char*)kmerIndex.getOffsets(),kmerIndex.getOffsetsSize()*sizeof(size_t), PrefilteringIndexReader::ENTRIESOFFSETS, 0);
dbw.alignToPageSize();
// write index
Debug(Debug::INFO) << "Write ENTRIESGRIDSIZE (" << PrefilteringIndexReader::ENTRIESGRIDSIZE << ")\n";
uint64_t gridResolution = static_cast<uint64_t >(kmerIndex.getGridResolution());
char *gridResolutionPtr = (char*) &gridResolution;
dbw.writeData(gridResolutionPtr, 1 * sizeof(uint64_t), PrefilteringIndexReader::ENTRIESGRIDSIZE, 0);
dbw.alignToPageSize();
// ENTRIESNUM
Debug(Debug::INFO) << "Write ENTRIESNUM (" << PrefilteringIndexReader::ENTRIESNUM << ")\n";
uint64_t entriesNum = kmerIndex.getTableEntriesNum();
char *entriesNumPtr = (char *) &entriesNum;
dbw.writeData(entriesNumPtr, 1 * sizeof(uint64_t), PrefilteringIndexReader::ENTRIESNUM, 0);
dbw.alignToPageSize();
}
template void LinsearchIndexReader::mergeAndWriteIndex<0>(DBWriter & dbw, std::vector<std::string> tmpFiles, int alphSize, int kmerSize);
template void LinsearchIndexReader::mergeAndWriteIndex<1>(DBWriter & dbw, std::vector<std::string> tmpFiles, int alphSize, int kmerSize);
template <int TYPE>
void LinsearchIndexReader::writeIndex(DBWriter & dbw,
KmerPosition<short> *hashSeqPair, size_t totalKmers,
int alphSize, int kmerSize) {
KmerIndex kmerIndex(alphSize - 1, kmerSize);
Debug(Debug::INFO) << "Write ENTRIES (" << PrefilteringIndexReader::ENTRIES << ")\n";
// write entries
dbw.writeStart(0);
for(size_t pos = 0; pos < totalKmers && hashSeqPair[pos].kmer != SIZE_T_MAX; pos++){
size_t kmer= hashSeqPair[pos].kmer;
bool isReverse = false;
if(TYPE == Parameters::DBTYPE_NUCLEOTIDES){
isReverse = (BIT_CHECK(hashSeqPair[pos].kmer, 63) == false);
kmer = BIT_CLEAR(kmer, 63);
}
if(kmerIndex.needsFlush(kmer) == true){
kmerIndex.flush(dbw);
}
kmerIndex.addElementSorted(kmer, hashSeqPair[pos].id, hashSeqPair[pos].pos, hashSeqPair[pos].seqLen, isReverse);
}
kmerIndex.flush(dbw);
dbw.writeEnd(PrefilteringIndexReader::ENTRIES, 0);
dbw.alignToPageSize();
// write index
Debug(Debug::INFO) << "Write ENTRIESOFFSETS (" << PrefilteringIndexReader::ENTRIESOFFSETS << ")\n";
kmerIndex.setupOffsetTable();
dbw.writeData((char*)kmerIndex.getOffsets(),kmerIndex.getOffsetsSize()*sizeof(size_t), PrefilteringIndexReader::ENTRIESOFFSETS, 0);
dbw.alignToPageSize();
// write index
Debug(Debug::INFO) << "Write ENTRIESGRIDSIZE (" << PrefilteringIndexReader::ENTRIESGRIDSIZE << ")\n";
uint64_t gridResolution = static_cast<uint64_t >(kmerIndex.getGridResolution());
char *gridResolutionPtr = (char*) &gridResolution;
dbw.writeData(gridResolutionPtr, 1 * sizeof(uint64_t), PrefilteringIndexReader::ENTRIESGRIDSIZE, 0);
dbw.alignToPageSize();
// ENTRIESNUM
Debug(Debug::INFO) << "Write ENTRIESNUM (" << PrefilteringIndexReader::ENTRIESNUM << ")\n";
uint64_t entriesNum = kmerIndex.getTableEntriesNum();
char *entriesNumPtr = (char *) &entriesNum;
dbw.writeData(entriesNumPtr, 1 * sizeof(uint64_t), PrefilteringIndexReader::ENTRIESNUM, 0);
dbw.alignToPageSize();
}
template void LinsearchIndexReader::writeIndex<0>(DBWriter & dbw,
KmerPosition<short> *hashSeqPair, size_t totalKmers,
int alphSize, int kmerSize);
template void LinsearchIndexReader::writeIndex<1>(DBWriter & dbw,
KmerPosition<short> *hashSeqPair, size_t totalKmers,
int alphSize, int kmerSize);
std::string LinsearchIndexReader::indexName(std::string baseName) {
std::string result(baseName);
result.append(".").append("linidx");
return result;
}
bool LinsearchIndexReader::checkIfIndexFile(DBReader<unsigned int> *pReader) {
char * version = pReader->getDataByDBKey(PrefilteringIndexReader::VERSION, 0);
if(version == NULL){
return false;
}
return (strncmp(version, PrefilteringIndexReader::CURRENT_VERSION, strlen(PrefilteringIndexReader::CURRENT_VERSION)) == 0 ) ? true : false;
}
void LinsearchIndexReader::writeKmerIndexToDisk(std::string fileName, KmerPosition<short> *kmers, size_t kmerCnt){
FILE* filePtr = fopen(fileName.c_str(), "wb");
if(filePtr == NULL) { perror(fileName.c_str()); EXIT(EXIT_FAILURE); }
fwrite(kmers, sizeof(KmerPosition<unsigned short>), kmerCnt, filePtr);
fclose(filePtr);
}
std::string LinsearchIndexReader::findIncompatibleParameter(DBReader<unsigned int> & index, Parameters &par, int dbtype) {
PrefilteringIndexData meta = PrefilteringIndexReader::getMetadata(&index);
if (meta.maxSeqLength != static_cast<int>(par.maxSeqLen))
return "maxSeqLen";
if (meta.seqType != dbtype)
return "seqType";
if (Parameters::isEqualDbtype(dbtype, Parameters::DBTYPE_NUCLEOTIDES) == false && meta.alphabetSize != par.alphabetSize)
return "alphabetSize";
if (meta.kmerSize != par.kmerSize)
return "kmerSize";
if (meta.mask != (par.maskMode > 0))
return "maskMode";
if (meta.spacedKmer != par.spacedKmer)
return "spacedKmer";
if (par.seedScoringMatrixFile != PrefilteringIndexReader::getSubstitutionMatrixName(&index))
return "seedScoringMatrixFile";
if (par.spacedKmerPattern != PrefilteringIndexReader::getSpacedPattern(&index))
return "spacedKmerPattern";
return "";
}
std::string LinsearchIndexReader::searchForIndex(const std::string& dbName) {
std::string outIndexName = dbName + ".linidx";
if (FileUtil::fileExists((outIndexName + ".dbtype").c_str()) == true) {
return outIndexName;
}
return "";
}
#undef SIZE_T_MAX