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(){
29 //initialize outputTypes
30 vector<string> tempOutNames;
31 outputTypes["fasta"] = tempOutNames;
32 outputTypes["filter"] = tempOutNames;
35 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
39 //**********************************************************************************************************************
40 vector<string> FilterSeqsCommand::getRequiredParameters(){
42 string Array[] = {"fasta"};
43 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
47 m->errorOut(e, "FilterSeqsCommand", "getRequiredParameters");
51 //**********************************************************************************************************************
52 vector<string> FilterSeqsCommand::getRequiredFiles(){
54 vector<string> myArray;
58 m->errorOut(e, "FilterSeqsCommand", "getRequiredFiles");
62 /**************************************************************************************/
63 FilterSeqsCommand::FilterSeqsCommand(string option) {
68 //allow user to run help
69 if(option == "help") { help(); abort = true; }
72 //valid paramters for this command
73 string Array[] = {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
74 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
76 OptionParser parser(option);
77 map<string,string> parameters = parser.getParameters();
79 ValidParameters validParameter("filter.seqs");
80 map<string,string>::iterator it;
82 //check to make sure all parameters are valid for command
83 for (it = parameters.begin(); it != parameters.end(); it++) {
84 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
87 //initialize outputTypes
88 vector<string> tempOutNames;
89 outputTypes["fasta"] = tempOutNames;
90 outputTypes["filter"] = tempOutNames;
92 //if the user changes the input directory command factory will send this info to us in the output parameter
93 string inputDir = validParameter.validFile(parameters, "inputdir", false);
94 if (inputDir == "not found"){ inputDir = ""; }
97 it = parameters.find("fasta");
98 //user has given a template file
99 if(it != parameters.end()){
100 path = m->hasPath(it->second);
101 //if the user has not given a path then, add inputdir. else leave path alone.
102 if (path == "") { parameters["fasta"] = inputDir + it->second; }
105 it = parameters.find("hard");
106 //user has given a template file
107 if(it != parameters.end()){
108 path = m->hasPath(it->second);
109 //if the user has not given a path then, add inputdir. else leave path alone.
110 if (path == "") { parameters["hard"] = inputDir + it->second; }
114 //check for required parameters
115 fasta = validParameter.validFile(parameters, "fasta", false);
116 if (fasta == "not found") { m->mothurOut("fasta is a required parameter for the filter.seqs command."); m->mothurOutEndLine(); abort = true; }
118 m->splitAtDash(fasta, fastafileNames);
120 //go through files and make sure they are good, if not, then disregard them
121 for (int i = 0; i < fastafileNames.size(); i++) {
122 if (inputDir != "") {
123 string path = m->hasPath(fastafileNames[i]);
124 //if the user has not given a path then, add inputdir. else leave path alone.
125 if (path == "") { fastafileNames[i] = inputDir + fastafileNames[i]; }
129 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
131 //if you can't open it, try default location
132 if (ableToOpen == 1) {
133 if (m->getDefaultPath() != "") { //default path is set
134 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
135 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
137 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
139 fastafileNames[i] = tryPath;
143 //if you can't open it, try default location
144 if (ableToOpen == 1) {
145 if (m->getOutputDir() != "") { //default path is set
146 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
147 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
149 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
151 fastafileNames[i] = tryPath;
157 if (ableToOpen == 1) {
158 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
159 //erase from file list
160 fastafileNames.erase(fastafileNames.begin()+i);
163 string simpleName = m->getSimpleName(fastafileNames[i]);
164 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
169 //make sure there is at least one valid file left
170 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
174 //if the user changes the output directory command factory will send this info to us in the output parameter
175 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
177 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
180 //check for optional parameter and set defaults
181 // ...at some point should added some additional type checking...
184 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
185 else if (hard == "not open") { hard = ""; abort = true; }
187 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
190 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
191 else { soft = (float)atoi(temp.c_str()) / 100.0; }
193 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
194 convert(temp, processors);
196 vertical = validParameter.validFile(parameters, "vertical", false);
197 if (vertical == "not found") {
198 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
199 else { vertical = "F"; }
206 catch(exception& e) {
207 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
212 //**********************************************************************************************************************
214 void FilterSeqsCommand::help(){
217 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
218 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n");
219 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");
220 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
221 m->mothurOut("The trump parameter .... The default is ...\n");
222 m->mothurOut("The soft parameter .... The default is ....\n");
223 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
224 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
225 m->mothurOut("The processors parameter allows you to specify the number of processors to use. The default is 1.\n");
226 m->mothurOut("The filter.seqs command should be in the following format: \n");
227 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
228 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
229 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
232 catch(exception& e) {
233 m->errorOut(e, "FilterSeqsCommand", "help");
238 /**************************************************************************************/
240 int FilterSeqsCommand::execute() {
243 if (abort == true) { return 0; }
246 m->openInputFile(fastafileNames[0], inFASTA);
248 Sequence testSeq(inFASTA);
249 alignmentLength = testSeq.getAlignLength();
252 ////////////create filter/////////////////
253 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
255 filter = createFilter();
257 m->mothurOutEndLine(); m->mothurOutEndLine();
259 if (m->control_pressed) { outputTypes.clear(); return 0; }
263 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
265 if (pid == 0) { //only one process should output the filter
270 //prevent giantic file name
272 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
273 else { filterFile = outputDir + filterFileName + ".filter"; }
275 m->openOutputFile(filterFile, outFilter);
276 outFilter << filter << endl;
278 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
284 ////////////run filter/////////////////
286 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
290 m->mothurOutEndLine(); m->mothurOutEndLine();
292 int filteredLength = 0;
293 for(int i=0;i<alignmentLength;i++){
294 if(filter[i] == '1'){ filteredLength++; }
297 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
300 m->mothurOutEndLine();
301 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
302 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
303 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
304 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
307 m->mothurOutEndLine();
308 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
309 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
310 m->mothurOutEndLine();
315 catch(exception& e) {
316 m->errorOut(e, "FilterSeqsCommand", "execute");
320 /**************************************************************************************/
321 int FilterSeqsCommand::filterSequences() {
326 for (int s = 0; s < fastafileNames.size(); s++) {
328 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
330 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
332 int pid, start, end, numSeqsPerProcessor, num;
334 vector<unsigned long int>MPIPos;
337 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
338 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
343 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
344 int inMode=MPI_MODE_RDONLY;
346 char outFilename[1024];
347 strcpy(outFilename, filteredFasta.c_str());
349 char inFileName[1024];
350 strcpy(inFileName, fastafileNames[s].c_str());
352 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
353 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
355 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
357 if (pid == 0) { //you are the root process
359 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
362 //send file positions to all processes
363 for(int i = 1; i < processors; i++) {
364 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
365 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
368 //figure out how many sequences you have to do
369 numSeqsPerProcessor = num / processors;
370 int startIndex = pid * numSeqsPerProcessor;
371 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
375 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
377 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
380 for(int i = 1; i < processors; i++) {
382 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
385 }else { //you are a child process
386 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
387 MPIPos.resize(num+1);
389 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
391 //figure out how many sequences you have to align
392 numSeqsPerProcessor = num / processors;
393 int startIndex = pid * numSeqsPerProcessor;
394 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
398 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
400 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
405 //tell parent you are done.
406 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
409 MPI_File_close(&outMPI);
410 MPI_File_close(&inMPI);
411 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
414 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
416 for (int i = 0; i < (positions.size()-1); i++) {
417 lines.push_back(new linePair(positions[i], positions[(i+1)]));
419 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
421 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
422 numSeqs += numFastaSeqs;
424 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
425 numSeqs += numFastaSeqs;
427 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
430 for(int i=1;i<processors;i++){
431 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
432 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
436 if (m->control_pressed) { return 1; }
438 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
439 numSeqs += numFastaSeqs;
441 if (m->control_pressed) { return 1; }
444 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
449 catch(exception& e) {
450 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
455 /**************************************************************************************/
456 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
458 string outputString = "";
462 for(int i=0;i<num;i++){
464 if (m->control_pressed) { return 0; }
467 int length = MPIPos[start+i+1] - MPIPos[start+i];
468 char* buf4 = new char[length];
469 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
471 string tempBuf = buf4;
472 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
473 istringstream iss (tempBuf,istringstream::in);
476 Sequence seq(iss); m->gobble(iss);
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];
489 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
491 if(count % 10 == 0){ //output to file
492 //send results to parent
493 int length = outputString.length();
494 char* buf = new char[length];
495 memcpy(buf, outputString.c_str(), length);
497 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
504 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
507 if(outputString != ""){ //output to file
508 //send results to parent
509 int length = outputString.length();
510 char* buf = new char[length];
511 memcpy(buf, outputString.c_str(), length);
513 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
518 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
522 catch(exception& e) {
523 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
528 /**************************************************************************************/
529 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
532 m->openOutputFile(outputFilename, out);
535 m->openInputFile(inputFilename, in);
537 in.seekg(filePos->start);
544 if (m->control_pressed) { in.close(); out.close(); return 0; }
546 Sequence seq(in); m->gobble(in);
547 if (seq.getName() != "") {
548 string align = seq.getAligned();
549 string filterSeq = "";
551 for(int j=0;j<alignmentLength;j++){
552 if(filter[j] == '1'){
553 filterSeq += align[j];
557 out << '>' << seq.getName() << endl << filterSeq << endl;
561 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
562 unsigned long int pos = in.tellg();
563 if ((pos == -1) || (pos >= filePos->end)) { break; }
565 if (in.eof()) { break; }
569 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
572 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
580 catch(exception& e) {
581 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
585 /**************************************************************************************************/
587 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
589 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
594 //loop through and create all the processes you want
595 while (process != processors) {
599 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
602 string filteredFasta = filename + toString(getpid()) + ".temp";
603 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
605 //pass numSeqs to parent
607 string tempFile = filename + toString(getpid()) + ".num.temp";
608 m->openOutputFile(tempFile, out);
613 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
616 //force parent to wait until all the processes are done
617 for (int i=0;i<processors;i++) {
618 int temp = processIDS[i];
622 for (int i = 0; i < processIDS.size(); i++) {
624 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
625 m->openInputFile(tempFile, in);
626 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
627 in.close(); remove(tempFile.c_str());
634 catch(exception& e) {
635 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
639 /**************************************************************************************/
640 string FilterSeqsCommand::createFilter() {
642 string filterString = "";
645 if (soft != 0) { F.setSoft(soft); }
646 if (trump != '*') { F.setTrump(trump); }
648 F.setLength(alignmentLength);
650 if(trump != '*' || m->isTrue(vertical) || soft != 0){
654 if(hard.compare("") != 0) { F.doHard(hard); }
655 else { F.setFilter(string(alignmentLength, '1')); }
658 if(trump != '*' || m->isTrue(vertical) || soft != 0){
659 for (int s = 0; s < fastafileNames.size(); s++) {
661 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
664 int pid, numSeqsPerProcessor, num;
666 vector<unsigned long int> MPIPos;
670 MPI_Comm_size(MPI_COMM_WORLD, &processors);
671 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
673 //char* tempFileName = new char(fastafileNames[s].length());
674 //tempFileName = &(fastafileNames[s][0]);
676 char tempFileName[1024];
677 strcpy(tempFileName, fastafileNames[s].c_str());
679 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
681 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
683 if (pid == 0) { //you are the root process
684 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
687 //send file positions to all processes
688 for(int i = 1; i < processors; i++) {
689 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
690 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
693 //figure out how many sequences you have to do
694 numSeqsPerProcessor = num / processors;
695 int startIndex = pid * numSeqsPerProcessor;
696 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
700 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
702 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
704 }else { //i am the child process
705 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
706 MPIPos.resize(num+1);
708 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
710 //figure out how many sequences you have to align
711 numSeqsPerProcessor = num / processors;
712 int startIndex = pid * numSeqsPerProcessor;
713 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
717 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
719 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
722 MPI_File_close(&inMPI);
723 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
726 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
728 for (int i = 0; i < (positions.size()-1); i++) {
729 lines.push_back(new linePair(positions[i], positions[(i+1)]));
731 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
733 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
734 numSeqs += numFastaSeqs;
736 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
737 numSeqs += numFastaSeqs;
740 if (m->control_pressed) { return filterString; }
742 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
743 numSeqs += numFastaSeqs;
744 if (m->control_pressed) { return filterString; }
754 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
757 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
759 if(trump != '*' || m->isTrue(vertical) || soft != 0){
761 if (pid == 0) { //only one process should output the filter
763 vector<int> temp; temp.resize(alignmentLength+1);
765 //get the frequencies from the child processes
766 for(int i = 1; i < processors; i++) {
768 for (int j = 0; j < 5; j++) {
770 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
771 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
773 if (receiveTag == Atag) { //you are recieveing the A frequencies
774 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
775 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
776 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
777 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
778 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
779 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
780 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
781 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
782 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
788 //send my fequency counts
790 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
792 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
794 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
796 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
797 F.gap.push_back(Gaptag);
798 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
803 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
805 if (pid == 0) { //only one process should output the filter
808 F.setNumSeqs(numSeqs);
809 if(m->isTrue(vertical) == 1) { F.doVertical(); }
810 if(soft != 0) { F.doSoft(); }
811 filterString = F.getFilter();
814 //send filter string to kids
815 //for(int i = 1; i < processors; i++) {
816 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
818 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
820 //recieve filterString
821 char* tempBuf = new char[alignmentLength];
822 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
823 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
825 filterString = tempBuf;
826 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
830 MPI_Barrier(MPI_COMM_WORLD);
835 catch(exception& e) {
836 m->errorOut(e, "FilterSeqsCommand", "createFilter");
840 /**************************************************************************************/
841 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
845 m->openInputFile(filename, in);
847 in.seekg(filePos->start);
854 if (m->control_pressed) { in.close(); return 1; }
856 Sequence seq(in); m->gobble(in);
857 if (seq.getName() != "") {
858 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
860 if(trump != '*') { F.doTrump(seq); }
861 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
866 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
867 unsigned long int pos = in.tellg();
868 if ((pos == -1) || (pos >= filePos->end)) { break; }
870 if (in.eof()) { break; }
874 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
877 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
882 catch(exception& e) {
883 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
888 /**************************************************************************************/
889 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
894 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
896 for(int i=0;i<num;i++){
898 if (m->control_pressed) { return 0; }
901 int length = MPIPos[start+i+1] - MPIPos[start+i];
903 char* buf4 = new char[length];
904 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
906 string tempBuf = buf4;
907 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
908 istringstream iss (tempBuf,istringstream::in);
913 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
915 if(trump != '*'){ F.doTrump(seq); }
916 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
920 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
924 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
928 catch(exception& e) {
929 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
934 /**************************************************************************************************/
936 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
938 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
943 //loop through and create all the processes you want
944 while (process != processors) {
948 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
951 //reset child's filter counts to 0;
952 F.a.clear(); F.a.resize(alignmentLength, 0);
953 F.t.clear(); F.t.resize(alignmentLength, 0);
954 F.g.clear(); F.g.resize(alignmentLength, 0);
955 F.c.clear(); F.c.resize(alignmentLength, 0);
956 F.gap.clear(); F.gap.resize(alignmentLength, 0);
958 num = driverCreateFilter(F, filename, lines[process]);
960 //write out filter counts to file
961 filename += toString(getpid()) + "filterValues.temp";
963 m->openOutputFile(filename, out);
966 out << F.getFilter() << endl;
967 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
968 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
969 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
970 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
971 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
973 //cout << F.getFilter() << endl;
977 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
980 //parent do your part
981 num = driverCreateFilter(F, filename, lines[0]);
983 //force parent to wait until all the processes are done
984 for (int i=0;i<(processors-1);i++) {
985 int temp = processIDS[i];
989 //parent reads in and combines Filter info
990 for (int i = 0; i < processIDS.size(); i++) {
991 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
993 m->openInputFile(tempFilename, in);
996 string tempFilterString;
998 in >> tempNum; m->gobble(in); num += tempNum;
1000 in >> tempFilterString;
1001 F.mergeFilter(tempFilterString);
1003 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1004 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1005 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1006 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1007 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1010 remove(tempFilename.c_str());
1016 catch(exception& e) {
1017 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1021 /**************************************************************************************/