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 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397
|
/* -*- mode: C++; c-basic-offset: 4; indent-tabs-mode: nil -*- */
#include <string>
#ifdef _OPENMP
#include <omp.h>
#else
#define omp_get_max_threads() (1)
#define omp_get_num_threads() (1)
#define omp_get_thread_num() (0)
#endif
#include "base/CommandLineParser.h"
#include "analysis/DNAVector.h"
#include "analysis/NonRedKmerTable.h"
//#include "analysis/CompMgr.h"
#include <queue>
#include <map>
#include <algorithm>
#define MAX_OPEN_FILES 1000
extern "C"
{
#include <stdio.h>
#include <unistd.h>
#include <sys/types.h>
#include <sys/stat.h>
#include <fcntl.h>
#include <string.h>
#include <errno.h>
}
class Assignment
{
public:
Assignment() {
m_comp = -1;
m_read = -1;
}
void SetComp(int c) {
m_comp = c;
}
void SetRead(int r) {
m_read = r;
}
int Read() const {return m_read;}
int Comp() const {return m_comp;}
bool operator < (const Assignment & a) const {
return (m_comp < a.m_comp);
}
private:
int m_comp;
int m_read;
};
int main(int argc,char** argv)
{
commandArg<string> aStringCmmd("-i","reads fasta");
commandArg<string> fStringCmmd("-f","fasta input file (concatenated flat components)");
commandArg<string> bStringCmmd("-o","output file");
commandArg<long> mmrCmmd("-max_mem_reads","Maximum number of reads to load into memory", -1);
commandArg<bool> strandCmmd("-strand","strand specific data", false);
commandArg<int> pctReadMapCmmd("-p", "percent of read kmers require mapping to component", 0);
commandArg<bool> verboseCmmd("-verbose", "prints more status info", false);
commandArg<int> numThreadsCmmd("-t", "number of threads (default: env OMP_NUM_THREADS)", 0);
commandLineParser P(argc,argv);
P.SetDescription("Assigns reads to graph components.");
P.registerArg(aStringCmmd);
P.registerArg(fStringCmmd);
P.registerArg(bStringCmmd);
P.registerArg(mmrCmmd);
P.registerArg(strandCmmd);
P.registerArg(pctReadMapCmmd);
P.registerArg(verboseCmmd);
P.registerArg(numThreadsCmmd);
P.parse();
cerr << "-------------------------------------------" << endl
<< "---- Chrysalis: ReadsToTranscripts --------" << endl
<< "-- (Place reads on Inchworm Bundles) ------" << endl
<< "-------------------------------------------" << endl << endl;
string aString = P.GetStringValueFor(aStringCmmd); // reads.fa
string readsToTranscriptsOutputFile = P.GetStringValueFor(bStringCmmd); // out directory
string fString = P.GetStringValueFor(fStringCmmd); // inchworm bundled.fasta
long max_mem_reads = P.GetLongValueFor(mmrCmmd);
bool bStrand = P.GetBoolValueFor(strandCmmd); // strand-specific data
int pctReadMapRequired = P.GetIntValueFor(pctReadMapCmmd);
int num_threads = P.GetIntValueFor(numThreadsCmmd);
bool VERBOSE = P.GetBoolValueFor(verboseCmmd);
vecDNAVector dna;
if(max_mem_reads > 0){
cerr << "Setting maximum number of reads to load in memory to " << max_mem_reads << endl;
} else {
max_mem_reads = 2147483647; // max int
}
if (num_threads > 0) {
cerr << "-setting num threads to: " << num_threads << endl;
omp_set_num_threads(num_threads);
}
cerr << "Reading bundled inchworm contigs... " << endl;
//ML: use 1M as read buffer size
dna.Read(fString, false, false, true, 1000000);
cerr << "done!" << endl;
//test.ReverseComplement();
int k = 25;
//ComponentFileMgr compMgr(bString);
multimap<int,int> componentReadMap;
unsigned long readCount = 0;
string readCountFile = readsToTranscriptsOutputFile + ".rcts.out";
//-------------------------------------------------//
//---- Build kmer index for Iworm Bundles ---------//
//-------------------------------------------------//
NonRedKmerTable kt(k);
kt.SetUp(dna, true);
kt.SetAllCounts(-1);
map<int,int> component_number_mapping;
cerr << "Assigning kmers to Iworm bundles ... ";
#pragma omp parallel for schedule(guided,100)
for (size_t i=0; i< dna.size(); i++) {
const DNAVector & d = dna[i];
string bundle_acc = dna.Name(i);
string component_no_string = bundle_acc.substr(3); // remove >s_ prefix
int component_no = atoi(component_no_string.c_str());
#pragma omp critical
{
component_number_mapping[i] = component_no;
//cerr << "\rBundle name: [" << i << "] = " << bundle_acc << ", so component_no = " << component_no << " ";
}
for (int j=0; j<=d.isize()-k; j++) {
kt.SetCount(d, j, i); // not really counting kmers here, but instead labeling kmers according to the iworm bundles they correspond to.
}
}
cerr << "done!" << endl;
// ------------------------------------------------//
// ------- Assign reads to Iworm Bundles ----------//
// ------------------------------------------------//
//ML: this implementation does not make use of class DNAVector for the reads
map<int,bool> seen_component; // open files on first write, append on later writes.
vector<int> read_pct_mapping_info;
read_pct_mapping_info.resize(100000, 0);
DNAStringStreamFast fastaReader;
fastaReader.ReadStream(aString);
//string readsToTranscriptsOutputFile = bString + "/readsToComponents.out"; // really mapping reads to components rather than transcripts, using spectral alignment
int outFile = open(readsToTranscriptsOutputFile.c_str(), O_WRONLY | O_CREAT | O_TRUNC , 0666);
cerr << "Processing reads:" << endl;
unsigned long total_reads_read = 0;
int reads_read = 0;
do {
vector<string> readSeqVector;
vector<string> readNameVector;
cerr << " reading another " << max_mem_reads << "... ";
for(int i=0;i<max_mem_reads;i++)
{
if(!fastaReader.NextToVector(readSeqVector, readNameVector)) break;
}
reads_read = readSeqVector.size();
if (reads_read > 0) {
cerr << "done. Read " << reads_read << " reads." << endl;
} else {
// reached EOF
cerr << "finished reading reads" << endl;
//MW: If nothing was read leave do loop right away
break;
}
componentReadMap.clear();
#pragma omp parallel for schedule(guided,100)
for (int i=0; i < (int)readSeqVector.size(); i++) {
string d = readSeqVector[i];
// ensure upper case
transform(d.begin(), d.end(), d.begin(), ::toupper);
// map kmers of read to Iworm bundles
svec<int> comp;
comp.reserve(4000);
int num_kmer_pos = d.size()-k + 1;
for (int j=0; j<=(int)d.size()-k; j++) {
int c = kt.GetCountReal(d, j); // the iworm bundle containing read kmer
if (c >= 0)
comp.push_back(c);
}
if (!bStrand) {
string dd = d;
DNAVector::ReverseComplement(dd);
for (int j=0; j<=(int)dd.size()-k; j++) {
int c = kt.GetCountReal(dd, j);
if (c >= 0)
comp.push_back(c);
}
}
// find the iworm bundle with the most kmer hits
Sort(comp);
int best = -1;
int max = 0;
int run = 0;
for (int j=1; j<comp.isize(); j++) {
if (comp[j] != comp[j-1] || j+1 == comp.isize()) {
if (run > max) {
max = run;
best = comp[j-1];
}
run = 0;
} else {
run++;
}
}
//cout << "Read " << i << " maps to " << best << " with " << max << " k-mers" << endl;
int pct_read_mapped = (int) ((float)max/num_kmer_pos*100 + 0.5);
//ML: the following calculates the same, but without floating point precision issues:
//int pct_read_mapped = ( 100*max + num_kmer_pos/2 ) / num_kmer_pos;
if(best != -1 && pct_read_mapped >= pctReadMapRequired) {
// Note: multiple threads shouldn't try to reorder the tree
// at the same time (which can happen during an insert)
#pragma omp critical
{
readCount++;
componentReadMap.insert(pair<int,int>(best,i)); // i = read_index, best = best_iworm_bundle
// store pct read mapping info for later reporting
if (i > (int) read_pct_mapping_info.size()-1) {
read_pct_mapping_info.resize(i+100000);
// cerr << "resizing vec to " << i << " + 100000 " << endl;
}
read_pct_mapping_info[i] = pct_read_mapped;
}
}
else {
// no mapping for read
if (VERBOSE)
cerr << "WARNING: No component mapping for read: " << readNameVector[i] << " : " << readSeqVector[i] << endl;
}
} // end of read assignment to components for this round of streaming
if (reads_read > 0) {
total_reads_read += reads_read;
cerr << "[" << total_reads_read << "] reads analyzed for mapping." << endl;
}
//-------------------------------------------------------
// Write reads to files based on components mapped to.
//cerr << "\nWriting to files... ";
// convert sorted map to vectors for threaded processing
vector<int> mapComponents; // list of component_ids
vector<vector<int> > componentReads; // list of read-vectors that correspond to the component_ids above (synched by index)
multimap<int,int>::iterator it;
int lastComponent = -1;
vector<int> tmpReadIDs; // temporary hold for the reads corresponding to the last component
for(it = componentReadMap.begin(); it != componentReadMap.end(); it++){
if(it->first != lastComponent){ // start new component, first store previous component info.
if(tmpReadIDs.size() > 0){
mapComponents.push_back(lastComponent);
componentReads.push_back(tmpReadIDs);
}
lastComponent = it->first;
tmpReadIDs.clear();
}
tmpReadIDs.push_back(it->second);
}
if(tmpReadIDs.size() > 0){ // get the last component info stored.
mapComponents.push_back(lastComponent);
componentReads.push_back(tmpReadIDs);
}
// write out to files in map order
#pragma omp parallel for schedule(static,1)
for(unsigned int i = 0; i < mapComponents.size(); i++){
// each thread accesses a unique component number, so no competition between threads in writing to component-based files.
int iworm_bundle_index = mapComponents[i];
int iworm_component_no = component_number_mapping[iworm_bundle_index];
/* changed to write to a single file - bhaas July 3, 2013
string name = compMgr.GetFileName(iworm_component_no, ".raw.fasta", false);
int outFile;
//ML: use systm calls open/write/close, buffer explicitly (faster than C++ buffered I/O)
if(seen_component.find(iworm_component_no) == seen_component.end()) {
// first write
if (VERBOSE)
cerr << "First write to component: " << iworm_component_no << endl;
outFile = open(name.c_str(), O_WRONLY | O_CREAT | O_TRUNC , 0666);
#pragma omp critical
seen_component[iworm_component_no] = true;
}
else {
// subsequent write as append
if (VERBOSE)
cerr << "Second write as append to component: " << iworm_component_no << endl;
outFile = open(name.c_str(), O_WRONLY | O_APPEND );
}
if(outFile<0) cerr << "cannot open " << name << endl;
*/
//ML: put everthing to write into buf
string tmpName;
char component_id_string[8];
sprintf(component_id_string, "%d", iworm_component_no);
string buf;
for(unsigned int j = 0; j < componentReads[i].size(); j++){
buf += component_id_string;
buf += "\t";
int read_index = componentReads[i][j];
DNAStringStreamFast::formatReadNameString(readNameVector[read_index], tmpName);
buf += tmpName;
char foo[8];
sprintf(foo, "%d%%\t", read_pct_mapping_info[read_index]);
buf += "\t";
buf += foo;
buf += readSeqVector[read_index];
buf += "\n";
}
// write mapping to file, avoid thread collisions
#pragma omp critical
if(write(outFile,buf.c_str(),buf.size())<0)
{
cerr << "error writing file " << readsToTranscriptsOutputFile << ": " << strerror(errno) << endl;
exit(0);
}
//close(outFile);
}
if(mapComponents.size()>0)
cerr << "[" << mapComponents.size() << "] components written." << endl;
//cerr << "done\n";
//clear out read component mappings
} while (reads_read > 0);
close(outFile);
cerr << "Done" << endl;
FILE * pReadCount = fopen(readCountFile.c_str(), "w");
fprintf(pReadCount, "%lu\n", readCount);
fclose(pReadCount);
return 0;
}
|