2 * filterseqscommand.cpp
5 * Created by Thomas Ryabin on 5/4/09.
6 * Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
10 #include "filterseqscommand.h"
11 #include "sequence.hpp"
14 //**********************************************************************************************************************
15 vector<string> FilterSeqsCommand::setParameters(){
17 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
18 CommandParameter phard("hard", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(phard);
19 CommandParameter ptrump("trump", "String", "", "*", "", "", "",false,false); parameters.push_back(ptrump);
20 CommandParameter psoft("soft", "Number", "", "0", "", "", "",false,false); parameters.push_back(psoft);
21 CommandParameter pvertical("vertical", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(pvertical);
22 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
23 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
24 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
26 vector<string> myArray;
27 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
31 m->errorOut(e, "FilterSeqsCommand", "setParameters");
35 //**********************************************************************************************************************
36 string FilterSeqsCommand::getHelpString(){
38 string helpString = "";
39 helpString += "The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n";
40 helpString += "The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n";
41 helpString += "The fasta parameter is required, unless you have a valid current fasta file. You may enter several fasta files to build the filter from and filter, by separating their names with -'s.\n";
42 helpString += "For example: fasta=abrecovery.fasta-amazon.fasta \n";
43 helpString += "The trump option will remove a column if the trump character is found at that position in any sequence of the alignment. Default=*, meaning no trump. \n";
44 helpString += "A soft mask removes any column where the dominant base (i.e. A, T, G, C, or U) does not occur in at least a designated percentage of sequences. Default=0.\n";
45 helpString += "The hard parameter allows you to enter a file containing the filter you want to use.\n";
46 helpString += "The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n";
47 helpString += "The processors parameter allows you to specify the number of processors to use. The default is 1.\n";
48 helpString += "The filter.seqs command should be in the following format: \n";
49 helpString += "filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n";
50 helpString += "Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n";
51 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n";
55 m->errorOut(e, "FilterSeqsCommand", "getHelpString");
59 //**********************************************************************************************************************
60 FilterSeqsCommand::FilterSeqsCommand(){
62 abort = true; calledHelp = true;
64 vector<string> tempOutNames;
65 outputTypes["fasta"] = tempOutNames;
66 outputTypes["filter"] = tempOutNames;
69 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
73 /**************************************************************************************/
74 FilterSeqsCommand::FilterSeqsCommand(string option) {
76 abort = false; calledHelp = false;
79 //allow user to run help
80 if(option == "help") { help(); abort = true; calledHelp = true; }
81 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
84 vector<string> myArray = setParameters();
86 OptionParser parser(option);
87 map<string,string> parameters = parser.getParameters();
89 ValidParameters validParameter("filter.seqs");
90 map<string,string>::iterator it;
92 //check to make sure all parameters are valid for command
93 for (it = parameters.begin(); it != parameters.end(); it++) {
94 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
97 //initialize outputTypes
98 vector<string> tempOutNames;
99 outputTypes["fasta"] = tempOutNames;
100 outputTypes["filter"] = tempOutNames;
102 //if the user changes the input directory command factory will send this info to us in the output parameter
103 string inputDir = validParameter.validFile(parameters, "inputdir", false);
104 if (inputDir == "not found"){ inputDir = ""; }
107 it = parameters.find("fasta");
108 //user has given a template file
109 if(it != parameters.end()){
110 path = m->hasPath(it->second);
111 //if the user has not given a path then, add inputdir. else leave path alone.
112 if (path == "") { parameters["fasta"] = inputDir + it->second; }
115 it = parameters.find("hard");
116 //user has given a template file
117 if(it != parameters.end()){
118 path = m->hasPath(it->second);
119 //if the user has not given a path then, add inputdir. else leave path alone.
120 if (path == "") { parameters["hard"] = inputDir + it->second; }
124 //check for required parameters
125 fasta = validParameter.validFile(parameters, "fasta", false);
126 if (fasta == "not found") {
127 fasta = m->getFastaFile();
128 if (fasta != "") { fastafileNames.push_back(fasta); m->mothurOut("Using " + fasta + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
129 else { m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
132 m->splitAtDash(fasta, fastafileNames);
134 //go through files and make sure they are good, if not, then disregard them
135 for (int i = 0; i < fastafileNames.size(); i++) {
138 if (fastafileNames[i] == "current") {
139 fastafileNames[i] = m->getFastaFile();
140 if (fastafileNames[i] != "") { m->mothurOut("Using " + fastafileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
142 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
143 //erase from file list
144 fastafileNames.erase(fastafileNames.begin()+i);
150 if (inputDir != "") {
151 string path = m->hasPath(fastafileNames[i]);
152 //if the user has not given a path then, add inputdir. else leave path alone.
153 if (path == "") { fastafileNames[i] = inputDir + fastafileNames[i]; }
157 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
159 //if you can't open it, try default location
160 if (ableToOpen == 1) {
161 if (m->getDefaultPath() != "") { //default path is set
162 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
163 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
165 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
167 fastafileNames[i] = tryPath;
171 //if you can't open it, try default location
172 if (ableToOpen == 1) {
173 if (m->getOutputDir() != "") { //default path is set
174 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
175 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
177 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
179 fastafileNames[i] = tryPath;
185 if (ableToOpen == 1) {
186 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
187 //erase from file list
188 fastafileNames.erase(fastafileNames.begin()+i);
191 string simpleName = m->getSimpleName(fastafileNames[i]);
192 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
193 m->setFastaFile(fastafileNames[i]);
199 //make sure there is at least one valid file left
200 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
204 //if the user changes the output directory command factory will send this info to us in the output parameter
205 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
207 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
210 //check for optional parameter and set defaults
211 // ...at some point should added some additional type checking...
214 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
215 else if (hard == "not open") { hard = ""; abort = true; }
217 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
220 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
221 else { soft = (float)atoi(temp.c_str()) / 100.0; }
223 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
224 m->setProcessors(temp);
225 m->mothurConvert(temp, processors);
227 vertical = validParameter.validFile(parameters, "vertical", false);
228 if (vertical == "not found") {
229 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
230 else { vertical = "F"; }
237 catch(exception& e) {
238 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
242 /**************************************************************************************/
244 int FilterSeqsCommand::execute() {
247 if (abort == true) { if (calledHelp) { return 0; } return 2; }
250 m->openInputFile(fastafileNames[0], inFASTA);
252 Sequence testSeq(inFASTA);
253 alignmentLength = testSeq.getAlignLength();
256 ////////////create filter/////////////////
257 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
259 filter = createFilter();
261 m->mothurOutEndLine(); m->mothurOutEndLine();
263 if (m->control_pressed) { outputTypes.clear(); return 0; }
267 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
269 if (pid == 0) { //only one process should output the filter
274 //prevent giantic file name
276 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
277 else { filterFile = outputDir + filterFileName + ".filter"; }
279 m->openOutputFile(filterFile, outFilter);
280 outFilter << filter << endl;
282 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
288 ////////////run filter/////////////////
290 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
294 m->mothurOutEndLine(); m->mothurOutEndLine();
296 int filteredLength = 0;
297 for(int i=0;i<alignmentLength;i++){
298 if(filter[i] == '1'){ filteredLength++; }
301 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); } return 0; }
304 m->mothurOutEndLine();
305 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
306 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
307 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
308 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
310 //set fasta file as new current fastafile
312 itTypes = outputTypes.find("fasta");
313 if (itTypes != outputTypes.end()) {
314 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
317 m->mothurOutEndLine();
318 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
319 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
320 m->mothurOutEndLine();
325 catch(exception& e) {
326 m->errorOut(e, "FilterSeqsCommand", "execute");
330 /**************************************************************************************/
331 int FilterSeqsCommand::filterSequences() {
336 for (int s = 0; s < fastafileNames.size(); s++) {
338 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
340 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
342 int pid, numSeqsPerProcessor, num;
344 vector<unsigned long long>MPIPos;
347 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
348 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
352 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
353 int inMode=MPI_MODE_RDONLY;
355 char outFilename[1024];
356 strcpy(outFilename, filteredFasta.c_str());
358 char inFileName[1024];
359 strcpy(inFileName, fastafileNames[s].c_str());
361 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
362 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
364 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
366 if (pid == 0) { //you are the root process
368 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
371 //send file positions to all processes
372 for(int i = 1; i < processors; i++) {
373 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
374 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
377 //figure out how many sequences you have to do
378 numSeqsPerProcessor = num / processors;
379 int startIndex = pid * numSeqsPerProcessor;
380 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
384 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
386 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
389 for(int i = 1; i < processors; i++) {
391 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
394 }else { //you are a child process
395 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
396 MPIPos.resize(num+1);
398 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
400 //figure out how many sequences you have to align
401 numSeqsPerProcessor = num / processors;
402 int startIndex = pid * numSeqsPerProcessor;
403 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
407 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
409 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
414 //tell parent you are done.
415 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
418 MPI_File_close(&outMPI);
419 MPI_File_close(&inMPI);
420 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
424 vector<unsigned long long> positions;
425 if (savedPositions.size() != 0) { positions = savedPositions[s]; }
427 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
428 positions = m->divideFile(fastafileNames[s], processors);
431 int numFastaSeqs = 0;
432 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
433 if (positions.size() < processors) { processors = positions.size(); }
437 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
438 //vector<unsigned long long> positions = m->divideFile(fastafileNames[s], processors);
440 for (int i = 0; i < (positions.size()-1); i++) {
441 lines.push_back(new linePair(positions[i], positions[(i+1)]));
445 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
446 numSeqs += numFastaSeqs;
448 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
449 numSeqs += numFastaSeqs;
452 if (m->control_pressed) { return 1; }
455 lines.push_back(new linePair(0, 1000));
456 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
457 numSeqs += numFastaSeqs;
459 int numFastaSeqs = positions.size()-1;
460 //positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
462 //figure out how many sequences you have to process
463 int numSeqsPerProcessor = numFastaSeqs / processors;
464 for (int i = 0; i < processors; i++) {
465 int startIndex = i * numSeqsPerProcessor;
466 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
467 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
470 numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
471 numSeqs += numFastaSeqs;
474 if (m->control_pressed) { return 1; }
477 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
482 catch(exception& e) {
483 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
488 /**************************************************************************************/
489 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long long>& MPIPos) {
491 string outputString = "";
495 for(int i=0;i<num;i++){
497 if (m->control_pressed) { return 0; }
500 int length = MPIPos[start+i+1] - MPIPos[start+i];
501 char* buf4 = new char[length];
502 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
504 string tempBuf = buf4;
505 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
506 istringstream iss (tempBuf,istringstream::in);
509 Sequence seq(iss); m->gobble(iss);
511 if (seq.getName() != "") {
512 string align = seq.getAligned();
513 string filterSeq = "";
515 for(int j=0;j<alignmentLength;j++){
516 if(filter[j] == '1'){
517 filterSeq += align[j];
522 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
524 if(count % 10 == 0){ //output to file
525 //send results to parent
526 int length = outputString.length();
527 char* buf = new char[length];
528 memcpy(buf, outputString.c_str(), length);
530 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
537 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
540 if(outputString != ""){ //output to file
541 //send results to parent
542 int length = outputString.length();
543 char* buf = new char[length];
544 memcpy(buf, outputString.c_str(), length);
546 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
551 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
555 catch(exception& e) {
556 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
561 /**************************************************************************************/
562 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
565 m->openOutputFile(outputFilename, out);
568 m->openInputFile(inputFilename, in);
570 in.seekg(filePos->start);
577 if (m->control_pressed) { in.close(); out.close(); return 0; }
579 Sequence seq(in); m->gobble(in);
580 if (seq.getName() != "") {
581 string align = seq.getAligned();
582 string filterSeq = "";
584 for(int j=0;j<alignmentLength;j++){
585 if(filter[j] == '1'){
586 filterSeq += align[j];
590 out << '>' << seq.getName() << endl << filterSeq << endl;
594 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
595 unsigned long long pos = in.tellg();
596 if ((pos == -1) || (pos >= filePos->end)) { break; }
598 if (in.eof()) { break; }
602 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
605 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
613 catch(exception& e) {
614 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
618 /**************************************************************************************************/
620 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename, string filteredFastaName) {
627 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
630 //loop through and create all the processes you want
631 while (process != processors) {
635 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
638 string filteredFasta = filename + toString(getpid()) + ".temp";
639 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
641 //pass numSeqs to parent
643 string tempFile = filename + toString(getpid()) + ".num.temp";
644 m->openOutputFile(tempFile, out);
650 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
651 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
656 num = driverRunFilter(F, filteredFastaName, filename, lines[0]);
658 //force parent to wait until all the processes are done
659 for (int i=0;i<processIDS.size();i++) {
660 int temp = processIDS[i];
664 for (int i = 0; i < processIDS.size(); i++) {
666 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
667 m->openInputFile(tempFile, in);
668 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
669 in.close(); m->mothurRemove(tempFile);
671 m->appendFiles((filename + toString(processIDS[i]) + ".temp"), filteredFastaName);
672 m->mothurRemove((filename + toString(processIDS[i]) + ".temp"));
677 //////////////////////////////////////////////////////////////////////////////////////////////////////
678 //Windows version shared memory, so be careful when passing variables through the filterData struct.
679 //Above fork() will clone, so memory is separate, but that's not the case with windows,
680 //Taking advantage of shared memory to allow both threads to add info to F.
681 //////////////////////////////////////////////////////////////////////////////////////////////////////
683 vector<filterRunData*> pDataArray;
684 DWORD dwThreadIdArray[processors-1];
685 HANDLE hThreadArray[processors-1];
687 //Create processor worker threads.
688 for( int i=0; i<processors-1; i++){
690 string extension = "";
691 if (i != 0) { extension = toString(i) + ".temp"; }
693 filterRunData* tempFilter = new filterRunData(filter, filename, (filteredFastaName + extension), m, lines[i]->start, lines[i]->end, alignmentLength, i);
694 pDataArray.push_back(tempFilter);
695 processIDS.push_back(i);
697 hThreadArray[i] = CreateThread(NULL, 0, MyRunFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
700 num = driverRunFilter(F, (filteredFastaName + toString(processors-1) + ".temp"), filename, lines[processors-1]);
702 //Wait until all threads have terminated.
703 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
705 //Close all thread handles and free memory allocations.
706 for(int i=0; i < pDataArray.size(); i++){
707 num += pDataArray[i]->count;
708 CloseHandle(hThreadArray[i]);
709 delete pDataArray[i];
712 for (int i = 1; i < processors; i++) {
713 m->appendFiles((filteredFastaName + toString(i) + ".temp"), filteredFastaName);
714 m->mothurRemove((filteredFastaName + toString(i) + ".temp"));
721 catch(exception& e) {
722 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
726 /**************************************************************************************/
727 string FilterSeqsCommand::createFilter() {
729 string filterString = "";
732 if (soft != 0) { F.setSoft(soft); }
733 if (trump != '*') { F.setTrump(trump); }
735 F.setLength(alignmentLength);
737 if(trump != '*' || m->isTrue(vertical) || soft != 0){
741 if(hard.compare("") != 0) { F.doHard(hard); }
742 else { F.setFilter(string(alignmentLength, '1')); }
745 if(trump != '*' || m->isTrue(vertical) || soft != 0){
746 for (int s = 0; s < fastafileNames.size(); s++) {
748 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
751 int pid, numSeqsPerProcessor, num;
753 vector<unsigned long long> MPIPos;
757 MPI_Comm_size(MPI_COMM_WORLD, &processors);
758 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
760 //char* tempFileName = new char(fastafileNames[s].length());
761 //tempFileName = &(fastafileNames[s][0]);
763 char tempFileName[1024];
764 strcpy(tempFileName, fastafileNames[s].c_str());
766 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
768 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
770 if (pid == 0) { //you are the root process
771 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
774 //send file positions to all processes
775 for(int i = 1; i < processors; i++) {
776 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
777 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
780 //figure out how many sequences you have to do
781 numSeqsPerProcessor = num / processors;
782 int startIndex = pid * numSeqsPerProcessor;
783 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
787 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
789 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
791 }else { //i am the child process
792 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
793 MPIPos.resize(num+1);
795 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
797 //figure out how many sequences you have to align
798 numSeqsPerProcessor = num / processors;
799 int startIndex = pid * numSeqsPerProcessor;
800 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
804 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
806 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
809 MPI_File_close(&inMPI);
810 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
814 vector<unsigned long long> positions;
815 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
816 positions = m->divideFile(fastafileNames[s], processors);
817 for (int i = 0; i < (positions.size()-1); i++) {
818 lines.push_back(new linePair(positions[i], positions[(i+1)]));
822 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
823 numSeqs += numFastaSeqs;
825 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
826 numSeqs += numFastaSeqs;
830 lines.push_back(new linePair(0, 1000));
831 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
832 numSeqs += numFastaSeqs;
834 int numFastaSeqs = 0;
835 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
836 if (positions.size() < processors) { processors = positions.size(); }
838 //figure out how many sequences you have to process
839 int numSeqsPerProcessor = numFastaSeqs / processors;
840 for (int i = 0; i < processors; i++) {
841 int startIndex = i * numSeqsPerProcessor;
842 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
843 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
846 numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
847 numSeqs += numFastaSeqs;
850 //save the file positions so we can reuse them in the runFilter function
851 savedPositions[s] = positions;
853 if (m->control_pressed) { return filterString; }
862 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
865 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
867 if(trump != '*' || m->isTrue(vertical) || soft != 0){
869 if (pid == 0) { //only one process should output the filter
871 vector<int> temp; temp.resize(alignmentLength+1);
873 //get the frequencies from the child processes
874 for(int i = 1; i < processors; i++) {
876 for (int j = 0; j < 5; j++) {
878 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
879 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
881 if (receiveTag == Atag) { //you are recieveing the A frequencies
882 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
883 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
884 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
885 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
886 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
887 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
888 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
889 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
890 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
896 //send my fequency counts
898 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
900 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
902 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
904 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
905 F.gap.push_back(Gaptag);
906 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
911 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
913 if (pid == 0) { //only one process should output the filter
916 F.setNumSeqs(numSeqs);
917 if(m->isTrue(vertical) == 1) { F.doVertical(); }
918 if(soft != 0) { F.doSoft(); }
919 filterString = F.getFilter();
922 //send filter string to kids
923 //for(int i = 1; i < processors; i++) {
924 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
926 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
928 //recieve filterString
929 char* tempBuf = new char[alignmentLength];
930 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
931 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
933 filterString = tempBuf;
934 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
938 MPI_Barrier(MPI_COMM_WORLD);
943 catch(exception& e) {
944 m->errorOut(e, "FilterSeqsCommand", "createFilter");
948 /**************************************************************************************/
949 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
953 m->openInputFile(filename, in);
955 in.seekg(filePos->start);
962 if (m->control_pressed) { in.close(); return 1; }
964 Sequence seq(in); m->gobble(in);
965 if (seq.getName() != "") {
966 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
968 if(trump != '*') { F.doTrump(seq); }
969 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
974 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
975 unsigned long long pos = in.tellg();
976 if ((pos == -1) || (pos >= filePos->end)) { break; }
978 if (in.eof()) { break; }
982 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
985 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
990 catch(exception& e) {
991 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
996 /**************************************************************************************/
997 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long long>& MPIPos) {
1002 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
1004 for(int i=0;i<num;i++){
1006 if (m->control_pressed) { return 0; }
1008 //read next sequence
1009 int length = MPIPos[start+i+1] - MPIPos[start+i];
1011 char* buf4 = new char[length];
1012 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
1014 string tempBuf = buf4;
1015 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
1016 istringstream iss (tempBuf,istringstream::in);
1021 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
1023 if(trump != '*'){ F.doTrump(seq); }
1024 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
1028 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
1032 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
1036 catch(exception& e) {
1037 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
1042 /**************************************************************************************************/
1044 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
1050 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1052 //loop through and create all the processes you want
1053 while (process != processors) {
1057 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
1059 }else if (pid == 0){
1060 //reset child's filter counts to 0;
1061 F.a.clear(); F.a.resize(alignmentLength, 0);
1062 F.t.clear(); F.t.resize(alignmentLength, 0);
1063 F.g.clear(); F.g.resize(alignmentLength, 0);
1064 F.c.clear(); F.c.resize(alignmentLength, 0);
1065 F.gap.clear(); F.gap.resize(alignmentLength, 0);
1067 num = driverCreateFilter(F, filename, lines[process]);
1069 //write out filter counts to file
1070 filename += toString(getpid()) + "filterValues.temp";
1072 m->openOutputFile(filename, out);
1075 out << F.getFilter() << endl;
1076 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
1077 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
1078 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
1079 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
1080 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
1082 //cout << F.getFilter() << endl;
1087 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
1088 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1093 //parent do your part
1094 num = driverCreateFilter(F, filename, lines[0]);
1096 //force parent to wait until all the processes are done
1097 for (int i=0;i<(processors-1);i++) {
1098 int temp = processIDS[i];
1102 //parent reads in and combines Filter info
1103 for (int i = 0; i < processIDS.size(); i++) {
1104 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1106 m->openInputFile(tempFilename, in);
1109 string tempFilterString;
1111 in >> tempNum; m->gobble(in); num += tempNum;
1113 in >> tempFilterString;
1114 F.mergeFilter(tempFilterString);
1116 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1117 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1118 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1119 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1120 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1123 m->mothurRemove(tempFilename);
1129 //////////////////////////////////////////////////////////////////////////////////////////////////////
1130 //Windows version shared memory, so be careful when passing variables through the filterData struct.
1131 //Above fork() will clone, so memory is separate, but that's not the case with windows,
1132 //Taking advantage of shared memory to allow both threads to add info to F.
1133 //////////////////////////////////////////////////////////////////////////////////////////////////////
1135 vector<filterData*> pDataArray;
1136 DWORD dwThreadIdArray[processors];
1137 HANDLE hThreadArray[processors];
1139 //Create processor worker threads.
1140 for( int i=0; i<processors; i++ ){
1142 filterData* tempFilter = new filterData(filename, m, lines[i]->start, lines[i]->end, alignmentLength, trump, vertical, soft, hard, i);
1143 pDataArray.push_back(tempFilter);
1144 processIDS.push_back(i);
1146 hThreadArray[i] = CreateThread(NULL, 0, MyCreateFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
1149 //Wait until all threads have terminated.
1150 WaitForMultipleObjects(processors, hThreadArray, TRUE, INFINITE);
1152 //Close all thread handles and free memory allocations.
1153 for(int i=0; i < pDataArray.size(); i++){
1154 num += pDataArray[i]->count;
1155 F.mergeFilter(pDataArray[i]->F.getFilter());
1157 for (int k = 0; k < alignmentLength; k++) { F.a[k] += pDataArray[i]->F.a[k]; }
1158 for (int k = 0; k < alignmentLength; k++) { F.t[k] += pDataArray[i]->F.t[k]; }
1159 for (int k = 0; k < alignmentLength; k++) { F.g[k] += pDataArray[i]->F.g[k]; }
1160 for (int k = 0; k < alignmentLength; k++) { F.c[k] += pDataArray[i]->F.c[k]; }
1161 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += pDataArray[i]->F.gap[k]; }
1163 CloseHandle(hThreadArray[i]);
1164 delete pDataArray[i];
1171 catch(exception& e) {
1172 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1176 /**************************************************************************************/