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
266 cout << pid << "is in create filter " << endl;
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 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
345 openInputFile(fastafileNames[s], inFASTA);
346 getNumSeqs(inFASTA, numFastaSeqs);
349 lines.push_back(new linePair(0, numFastaSeqs));
351 numSeqs += numFastaSeqs;
353 driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
355 setLines(fastafileNames[s]);
356 createProcessesRunFilter(filter, fastafileNames[s]);
358 rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
361 for(int i=1;i<processors;i++){
362 appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
363 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
367 if (m->control_pressed) { return 1; }
371 openInputFile(fastafileNames[s], inFASTA);
372 getNumSeqs(inFASTA, numFastaSeqs);
375 lines.push_back(new linePair(0, numFastaSeqs));
377 numSeqs += numFastaSeqs;
379 driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
381 if (m->control_pressed) { return 1; }
384 outputNames.push_back(filteredFasta);
389 catch(exception& e) {
390 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
395 /**************************************************************************************/
396 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
398 string outputString = "";
402 for(int i=0;i<num;i++){
404 if (m->control_pressed) { return 0; }
407 int length = MPIPos[start+i+1] - MPIPos[start+i];
408 char* buf4 = new char[length];
409 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
411 string tempBuf = buf4;
412 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
413 istringstream iss (tempBuf,istringstream::in);
416 Sequence seq(iss); gobble(iss);
418 if (seq.getName() != "") {
419 string align = seq.getAligned();
420 string filterSeq = "";
422 for(int j=0;j<alignmentLength;j++){
423 if(filter[j] == '1'){
424 filterSeq += align[j];
429 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
431 if(count % 10 == 0){ //output to file
432 //send results to parent
433 int length = outputString.length();
434 char* buf = new char[length];
435 memcpy(buf, outputString.c_str(), length);
437 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
444 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
447 if(outputString != ""){ //output to file
448 //send results to parent
449 int length = outputString.length();
450 char* buf = new char[length];
451 memcpy(buf, outputString.c_str(), length);
453 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
458 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
462 catch(exception& e) {
463 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
468 /**************************************************************************************/
469 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* line) {
472 openOutputFile(outputFilename, out);
475 openInputFile(inputFilename, in);
477 in.seekg(line->start);
479 for(int i=0;i<line->num;i++){
481 if (m->control_pressed) { in.close(); out.close(); return 0; }
484 if (seq.getName() != "") {
485 string align = seq.getAligned();
486 string filterSeq = "";
488 for(int j=0;j<alignmentLength;j++){
489 if(filter[j] == '1'){
490 filterSeq += align[j];
494 out << '>' << seq.getName() << endl << filterSeq << endl;
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 driverRunFilter(F, filteredFasta, filename, lines[process]);
528 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
531 //force parent to wait until all the processes are done
532 for (int i=0;i<processors;i++) {
533 int temp = processIDS[i];
540 catch(exception& e) {
541 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
545 /**************************************************************************************/
546 string FilterSeqsCommand::createFilter() {
548 string filterString = "";
551 if (soft != 0) { F.setSoft(soft); }
552 if (trump != '*') { F.setTrump(trump); }
554 F.setLength(alignmentLength);
556 if(trump != '*' || isTrue(vertical) || soft != 0){
560 if(hard.compare("") != 0) { F.doHard(hard); }
561 else { F.setFilter(string(alignmentLength, '1')); }
564 if(trump != '*' || isTrue(vertical) || soft != 0){
565 for (int s = 0; s < fastafileNames.size(); s++) {
567 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
570 int pid, numSeqsPerProcessor, num;
572 vector<unsigned long int> MPIPos;
576 MPI_Comm_size(MPI_COMM_WORLD, &processors);
577 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
579 //char* tempFileName = new char(fastafileNames[s].length());
580 //tempFileName = &(fastafileNames[s][0]);
582 char tempFileName[1024];
583 strcpy(tempFileName, fastafileNames[s].c_str());
585 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
587 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
589 if (pid == 0) { //you are the root process
590 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
593 //send file positions to all processes
594 for(int i = 1; i < processors; i++) {
595 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
596 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
599 //figure out how many sequences you have to do
600 numSeqsPerProcessor = num / processors;
601 int startIndex = pid * numSeqsPerProcessor;
602 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
606 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
608 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
610 }else { //i am the child process
611 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
612 MPIPos.resize(num+1);
614 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
616 //figure out how many sequences you have to align
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; }
628 MPI_File_close(&inMPI);
629 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
632 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
636 openInputFile(fastafileNames[s], inFASTA);
637 getNumSeqs(inFASTA, numFastaSeqs);
640 numSeqs += numFastaSeqs;
642 lines.push_back(new linePair(0, numFastaSeqs));
644 driverCreateFilter(F, fastafileNames[s], lines[0]);
646 setLines(fastafileNames[s]);
647 createProcessesCreateFilter(F, fastafileNames[s]);
650 if (m->control_pressed) { return filterString; }
654 openInputFile(fastafileNames[s], inFASTA);
655 getNumSeqs(inFASTA, numFastaSeqs);
658 numSeqs += numFastaSeqs;
660 lines.push_back(new linePair(0, numFastaSeqs));
662 driverCreateFilter(F, fastafileNames[s], lines[0]);
663 if (m->control_pressed) { return filterString; }
673 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
676 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
678 if(trump != '*' || isTrue(vertical) || soft != 0){
680 if (pid == 0) { //only one process should output the filter
682 vector<int> temp; temp.resize(alignmentLength+1);
684 //get the frequencies from the child processes
685 for(int i = 1; i < processors; i++) {
687 for (int j = 0; j < 5; j++) {
689 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
690 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
692 if (receiveTag == Atag) { //you are recieveing the A frequencies
693 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
694 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
695 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
696 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
697 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
698 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
699 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
700 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
701 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
707 //send my fequency counts
709 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
711 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
713 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
715 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
716 F.gap.push_back(Gaptag);
717 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
722 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
724 if (pid == 0) { //only one process should output the filter
726 F.setNumSeqs(numSeqs);
728 if(isTrue(vertical) == 1) { F.doVertical(); }
729 if(soft != 0) { F.doSoft(); }
731 filterString = F.getFilter();
734 //send filter string to kids
735 //for(int i = 1; i < processors; i++) {
736 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
738 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
740 //recieve filterString
741 char* tempBuf = new char[alignmentLength];
742 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
743 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
745 filterString = tempBuf;
746 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
750 MPI_Barrier(MPI_COMM_WORLD);
755 catch(exception& e) {
756 m->errorOut(e, "FilterSeqsCommand", "createFilter");
760 /**************************************************************************************/
761 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* line) {
765 openInputFile(filename, in);
767 in.seekg(line->start);
769 for(int i=0;i<line->num;i++){
771 if (m->control_pressed) { in.close(); return 1; }
774 if (seq.getName() != "") {
775 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
777 if(trump != '*'){ F.doTrump(seq); }
778 if(isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
783 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
787 if((line->num) % 100 != 0){ m->mothurOut(toString(line->num)); m->mothurOutEndLine(); }
793 catch(exception& e) {
794 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
799 /**************************************************************************************/
800 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {
805 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
807 for(int i=0;i<num;i++){
809 if (m->control_pressed) { return 0; }
812 int length = MPIPos[start+i+1] - MPIPos[start+i];
814 char* buf4 = new char[length];
815 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
817 string tempBuf = buf4;
818 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
819 istringstream iss (tempBuf,istringstream::in);
824 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
826 if(trump != '*'){ F.doTrump(seq); }
827 if(isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
831 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
835 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
839 catch(exception& e) {
840 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
845 /**************************************************************************************************/
847 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
849 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
854 //loop through and create all the processes you want
855 while (process != processors) {
859 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
862 driverCreateFilter(F, filename, lines[process]);
864 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
867 //force parent to wait until all the processes are done
868 for (int i=0;i<processors;i++) {
869 int temp = processIDS[i];
876 catch(exception& e) {
877 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
881 /**************************************************************************************************/
883 int FilterSeqsCommand::setLines(string filename) {
886 vector<unsigned long int> positions;
890 openInputFile(filename, inFASTA);
893 while(!inFASTA.eof()){
894 input = getline(inFASTA);
896 if (input.length() != 0) {
897 if(input[0] == '>'){ unsigned long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
902 int numFastaSeqs = positions.size();
905 unsigned long int size;
907 //get num bytes in file
908 pFile = fopen (filename.c_str(),"rb");
909 if (pFile==NULL) perror ("Error opening file");
911 fseek (pFile, 0, SEEK_END);
916 numSeqs += numFastaSeqs;
918 int numSeqsPerProcessor = numFastaSeqs / processors;
920 for (int i = 0; i < processors; i++) {
922 unsigned long int startPos = positions[ i * numSeqsPerProcessor ];
923 if(i == processors - 1){
924 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
925 bufferSizes.push_back(size - startPos);
927 unsigned long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
928 bufferSizes.push_back(myEnd-startPos);
930 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
935 catch(exception& e) {
936 m->errorOut(e, "FilterSeqsCommand", "setLines");
940 /**************************************************************************************/