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 = 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 = 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 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 = 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 = 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() + getSimpleName(fastafileNames[i]);
82 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
83 ableToOpen = 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 = 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 += 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 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 string filterFile = outputDir + filterFileName + ".filter";
202 openOutputFile(filterFile, outFilter);
203 outFilter << filter << endl;
205 outputNames.push_back(filterFile);
211 ////////////run filter/////////////////
213 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
217 m->mothurOutEndLine(); m->mothurOutEndLine();
219 int filteredLength = 0;
220 for(int i=0;i<alignmentLength;i++){
221 if(filter[i] == '1'){ filteredLength++; }
224 if (m->control_pressed) { for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
227 m->mothurOutEndLine();
228 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
229 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
230 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
231 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
234 m->mothurOutEndLine();
235 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
236 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
237 m->mothurOutEndLine();
242 catch(exception& e) {
243 m->errorOut(e, "FilterSeqsCommand", "execute");
247 /**************************************************************************************/
248 int FilterSeqsCommand::filterSequences() {
253 for (int s = 0; s < fastafileNames.size(); s++) {
255 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
257 string filteredFasta = outputDir + getRootName(getSimpleName(fastafileNames[s])) + "filter.fasta";
259 int pid, start, end, numSeqsPerProcessor, num;
261 vector<unsigned long int>MPIPos;
264 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
265 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
270 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
271 int inMode=MPI_MODE_RDONLY;
273 char outFilename[1024];
274 strcpy(outFilename, filteredFasta.c_str());
276 char inFileName[1024];
277 strcpy(inFileName, fastafileNames[s].c_str());
279 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
280 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
282 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
284 if (pid == 0) { //you are the root process
286 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
289 //send file positions to all processes
290 for(int i = 1; i < processors; i++) {
291 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
292 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
295 //figure out how many sequences you have to do
296 numSeqsPerProcessor = num / processors;
297 int startIndex = pid * numSeqsPerProcessor;
298 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
302 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
304 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
307 for(int i = 1; i < processors; i++) {
309 MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
312 }else { //you are a child process
313 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
314 MPIPos.resize(num+1);
316 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
318 //figure out how many sequences you have to align
319 numSeqsPerProcessor = num / processors;
320 int startIndex = pid * numSeqsPerProcessor;
321 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
325 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
327 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
332 //tell parent you are done.
333 MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
336 MPI_File_close(&outMPI);
337 MPI_File_close(&inMPI);
338 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
341 vector<unsigned long int> positions = divideFile(fastafileNames[s], processors);
343 for (int i = 0; i < (positions.size()-1); i++) {
344 lines.push_back(new linePair(positions[i], positions[(i+1)]));
346 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
348 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
349 numSeqs += numFastaSeqs;
351 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]);
352 numSeqs += numFastaSeqs;
354 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
357 for(int i=1;i<processors;i++){
358 appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
359 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
363 if (m->control_pressed) { return 1; }
365 numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
366 numSeqs += numFastaSeqs;
368 if (m->control_pressed) { return 1; }
371 outputNames.push_back(filteredFasta);
376 catch(exception& e) {
377 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
382 /**************************************************************************************/
383 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
385 string outputString = "";
389 for(int i=0;i<num;i++){
391 if (m->control_pressed) { return 0; }
394 int length = MPIPos[start+i+1] - MPIPos[start+i];
395 char* buf4 = new char[length];
396 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
398 string tempBuf = buf4;
399 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
400 istringstream iss (tempBuf,istringstream::in);
403 Sequence seq(iss); gobble(iss);
405 if (seq.getName() != "") {
406 string align = seq.getAligned();
407 string filterSeq = "";
409 for(int j=0;j<alignmentLength;j++){
410 if(filter[j] == '1'){
411 filterSeq += align[j];
416 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
418 if(count % 10 == 0){ //output to file
419 //send results to parent
420 int length = outputString.length();
421 char* buf = new char[length];
422 memcpy(buf, outputString.c_str(), length);
424 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
431 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
434 if(outputString != ""){ //output to file
435 //send results to parent
436 int length = outputString.length();
437 char* buf = new char[length];
438 memcpy(buf, outputString.c_str(), length);
440 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
445 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
449 catch(exception& e) {
450 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
455 /**************************************************************************************/
456 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
459 openOutputFile(outputFilename, out);
462 openInputFile(inputFilename, in);
464 in.seekg(filePos->start);
471 if (m->control_pressed) { in.close(); out.close(); return 0; }
473 Sequence seq(in); gobble(in);
474 if (seq.getName() != "") {
475 string align = seq.getAligned();
476 string filterSeq = "";
478 for(int j=0;j<alignmentLength;j++){
479 if(filter[j] == '1'){
480 filterSeq += align[j];
484 out << '>' << seq.getName() << endl << filterSeq << endl;
488 unsigned long int pos = in.tellg();
489 if ((pos == -1) || (pos >= filePos->end)) { break; }
492 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
495 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
503 catch(exception& e) {
504 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
508 /**************************************************************************************************/
510 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
512 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
517 //loop through and create all the processes you want
518 while (process != processors) {
522 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
525 string filteredFasta = filename + toString(getpid()) + ".temp";
526 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
528 //pass numSeqs to parent
530 string tempFile = filename + toString(getpid()) + ".num.temp";
531 openOutputFile(tempFile, out);
536 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
539 //force parent to wait until all the processes are done
540 for (int i=0;i<processors;i++) {
541 int temp = processIDS[i];
545 for (int i = 0; i < processIDS.size(); i++) {
547 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
548 openInputFile(tempFile, in);
549 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
550 in.close(); remove(tempFile.c_str());
557 catch(exception& e) {
558 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
562 /**************************************************************************************/
563 string FilterSeqsCommand::createFilter() {
565 string filterString = "";
568 if (soft != 0) { F.setSoft(soft); }
569 if (trump != '*') { F.setTrump(trump); }
571 F.setLength(alignmentLength);
573 if(trump != '*' || isTrue(vertical) || soft != 0){
577 if(hard.compare("") != 0) { F.doHard(hard); }
578 else { F.setFilter(string(alignmentLength, '1')); }
581 if(trump != '*' || isTrue(vertical) || soft != 0){
582 for (int s = 0; s < fastafileNames.size(); s++) {
584 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
587 int pid, numSeqsPerProcessor, num;
589 vector<unsigned long int> MPIPos;
593 MPI_Comm_size(MPI_COMM_WORLD, &processors);
594 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
596 //char* tempFileName = new char(fastafileNames[s].length());
597 //tempFileName = &(fastafileNames[s][0]);
599 char tempFileName[1024];
600 strcpy(tempFileName, fastafileNames[s].c_str());
602 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
604 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
606 if (pid == 0) { //you are the root process
607 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
610 //send file positions to all processes
611 for(int i = 1; i < processors; i++) {
612 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
613 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
616 //figure out how many sequences you have to do
617 numSeqsPerProcessor = num / processors;
618 int startIndex = pid * numSeqsPerProcessor;
619 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
623 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
625 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
627 }else { //i am the child process
628 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
629 MPIPos.resize(num+1);
631 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
633 //figure out how many sequences you have to align
634 numSeqsPerProcessor = num / processors;
635 int startIndex = pid * numSeqsPerProcessor;
636 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
640 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
642 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
645 MPI_File_close(&inMPI);
646 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
649 vector<unsigned long int> positions = divideFile(fastafileNames[s], processors);
651 for (int i = 0; i < (positions.size()-1); i++) {
652 lines.push_back(new linePair(positions[i], positions[(i+1)]));
654 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
656 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
657 numSeqs += numFastaSeqs;
659 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
660 numSeqs += numFastaSeqs;
663 if (m->control_pressed) { return filterString; }
665 numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
666 numSeqs += numFastaSeqs;
667 if (m->control_pressed) { return filterString; }
677 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
680 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
682 if(trump != '*' || isTrue(vertical) || soft != 0){
684 if (pid == 0) { //only one process should output the filter
686 vector<int> temp; temp.resize(alignmentLength+1);
688 //get the frequencies from the child processes
689 for(int i = 1; i < processors; i++) {
691 for (int j = 0; j < 5; j++) {
693 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
694 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
696 if (receiveTag == Atag) { //you are recieveing the A frequencies
697 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
698 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
699 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
700 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
701 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
702 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
703 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
704 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
705 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
711 //send my fequency counts
713 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
715 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
717 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
719 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
720 F.gap.push_back(Gaptag);
721 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
726 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
728 if (pid == 0) { //only one process should output the filter
730 F.setNumSeqs(numSeqs);
731 if(isTrue(vertical) == 1) { F.doVertical(); }
732 if(soft != 0) { F.doSoft(); }
733 filterString = F.getFilter();
736 //send filter string to kids
737 //for(int i = 1; i < processors; i++) {
738 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
740 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
742 //recieve filterString
743 char* tempBuf = new char[alignmentLength];
744 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
745 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
747 filterString = tempBuf;
748 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
752 MPI_Barrier(MPI_COMM_WORLD);
757 catch(exception& e) {
758 m->errorOut(e, "FilterSeqsCommand", "createFilter");
762 /**************************************************************************************/
763 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
767 openInputFile(filename, in);
769 in.seekg(filePos->start);
776 if (m->control_pressed) { in.close(); return 1; }
778 Sequence seq(in); gobble(in);
779 if (seq.getName() != "") {
780 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
782 if(trump != '*') { F.doTrump(seq); }
783 if(isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
788 unsigned long int pos = in.tellg();
789 if ((pos == -1) || (pos >= filePos->end)) { break; }
792 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
795 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
800 catch(exception& e) {
801 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
806 /**************************************************************************************/
807 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
812 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
814 for(int i=0;i<num;i++){
816 if (m->control_pressed) { return 0; }
819 int length = MPIPos[start+i+1] - MPIPos[start+i];
821 char* buf4 = new char[length];
822 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
824 string tempBuf = buf4;
825 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
826 istringstream iss (tempBuf,istringstream::in);
831 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
833 if(trump != '*'){ F.doTrump(seq); }
834 if(isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
838 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
842 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
846 catch(exception& e) {
847 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
852 /**************************************************************************************************/
854 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
856 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
861 //loop through and create all the processes you want
862 while (process != processors) {
866 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
869 //reset child's filter counts to 0;
870 F.a.clear(); F.a.resize(alignmentLength, 0);
871 F.t.clear(); F.t.resize(alignmentLength, 0);
872 F.g.clear(); F.g.resize(alignmentLength, 0);
873 F.c.clear(); F.c.resize(alignmentLength, 0);
874 F.gap.clear(); F.gap.resize(alignmentLength, 0);
876 num = driverCreateFilter(F, filename, lines[process]);
878 //write out filter counts to file
879 filename += toString(getpid()) + "filterValues.temp";
881 openOutputFile(filename, out);
884 out << F.getFilter() << endl;
885 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
886 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
887 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
888 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
889 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
891 //cout << F.getFilter() << endl;
895 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
898 //force parent to wait until all the processes are done
899 for (int i=0;i<processors;i++) {
900 int temp = processIDS[i];
904 //parent reads in and combine Filter info
905 for (int i = 0; i < processIDS.size(); i++) {
906 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
908 openInputFile(tempFilename, in);
911 string tempFilterString;
913 in >> tempNum; gobble(in); num += tempNum;
915 in >> tempFilterString;
916 F.mergeFilter(tempFilterString);
918 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } gobble(in);
919 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } gobble(in);
920 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } gobble(in);
921 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } gobble(in);
922 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } gobble(in);
925 remove(tempFilename.c_str());
931 catch(exception& e) {
932 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
936 /**************************************************************************************/