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"
13 //**********************************************************************************************************************
14 vector<string> FilterSeqsCommand::getValidParameters(){
16 string Array[] = {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
17 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
21 m->errorOut(e, "FilterSeqsCommand", "getValidParameters");
25 //**********************************************************************************************************************
26 FilterSeqsCommand::FilterSeqsCommand(){
28 abort = true; calledHelp = true;
29 vector<string> tempOutNames;
30 outputTypes["fasta"] = tempOutNames;
31 outputTypes["filter"] = tempOutNames;
34 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
38 //**********************************************************************************************************************
39 vector<string> FilterSeqsCommand::getRequiredParameters(){
41 string Array[] = {"fasta"};
42 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
46 m->errorOut(e, "FilterSeqsCommand", "getRequiredParameters");
50 //**********************************************************************************************************************
51 vector<string> FilterSeqsCommand::getRequiredFiles(){
53 vector<string> myArray;
57 m->errorOut(e, "FilterSeqsCommand", "getRequiredFiles");
61 /**************************************************************************************/
62 FilterSeqsCommand::FilterSeqsCommand(string option) {
64 abort = false; calledHelp = false;
67 //allow user to run help
68 if(option == "help") { help(); abort = true; calledHelp = true; }
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)));
75 OptionParser parser(option);
76 map<string,string> parameters = parser.getParameters();
78 ValidParameters validParameter("filter.seqs");
79 map<string,string>::iterator it;
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; }
86 //initialize outputTypes
87 vector<string> tempOutNames;
88 outputTypes["fasta"] = tempOutNames;
89 outputTypes["filter"] = tempOutNames;
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 = ""; }
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; }
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; }
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; }
117 m->splitAtDash(fasta, fastafileNames);
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]; }
128 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
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();
136 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
138 fastafileNames[i] = tryPath;
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();
148 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
150 fastafileNames[i] = tryPath;
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);
162 string simpleName = m->getSimpleName(fastafileNames[i]);
163 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
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; }
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"){
176 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
179 //check for optional parameter and set defaults
180 // ...at some point should added some additional type checking...
183 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
184 else if (hard == "not open") { hard = ""; abort = true; }
186 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
189 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
190 else { soft = (float)atoi(temp.c_str()) / 100.0; }
192 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
193 convert(temp, processors);
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"; }
205 catch(exception& e) {
206 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
211 //**********************************************************************************************************************
213 void FilterSeqsCommand::help(){
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");
231 catch(exception& e) {
232 m->errorOut(e, "FilterSeqsCommand", "help");
237 /**************************************************************************************/
239 int FilterSeqsCommand::execute() {
242 if (abort == true) { if (calledHelp) { return 0; } return 2; }
245 m->openInputFile(fastafileNames[0], inFASTA);
247 Sequence testSeq(inFASTA);
248 alignmentLength = testSeq.getAlignLength();
251 ////////////create filter/////////////////
252 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
254 filter = createFilter();
256 m->mothurOutEndLine(); m->mothurOutEndLine();
258 if (m->control_pressed) { outputTypes.clear(); return 0; }
262 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
264 if (pid == 0) { //only one process should output the filter
269 //prevent giantic file name
271 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
272 else { filterFile = outputDir + filterFileName + ".filter"; }
274 m->openOutputFile(filterFile, outFilter);
275 outFilter << filter << endl;
277 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
283 ////////////run filter/////////////////
285 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
289 m->mothurOutEndLine(); m->mothurOutEndLine();
291 int filteredLength = 0;
292 for(int i=0;i<alignmentLength;i++){
293 if(filter[i] == '1'){ filteredLength++; }
296 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
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();
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();
314 catch(exception& e) {
315 m->errorOut(e, "FilterSeqsCommand", "execute");
319 /**************************************************************************************/
320 int FilterSeqsCommand::filterSequences() {
325 for (int s = 0; s < fastafileNames.size(); s++) {
327 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
329 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
331 int pid, numSeqsPerProcessor, num;
333 vector<unsigned long int>MPIPos;
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
341 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
342 int inMode=MPI_MODE_RDONLY;
344 char outFilename[1024];
345 strcpy(outFilename, filteredFasta.c_str());
347 char inFileName[1024];
348 strcpy(inFileName, fastafileNames[s].c_str());
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);
353 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
355 if (pid == 0) { //you are the root process
357 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
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);
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; }
373 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
375 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
378 for(int i = 1; i < processors; i++) {
380 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
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);
387 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
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; }
396 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
398 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
403 //tell parent you are done.
404 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
407 MPI_File_close(&outMPI);
408 MPI_File_close(&inMPI);
409 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
412 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
414 for (int i = 0; i < (positions.size()-1); i++) {
415 lines.push_back(new linePair(positions[i], positions[(i+1)]));
417 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
419 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
420 numSeqs += numFastaSeqs;
422 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
423 numSeqs += numFastaSeqs;
425 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
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());
434 if (m->control_pressed) { return 1; }
436 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
437 numSeqs += numFastaSeqs;
439 if (m->control_pressed) { return 1; }
442 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
447 catch(exception& e) {
448 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
453 /**************************************************************************************/
454 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
456 string outputString = "";
460 for(int i=0;i<num;i++){
462 if (m->control_pressed) { return 0; }
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);
469 string tempBuf = buf4;
470 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
471 istringstream iss (tempBuf,istringstream::in);
474 Sequence seq(iss); m->gobble(iss);
476 if (seq.getName() != "") {
477 string align = seq.getAligned();
478 string filterSeq = "";
480 for(int j=0;j<alignmentLength;j++){
481 if(filter[j] == '1'){
482 filterSeq += align[j];
487 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
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);
495 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
502 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
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);
511 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
516 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
520 catch(exception& e) {
521 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
526 /**************************************************************************************/
527 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
530 m->openOutputFile(outputFilename, out);
533 m->openInputFile(inputFilename, in);
535 in.seekg(filePos->start);
542 if (m->control_pressed) { in.close(); out.close(); return 0; }
544 Sequence seq(in); m->gobble(in);
545 if (seq.getName() != "") {
546 string align = seq.getAligned();
547 string filterSeq = "";
549 for(int j=0;j<alignmentLength;j++){
550 if(filter[j] == '1'){
551 filterSeq += align[j];
555 out << '>' << seq.getName() << endl << filterSeq << endl;
559 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
560 unsigned long int pos = in.tellg();
561 if ((pos == -1) || (pos >= filePos->end)) { break; }
563 if (in.eof()) { break; }
567 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
570 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
578 catch(exception& e) {
579 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
583 /**************************************************************************************************/
585 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
587 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
592 //loop through and create all the processes you want
593 while (process != processors) {
597 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
600 string filteredFasta = filename + toString(getpid()) + ".temp";
601 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
603 //pass numSeqs to parent
605 string tempFile = filename + toString(getpid()) + ".num.temp";
606 m->openOutputFile(tempFile, out);
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); }
618 //force parent to wait until all the processes are done
619 for (int i=0;i<processors;i++) {
620 int temp = processIDS[i];
624 for (int i = 0; i < processIDS.size(); i++) {
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());
636 catch(exception& e) {
637 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
641 /**************************************************************************************/
642 string FilterSeqsCommand::createFilter() {
644 string filterString = "";
647 if (soft != 0) { F.setSoft(soft); }
648 if (trump != '*') { F.setTrump(trump); }
650 F.setLength(alignmentLength);
652 if(trump != '*' || m->isTrue(vertical) || soft != 0){
656 if(hard.compare("") != 0) { F.doHard(hard); }
657 else { F.setFilter(string(alignmentLength, '1')); }
660 if(trump != '*' || m->isTrue(vertical) || soft != 0){
661 for (int s = 0; s < fastafileNames.size(); s++) {
663 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
666 int pid, numSeqsPerProcessor, num;
668 vector<unsigned long int> MPIPos;
672 MPI_Comm_size(MPI_COMM_WORLD, &processors);
673 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
675 //char* tempFileName = new char(fastafileNames[s].length());
676 //tempFileName = &(fastafileNames[s][0]);
678 char tempFileName[1024];
679 strcpy(tempFileName, fastafileNames[s].c_str());
681 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
683 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
685 if (pid == 0) { //you are the root process
686 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
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);
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; }
702 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
704 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
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);
710 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
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; }
719 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
721 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
724 MPI_File_close(&inMPI);
725 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
728 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
730 for (int i = 0; i < (positions.size()-1); i++) {
731 lines.push_back(new linePair(positions[i], positions[(i+1)]));
733 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
735 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
736 numSeqs += numFastaSeqs;
738 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
739 numSeqs += numFastaSeqs;
742 if (m->control_pressed) { return filterString; }
744 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
745 numSeqs += numFastaSeqs;
746 if (m->control_pressed) { return filterString; }
756 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
759 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
761 if(trump != '*' || m->isTrue(vertical) || soft != 0){
763 if (pid == 0) { //only one process should output the filter
765 vector<int> temp; temp.resize(alignmentLength+1);
767 //get the frequencies from the child processes
768 for(int i = 1; i < processors; i++) {
770 for (int j = 0; j < 5; j++) {
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
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]; }
790 //send my fequency counts
792 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
794 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
796 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
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);
805 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
807 if (pid == 0) { //only one process should output the filter
810 F.setNumSeqs(numSeqs);
811 if(m->isTrue(vertical) == 1) { F.doVertical(); }
812 if(soft != 0) { F.doSoft(); }
813 filterString = F.getFilter();
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);
820 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
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);
827 filterString = tempBuf;
828 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
832 MPI_Barrier(MPI_COMM_WORLD);
837 catch(exception& e) {
838 m->errorOut(e, "FilterSeqsCommand", "createFilter");
842 /**************************************************************************************/
843 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
847 m->openInputFile(filename, in);
849 in.seekg(filePos->start);
856 if (m->control_pressed) { in.close(); return 1; }
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; }
862 if(trump != '*') { F.doTrump(seq); }
863 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
868 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
869 unsigned long int pos = in.tellg();
870 if ((pos == -1) || (pos >= filePos->end)) { break; }
872 if (in.eof()) { break; }
876 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
879 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
884 catch(exception& e) {
885 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
890 /**************************************************************************************/
891 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
896 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
898 for(int i=0;i<num;i++){
900 if (m->control_pressed) { return 0; }
903 int length = MPIPos[start+i+1] - MPIPos[start+i];
905 char* buf4 = new char[length];
906 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
908 string tempBuf = buf4;
909 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
910 istringstream iss (tempBuf,istringstream::in);
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); }
917 if(trump != '*'){ F.doTrump(seq); }
918 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
922 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
926 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
930 catch(exception& e) {
931 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
936 /**************************************************************************************************/
938 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
940 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
945 //loop through and create all the processes you want
946 while (process != processors) {
950 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
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);
960 num = driverCreateFilter(F, filename, lines[process]);
962 //write out filter counts to file
963 filename += toString(getpid()) + "filterValues.temp";
965 m->openOutputFile(filename, out);
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;
975 //cout << F.getFilter() << endl;
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); }
986 //parent do your part
987 num = driverCreateFilter(F, filename, lines[0]);
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];
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";
999 m->openInputFile(tempFilename, in);
1002 string tempFilterString;
1004 in >> tempNum; m->gobble(in); num += tempNum;
1006 in >> tempFilterString;
1007 F.mergeFilter(tempFilterString);
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);
1016 remove(tempFilename.c_str());
1022 catch(exception& e) {
1023 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1027 /**************************************************************************************/