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;
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]; }
77 ableToOpen = openInputFile(fastafileNames[i], in);
78 if (ableToOpen == 1) {
79 m->mothurOut(fastafileNames[i] + " will be disregarded."); m->mothurOutEndLine();
80 //erase from file list
81 fastafileNames.erase(fastafileNames.begin()+i);
84 string simpleName = getSimpleName(fastafileNames[i]);
85 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
90 //make sure there is at least one valid file left
91 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
95 //if the user changes the output directory command factory will send this info to us in the output parameter
96 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
98 outputDir += hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
101 //check for optional parameter and set defaults
102 // ...at some point should added some additional type checking...
105 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
106 else if (hard == "not open") { hard = ""; abort = true; }
108 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
111 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
112 else { soft = (float)atoi(temp.c_str()) / 100.0; }
114 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
115 convert(temp, processors);
117 vertical = validParameter.validFile(parameters, "vertical", false);
118 if (vertical == "not found") {
119 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
120 else { vertical = "F"; }
127 catch(exception& e) {
128 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
133 //**********************************************************************************************************************
135 void FilterSeqsCommand::help(){
138 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
139 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard and vertical. \n");
140 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");
141 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
142 m->mothurOut("The trump parameter .... The default is ...\n");
143 m->mothurOut("The soft parameter .... The default is ....\n");
144 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
145 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
146 m->mothurOut("The filter.seqs command should be in the following format: \n");
147 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
148 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
149 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
152 catch(exception& e) {
153 m->errorOut(e, "FilterSeqsCommand", "help");
158 /**************************************************************************************/
160 int FilterSeqsCommand::execute() {
163 if (abort == true) { return 0; }
166 openInputFile(fastafileNames[0], inFASTA);
168 Sequence testSeq(inFASTA);
169 alignmentLength = testSeq.getAlignLength();
172 ////////////create filter/////////////////
173 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
175 filter = createFilter();
177 m->mothurOutEndLine(); m->mothurOutEndLine();
179 if (m->control_pressed) { return 0; }
183 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
185 if (pid == 0) { //only one process should output the filter
190 string filterFile = outputDir + filterFileName + ".filter";
191 openOutputFile(filterFile, outFilter);
192 outFilter << filter << endl;
194 outputNames.push_back(filterFile);
200 ////////////run filter/////////////////
202 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
206 m->mothurOutEndLine(); m->mothurOutEndLine();
208 int filteredLength = 0;
209 for(int i=0;i<alignmentLength;i++){
210 if(filter[i] == '1'){ filteredLength++; }
213 if (m->control_pressed) { for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } return 0; }
216 m->mothurOutEndLine();
217 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
218 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
219 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
220 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
223 m->mothurOutEndLine();
224 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
225 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
226 m->mothurOutEndLine();
231 catch(exception& e) {
232 m->errorOut(e, "FilterSeqsCommand", "execute");
236 /**************************************************************************************/
237 int FilterSeqsCommand::filterSequences() {
242 for (int s = 0; s < fastafileNames.size(); s++) {
244 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
246 string filteredFasta = outputDir + getRootName(getSimpleName(fastafileNames[s])) + "filter.fasta";
248 int pid, start, end, numSeqsPerProcessor, num;
253 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
254 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
259 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
260 int inMode=MPI_MODE_RDONLY;
262 char outFilename[1024];
263 strcpy(outFilename, filteredFasta.c_str());
265 char inFileName[1024];
266 strcpy(inFileName, fastafileNames[s].c_str());
268 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
269 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
271 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
273 if (pid == 0) { //you are the root process
275 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
278 //send file positions to all processes
279 for(int i = 1; i < processors; i++) {
280 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
281 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
284 //figure out how many sequences you have to do
285 numSeqsPerProcessor = num / processors;
286 int startIndex = pid * numSeqsPerProcessor;
287 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
291 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
293 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
296 for(int i = 1; i < processors; i++) {
298 MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
301 }else { //you are a child process
302 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
303 MPIPos.resize(num+1);
305 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
307 //figure out how many sequences you have to align
308 numSeqsPerProcessor = num / processors;
309 int startIndex = pid * numSeqsPerProcessor;
310 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
314 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
316 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
321 //tell parent you are done.
322 MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
325 MPI_File_close(&outMPI);
326 MPI_File_close(&inMPI);
327 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
330 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
334 openInputFile(fastafileNames[s], inFASTA);
335 getNumSeqs(inFASTA, numFastaSeqs);
338 lines.push_back(new linePair(0, numFastaSeqs));
340 numSeqs += numFastaSeqs;
342 driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
344 setLines(fastafileNames[s]);
345 createProcessesRunFilter(filter, fastafileNames[s]);
347 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
350 for(int i=1;i<processors;i++){
351 appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
352 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
356 if (m->control_pressed) { return 1; }
360 openInputFile(fastafileNames[s], inFASTA);
361 getNumSeqs(inFASTA, numFastaSeqs);
364 lines.push_back(new linePair(0, numFastaSeqs));
366 numSeqs += numFastaSeqs;
368 driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
370 if (m->control_pressed) { return 1; }
373 outputNames.push_back(filteredFasta);
378 catch(exception& e) {
379 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
384 /**************************************************************************************/
385 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<long>& MPIPos) {
387 string outputString = "";
391 for(int i=0;i<num;i++){
393 if (m->control_pressed) { return 0; }
396 int length = MPIPos[start+i+1] - MPIPos[start+i];
397 char* buf4 = new char[length];
398 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
400 string tempBuf = buf4;
401 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
402 istringstream iss (tempBuf,istringstream::in);
405 Sequence seq(iss); gobble(iss);
407 if (seq.getName() != "") {
408 string align = seq.getAligned();
409 string filterSeq = "";
411 for(int j=0;j<alignmentLength;j++){
412 if(filter[j] == '1'){
413 filterSeq += align[j];
418 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
420 if(count % 10 == 0){ //output to file
421 //send results to parent
422 int length = outputString.length();
423 char* buf = new char[length];
424 memcpy(buf, outputString.c_str(), length);
426 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
433 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
436 if(outputString != ""){ //output to file
437 //send results to parent
438 int length = outputString.length();
439 char* buf = new char[length];
440 memcpy(buf, outputString.c_str(), length);
442 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
447 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
451 catch(exception& e) {
452 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
457 /**************************************************************************************/
458 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* line) {
461 openOutputFile(outputFilename, out);
464 openInputFile(inputFilename, in);
466 in.seekg(line->start);
468 for(int i=0;i<line->num;i++){
470 if (m->control_pressed) { in.close(); out.close(); return 0; }
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];
483 out << '>' << seq.getName() << endl << filterSeq << endl;
492 catch(exception& e) {
493 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
497 /**************************************************************************************************/
499 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
501 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
506 //loop through and create all the processes you want
507 while (process != processors) {
511 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
514 string filteredFasta = filename + toString(getpid()) + ".temp";
515 driverRunFilter(F, filteredFasta, filename, lines[process]);
517 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
520 //force parent to wait until all the processes are done
521 for (int i=0;i<processors;i++) {
522 int temp = processIDS[i];
529 catch(exception& e) {
530 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
534 /**************************************************************************************/
535 string FilterSeqsCommand::createFilter() {
537 string filterString = "";
540 if (soft != 0) { F.setSoft(soft); }
541 if (trump != '*') { F.setTrump(trump); }
543 F.setLength(alignmentLength);
545 if(trump != '*' || isTrue(vertical) || soft != 0){
549 if(hard.compare("") != 0) { F.doHard(hard); }
550 else { F.setFilter(string(alignmentLength, '1')); }
553 if(trump != '*' || isTrue(vertical) || soft != 0){
554 for (int s = 0; s < fastafileNames.size(); s++) {
556 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
559 int pid, numSeqsPerProcessor, num;
565 MPI_Comm_size(MPI_COMM_WORLD, &processors);
566 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
568 //char* tempFileName = new char(fastafileNames[s].length());
569 //tempFileName = &(fastafileNames[s][0]);
571 char tempFileName[1024];
572 strcpy(tempFileName, fastafileNames[s].c_str());
574 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
576 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
578 if (pid == 0) { //you are the root process
579 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
582 //send file positions to all processes
583 for(int i = 1; i < processors; i++) {
584 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
585 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
588 //figure out how many sequences you have to do
589 numSeqsPerProcessor = num / processors;
590 int startIndex = pid * numSeqsPerProcessor;
591 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
595 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
597 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
599 }else { //i am the child process
600 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
601 MPIPos.resize(num+1);
603 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
605 //figure out how many sequences you have to align
606 numSeqsPerProcessor = num / processors;
607 int startIndex = pid * numSeqsPerProcessor;
608 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
612 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
614 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
617 MPI_File_close(&inMPI);
618 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
621 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
625 openInputFile(fastafileNames[s], inFASTA);
626 getNumSeqs(inFASTA, numFastaSeqs);
629 numSeqs += numFastaSeqs;
631 lines.push_back(new linePair(0, numFastaSeqs));
633 driverCreateFilter(F, fastafileNames[s], lines[0]);
635 setLines(fastafileNames[s]);
636 createProcessesCreateFilter(F, fastafileNames[s]);
639 if (m->control_pressed) { return filterString; }
643 openInputFile(fastafileNames[s], inFASTA);
644 getNumSeqs(inFASTA, numFastaSeqs);
647 numSeqs += numFastaSeqs;
649 lines.push_back(new linePair(0, numFastaSeqs));
651 driverCreateFilter(F, fastafileNames[s], lines[0]);
652 if (m->control_pressed) { return filterString; }
662 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
665 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
667 if(trump != '*' || isTrue(vertical) || soft != 0){
669 if (pid == 0) { //only one process should output the filter
671 vector<int> temp; temp.resize(alignmentLength+1);
673 //get the frequencies from the child processes
674 for(int i = 1; i < processors; i++) {
676 for (int j = 0; j < 5; j++) {
678 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
679 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
681 if (receiveTag == Atag) { //you are recieveing the A frequencies
682 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
683 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
684 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
685 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
686 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
687 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
688 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
689 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
690 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
696 //send my fequency counts
698 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
700 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
702 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
704 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
705 F.gap.push_back(Gaptag);
706 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
711 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
713 if (pid == 0) { //only one process should output the filter
715 F.setNumSeqs(numSeqs);
717 if(isTrue(vertical) == 1) { F.doVertical(); }
718 if(soft != 0) { F.doSoft(); }
720 filterString = F.getFilter();
723 //send filter string to kids
724 //for(int i = 1; i < processors; i++) {
725 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
727 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
729 //recieve filterString
730 char* tempBuf = new char[alignmentLength];
731 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
732 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
734 filterString = tempBuf;
735 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
739 MPI_Barrier(MPI_COMM_WORLD);
744 catch(exception& e) {
745 m->errorOut(e, "FilterSeqsCommand", "createFilter");
749 /**************************************************************************************/
750 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* line) {
754 openInputFile(filename, in);
756 in.seekg(line->start);
758 for(int i=0;i<line->num;i++){
760 if (m->control_pressed) { in.close(); return 1; }
763 if (seq.getName() != "") {
764 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
766 if(trump != '*'){ F.doTrump(seq); }
767 if(isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
772 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
776 if((line->num) % 100 != 0){ m->mothurOut(toString(line->num)); m->mothurOutEndLine(); }
782 catch(exception& e) {
783 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
788 /**************************************************************************************/
789 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<long>& MPIPos) {
794 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
796 for(int i=0;i<num;i++){
798 if (m->control_pressed) { return 0; }
801 int length = MPIPos[start+i+1] - MPIPos[start+i];
803 char* buf4 = new char[length];
804 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
806 string tempBuf = buf4;
807 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
808 istringstream iss (tempBuf,istringstream::in);
813 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
815 if(trump != '*'){ F.doTrump(seq); }
816 if(isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
820 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
824 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
828 catch(exception& e) {
829 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
834 /**************************************************************************************************/
836 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
838 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
843 //loop through and create all the processes you want
844 while (process != processors) {
848 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
851 driverCreateFilter(F, filename, lines[process]);
853 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
856 //force parent to wait until all the processes are done
857 for (int i=0;i<processors;i++) {
858 int temp = processIDS[i];
865 catch(exception& e) {
866 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
870 /**************************************************************************************************/
872 int FilterSeqsCommand::setLines(string filename) {
875 vector<unsigned long int> positions;
879 openInputFile(filename, inFASTA);
882 while(!inFASTA.eof()){
883 input = getline(inFASTA);
885 if (input.length() != 0) {
886 if(input[0] == '>'){ unsigned long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
891 int numFastaSeqs = positions.size();
894 unsigned long int size;
896 //get num bytes in file
897 pFile = fopen (filename.c_str(),"rb");
898 if (pFile==NULL) perror ("Error opening file");
900 fseek (pFile, 0, SEEK_END);
905 numSeqs += numFastaSeqs;
907 int numSeqsPerProcessor = numFastaSeqs / processors;
909 for (int i = 0; i < processors; i++) {
911 unsigned long int startPos = positions[ i * numSeqsPerProcessor ];
912 if(i == processors - 1){
913 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
914 bufferSizes.push_back(size - startPos);
916 unsigned long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
917 bufferSizes.push_back(myEnd-startPos);
919 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
924 catch(exception& e) {
925 m->errorOut(e, "FilterSeqsCommand", "setLines");
929 /**************************************************************************************/