]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
modified mpi code to save ram by writing out every 10 seqs.
[mothur.git] / filterseqscommand.cpp
1 /*
2  *  filterseqscommand.cpp
3  *  Mothur
4  *
5  *  Created by Thomas Ryabin on 5/4/09.
6  *  Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
7  *
8  */
9
10 #include "filterseqscommand.h"
11 #include "sequence.hpp"
12
13 /**************************************************************************************/
14
15 FilterSeqsCommand::FilterSeqsCommand(string option)  {
16         try {
17                 abort = false;
18                 filterFileName = "";
19                 
20                 //allow user to run help
21                 if(option == "help") { help(); abort = true; }
22                 
23                 else {
24                         //valid paramters for this command
25                         string Array[] =  {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
26                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
27                         
28                         OptionParser parser(option);
29                         map<string,string> parameters = parser.getParameters();
30                         
31                         ValidParameters validParameter;
32                         map<string,string>::iterator it;
33                         
34                         //check to make sure all parameters are valid for command
35                         for (it = parameters.begin(); it != parameters.end(); it++) { 
36                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
37                         }
38                         
39                         //if the user changes the input directory command factory will send this info to us in the output parameter 
40                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
41                         if (inputDir == "not found"){   inputDir = "";          }
42                         else {
43                                 string path;
44                                 it = parameters.find("fasta");
45                                 //user has given a template file
46                                 if(it != parameters.end()){ 
47                                         path = hasPath(it->second);
48                                         //if the user has not given a path then, add inputdir. else leave path alone.
49                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
50                                 }
51                                 
52                                 it = parameters.find("hard");
53                                 //user has given a template file
54                                 if(it != parameters.end()){ 
55                                         path = hasPath(it->second);
56                                         //if the user has not given a path then, add inputdir. else leave path alone.
57                                         if (path == "") {       parameters["hard"] = inputDir + it->second;             }
58                                 }
59                         }
60                         
61                         //check for required parameters
62                         fasta = validParameter.validFile(parameters, "fasta", false);
63                         if (fasta == "not found") { m->mothurOut("fasta is a required parameter for the filter.seqs command."); m->mothurOutEndLine(); abort = true;  }
64                         else { 
65                                 splitAtDash(fasta, fastafileNames);
66                                 
67                                 //go through files and make sure they are good, if not, then disregard them
68                                 for (int i = 0; i < fastafileNames.size(); i++) {
69                                         if (inputDir != "") {
70                                                 string path = hasPath(fastafileNames[i]);
71                                                 //if the user has not given a path then, add inputdir. else leave path alone.
72                                                 if (path == "") {       fastafileNames[i] = inputDir + fastafileNames[i];               }
73                                         }
74
75                                         int ableToOpen;
76                                         ifstream in;
77                                         ableToOpen = openInputFile(fastafileNames[i], in);
78                                         if (ableToOpen == 1) { 
79                                                 m->mothurOut(fastafileNames[i] + " will be disregarded."); m->mothurOutEndLine(); 
80                                                 //erase from file list
81                                                 fastafileNames.erase(fastafileNames.begin()+i);
82                                                 i--;
83                                         }else{  
84                                                 string simpleName = getSimpleName(fastafileNames[i]);
85                                                 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
86                                         }
87                                         in.close();
88                                 }
89                                 
90                                 //make sure there is at least one valid file left
91                                 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
92                         }
93                         
94                         if (!abort) {
95                                 //if the user changes the output directory command factory will send this info to us in the output parameter 
96                                 outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
97                                         outputDir = ""; 
98                                         outputDir += hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it  
99                                 }
100                         }
101                         //check for optional parameter and set defaults
102                         // ...at some point should added some additional type checking...
103                         
104                         string temp;
105                         temp = validParameter.validFile(parameters, "trump", false);                    if (temp == "not found") { temp = "*"; }
106                         trump = temp[0];
107                         
108                         temp = validParameter.validFile(parameters, "soft", false);                             if (temp == "not found") { soft = 0; }
109                         else {  soft = (float)atoi(temp.c_str()) / 100.0;  }
110                         
111                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
112                         convert(temp, processors); 
113                         
114                         hard = validParameter.validFile(parameters, "hard", true);                              if (hard == "not found") { hard = ""; }
115                         else if (hard == "not open") { abort = true; }  
116                         
117                         vertical = validParameter.validFile(parameters, "vertical", false);             if (vertical == "not found") { vertical = "T"; }
118                         
119                         numSeqs = 0;
120                         
121                 }
122                 
123         }
124         catch(exception& e) {
125                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
126                 exit(1);
127         }
128 }
129
130 //**********************************************************************************************************************
131
132 void FilterSeqsCommand::help(){
133         try {
134                                 
135                 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
136                 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard and vertical. \n");
137                 m->mothurOut("The fasta parameter is required. You may enter several fasta files to build the filter from and filter, by separating their names with -'s.\n");
138                 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
139                 m->mothurOut("The trump parameter .... The default is ...\n");
140                 m->mothurOut("The soft parameter .... The default is ....\n");
141                 m->mothurOut("The hard parameter .... The default is ....\n");
142                 m->mothurOut("The vertical parameter .... The default is T.\n");
143                 m->mothurOut("The filter.seqs command should be in the following format: \n");
144                 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump, soft=yourSoft, hard=yourHard, vertical=yourVertical) \n");
145                 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=..., soft=..., hard=..., vertical=T).\n");
146                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
147                 
148         }
149         catch(exception& e) {
150                 m->errorOut(e, "FilterSeqsCommand", "help");
151                 exit(1);
152         }
153 }
154
155 /**************************************************************************************/
156
157 int FilterSeqsCommand::execute() {      
158         try {
159         
160                 if (abort == true) { return 0; }
161                 vector<string> outputNames;
162                 
163                 ifstream inFASTA;
164                 openInputFile(fastafileNames[0], inFASTA);
165                 
166                 Sequence testSeq(inFASTA);
167                 alignmentLength = testSeq.getAlignLength();
168                 inFASTA.close();
169                 
170                 ////////////create filter/////////////////
171                 
172                 filter = createFilter();
173                 
174                 ofstream outFilter;
175                 
176                 string filterFile = outputDir + filterFileName + ".filter";
177                 openOutputFile(filterFile, outFilter);
178                 outFilter << filter << endl;
179                 outFilter.close();
180                 outputNames.push_back(filterFile);
181                 
182                 
183                 ////////////run filter/////////////////
184                 
185                 numSeqs = 0;
186                 for (int i = 0; i < fastafileNames.size(); i++) {
187                         ifstream in;
188                         openInputFile(fastafileNames[i], in);
189                         string filteredFasta = outputDir + getRootName(getSimpleName(fastafileNames[i])) + "filter.fasta";
190                         ofstream outFASTA;
191                         openOutputFile(filteredFasta, outFASTA);
192                         outputNames.push_back(filteredFasta);
193                         
194                         
195                         while(!in.eof()){
196                                 if (m->control_pressed) { in.close(); outFASTA.close(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }  return 0; }
197                                 
198                                 Sequence seq(in);
199                                 if (seq.getName() != "") {
200                                         string align = seq.getAligned();
201                                         string filterSeq = "";
202                                         
203                                         for(int j=0;j<alignmentLength;j++){
204                                                 if(filter[j] == '1'){
205                                                         filterSeq += align[j];
206                                                 }
207                                         }
208                                         
209                                         outFASTA << '>' << seq.getName() << endl << filterSeq << endl;
210                                         numSeqs++;
211                                 }
212                                 gobble(in);
213                         }
214                         outFASTA.close();
215                         in.close();
216                 }
217                 
218                 int filteredLength = 0;
219                 for(int i=0;i<alignmentLength;i++){
220                         if(filter[i] == '1'){   filteredLength++;       }
221                 }
222                 
223                 if (m->control_pressed) {  for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }  return 0; }
224
225                 
226                 m->mothurOutEndLine();
227                 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
228                 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
229                 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
230                 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
231                 
232                 
233                 m->mothurOutEndLine();
234                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
235                 for(int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();      }
236                 m->mothurOutEndLine();
237                 
238                 return 0;
239                 
240         }
241         catch(exception& e) {
242                 m->errorOut(e, "FilterSeqsCommand", "execute");
243                 exit(1);
244         }
245 }
246 /**************************************************************************************/
247 string FilterSeqsCommand::createFilter() {      
248         try {
249                 string filterString = "";
250                 
251                 Filters F;
252                 
253                 if (soft != 0)                  {  F.setSoft(soft);             }
254                 if (trump != '*')               {  F.setTrump(trump);   }
255                 
256                 F.setLength(alignmentLength);
257                 
258                 if(soft != 0 || isTrue(vertical)){
259                         F.initialize();
260                 }
261                 
262                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
263                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
264                 
265                 numSeqs = 0;
266                 
267                 if(trump != '*' || isTrue(vertical) || soft != 0){
268                         for (int s = 0; s < fastafileNames.size(); s++) {
269                         
270                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
271                         
272 #ifdef USE_MPI  
273                                 int pid, rc, ierr; 
274                                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
275                                 int tag = 2001;
276                                 
277                                 MPI_Status status; 
278                                 MPI_File in; 
279                                 rc = MPI_Comm_size(MPI_COMM_WORLD, &processors);
280                                 rc = MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
281                                                 
282                                                         
283                                 char* tempFileName = new char(fastafileNames[s].length());
284                                 tempFileName = &(fastafileNames[s][0]);
285                                 
286                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &in);  //comm, filename, mode, info, filepointer
287                                 
288                                 if (pid == 0) { //you are the root process
289                                                 setLines(fastafileNames[s]);
290                                                 
291                                                 for (int j = 0; j < lines.size(); j++) { //each process
292                                                         if (j != 0) { //don't send to yourself
293                                                                 MPI_Send(&lines[j]->start, 1, MPI_INT, j, tag, MPI_COMM_WORLD); //start position in file
294                                                                 MPI_Send(&lines[j]->numSeqs, 1, MPI_INT, j, tag, MPI_COMM_WORLD); //how many sequences we are sending
295                                                                 MPI_Send(&bufferSizes[j], 1, MPI_INT, j, tag, MPI_COMM_WORLD); //how bytes for the read
296                                                         }
297                                                 }
298                                 
299                                                 char buf[bufferSizes[0]];
300                                                 MPI_File_read_at(in, 0, buf, bufferSizes[0], MPI_CHAR, &status);
301                         
302                                                 MPICreateFilter(F, buf);
303                                                                                                 
304                                                 vector<int> temp; temp.resize(alignmentLength+1);
305                                                 
306                                                 //get the frequencies from the child processes
307                                                 for(int i = 0; i < ((processors-1)*5); i++) { 
308                                 cout << "i = " << i << endl;
309                                 //vector<int> trial; trial.resize(10);
310                                 //cout << "trials size = " << trial.size() << endl;
311                                 //int ierr = MPI_Recv(&trial[0], 10, MPI_INT, MPI_ANY_SOURCE, tag, MPI_COMM_WORLD, &status);
312                                                         int ierr = MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, MPI_ANY_SOURCE, tag, MPI_COMM_WORLD, &status); 
313                                         cout << "recieved something" << endl;           
314                         //for (int g = 0; g < trial.size(); g++) {  cout << trial[g] << '\t';  } cout << endl;
315                                                         int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
316                                         cout << "reciveve tag = " << receiveTag << endl;
317                                         for (int k = 0; k < alignmentLength; k++) {             cout << k << '\t' << temp[k] << endl;   }
318                                         cout << "done " << endl << endl;
319                                                         
320                                                         int sender = status.MPI_SOURCE; 
321                                                         
322                                                         if (receiveTag == Atag) { //you are recieveing the A frequencies
323                                                                 for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
324                                                         }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
325                                                                 for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
326                                                         }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
327                                                                 for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
328                                                         }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
329                                                                 for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
330                                                         }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
331                                                                 for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
332                                                         }
333                                                         
334                                                         m->mothurOut("receive tag = " + toString(receiveTag) + " " + toString(sender) + " is complete."); m->mothurOutEndLine();
335                                                 } 
336
337                                                 
338                                 }else { //i am the child process
339                                         int startPos, numLines, bufferSize;
340                                         ierr = MPI_Recv(&startPos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
341                                         ierr = MPI_Recv(&numLines, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
342                                         ierr = MPI_Recv(&bufferSize, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
343                                                                         
344                                         //send freqs
345                                         char buf2[bufferSize];
346                                         MPI_File_read_at( in, startPos, buf2, bufferSize, MPI_CHAR, &status);
347                                                 
348                                         MPICreateFilter(F, buf2);
349                                         
350                                         //send my fequency counts
351                                         F.a.push_back(Atag);
352                                         int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, tag, MPI_COMM_WORLD);
353                                         F.t.push_back(Ttag);
354                                         ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, tag, MPI_COMM_WORLD);
355                                         F.c.push_back(Ctag);
356                                         ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, tag, MPI_COMM_WORLD);
357                                         F.g.push_back(Gtag);
358                                         ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, tag, MPI_COMM_WORLD);
359                                         F.gap.push_back(Gaptag);
360                                         ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, tag, MPI_COMM_WORLD);
361                                 }
362                                 
363                                 MPI_Barrier(MPI_COMM_WORLD);
364                                 
365 #else
366                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
367                                 if(processors == 1){
368                                         ifstream inFASTA;
369                                         openInputFile(fastafileNames[s], inFASTA);
370                                         int numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
371                                         inFASTA.close();
372                                         
373                                         numSeqs += numFastaSeqs;
374                                         
375                                         lines.push_back(new linePair(0, numFastaSeqs));
376                                         
377                                         driverCreateFilter(F, fastafileNames[s], lines[0]);
378                                 }else{
379                                         
380                                         setLines(fastafileNames[s]);                                    
381                                         createProcessesCreateFilter(F, fastafileNames[s]); 
382                                 }
383                 #else
384                                 ifstream inFASTA;
385                                 openInputFile(fastafileNames[s], inFASTA);
386                                 int numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
387                                 inFASTA.close();
388                                 
389                                 numSeqs += numFastaSeqs;
390                                 
391                                 lines.push_back(new linePair(0, numFastaSeqs));
392                                 
393                                 driverCreateFilter(F, lines[0], fastafileNames[s]);
394                 #endif
395 #endif
396                         
397                         }
398                 }
399
400
401         cout << "made it here, numSeqs = " << numSeqs << endl;  
402                 F.setNumSeqs(numSeqs);
403                                 
404                 if(isTrue(vertical) == 1)       {       F.doVertical(); }
405                 if(soft != 0)                           {       F.doSoft();             }
406 //cout << "Filter String = " << F.getFilter() << endl;                  
407                 filterString = F.getFilter();
408
409                 return filterString;
410         }
411         catch(exception& e) {
412                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
413                 exit(1);
414         }
415 }
416 /**************************************************************************************/
417 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* line) {        
418         try {
419                 
420                 ifstream in;
421                 openInputFile(filename, in);
422                                 
423                 in.seekg(line->start);
424                 
425                 for(int i=0;i<line->numSeqs;i++){
426                                 
427                         if (m->control_pressed) { in.close(); return 1; }
428                                         
429                         Sequence seq(in);
430                         if (seq.getName() != "") {
431                                         if(trump != '*'){       F.doTrump(seq); }
432                                         if(isTrue(vertical) || soft != 0){      F.getFreqs(seq);        }
433                                         cout.flush();
434                         }
435                         
436                         //report progress
437                         if((i+1) % 100 == 0){   m->mothurOut(toString(i+1)); m->mothurOutEndLine();             }
438                 }
439                 
440                 //report progress
441                 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine();           }
442                 
443                 in.close();
444                 
445                 return 0;
446         }
447         catch(exception& e) {
448                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
449                 exit(1);
450         }
451 }
452 /**************************************************************************************/
453 int FilterSeqsCommand::MPICreateFilter(Filters& F, string input) {      
454         try {
455                 
456                 vector<string> seqStrings;
457                 parseBuffer(input, seqStrings);
458                 
459                 for(int i=0;i<seqStrings.size();i++){
460                                 
461                         if (m->control_pressed) { return 1; }
462                         
463                         Sequence seq("", seqStrings[i]);
464                         
465                         if(trump != '*'){       F.doTrump(seq); }
466                         if(isTrue(vertical) || soft != 0){      F.getFreqs(seq);        }
467                         cout.flush();
468                                                 
469                         //report progress
470                         if((i+1) % 100 == 0){   m->mothurOut(toString(i+1)); m->mothurOutEndLine();             }
471                 }
472                 
473                 //report progress
474                 if((seqStrings.size()) % 100 != 0){     m->mothurOut(toString(seqStrings.size())); m->mothurOutEndLine();               }
475                 
476                 return 0;
477         }
478         catch(exception& e) {
479                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
480                 exit(1);
481         }
482 }
483
484 /**************************************************************************************************/
485
486 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
487         try {
488 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
489                 int process = 0;
490                 int exitCommand = 1;
491                 vector<int> processIDS;
492                 
493                 //loop through and create all the processes you want
494                 while (process != processors) {
495                         int pid = vfork();
496                         
497                         if (pid > 0) {
498                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
499                                 process++;
500                         }else if (pid == 0){
501                                 driverCreateFilter(F, filename, lines[process]);
502                                 exit(0);
503                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
504                 }
505                 
506                 //force parent to wait until all the processes are done
507                 for (int i=0;i<processors;i++) { 
508                         int temp = processIDS[i];
509                         wait(&temp);
510                 }
511                 
512                 return exitCommand;
513 #endif          
514         }
515         catch(exception& e) {
516                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
517                 exit(1);
518         }
519 }
520 /**************************************************************************************************/
521
522 int FilterSeqsCommand::setLines(string filename) {
523         try {
524                 vector<int> positions;
525                 map<int, int> buf;
526                 bufferSizes.clear();
527                 
528                 int pid;
529                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
530         
531                 ifstream inFASTA;
532                 openInputFile(filename, inFASTA);
533                         
534                 string input;
535                 int numbuf = 0;
536                 while(!inFASTA.eof()){  
537                         input = getline(inFASTA);
538
539                         if (input.length() != 0) {
540                                 numbuf += input.length();
541                                 if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  buf[(pos - input.length() - 1)] = numbuf; }
542                         }
543                 }
544
545                 inFASTA.close();
546                 
547                 int numFastaSeqs = positions.size();
548                 
549                 numSeqs += numFastaSeqs;
550                 
551                 int numSeqsPerProcessor = numFastaSeqs / processors;
552                 
553                 for (int i = 0; i < processors; i++) {
554
555                         long int startPos = positions[ i * numSeqsPerProcessor ];
556                         if(i == processors - 1){
557                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
558                                 bufferSizes.push_back(numbuf-buf[startPos]);
559                         }else{  
560                                 int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
561                                 bufferSizes.push_back(buf[myEnd]-buf[startPos]);
562                         }
563                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
564                 }
565                 
566                 return 0;
567         }
568         catch(exception& e) {
569                 m->errorOut(e, "FilterSeqsCommand", "setLines");
570                 exit(1);
571         }
572 }
573 /**************************************************************************************************/
574 int FilterSeqsCommand::parseBuffer(string file, vector<string>& seqs) {
575         try {
576                 
577                 istringstream iss (file,istringstream::in);
578                 string name, seqstring;
579         
580                 while (iss) {
581                         
582                         if (m->control_pressed) { return 0; }
583                                 
584                         Sequence seq(iss); gobble(iss);
585         
586                         if (seq.getName() != "") {
587                                 seqs.push_back(seq.getAligned());       
588                         }
589                 }
590                 
591                 return 0;
592         }
593         catch(exception& e) {
594                 m->errorOut(e, "FilterSeqsCommand", "parseBuffer");
595                 exit(1);
596         }
597 }
598 /**************************************************************************************/