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('.'));
198 //make sure there is at least one valid file left
199 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
203 //if the user changes the output directory command factory will send this info to us in the output parameter
204 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
206 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
209 //check for optional parameter and set defaults
210 // ...at some point should added some additional type checking...
213 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
214 else if (hard == "not open") { hard = ""; abort = true; }
216 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
219 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
220 else { soft = (float)atoi(temp.c_str()) / 100.0; }
222 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
223 m->setProcessors(temp);
224 convert(temp, processors);
226 vertical = validParameter.validFile(parameters, "vertical", false);
227 if (vertical == "not found") {
228 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
229 else { vertical = "F"; }
236 catch(exception& e) {
237 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
241 /**************************************************************************************/
243 int FilterSeqsCommand::execute() {
246 if (abort == true) { if (calledHelp) { return 0; } return 2; }
249 m->openInputFile(fastafileNames[0], inFASTA);
251 Sequence testSeq(inFASTA);
252 alignmentLength = testSeq.getAlignLength();
255 ////////////create filter/////////////////
256 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
258 filter = createFilter();
260 m->mothurOutEndLine(); m->mothurOutEndLine();
262 if (m->control_pressed) { outputTypes.clear(); return 0; }
266 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
268 if (pid == 0) { //only one process should output the filter
273 //prevent giantic file name
275 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
276 else { filterFile = outputDir + filterFileName + ".filter"; }
278 m->openOutputFile(filterFile, outFilter);
279 outFilter << filter << endl;
281 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
287 ////////////run filter/////////////////
289 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
293 m->mothurOutEndLine(); m->mothurOutEndLine();
295 int filteredLength = 0;
296 for(int i=0;i<alignmentLength;i++){
297 if(filter[i] == '1'){ filteredLength++; }
300 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
303 m->mothurOutEndLine();
304 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
305 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
306 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
307 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
309 //set fasta file as new current fastafile
311 itTypes = outputTypes.find("fasta");
312 if (itTypes != outputTypes.end()) {
313 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
316 m->mothurOutEndLine();
317 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
318 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
319 m->mothurOutEndLine();
324 catch(exception& e) {
325 m->errorOut(e, "FilterSeqsCommand", "execute");
329 /**************************************************************************************/
330 int FilterSeqsCommand::filterSequences() {
335 for (int s = 0; s < fastafileNames.size(); s++) {
337 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
339 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
341 int pid, numSeqsPerProcessor, num;
343 vector<unsigned long int>MPIPos;
346 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
347 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
351 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
352 int inMode=MPI_MODE_RDONLY;
354 char outFilename[1024];
355 strcpy(outFilename, filteredFasta.c_str());
357 char inFileName[1024];
358 strcpy(inFileName, fastafileNames[s].c_str());
360 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
361 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
363 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
365 if (pid == 0) { //you are the root process
367 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
370 //send file positions to all processes
371 for(int i = 1; i < processors; i++) {
372 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
373 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
376 //figure out how many sequences you have to do
377 numSeqsPerProcessor = num / processors;
378 int startIndex = pid * numSeqsPerProcessor;
379 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
383 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
385 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
388 for(int i = 1; i < processors; i++) {
390 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
393 }else { //you are a child process
394 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
395 MPIPos.resize(num+1);
397 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
399 //figure out how many sequences you have to align
400 numSeqsPerProcessor = num / processors;
401 int startIndex = pid * numSeqsPerProcessor;
402 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
406 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
408 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
413 //tell parent you are done.
414 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
417 MPI_File_close(&outMPI);
418 MPI_File_close(&inMPI);
419 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
422 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
424 for (int i = 0; i < (positions.size()-1); i++) {
425 lines.push_back(new linePair(positions[i], positions[(i+1)]));
427 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
429 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
430 numSeqs += numFastaSeqs;
432 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
433 numSeqs += numFastaSeqs;
435 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
438 for(int i=1;i<processors;i++){
439 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
440 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
444 if (m->control_pressed) { return 1; }
446 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
447 numSeqs += numFastaSeqs;
449 if (m->control_pressed) { return 1; }
452 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
457 catch(exception& e) {
458 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
463 /**************************************************************************************/
464 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
466 string outputString = "";
470 for(int i=0;i<num;i++){
472 if (m->control_pressed) { return 0; }
475 int length = MPIPos[start+i+1] - MPIPos[start+i];
476 char* buf4 = new char[length];
477 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
479 string tempBuf = buf4;
480 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
481 istringstream iss (tempBuf,istringstream::in);
484 Sequence seq(iss); m->gobble(iss);
486 if (seq.getName() != "") {
487 string align = seq.getAligned();
488 string filterSeq = "";
490 for(int j=0;j<alignmentLength;j++){
491 if(filter[j] == '1'){
492 filterSeq += align[j];
497 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
499 if(count % 10 == 0){ //output to file
500 //send results to parent
501 int length = outputString.length();
502 char* buf = new char[length];
503 memcpy(buf, outputString.c_str(), length);
505 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
512 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
515 if(outputString != ""){ //output to file
516 //send results to parent
517 int length = outputString.length();
518 char* buf = new char[length];
519 memcpy(buf, outputString.c_str(), length);
521 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
526 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
530 catch(exception& e) {
531 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
536 /**************************************************************************************/
537 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
540 m->openOutputFile(outputFilename, out);
543 m->openInputFile(inputFilename, in);
545 in.seekg(filePos->start);
552 if (m->control_pressed) { in.close(); out.close(); return 0; }
554 Sequence seq(in); m->gobble(in);
555 if (seq.getName() != "") {
556 string align = seq.getAligned();
557 string filterSeq = "";
559 for(int j=0;j<alignmentLength;j++){
560 if(filter[j] == '1'){
561 filterSeq += align[j];
565 out << '>' << seq.getName() << endl << filterSeq << endl;
569 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
570 unsigned long int pos = in.tellg();
571 if ((pos == -1) || (pos >= filePos->end)) { break; }
573 if (in.eof()) { break; }
577 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
580 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
588 catch(exception& e) {
589 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
593 /**************************************************************************************************/
595 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
597 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
602 //loop through and create all the processes you want
603 while (process != processors) {
607 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
610 string filteredFasta = filename + toString(getpid()) + ".temp";
611 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
613 //pass numSeqs to parent
615 string tempFile = filename + toString(getpid()) + ".num.temp";
616 m->openOutputFile(tempFile, out);
622 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
623 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
628 //force parent to wait until all the processes are done
629 for (int i=0;i<processors;i++) {
630 int temp = processIDS[i];
634 for (int i = 0; i < processIDS.size(); i++) {
636 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
637 m->openInputFile(tempFile, in);
638 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
639 in.close(); remove(tempFile.c_str());
646 catch(exception& e) {
647 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
651 /**************************************************************************************/
652 string FilterSeqsCommand::createFilter() {
654 string filterString = "";
657 if (soft != 0) { F.setSoft(soft); }
658 if (trump != '*') { F.setTrump(trump); }
660 F.setLength(alignmentLength);
662 if(trump != '*' || m->isTrue(vertical) || soft != 0){
666 if(hard.compare("") != 0) { F.doHard(hard); }
667 else { F.setFilter(string(alignmentLength, '1')); }
670 if(trump != '*' || m->isTrue(vertical) || soft != 0){
671 for (int s = 0; s < fastafileNames.size(); s++) {
673 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
676 int pid, numSeqsPerProcessor, num;
678 vector<unsigned long int> MPIPos;
682 MPI_Comm_size(MPI_COMM_WORLD, &processors);
683 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
685 //char* tempFileName = new char(fastafileNames[s].length());
686 //tempFileName = &(fastafileNames[s][0]);
688 char tempFileName[1024];
689 strcpy(tempFileName, fastafileNames[s].c_str());
691 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
693 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
695 if (pid == 0) { //you are the root process
696 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
699 //send file positions to all processes
700 for(int i = 1; i < processors; i++) {
701 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
702 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
705 //figure out how many sequences you have to do
706 numSeqsPerProcessor = num / processors;
707 int startIndex = pid * numSeqsPerProcessor;
708 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
712 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
714 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
716 }else { //i am the child process
717 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
718 MPIPos.resize(num+1);
720 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
722 //figure out how many sequences you have to align
723 numSeqsPerProcessor = num / processors;
724 int startIndex = pid * numSeqsPerProcessor;
725 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
729 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
731 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
734 MPI_File_close(&inMPI);
735 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
738 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
740 for (int i = 0; i < (positions.size()-1); i++) {
741 lines.push_back(new linePair(positions[i], positions[(i+1)]));
743 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
745 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
746 numSeqs += numFastaSeqs;
748 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
749 numSeqs += numFastaSeqs;
752 if (m->control_pressed) { return filterString; }
754 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
755 numSeqs += numFastaSeqs;
756 if (m->control_pressed) { return filterString; }
766 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
769 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
771 if(trump != '*' || m->isTrue(vertical) || soft != 0){
773 if (pid == 0) { //only one process should output the filter
775 vector<int> temp; temp.resize(alignmentLength+1);
777 //get the frequencies from the child processes
778 for(int i = 1; i < processors; i++) {
780 for (int j = 0; j < 5; j++) {
782 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
783 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
785 if (receiveTag == Atag) { //you are recieveing the A frequencies
786 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
787 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
788 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
789 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
790 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
791 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
792 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
793 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
794 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
800 //send my fequency counts
802 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
804 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
806 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
808 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
809 F.gap.push_back(Gaptag);
810 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
815 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
817 if (pid == 0) { //only one process should output the filter
820 F.setNumSeqs(numSeqs);
821 if(m->isTrue(vertical) == 1) { F.doVertical(); }
822 if(soft != 0) { F.doSoft(); }
823 filterString = F.getFilter();
826 //send filter string to kids
827 //for(int i = 1; i < processors; i++) {
828 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
830 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
832 //recieve filterString
833 char* tempBuf = new char[alignmentLength];
834 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
835 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
837 filterString = tempBuf;
838 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
842 MPI_Barrier(MPI_COMM_WORLD);
847 catch(exception& e) {
848 m->errorOut(e, "FilterSeqsCommand", "createFilter");
852 /**************************************************************************************/
853 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
857 m->openInputFile(filename, in);
859 in.seekg(filePos->start);
866 if (m->control_pressed) { in.close(); return 1; }
868 Sequence seq(in); m->gobble(in);
869 if (seq.getName() != "") {
870 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
872 if(trump != '*') { F.doTrump(seq); }
873 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
878 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
879 unsigned long int pos = in.tellg();
880 if ((pos == -1) || (pos >= filePos->end)) { break; }
882 if (in.eof()) { break; }
886 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
889 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
894 catch(exception& e) {
895 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
900 /**************************************************************************************/
901 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
906 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
908 for(int i=0;i<num;i++){
910 if (m->control_pressed) { return 0; }
913 int length = MPIPos[start+i+1] - MPIPos[start+i];
915 char* buf4 = new char[length];
916 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
918 string tempBuf = buf4;
919 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
920 istringstream iss (tempBuf,istringstream::in);
925 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
927 if(trump != '*'){ F.doTrump(seq); }
928 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
932 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
936 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
940 catch(exception& e) {
941 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
946 /**************************************************************************************************/
948 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
950 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
955 //loop through and create all the processes you want
956 while (process != processors) {
960 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
963 //reset child's filter counts to 0;
964 F.a.clear(); F.a.resize(alignmentLength, 0);
965 F.t.clear(); F.t.resize(alignmentLength, 0);
966 F.g.clear(); F.g.resize(alignmentLength, 0);
967 F.c.clear(); F.c.resize(alignmentLength, 0);
968 F.gap.clear(); F.gap.resize(alignmentLength, 0);
970 num = driverCreateFilter(F, filename, lines[process]);
972 //write out filter counts to file
973 filename += toString(getpid()) + "filterValues.temp";
975 m->openOutputFile(filename, out);
978 out << F.getFilter() << endl;
979 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
980 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
981 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
982 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
983 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
985 //cout << F.getFilter() << endl;
990 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
991 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
996 //parent do your part
997 num = driverCreateFilter(F, filename, lines[0]);
999 //force parent to wait until all the processes are done
1000 for (int i=0;i<(processors-1);i++) {
1001 int temp = processIDS[i];
1005 //parent reads in and combines Filter info
1006 for (int i = 0; i < processIDS.size(); i++) {
1007 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1009 m->openInputFile(tempFilename, in);
1012 string tempFilterString;
1014 in >> tempNum; m->gobble(in); num += tempNum;
1016 in >> tempFilterString;
1017 F.mergeFilter(tempFilterString);
1019 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1020 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1021 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1022 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1023 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1026 remove(tempFilename.c_str());
1032 catch(exception& e) {
1033 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1037 /**************************************************************************************/