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 /**************************************************************************************/
15 FilterSeqsCommand::FilterSeqsCommand(string option) {
20 //allow user to run help
21 if(option == "help") { help(); abort = true; }
24 //valid paramters for this command
25 string Array[] = {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
26 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
28 OptionParser parser(option);
29 map<string,string> parameters = parser.getParameters();
31 ValidParameters validParameter("filter.seqs");
32 map<string,string>::iterator it;
34 //check to make sure all parameters are valid for command
35 for (it = parameters.begin(); it != parameters.end(); it++) {
36 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
39 //if the user changes the input directory command factory will send this info to us in the output parameter
40 string inputDir = validParameter.validFile(parameters, "inputdir", false);
41 if (inputDir == "not found"){ inputDir = ""; }
44 it = parameters.find("fasta");
45 //user has given a template file
46 if(it != parameters.end()){
47 path = m->hasPath(it->second);
48 //if the user has not given a path then, add inputdir. else leave path alone.
49 if (path == "") { parameters["fasta"] = inputDir + it->second; }
52 it = parameters.find("hard");
53 //user has given a template file
54 if(it != parameters.end()){
55 path = m->hasPath(it->second);
56 //if the user has not given a path then, add inputdir. else leave path alone.
57 if (path == "") { parameters["hard"] = inputDir + it->second; }
61 //check for required parameters
62 fasta = validParameter.validFile(parameters, "fasta", false);
63 if (fasta == "not found") { m->mothurOut("fasta is a required parameter for the filter.seqs command."); m->mothurOutEndLine(); abort = true; }
65 m->splitAtDash(fasta, fastafileNames);
67 //go through files and make sure they are good, if not, then disregard them
68 for (int i = 0; i < fastafileNames.size(); i++) {
70 string path = m->hasPath(fastafileNames[i]);
71 //if the user has not given a path then, add inputdir. else leave path alone.
72 if (path == "") { fastafileNames[i] = inputDir + fastafileNames[i]; }
76 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
78 //if you can't open it, try default location
79 if (ableToOpen == 1) {
80 if (m->getDefaultPath() != "") { //default path is set
81 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
82 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
83 ableToOpen = m->openInputFile(tryPath, in, "noerror");
84 fastafileNames[i] = tryPath;
89 if (ableToOpen == 1) {
90 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
91 //erase from file list
92 fastafileNames.erase(fastafileNames.begin()+i);
95 string simpleName = m->getSimpleName(fastafileNames[i]);
96 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
101 //make sure there is at least one valid file left
102 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
106 //if the user changes the output directory command factory will send this info to us in the output parameter
107 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
109 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
112 //check for optional parameter and set defaults
113 // ...at some point should added some additional type checking...
116 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
117 else if (hard == "not open") { hard = ""; abort = true; }
119 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
122 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
123 else { soft = (float)atoi(temp.c_str()) / 100.0; }
125 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
126 convert(temp, processors);
128 vertical = validParameter.validFile(parameters, "vertical", false);
129 if (vertical == "not found") {
130 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
131 else { vertical = "F"; }
138 catch(exception& e) {
139 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
144 //**********************************************************************************************************************
146 void FilterSeqsCommand::help(){
149 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
150 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard and vertical. \n");
151 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");
152 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
153 m->mothurOut("The trump parameter .... The default is ...\n");
154 m->mothurOut("The soft parameter .... The default is ....\n");
155 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
156 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
157 m->mothurOut("The filter.seqs command should be in the following format: \n");
158 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
159 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
160 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
163 catch(exception& e) {
164 m->errorOut(e, "FilterSeqsCommand", "help");
169 /**************************************************************************************/
171 int FilterSeqsCommand::execute() {
174 if (abort == true) { return 0; }
177 m->openInputFile(fastafileNames[0], inFASTA);
179 Sequence testSeq(inFASTA);
180 alignmentLength = testSeq.getAlignLength();
183 ////////////create filter/////////////////
184 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
186 filter = createFilter();
188 m->mothurOutEndLine(); m->mothurOutEndLine();
190 if (m->control_pressed) { return 0; }
194 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
196 if (pid == 0) { //only one process should output the filter
201 //prevent giantic file name
203 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
204 else { filterFile = outputDir + filterFileName + ".filter"; }
206 m->openOutputFile(filterFile, outFilter);
207 outFilter << filter << endl;
209 outputNames.push_back(filterFile);
215 ////////////run filter/////////////////
217 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
221 m->mothurOutEndLine(); m->mothurOutEndLine();
223 int filteredLength = 0;
224 for(int i=0;i<alignmentLength;i++){
225 if(filter[i] == '1'){ filteredLength++; }
228 if (m->control_pressed) { for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
231 m->mothurOutEndLine();
232 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
233 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
234 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
235 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
238 m->mothurOutEndLine();
239 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
240 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
241 m->mothurOutEndLine();
246 catch(exception& e) {
247 m->errorOut(e, "FilterSeqsCommand", "execute");
251 /**************************************************************************************/
252 int FilterSeqsCommand::filterSequences() {
257 for (int s = 0; s < fastafileNames.size(); s++) {
259 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
261 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
263 int pid, start, end, numSeqsPerProcessor, num;
265 vector<unsigned long int>MPIPos;
268 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
269 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
274 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
275 int inMode=MPI_MODE_RDONLY;
277 char outFilename[1024];
278 strcpy(outFilename, filteredFasta.c_str());
280 char inFileName[1024];
281 strcpy(inFileName, fastafileNames[s].c_str());
283 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
284 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
286 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
288 if (pid == 0) { //you are the root process
290 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
293 //send file positions to all processes
294 for(int i = 1; i < processors; i++) {
295 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
296 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
299 //figure out how many sequences you have to do
300 numSeqsPerProcessor = num / processors;
301 int startIndex = pid * numSeqsPerProcessor;
302 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
306 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
308 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
311 for(int i = 1; i < processors; i++) {
313 MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
316 }else { //you are a child process
317 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
318 MPIPos.resize(num+1);
320 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
322 //figure out how many sequences you have to align
323 numSeqsPerProcessor = num / processors;
324 int startIndex = pid * numSeqsPerProcessor;
325 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
329 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
331 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
336 //tell parent you are done.
337 MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
340 MPI_File_close(&outMPI);
341 MPI_File_close(&inMPI);
342 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
345 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
347 for (int i = 0; i < (positions.size()-1); i++) {
348 lines.push_back(new linePair(positions[i], positions[(i+1)]));
350 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
352 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
353 numSeqs += numFastaSeqs;
355 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
356 numSeqs += numFastaSeqs;
358 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
361 for(int i=1;i<processors;i++){
362 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
363 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
367 if (m->control_pressed) { return 1; }
369 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
370 numSeqs += numFastaSeqs;
372 if (m->control_pressed) { return 1; }
375 outputNames.push_back(filteredFasta);
380 catch(exception& e) {
381 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
386 /**************************************************************************************/
387 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
389 string outputString = "";
393 for(int i=0;i<num;i++){
395 if (m->control_pressed) { return 0; }
398 int length = MPIPos[start+i+1] - MPIPos[start+i];
399 char* buf4 = new char[length];
400 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
402 string tempBuf = buf4;
403 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
404 istringstream iss (tempBuf,istringstream::in);
407 Sequence seq(iss); m->gobble(iss);
409 if (seq.getName() != "") {
410 string align = seq.getAligned();
411 string filterSeq = "";
413 for(int j=0;j<alignmentLength;j++){
414 if(filter[j] == '1'){
415 filterSeq += align[j];
420 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
422 if(count % 10 == 0){ //output to file
423 //send results to parent
424 int length = outputString.length();
425 char* buf = new char[length];
426 memcpy(buf, outputString.c_str(), length);
428 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
435 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
438 if(outputString != ""){ //output to file
439 //send results to parent
440 int length = outputString.length();
441 char* buf = new char[length];
442 memcpy(buf, outputString.c_str(), length);
444 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
449 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
453 catch(exception& e) {
454 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
459 /**************************************************************************************/
460 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
463 m->openOutputFile(outputFilename, out);
466 m->openInputFile(inputFilename, in);
468 in.seekg(filePos->start);
475 if (m->control_pressed) { in.close(); out.close(); return 0; }
477 Sequence seq(in); m->gobble(in);
478 if (seq.getName() != "") {
479 string align = seq.getAligned();
480 string filterSeq = "";
482 for(int j=0;j<alignmentLength;j++){
483 if(filter[j] == '1'){
484 filterSeq += align[j];
488 out << '>' << seq.getName() << endl << filterSeq << endl;
492 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
493 unsigned long int pos = in.tellg();
494 if ((pos == -1) || (pos >= filePos->end)) { break; }
496 if (in.eof()) { break; }
500 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
503 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
511 catch(exception& e) {
512 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
516 /**************************************************************************************************/
518 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
520 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
525 //loop through and create all the processes you want
526 while (process != processors) {
530 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
533 string filteredFasta = filename + toString(getpid()) + ".temp";
534 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
536 //pass numSeqs to parent
538 string tempFile = filename + toString(getpid()) + ".num.temp";
539 m->openOutputFile(tempFile, out);
544 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
547 //force parent to wait until all the processes are done
548 for (int i=0;i<processors;i++) {
549 int temp = processIDS[i];
553 for (int i = 0; i < processIDS.size(); i++) {
555 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
556 m->openInputFile(tempFile, in);
557 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
558 in.close(); remove(tempFile.c_str());
565 catch(exception& e) {
566 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
570 /**************************************************************************************/
571 string FilterSeqsCommand::createFilter() {
573 string filterString = "";
576 if (soft != 0) { F.setSoft(soft); }
577 if (trump != '*') { F.setTrump(trump); }
579 F.setLength(alignmentLength);
581 if(trump != '*' || m->isTrue(vertical) || soft != 0){
585 if(hard.compare("") != 0) { F.doHard(hard); }
586 else { F.setFilter(string(alignmentLength, '1')); }
589 if(trump != '*' || m->isTrue(vertical) || soft != 0){
590 for (int s = 0; s < fastafileNames.size(); s++) {
592 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
595 int pid, numSeqsPerProcessor, num;
597 vector<unsigned long int> MPIPos;
601 MPI_Comm_size(MPI_COMM_WORLD, &processors);
602 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
604 //char* tempFileName = new char(fastafileNames[s].length());
605 //tempFileName = &(fastafileNames[s][0]);
607 char tempFileName[1024];
608 strcpy(tempFileName, fastafileNames[s].c_str());
610 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
612 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
614 if (pid == 0) { //you are the root process
615 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
618 //send file positions to all processes
619 for(int i = 1; i < processors; i++) {
620 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
621 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
624 //figure out how many sequences you have to do
625 numSeqsPerProcessor = num / processors;
626 int startIndex = pid * numSeqsPerProcessor;
627 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
631 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
633 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
635 }else { //i am the child process
636 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
637 MPIPos.resize(num+1);
639 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
641 //figure out how many sequences you have to align
642 numSeqsPerProcessor = num / processors;
643 int startIndex = pid * numSeqsPerProcessor;
644 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
648 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
650 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
653 MPI_File_close(&inMPI);
654 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
657 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
659 for (int i = 0; i < (positions.size()-1); i++) {
660 lines.push_back(new linePair(positions[i], positions[(i+1)]));
662 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
664 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
665 numSeqs += numFastaSeqs;
667 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
668 numSeqs += numFastaSeqs;
671 if (m->control_pressed) { return filterString; }
673 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
674 numSeqs += numFastaSeqs;
675 if (m->control_pressed) { return filterString; }
685 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
688 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
690 if(trump != '*' || m->isTrue(vertical) || soft != 0){
692 if (pid == 0) { //only one process should output the filter
694 vector<int> temp; temp.resize(alignmentLength+1);
696 //get the frequencies from the child processes
697 for(int i = 1; i < processors; i++) {
699 for (int j = 0; j < 5; j++) {
701 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
702 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
704 if (receiveTag == Atag) { //you are recieveing the A frequencies
705 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
706 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
707 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
708 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
709 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
710 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
711 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
712 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
713 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
719 //send my fequency counts
721 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
723 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
725 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
727 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
728 F.gap.push_back(Gaptag);
729 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
734 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
736 if (pid == 0) { //only one process should output the filter
739 F.setNumSeqs(numSeqs);
740 if(m->isTrue(vertical) == 1) { F.doVertical(); }
741 if(soft != 0) { F.doSoft(); }
742 filterString = F.getFilter();
745 //send filter string to kids
746 //for(int i = 1; i < processors; i++) {
747 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
749 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
751 //recieve filterString
752 char* tempBuf = new char[alignmentLength];
753 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
754 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
756 filterString = tempBuf;
757 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
761 MPI_Barrier(MPI_COMM_WORLD);
766 catch(exception& e) {
767 m->errorOut(e, "FilterSeqsCommand", "createFilter");
771 /**************************************************************************************/
772 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
776 m->openInputFile(filename, in);
778 in.seekg(filePos->start);
785 if (m->control_pressed) { in.close(); return 1; }
787 Sequence seq(in); m->gobble(in);
788 if (seq.getName() != "") {
789 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
791 if(trump != '*') { F.doTrump(seq); }
792 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
797 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
798 unsigned long int pos = in.tellg();
799 if ((pos == -1) || (pos >= filePos->end)) { break; }
801 if (in.eof()) { break; }
805 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
808 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
813 catch(exception& e) {
814 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
819 /**************************************************************************************/
820 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
825 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
827 for(int i=0;i<num;i++){
829 if (m->control_pressed) { return 0; }
832 int length = MPIPos[start+i+1] - MPIPos[start+i];
834 char* buf4 = new char[length];
835 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
837 string tempBuf = buf4;
838 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
839 istringstream iss (tempBuf,istringstream::in);
844 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
846 if(trump != '*'){ F.doTrump(seq); }
847 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
851 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
855 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
859 catch(exception& e) {
860 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
865 /**************************************************************************************************/
867 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
869 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
874 //loop through and create all the processes you want
875 while (process != processors) {
879 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
882 //reset child's filter counts to 0;
883 F.a.clear(); F.a.resize(alignmentLength, 0);
884 F.t.clear(); F.t.resize(alignmentLength, 0);
885 F.g.clear(); F.g.resize(alignmentLength, 0);
886 F.c.clear(); F.c.resize(alignmentLength, 0);
887 F.gap.clear(); F.gap.resize(alignmentLength, 0);
889 num = driverCreateFilter(F, filename, lines[process]);
891 //write out filter counts to file
892 filename += toString(getpid()) + "filterValues.temp";
894 m->openOutputFile(filename, out);
897 out << F.getFilter() << endl;
898 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
899 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
900 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
901 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
902 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
904 //cout << F.getFilter() << endl;
908 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
911 //force parent to wait until all the processes are done
912 for (int i=0;i<processors;i++) {
913 int temp = processIDS[i];
917 //parent reads in and combine Filter info
918 for (int i = 0; i < processIDS.size(); i++) {
919 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
921 m->openInputFile(tempFilename, in);
924 string tempFilterString;
926 in >> tempNum; m->gobble(in); num += tempNum;
928 in >> tempFilterString;
929 F.mergeFilter(tempFilterString);
931 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
932 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
933 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
934 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
935 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
938 remove(tempFilename.c_str());
944 catch(exception& e) {
945 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
949 /**************************************************************************************/