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 //initialize outputTypes
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) {
67 //allow user to run help
68 if(option == "help") { help(); abort = 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) { return 0; }
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, start, end, 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
342 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
343 int inMode=MPI_MODE_RDONLY;
345 char outFilename[1024];
346 strcpy(outFilename, filteredFasta.c_str());
348 char inFileName[1024];
349 strcpy(inFileName, fastafileNames[s].c_str());
351 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
352 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
354 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
356 if (pid == 0) { //you are the root process
358 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
361 //send file positions to all processes
362 for(int i = 1; i < processors; i++) {
363 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
364 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
367 //figure out how many sequences you have to do
368 numSeqsPerProcessor = num / processors;
369 int startIndex = pid * numSeqsPerProcessor;
370 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
374 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
376 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
379 for(int i = 1; i < processors; i++) {
381 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
384 }else { //you are a child process
385 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
386 MPIPos.resize(num+1);
388 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
390 //figure out how many sequences you have to align
391 numSeqsPerProcessor = num / processors;
392 int startIndex = pid * numSeqsPerProcessor;
393 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
397 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
399 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
404 //tell parent you are done.
405 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
408 MPI_File_close(&outMPI);
409 MPI_File_close(&inMPI);
410 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
413 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
415 for (int i = 0; i < (positions.size()-1); i++) {
416 lines.push_back(new linePair(positions[i], positions[(i+1)]));
418 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
420 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
421 numSeqs += numFastaSeqs;
423 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
424 numSeqs += numFastaSeqs;
426 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
429 for(int i=1;i<processors;i++){
430 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
431 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
435 if (m->control_pressed) { return 1; }
437 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
438 numSeqs += numFastaSeqs;
440 if (m->control_pressed) { return 1; }
443 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
448 catch(exception& e) {
449 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
454 /**************************************************************************************/
455 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
457 string outputString = "";
461 for(int i=0;i<num;i++){
463 if (m->control_pressed) { return 0; }
466 int length = MPIPos[start+i+1] - MPIPos[start+i];
467 char* buf4 = new char[length];
468 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
470 string tempBuf = buf4;
471 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
472 istringstream iss (tempBuf,istringstream::in);
475 Sequence seq(iss); m->gobble(iss);
477 if (seq.getName() != "") {
478 string align = seq.getAligned();
479 string filterSeq = "";
481 for(int j=0;j<alignmentLength;j++){
482 if(filter[j] == '1'){
483 filterSeq += align[j];
488 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
490 if(count % 10 == 0){ //output to file
491 //send results to parent
492 int length = outputString.length();
493 char* buf = new char[length];
494 memcpy(buf, outputString.c_str(), length);
496 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
503 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
506 if(outputString != ""){ //output to file
507 //send results to parent
508 int length = outputString.length();
509 char* buf = new char[length];
510 memcpy(buf, outputString.c_str(), length);
512 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
517 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
521 catch(exception& e) {
522 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
527 /**************************************************************************************/
528 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
531 m->openOutputFile(outputFilename, out);
534 m->openInputFile(inputFilename, in);
536 in.seekg(filePos->start);
543 if (m->control_pressed) { in.close(); out.close(); return 0; }
545 Sequence seq(in); m->gobble(in);
546 if (seq.getName() != "") {
547 string align = seq.getAligned();
548 string filterSeq = "";
550 for(int j=0;j<alignmentLength;j++){
551 if(filter[j] == '1'){
552 filterSeq += align[j];
556 out << '>' << seq.getName() << endl << filterSeq << endl;
560 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
561 unsigned long int pos = in.tellg();
562 if ((pos == -1) || (pos >= filePos->end)) { break; }
564 if (in.eof()) { break; }
568 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
571 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
579 catch(exception& e) {
580 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
584 /**************************************************************************************************/
586 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
588 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
593 //loop through and create all the processes you want
594 while (process != processors) {
598 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
601 string filteredFasta = filename + toString(getpid()) + ".temp";
602 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
604 //pass numSeqs to parent
606 string tempFile = filename + toString(getpid()) + ".num.temp";
607 m->openOutputFile(tempFile, out);
612 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
615 //force parent to wait until all the processes are done
616 for (int i=0;i<processors;i++) {
617 int temp = processIDS[i];
621 for (int i = 0; i < processIDS.size(); i++) {
623 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
624 m->openInputFile(tempFile, in);
625 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
626 in.close(); remove(tempFile.c_str());
633 catch(exception& e) {
634 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
638 /**************************************************************************************/
639 string FilterSeqsCommand::createFilter() {
641 string filterString = "";
644 if (soft != 0) { F.setSoft(soft); }
645 if (trump != '*') { F.setTrump(trump); }
647 F.setLength(alignmentLength);
649 if(trump != '*' || m->isTrue(vertical) || soft != 0){
653 if(hard.compare("") != 0) { F.doHard(hard); }
654 else { F.setFilter(string(alignmentLength, '1')); }
657 if(trump != '*' || m->isTrue(vertical) || soft != 0){
658 for (int s = 0; s < fastafileNames.size(); s++) {
660 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
663 int pid, numSeqsPerProcessor, num;
665 vector<unsigned long int> MPIPos;
669 MPI_Comm_size(MPI_COMM_WORLD, &processors);
670 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
672 //char* tempFileName = new char(fastafileNames[s].length());
673 //tempFileName = &(fastafileNames[s][0]);
675 char tempFileName[1024];
676 strcpy(tempFileName, fastafileNames[s].c_str());
678 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
680 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
682 if (pid == 0) { //you are the root process
683 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
686 //send file positions to all processes
687 for(int i = 1; i < processors; i++) {
688 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
689 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
692 //figure out how many sequences you have to do
693 numSeqsPerProcessor = num / processors;
694 int startIndex = pid * numSeqsPerProcessor;
695 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
699 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
701 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
703 }else { //i am the child process
704 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
705 MPIPos.resize(num+1);
707 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
709 //figure out how many sequences you have to align
710 numSeqsPerProcessor = num / processors;
711 int startIndex = pid * numSeqsPerProcessor;
712 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
716 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
718 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
721 MPI_File_close(&inMPI);
722 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
725 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
727 for (int i = 0; i < (positions.size()-1); i++) {
728 lines.push_back(new linePair(positions[i], positions[(i+1)]));
730 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
732 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
733 numSeqs += numFastaSeqs;
735 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
736 numSeqs += numFastaSeqs;
739 if (m->control_pressed) { return filterString; }
741 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
742 numSeqs += numFastaSeqs;
743 if (m->control_pressed) { return filterString; }
753 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
756 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
758 if(trump != '*' || m->isTrue(vertical) || soft != 0){
760 if (pid == 0) { //only one process should output the filter
762 vector<int> temp; temp.resize(alignmentLength+1);
764 //get the frequencies from the child processes
765 for(int i = 1; i < processors; i++) {
767 for (int j = 0; j < 5; j++) {
769 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
770 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
772 if (receiveTag == Atag) { //you are recieveing the A frequencies
773 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
774 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
775 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
776 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
777 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
778 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
779 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
780 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
781 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
787 //send my fequency counts
789 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
791 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
793 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
795 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
796 F.gap.push_back(Gaptag);
797 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
802 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
804 if (pid == 0) { //only one process should output the filter
807 F.setNumSeqs(numSeqs);
808 if(m->isTrue(vertical) == 1) { F.doVertical(); }
809 if(soft != 0) { F.doSoft(); }
810 filterString = F.getFilter();
813 //send filter string to kids
814 //for(int i = 1; i < processors; i++) {
815 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
817 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
819 //recieve filterString
820 char* tempBuf = new char[alignmentLength];
821 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
822 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
824 filterString = tempBuf;
825 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
829 MPI_Barrier(MPI_COMM_WORLD);
834 catch(exception& e) {
835 m->errorOut(e, "FilterSeqsCommand", "createFilter");
839 /**************************************************************************************/
840 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
844 m->openInputFile(filename, in);
846 in.seekg(filePos->start);
853 if (m->control_pressed) { in.close(); return 1; }
855 Sequence seq(in); m->gobble(in);
856 if (seq.getName() != "") {
857 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
859 if(trump != '*') { F.doTrump(seq); }
860 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
865 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
866 unsigned long int pos = in.tellg();
867 if ((pos == -1) || (pos >= filePos->end)) { break; }
869 if (in.eof()) { break; }
873 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
876 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
881 catch(exception& e) {
882 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
887 /**************************************************************************************/
888 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
893 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
895 for(int i=0;i<num;i++){
897 if (m->control_pressed) { return 0; }
900 int length = MPIPos[start+i+1] - MPIPos[start+i];
902 char* buf4 = new char[length];
903 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
905 string tempBuf = buf4;
906 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
907 istringstream iss (tempBuf,istringstream::in);
912 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
914 if(trump != '*'){ F.doTrump(seq); }
915 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
919 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
923 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
927 catch(exception& e) {
928 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
933 /**************************************************************************************************/
935 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
937 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
942 //loop through and create all the processes you want
943 while (process != processors) {
947 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
950 //reset child's filter counts to 0;
951 F.a.clear(); F.a.resize(alignmentLength, 0);
952 F.t.clear(); F.t.resize(alignmentLength, 0);
953 F.g.clear(); F.g.resize(alignmentLength, 0);
954 F.c.clear(); F.c.resize(alignmentLength, 0);
955 F.gap.clear(); F.gap.resize(alignmentLength, 0);
957 num = driverCreateFilter(F, filename, lines[process]);
959 //write out filter counts to file
960 filename += toString(getpid()) + "filterValues.temp";
962 m->openOutputFile(filename, out);
965 out << F.getFilter() << endl;
966 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
967 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
968 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
969 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
970 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
972 //cout << F.getFilter() << endl;
976 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
979 //parent do your part
980 num = driverCreateFilter(F, filename, lines[0]);
982 //force parent to wait until all the processes are done
983 for (int i=0;i<(processors-1);i++) {
984 int temp = processIDS[i];
988 //parent reads in and combines Filter info
989 for (int i = 0; i < processIDS.size(); i++) {
990 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
992 m->openInputFile(tempFilename, in);
995 string tempFilterString;
997 in >> tempNum; m->gobble(in); num += tempNum;
999 in >> tempFilterString;
1000 F.mergeFilter(tempFilterString);
1002 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1003 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1004 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1005 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1006 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1009 remove(tempFilename.c_str());
1015 catch(exception& e) {
1016 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1020 /**************************************************************************************/