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();
135 ableToOpen = m->openInputFile(tryPath, in, "noerror");
136 fastafileNames[i] = tryPath;
140 //if you can't open it, try default location
141 if (ableToOpen == 1) {
142 if (m->getOutputDir() != "") { //default path is set
143 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
144 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
145 ableToOpen = m->openInputFile(tryPath, in, "noerror");
146 fastafileNames[i] = tryPath;
152 if (ableToOpen == 1) {
153 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
154 //erase from file list
155 fastafileNames.erase(fastafileNames.begin()+i);
158 string simpleName = m->getSimpleName(fastafileNames[i]);
159 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
164 //make sure there is at least one valid file left
165 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
169 //if the user changes the output directory command factory will send this info to us in the output parameter
170 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
172 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
175 //check for optional parameter and set defaults
176 // ...at some point should added some additional type checking...
179 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
180 else if (hard == "not open") { hard = ""; abort = true; }
182 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
185 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
186 else { soft = (float)atoi(temp.c_str()) / 100.0; }
188 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
189 convert(temp, processors);
191 vertical = validParameter.validFile(parameters, "vertical", false);
192 if (vertical == "not found") {
193 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
194 else { vertical = "F"; }
201 catch(exception& e) {
202 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
207 //**********************************************************************************************************************
209 void FilterSeqsCommand::help(){
212 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
213 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n");
214 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");
215 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
216 m->mothurOut("The trump parameter .... The default is ...\n");
217 m->mothurOut("The soft parameter .... The default is ....\n");
218 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
219 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
220 m->mothurOut("The processors parameter allows you to specify the number of processors to use. The default is 1.\n");
221 m->mothurOut("The filter.seqs command should be in the following format: \n");
222 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
223 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
224 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
227 catch(exception& e) {
228 m->errorOut(e, "FilterSeqsCommand", "help");
233 /**************************************************************************************/
235 int FilterSeqsCommand::execute() {
238 if (abort == true) { return 0; }
241 m->openInputFile(fastafileNames[0], inFASTA);
243 Sequence testSeq(inFASTA);
244 alignmentLength = testSeq.getAlignLength();
247 ////////////create filter/////////////////
248 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
250 filter = createFilter();
252 m->mothurOutEndLine(); m->mothurOutEndLine();
254 if (m->control_pressed) { outputTypes.clear(); return 0; }
258 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
260 if (pid == 0) { //only one process should output the filter
265 //prevent giantic file name
267 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
268 else { filterFile = outputDir + filterFileName + ".filter"; }
270 m->openOutputFile(filterFile, outFilter);
271 outFilter << filter << endl;
273 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
279 ////////////run filter/////////////////
281 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
285 m->mothurOutEndLine(); m->mothurOutEndLine();
287 int filteredLength = 0;
288 for(int i=0;i<alignmentLength;i++){
289 if(filter[i] == '1'){ filteredLength++; }
292 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
295 m->mothurOutEndLine();
296 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
297 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
298 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
299 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
302 m->mothurOutEndLine();
303 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
304 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
305 m->mothurOutEndLine();
310 catch(exception& e) {
311 m->errorOut(e, "FilterSeqsCommand", "execute");
315 /**************************************************************************************/
316 int FilterSeqsCommand::filterSequences() {
321 for (int s = 0; s < fastafileNames.size(); s++) {
323 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
325 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
327 int pid, start, end, numSeqsPerProcessor, num;
329 vector<unsigned long int>MPIPos;
332 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
333 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
338 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
339 int inMode=MPI_MODE_RDONLY;
341 char outFilename[1024];
342 strcpy(outFilename, filteredFasta.c_str());
344 char inFileName[1024];
345 strcpy(inFileName, fastafileNames[s].c_str());
347 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
348 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
350 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
352 if (pid == 0) { //you are the root process
354 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
357 //send file positions to all processes
358 for(int i = 1; i < processors; i++) {
359 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
360 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
363 //figure out how many sequences you have to do
364 numSeqsPerProcessor = num / processors;
365 int startIndex = pid * numSeqsPerProcessor;
366 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
370 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
372 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
375 for(int i = 1; i < processors; i++) {
377 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
380 }else { //you are a child process
381 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
382 MPIPos.resize(num+1);
384 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
386 //figure out how many sequences you have to align
387 numSeqsPerProcessor = num / processors;
388 int startIndex = pid * numSeqsPerProcessor;
389 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
393 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
395 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
400 //tell parent you are done.
401 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
404 MPI_File_close(&outMPI);
405 MPI_File_close(&inMPI);
406 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
409 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
411 for (int i = 0; i < (positions.size()-1); i++) {
412 lines.push_back(new linePair(positions[i], positions[(i+1)]));
414 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
416 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
417 numSeqs += numFastaSeqs;
419 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
420 numSeqs += numFastaSeqs;
422 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
425 for(int i=1;i<processors;i++){
426 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
427 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
431 if (m->control_pressed) { return 1; }
433 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
434 numSeqs += numFastaSeqs;
436 if (m->control_pressed) { return 1; }
439 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
444 catch(exception& e) {
445 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
450 /**************************************************************************************/
451 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
453 string outputString = "";
457 for(int i=0;i<num;i++){
459 if (m->control_pressed) { return 0; }
462 int length = MPIPos[start+i+1] - MPIPos[start+i];
463 char* buf4 = new char[length];
464 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
466 string tempBuf = buf4;
467 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
468 istringstream iss (tempBuf,istringstream::in);
471 Sequence seq(iss); m->gobble(iss);
473 if (seq.getName() != "") {
474 string align = seq.getAligned();
475 string filterSeq = "";
477 for(int j=0;j<alignmentLength;j++){
478 if(filter[j] == '1'){
479 filterSeq += align[j];
484 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
486 if(count % 10 == 0){ //output to file
487 //send results to parent
488 int length = outputString.length();
489 char* buf = new char[length];
490 memcpy(buf, outputString.c_str(), length);
492 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
499 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
502 if(outputString != ""){ //output to file
503 //send results to parent
504 int length = outputString.length();
505 char* buf = new char[length];
506 memcpy(buf, outputString.c_str(), length);
508 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
513 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
517 catch(exception& e) {
518 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
523 /**************************************************************************************/
524 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
527 m->openOutputFile(outputFilename, out);
530 m->openInputFile(inputFilename, in);
532 in.seekg(filePos->start);
539 if (m->control_pressed) { in.close(); out.close(); return 0; }
541 Sequence seq(in); m->gobble(in);
542 if (seq.getName() != "") {
543 string align = seq.getAligned();
544 string filterSeq = "";
546 for(int j=0;j<alignmentLength;j++){
547 if(filter[j] == '1'){
548 filterSeq += align[j];
552 out << '>' << seq.getName() << endl << filterSeq << endl;
556 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
557 unsigned long int pos = in.tellg();
558 if ((pos == -1) || (pos >= filePos->end)) { break; }
560 if (in.eof()) { break; }
564 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
567 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
575 catch(exception& e) {
576 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
580 /**************************************************************************************************/
582 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
584 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
589 //loop through and create all the processes you want
590 while (process != processors) {
594 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
597 string filteredFasta = filename + toString(getpid()) + ".temp";
598 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
600 //pass numSeqs to parent
602 string tempFile = filename + toString(getpid()) + ".num.temp";
603 m->openOutputFile(tempFile, out);
608 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
611 //force parent to wait until all the processes are done
612 for (int i=0;i<processors;i++) {
613 int temp = processIDS[i];
617 for (int i = 0; i < processIDS.size(); i++) {
619 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
620 m->openInputFile(tempFile, in);
621 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
622 in.close(); remove(tempFile.c_str());
629 catch(exception& e) {
630 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
634 /**************************************************************************************/
635 string FilterSeqsCommand::createFilter() {
637 string filterString = "";
640 if (soft != 0) { F.setSoft(soft); }
641 if (trump != '*') { F.setTrump(trump); }
643 F.setLength(alignmentLength);
645 if(trump != '*' || m->isTrue(vertical) || soft != 0){
649 if(hard.compare("") != 0) { F.doHard(hard); }
650 else { F.setFilter(string(alignmentLength, '1')); }
653 if(trump != '*' || m->isTrue(vertical) || soft != 0){
654 for (int s = 0; s < fastafileNames.size(); s++) {
656 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
659 int pid, numSeqsPerProcessor, num;
661 vector<unsigned long int> MPIPos;
665 MPI_Comm_size(MPI_COMM_WORLD, &processors);
666 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
668 //char* tempFileName = new char(fastafileNames[s].length());
669 //tempFileName = &(fastafileNames[s][0]);
671 char tempFileName[1024];
672 strcpy(tempFileName, fastafileNames[s].c_str());
674 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
676 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
678 if (pid == 0) { //you are the root process
679 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
682 //send file positions to all processes
683 for(int i = 1; i < processors; i++) {
684 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
685 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
688 //figure out how many sequences you have to do
689 numSeqsPerProcessor = num / processors;
690 int startIndex = pid * numSeqsPerProcessor;
691 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
695 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
697 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
699 }else { //i am the child process
700 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
701 MPIPos.resize(num+1);
703 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
705 //figure out how many sequences you have to align
706 numSeqsPerProcessor = num / processors;
707 int startIndex = pid * numSeqsPerProcessor;
708 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
712 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
714 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
717 MPI_File_close(&inMPI);
718 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
721 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
723 for (int i = 0; i < (positions.size()-1); i++) {
724 lines.push_back(new linePair(positions[i], positions[(i+1)]));
726 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
728 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
729 numSeqs += numFastaSeqs;
731 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
732 numSeqs += numFastaSeqs;
735 if (m->control_pressed) { return filterString; }
737 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
738 numSeqs += numFastaSeqs;
739 if (m->control_pressed) { return filterString; }
749 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
752 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
754 if(trump != '*' || m->isTrue(vertical) || soft != 0){
756 if (pid == 0) { //only one process should output the filter
758 vector<int> temp; temp.resize(alignmentLength+1);
760 //get the frequencies from the child processes
761 for(int i = 1; i < processors; i++) {
763 for (int j = 0; j < 5; j++) {
765 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
766 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
768 if (receiveTag == Atag) { //you are recieveing the A frequencies
769 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
770 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
771 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
772 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
773 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
774 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
775 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
776 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
777 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
783 //send my fequency counts
785 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
787 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
789 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
791 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
792 F.gap.push_back(Gaptag);
793 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
798 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
800 if (pid == 0) { //only one process should output the filter
803 F.setNumSeqs(numSeqs);
804 if(m->isTrue(vertical) == 1) { F.doVertical(); }
805 if(soft != 0) { F.doSoft(); }
806 filterString = F.getFilter();
809 //send filter string to kids
810 //for(int i = 1; i < processors; i++) {
811 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
813 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
815 //recieve filterString
816 char* tempBuf = new char[alignmentLength];
817 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
818 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
820 filterString = tempBuf;
821 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
825 MPI_Barrier(MPI_COMM_WORLD);
830 catch(exception& e) {
831 m->errorOut(e, "FilterSeqsCommand", "createFilter");
835 /**************************************************************************************/
836 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
840 m->openInputFile(filename, in);
842 in.seekg(filePos->start);
849 if (m->control_pressed) { in.close(); return 1; }
851 Sequence seq(in); m->gobble(in);
852 if (seq.getName() != "") {
853 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
855 if(trump != '*') { F.doTrump(seq); }
856 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
861 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
862 unsigned long int pos = in.tellg();
863 if ((pos == -1) || (pos >= filePos->end)) { break; }
865 if (in.eof()) { break; }
869 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
872 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
877 catch(exception& e) {
878 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
883 /**************************************************************************************/
884 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
889 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
891 for(int i=0;i<num;i++){
893 if (m->control_pressed) { return 0; }
896 int length = MPIPos[start+i+1] - MPIPos[start+i];
898 char* buf4 = new char[length];
899 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
901 string tempBuf = buf4;
902 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
903 istringstream iss (tempBuf,istringstream::in);
908 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
910 if(trump != '*'){ F.doTrump(seq); }
911 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
915 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
919 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
923 catch(exception& e) {
924 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
929 /**************************************************************************************************/
931 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
933 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
938 //loop through and create all the processes you want
939 while (process != processors) {
943 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
946 //reset child's filter counts to 0;
947 F.a.clear(); F.a.resize(alignmentLength, 0);
948 F.t.clear(); F.t.resize(alignmentLength, 0);
949 F.g.clear(); F.g.resize(alignmentLength, 0);
950 F.c.clear(); F.c.resize(alignmentLength, 0);
951 F.gap.clear(); F.gap.resize(alignmentLength, 0);
953 num = driverCreateFilter(F, filename, lines[process]);
955 //write out filter counts to file
956 filename += toString(getpid()) + "filterValues.temp";
958 m->openOutputFile(filename, out);
961 out << F.getFilter() << endl;
962 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
963 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
964 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
965 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
966 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
968 //cout << F.getFilter() << endl;
972 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
975 //parent do your part
976 num = driverCreateFilter(F, filename, lines[0]);
978 //force parent to wait until all the processes are done
979 for (int i=0;i<(processors-1);i++) {
980 int temp = processIDS[i];
984 //parent reads in and combines Filter info
985 for (int i = 0; i < processIDS.size(); i++) {
986 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
988 m->openInputFile(tempFilename, in);
991 string tempFilterString;
993 in >> tempNum; m->gobble(in); num += tempNum;
995 in >> tempFilterString;
996 F.mergeFilter(tempFilterString);
998 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
999 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1000 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1001 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1002 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1005 remove(tempFilename.c_str());
1011 catch(exception& e) {
1012 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1016 /**************************************************************************************/