]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
added [ERROR] flag if command aborts
[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 vector<string> FilterSeqsCommand::getValidParameters(){ 
15         try {
16                 string Array[] =  {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
17                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
18                 return myArray;
19         }
20         catch(exception& e) {
21                 m->errorOut(e, "FilterSeqsCommand", "getValidParameters");
22                 exit(1);
23         }
24 }
25 //**********************************************************************************************************************
26 FilterSeqsCommand::FilterSeqsCommand(){ 
27         try {
28                 abort = true; calledHelp = true; 
29                 vector<string> tempOutNames;
30                 outputTypes["fasta"] = tempOutNames;
31                 outputTypes["filter"] = tempOutNames;
32         }
33         catch(exception& e) {
34                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
35                 exit(1);
36         }
37 }
38 //**********************************************************************************************************************
39 vector<string> FilterSeqsCommand::getRequiredParameters(){      
40         try {
41                 string Array[] =  {"fasta"};
42                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
43                 return myArray;
44         }
45         catch(exception& e) {
46                 m->errorOut(e, "FilterSeqsCommand", "getRequiredParameters");
47                 exit(1);
48         }
49 }
50 //**********************************************************************************************************************
51 vector<string> FilterSeqsCommand::getRequiredFiles(){   
52         try {
53                 vector<string> myArray;
54                 return myArray;
55         }
56         catch(exception& e) {
57                 m->errorOut(e, "FilterSeqsCommand", "getRequiredFiles");
58                 exit(1);
59         }
60 }
61 /**************************************************************************************/
62 FilterSeqsCommand::FilterSeqsCommand(string option)  {
63         try {
64                 abort = false; calledHelp = false;   
65                 filterFileName = "";
66                 
67                 //allow user to run help
68                 if(option == "help") { help(); abort = true; calledHelp = true; }
69                 
70                 else {
71                         //valid paramters for this command
72                         string Array[] =  {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
73                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
74                         
75                         OptionParser parser(option);
76                         map<string,string> parameters = parser.getParameters();
77                         
78                         ValidParameters validParameter("filter.seqs");
79                         map<string,string>::iterator it;
80                         
81                         //check to make sure all parameters are valid for command
82                         for (it = parameters.begin(); it != parameters.end(); it++) { 
83                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
84                         }
85                         
86                         //initialize outputTypes
87                         vector<string> tempOutNames;
88                         outputTypes["fasta"] = tempOutNames;
89                         outputTypes["filter"] = tempOutNames;
90                 
91                         //if the user changes the input directory command factory will send this info to us in the output parameter 
92                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
93                         if (inputDir == "not found"){   inputDir = "";          }
94                         else {
95                                 string path;
96                                 it = parameters.find("fasta");
97                                 //user has given a template file
98                                 if(it != parameters.end()){ 
99                                         path = m->hasPath(it->second);
100                                         //if the user has not given a path then, add inputdir. else leave path alone.
101                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
102                                 }
103                                 
104                                 it = parameters.find("hard");
105                                 //user has given a template file
106                                 if(it != parameters.end()){ 
107                                         path = m->hasPath(it->second);
108                                         //if the user has not given a path then, add inputdir. else leave path alone.
109                                         if (path == "") {       parameters["hard"] = inputDir + it->second;             }
110                                 }
111                         }
112                         
113                         //check for required parameters
114                         fasta = validParameter.validFile(parameters, "fasta", false);
115                         if (fasta == "not found") { m->mothurOut("fasta is a required parameter for the filter.seqs command."); m->mothurOutEndLine(); abort = true;  }
116                         else { 
117                                 m->splitAtDash(fasta, fastafileNames);
118                                 
119                                 //go through files and make sure they are good, if not, then disregard them
120                                 for (int i = 0; i < fastafileNames.size(); i++) {
121                                         if (inputDir != "") {
122                                                 string path = m->hasPath(fastafileNames[i]);
123                                                 //if the user has not given a path then, add inputdir. else leave path alone.
124                                                 if (path == "") {       fastafileNames[i] = inputDir + fastafileNames[i];               }
125                                         }
126
127                                         ifstream in;
128                                         int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
129                                 
130                                         //if you can't open it, try default location
131                                         if (ableToOpen == 1) {
132                                                 if (m->getDefaultPath() != "") { //default path is set
133                                                         string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
134                                                         m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
135                                                         ifstream in2;
136                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
137                                                         in2.close();
138                                                         fastafileNames[i] = tryPath;
139                                                 }
140                                         }
141                                         
142                                         //if you can't open it, try default location
143                                         if (ableToOpen == 1) {
144                                                 if (m->getOutputDir() != "") { //default path is set
145                                                         string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
146                                                         m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
147                                                         ifstream in2;
148                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
149                                                         in2.close();
150                                                         fastafileNames[i] = tryPath;
151                                                 }
152                                         }
153                                         
154                                         in.close();
155                                         
156                                         if (ableToOpen == 1) { 
157                                                 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
158                                                 //erase from file list
159                                                 fastafileNames.erase(fastafileNames.begin()+i);
160                                                 i--;
161                                         }else{  
162                                                 string simpleName = m->getSimpleName(fastafileNames[i]);
163                                                 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
164                                         }
165                                         in.close();
166                                 }
167                                 
168                                 //make sure there is at least one valid file left
169                                 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
170                         }
171                         
172                         if (!abort) {
173                                 //if the user changes the output directory command factory will send this info to us in the output parameter 
174                                 outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
175                                         outputDir = ""; 
176                                         outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it       
177                                 }
178                         }
179                         //check for optional parameter and set defaults
180                         // ...at some point should added some additional type checking...
181                         
182                         string temp;
183                         hard = validParameter.validFile(parameters, "hard", true);                              if (hard == "not found") { hard = ""; }
184                         else if (hard == "not open") { hard = ""; abort = true; }       
185
186                         temp = validParameter.validFile(parameters, "trump", false);                    if (temp == "not found") { temp = "*"; }
187                         trump = temp[0];
188                         
189                         temp = validParameter.validFile(parameters, "soft", false);                             if (temp == "not found") { soft = 0; }
190                         else {  soft = (float)atoi(temp.c_str()) / 100.0;  }
191                         
192                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
193                         convert(temp, processors); 
194                         
195                         vertical = validParameter.validFile(parameters, "vertical", false);             
196                         if (vertical == "not found") { 
197                                 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
198                                 else { vertical = "F";  }
199                         }
200                         
201                         numSeqs = 0;
202                 }
203                 
204         }
205         catch(exception& e) {
206                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
207                 exit(1);
208         }
209 }
210
211 //**********************************************************************************************************************
212
213 void FilterSeqsCommand::help(){
214         try {
215                                 
216                 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
217                 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n");
218                 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");
219                 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
220                 m->mothurOut("The trump parameter .... The default is ...\n");
221                 m->mothurOut("The soft parameter .... The default is ....\n");
222                 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
223                 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
224                 m->mothurOut("The processors parameter allows you to specify the number of processors to use. The default is 1.\n");
225                 m->mothurOut("The filter.seqs command should be in the following format: \n");
226                 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
227                 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
228                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
229                 
230         }
231         catch(exception& e) {
232                 m->errorOut(e, "FilterSeqsCommand", "help");
233                 exit(1);
234         }
235 }
236
237 /**************************************************************************************/
238
239 int FilterSeqsCommand::execute() {      
240         try {
241         
242                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
243                 
244                 ifstream inFASTA;
245                 m->openInputFile(fastafileNames[0], inFASTA);
246                 
247                 Sequence testSeq(inFASTA);
248                 alignmentLength = testSeq.getAlignLength();
249                 inFASTA.close();
250                 
251                 ////////////create filter/////////////////
252                 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
253                 
254                 filter = createFilter();
255                 
256                 m->mothurOutEndLine();  m->mothurOutEndLine();
257                 
258                 if (m->control_pressed) { outputTypes.clear(); return 0; }
259                 
260                 #ifdef USE_MPI
261                         int pid;
262                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
263                                         
264                         if (pid == 0) { //only one process should output the filter
265                 #endif
266                 
267                 ofstream outFilter;
268                 
269                 //prevent giantic file name
270                 string filterFile;
271                 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
272                 else {  filterFile = outputDir + filterFileName + ".filter";  }
273                 
274                 m->openOutputFile(filterFile, outFilter);
275                 outFilter << filter << endl;
276                 outFilter.close();
277                 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
278                 
279                 #ifdef USE_MPI
280                         }
281                 #endif
282                 
283                 ////////////run filter/////////////////
284                 
285                 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
286                 
287                 filterSequences();
288                 
289                 m->mothurOutEndLine();  m->mothurOutEndLine();
290                                         
291                 int filteredLength = 0;
292                 for(int i=0;i<alignmentLength;i++){
293                         if(filter[i] == '1'){   filteredLength++;       }
294                 }
295                 
296                 if (m->control_pressed) {  outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }  return 0; }
297
298                 
299                 m->mothurOutEndLine();
300                 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
301                 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
302                 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
303                 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
304                 
305                 
306                 m->mothurOutEndLine();
307                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
308                 for(int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();      }
309                 m->mothurOutEndLine();
310                 
311                 return 0;
312                 
313         }
314         catch(exception& e) {
315                 m->errorOut(e, "FilterSeqsCommand", "execute");
316                 exit(1);
317         }
318 }
319 /**************************************************************************************/
320 int FilterSeqsCommand::filterSequences() {      
321         try {
322                 
323                 numSeqs = 0;
324                 
325                 for (int s = 0; s < fastafileNames.size(); s++) {
326                         
327                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
328                                 
329                                 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
330 #ifdef USE_MPI  
331                                 int pid, numSeqsPerProcessor, num; 
332                                 int tag = 2001;
333                                 vector<unsigned long int>MPIPos;
334                                                 
335                                 MPI_Status status; 
336                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
337                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
338                                 
339                                 MPI_File outMPI;
340                                 MPI_File inMPI;
341                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
342                                 int inMode=MPI_MODE_RDONLY; 
343                                 
344                                 char outFilename[1024];
345                                 strcpy(outFilename, filteredFasta.c_str());
346                         
347                                 char inFileName[1024];
348                                 strcpy(inFileName, fastafileNames[s].c_str());
349                                 
350                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
351                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
352
353                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
354
355                                 if (pid == 0) { //you are the root process 
356                                         
357                                         MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
358                                         numSeqs += num;
359                                         
360                                         //send file positions to all processes
361                                         for(int i = 1; i < processors; i++) { 
362                                                 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
363                                                 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
364                                         }
365                                         
366                                         //figure out how many sequences you have to do
367                                         numSeqsPerProcessor = num / processors;
368                                         int startIndex =  pid * numSeqsPerProcessor;
369                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
370                                         
371                                 
372                                         //do your part
373                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
374                                         
375                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
376                                         
377                                         //wait on chidren
378                                         for(int i = 1; i < processors; i++) { 
379                                                 char buf[5];
380                                                 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
381                                         }
382                                         
383                                 }else { //you are a child process
384                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
385                                         MPIPos.resize(num+1);
386                                         numSeqs += num;
387                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
388                                         
389                                         //figure out how many sequences you have to align
390                                         numSeqsPerProcessor = num / processors;
391                                         int startIndex =  pid * numSeqsPerProcessor;
392                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
393                                         
394                                         
395                                         //align your part
396                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);           
397                                         
398                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
399                                         
400                                         char buf[5];
401                                         strcpy(buf, "done"); 
402                                         
403                                         //tell parent you are done.
404                                         MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
405                                 }
406                                 
407                                 MPI_File_close(&outMPI);
408                                 MPI_File_close(&inMPI);
409                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
410                                 
411 #else
412                         vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
413                                 
414                         for (int i = 0; i < (positions.size()-1); i++) {
415                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
416                         }       
417                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
418                                 if(processors == 1){
419                                         int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
420                                         numSeqs += numFastaSeqs;
421                                 }else{
422                                         int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]); 
423                                         numSeqs += numFastaSeqs;
424                                 
425                                         rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
426                                 
427                                         //append fasta files
428                                         for(int i=1;i<processors;i++){
429                                                 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
430                                                 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
431                                         }
432                                 }
433                                 
434                                 if (m->control_pressed) {  return 1; }
435                 #else
436                                 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
437                                 numSeqs += numFastaSeqs;
438
439                                 if (m->control_pressed) {  return 1; }
440                 #endif
441 #endif
442                         outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
443                 }
444
445                 return 0;
446         }
447         catch(exception& e) {
448                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
449                 exit(1);
450         }
451 }
452 #ifdef USE_MPI
453 /**************************************************************************************/
454 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) { 
455         try {
456                 string outputString = "";
457                 int count = 0;
458                 MPI_Status status; 
459                 
460                 for(int i=0;i<num;i++){
461                 
462                         if (m->control_pressed) { return 0; }
463                 
464                         //read next sequence
465                         int length = MPIPos[start+i+1] - MPIPos[start+i];
466                         char* buf4 = new char[length];
467                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
468                         
469                         string tempBuf = buf4;
470                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
471                         istringstream iss (tempBuf,istringstream::in);
472                         delete buf4;
473         
474                         Sequence seq(iss);  m->gobble(iss);
475                         
476                         if (seq.getName() != "") {
477                                 string align = seq.getAligned();
478                                 string filterSeq = "";
479                                         
480                                 for(int j=0;j<alignmentLength;j++){
481                                         if(filter[j] == '1'){
482                                                 filterSeq += align[j];
483                                         }
484                                 }
485                                 
486                                 count++;
487                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
488                                 
489                                 if(count % 10 == 0){ //output to file 
490                                         //send results to parent
491                                         int length = outputString.length();
492                                         char* buf = new char[length];
493                                         memcpy(buf, outputString.c_str(), length);
494                                 
495                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
496                                         outputString = "";
497                                         delete buf;
498                                 }
499
500                         }
501                         
502                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
503                 }
504                 
505                 if(outputString != ""){ //output to file 
506                         //send results to parent
507                         int length = outputString.length();
508                         char* buf = new char[length];
509                         memcpy(buf, outputString.c_str(), length);
510                         
511                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
512                         outputString = "";
513                         delete buf;
514                 }
515                 
516                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
517                         
518                 return 0;
519         }
520         catch(exception& e) {
521                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
522                 exit(1);
523         }
524 }
525 #endif
526 /**************************************************************************************/
527 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {      
528         try {
529                 ofstream out;
530                 m->openOutputFile(outputFilename, out);
531                 
532                 ifstream in;
533                 m->openInputFile(inputFilename, in);
534                                 
535                 in.seekg(filePos->start);
536
537                 bool done = false;
538                 int count = 0;
539         
540                 while (!done) {
541                                 
542                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
543                                 
544                                 Sequence seq(in); m->gobble(in);
545                                 if (seq.getName() != "") {
546                                         string align = seq.getAligned();
547                                         string filterSeq = "";
548                                         
549                                         for(int j=0;j<alignmentLength;j++){
550                                                 if(filter[j] == '1'){
551                                                         filterSeq += align[j];
552                                                 }
553                                         }
554                                         
555                                         out << '>' << seq.getName() << endl << filterSeq << endl;
556                                 count++;
557                         }
558                         
559                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
560                                 unsigned long int pos = in.tellg();
561                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
562                         #else
563                                 if (in.eof()) { break; }
564                         #endif
565                         
566                         //report progress
567                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
568                 }
569                 //report progress
570                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
571                 
572                 
573                 out.close();
574                 in.close();
575                 
576                 return count;
577         }
578         catch(exception& e) {
579                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
580                 exit(1);
581         }
582 }
583 /**************************************************************************************************/
584
585 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
586         try {
587 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
588                 int process = 0;
589                 int num = 0;
590                 processIDS.clear();
591                 
592                 //loop through and create all the processes you want
593                 while (process != processors) {
594                         int pid = fork();
595                         
596                         if (pid > 0) {
597                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
598                                 process++;
599                         }else if (pid == 0){
600                                 string filteredFasta = filename + toString(getpid()) + ".temp";
601                                 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
602                                 
603                                 //pass numSeqs to parent
604                                 ofstream out;
605                                 string tempFile = filename +  toString(getpid()) + ".num.temp";
606                                 m->openOutputFile(tempFile, out);
607                                 out << num << endl;
608                                 out.close();
609                                 
610                                 exit(0);
611                         }else { 
612                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
613                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
614                                 exit(0);
615                         }
616                 }
617                 
618                 //force parent to wait until all the processes are done
619                 for (int i=0;i<processors;i++) { 
620                         int temp = processIDS[i];
621                         wait(&temp);
622                 }       
623                                         
624                 for (int i = 0; i < processIDS.size(); i++) {
625                         ifstream in;
626                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
627                         m->openInputFile(tempFile, in);
628                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
629                         in.close(); remove(tempFile.c_str());
630                 }
631
632                 
633                 return num;
634 #endif          
635         }
636         catch(exception& e) {
637                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
638                 exit(1);
639         }
640 }
641 /**************************************************************************************/
642 string FilterSeqsCommand::createFilter() {      
643         try {
644                 string filterString = "";                       
645                 Filters F;
646                 
647                 if (soft != 0)                  {  F.setSoft(soft);             }
648                 if (trump != '*')               {  F.setTrump(trump);   }
649                 
650                 F.setLength(alignmentLength);
651                 
652                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
653                         F.initialize();
654                 }
655                 
656                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
657                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
658                 
659                 numSeqs = 0;
660                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
661                         for (int s = 0; s < fastafileNames.size(); s++) {
662                         
663                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
664                         
665 #ifdef USE_MPI  
666                                 int pid, numSeqsPerProcessor, num; 
667                                 int tag = 2001;
668                                 vector<unsigned long int> MPIPos;
669                                 
670                                 MPI_Status status; 
671                                 MPI_File inMPI; 
672                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
673                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
674                                                         
675                                 //char* tempFileName = new char(fastafileNames[s].length());
676                                 //tempFileName = &(fastafileNames[s][0]);
677                                 
678                                 char tempFileName[1024];
679                                 strcpy(tempFileName, fastafileNames[s].c_str());
680                 
681                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
682                                 
683                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
684                                 
685                                 if (pid == 0) { //you are the root process
686                                                 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
687                                                 numSeqs += num;
688                                                 
689                                                 //send file positions to all processes
690                                                 for(int i = 1; i < processors; i++) { 
691                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
692                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
693                                                 }
694                                                                 
695                                                 //figure out how many sequences you have to do
696                                                 numSeqsPerProcessor = num / processors;
697                                                 int startIndex =  pid * numSeqsPerProcessor;
698                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
699                                                 
700                                 
701                                                 //do your part
702                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
703                                                 
704                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
705                                                                                                 
706                                 }else { //i am the child process
707                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
708                                         MPIPos.resize(num+1);
709                                         numSeqs += num;
710                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
711                                         
712                                         //figure out how many sequences you have to align
713                                         numSeqsPerProcessor = num / processors;
714                                         int startIndex =  pid * numSeqsPerProcessor;
715                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
716                                         
717                                         
718                                         //do your part
719                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
720                                         
721                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
722                                 }
723                                 
724                                 MPI_File_close(&inMPI);
725                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
726                                 
727 #else
728                 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
729                                 
730                 for (int i = 0; i < (positions.size()-1); i++) {
731                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
732                 }       
733                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
734                                 if(processors == 1){
735                                         int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
736                                         numSeqs += numFastaSeqs;
737                                 }else{
738                                         int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
739                                         numSeqs += numFastaSeqs;
740                                 }
741                                 
742                                 if (m->control_pressed) {  return filterString; }
743                 #else
744                                 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
745                                 numSeqs += numFastaSeqs;
746                                 if (m->control_pressed) {  return filterString; }
747                 #endif
748 #endif
749                         
750                         }
751                 }
752
753
754 #ifdef USE_MPI  
755                 int pid;
756                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
757                 MPI_Status status;
758                 
759                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
760                 
761                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
762                         
763                         if (pid == 0) { //only one process should output the filter
764                         
765                                 vector<int> temp; temp.resize(alignmentLength+1);
766                                                                 
767                                 //get the frequencies from the child processes
768                                 for(int i = 1; i < processors; i++) { 
769                                 
770                                         for (int j = 0; j < 5; j++) {
771                                         
772                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
773                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
774                                                 
775                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
776                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
777                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
778                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
779                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
780                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
781                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
782                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
783                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
784                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
785                                                 }
786                                         }
787                                 } 
788                         }else{
789                         
790                                 //send my fequency counts
791                                 F.a.push_back(Atag);
792                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
793                                 F.t.push_back(Ttag);
794                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
795                                 F.c.push_back(Ctag);
796                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
797                                 F.g.push_back(Gtag);
798                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
799                                 F.gap.push_back(Gaptag);
800                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
801                         }
802                         
803                 }
804                 
805                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
806                 
807                 if (pid == 0) { //only one process should output the filter
808 #endif
809
810                 F.setNumSeqs(numSeqs);
811                 if(m->isTrue(vertical) == 1)    {       F.doVertical(); }
812                 if(soft != 0)                           {       F.doSoft();             }
813                 filterString = F.getFilter();
814                 
815 #ifdef USE_MPI
816                 //send filter string to kids
817                 //for(int i = 1; i < processors; i++) { 
818                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
819                 //}
820                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
821         }else{
822                 //recieve filterString
823                 char* tempBuf = new char[alignmentLength];
824                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
825                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
826                 
827                 filterString = tempBuf;
828                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
829                 delete tempBuf; 
830         }
831         
832         MPI_Barrier(MPI_COMM_WORLD);
833 #endif
834                                 
835                 return filterString;
836         }
837         catch(exception& e) {
838                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
839                 exit(1);
840         }
841 }
842 /**************************************************************************************/
843 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {     
844         try {
845                 
846                 ifstream in;
847                 m->openInputFile(filename, in);
848                                 
849                 in.seekg(filePos->start);
850
851                 bool done = false;
852                 int count = 0;
853         
854                 while (!done) {
855                                 
856                         if (m->control_pressed) { in.close(); return 1; }
857                                         
858                         Sequence seq(in); m->gobble(in);
859                         if (seq.getName() != "") {
860                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
861                                         
862                                         if(trump != '*')                        {       F.doTrump(seq);         }
863                                         if(m->isTrue(vertical) || soft != 0)    {       F.getFreqs(seq);        }
864                                         cout.flush();
865                                         count++;
866                         }
867                         
868                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
869                                 unsigned long int pos = in.tellg();
870                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
871                         #else
872                                 if (in.eof()) { break; }
873                         #endif
874                         
875                         //report progress
876                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
877                 }
878                 //report progress
879                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
880                 in.close();
881                 
882                 return count;
883         }
884         catch(exception& e) {
885                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
886                 exit(1);
887         }
888 }
889 #ifdef USE_MPI
890 /**************************************************************************************/
891 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {    
892         try {
893                 
894                 MPI_Status status; 
895                 int pid;
896                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
897                 
898                 for(int i=0;i<num;i++){
899                         
900                         if (m->control_pressed) { return 0; }
901                         
902                         //read next sequence
903                         int length = MPIPos[start+i+1] - MPIPos[start+i];
904         
905                         char* buf4 = new char[length];
906                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
907                         
908                         string tempBuf = buf4;
909                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
910                         istringstream iss (tempBuf,istringstream::in);
911                         delete buf4;
912
913                         Sequence seq(iss);  
914
915                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
916                         
917                         if(trump != '*'){       F.doTrump(seq); }
918                         if(m->isTrue(vertical) || soft != 0){   F.getFreqs(seq);        }
919                         cout.flush();
920                                                 
921                         //report progress
922                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
923                 }
924                 
925                 //report progress
926                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
927                 
928                 return 0;
929         }
930         catch(exception& e) {
931                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
932                 exit(1);
933         }
934 }
935 #endif
936 /**************************************************************************************************/
937
938 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
939         try {
940 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
941                 int process = 1;
942                 int num = 0;
943                 processIDS.clear();
944                 
945                 //loop through and create all the processes you want
946                 while (process != processors) {
947                         int pid = fork();
948                         
949                         if (pid > 0) {
950                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
951                                 process++;
952                         }else if (pid == 0){
953                                 //reset child's filter counts to 0;
954                                 F.a.clear(); F.a.resize(alignmentLength, 0);
955                                 F.t.clear(); F.t.resize(alignmentLength, 0);
956                                 F.g.clear(); F.g.resize(alignmentLength, 0);
957                                 F.c.clear(); F.c.resize(alignmentLength, 0);
958                                 F.gap.clear(); F.gap.resize(alignmentLength, 0);
959                                 
960                                 num = driverCreateFilter(F, filename, lines[process]);
961                                 
962                                 //write out filter counts to file
963                                 filename += toString(getpid()) + "filterValues.temp";
964                                 ofstream out;
965                                 m->openOutputFile(filename, out);
966                                 
967                                 out << num << endl;
968                                 out << F.getFilter() << endl;
969                                 for (int k = 0; k < alignmentLength; k++) {             out << F.a[k] << '\t'; }  out << endl;
970                                 for (int k = 0; k < alignmentLength; k++) {             out << F.t[k] << '\t'; }  out << endl;
971                                 for (int k = 0; k < alignmentLength; k++) {             out << F.g[k] << '\t'; }  out << endl;
972                                 for (int k = 0; k < alignmentLength; k++) {             out << F.c[k] << '\t'; }  out << endl;
973                                 for (int k = 0; k < alignmentLength; k++) {             out << F.gap[k] << '\t'; }  out << endl;
974
975                                 //cout << F.getFilter() << endl;
976                                 out.close();
977                                 
978                                 exit(0);
979                         }else { 
980                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
981                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
982                                 exit(0);
983                         }
984                 }
985                 
986                 //parent do your part
987                 num = driverCreateFilter(F, filename, lines[0]);
988                 
989                 //force parent to wait until all the processes are done
990                 for (int i=0;i<(processors-1);i++) { 
991                         int temp = processIDS[i];
992                         wait(&temp);
993                 }
994                 
995                 //parent reads in and combines Filter info
996                 for (int i = 0; i < processIDS.size(); i++) {
997                         string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
998                         ifstream in;
999                         m->openInputFile(tempFilename, in);
1000                         
1001                         int temp, tempNum;
1002                         string tempFilterString;
1003
1004                         in >> tempNum; m->gobble(in); num += tempNum;
1005
1006                         in >> tempFilterString;
1007                         F.mergeFilter(tempFilterString);
1008
1009                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.a[k] += temp; }           m->gobble(in);
1010                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.t[k] += temp; }           m->gobble(in);
1011                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.g[k] += temp; }           m->gobble(in);
1012                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.c[k] += temp; }           m->gobble(in);
1013                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.gap[k] += temp; } m->gobble(in);
1014                                 
1015                         in.close();
1016                         remove(tempFilename.c_str());
1017                 }
1018                 
1019                 return num;
1020 #endif          
1021         }
1022         catch(exception& e) {
1023                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1024                 exit(1);
1025         }
1026 }
1027 /**************************************************************************************/