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"
14 //**********************************************************************************************************************
15 vector<string> FilterSeqsCommand::setParameters(){
17 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
18 CommandParameter phard("hard", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(phard);
19 CommandParameter ptrump("trump", "String", "", "*", "", "", "",false,false); parameters.push_back(ptrump);
20 CommandParameter psoft("soft", "Number", "", "0", "", "", "",false,false); parameters.push_back(psoft);
21 CommandParameter pvertical("vertical", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(pvertical);
22 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
23 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
24 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
26 vector<string> myArray;
27 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
31 m->errorOut(e, "FilterSeqsCommand", "setParameters");
35 //**********************************************************************************************************************
36 string FilterSeqsCommand::getHelpString(){
38 string helpString = "";
39 helpString += "The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n";
40 helpString += "The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n";
41 helpString += "The fasta parameter is required, unless you have a valid current fasta file. You may enter several fasta files to build the filter from and filter, by separating their names with -'s.\n";
42 helpString += "For example: fasta=abrecovery.fasta-amazon.fasta \n";
43 helpString += "The trump option will remove a column if the trump character is found at that position in any sequence of the alignment. Default=*, meaning no trump. \n";
44 helpString += "A soft mask removes any column where the dominant base (i.e. A, T, G, C, or U) does not occur in at least a designated percentage of sequences. Default=0.\n";
45 helpString += "The hard parameter allows you to enter a file containing the filter you want to use.\n";
46 helpString += "The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n";
47 helpString += "The processors parameter allows you to specify the number of processors to use. The default is 1.\n";
48 helpString += "The filter.seqs command should be in the following format: \n";
49 helpString += "filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n";
50 helpString += "Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n";
51 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n";
55 m->errorOut(e, "FilterSeqsCommand", "getHelpString");
59 //**********************************************************************************************************************
60 string FilterSeqsCommand::getOutputFileNameTag(string type, string inputName=""){
62 string outputFileName = "";
63 map<string, vector<string> >::iterator it;
65 //is this a type this command creates
66 it = outputTypes.find(type);
67 if (it == outputTypes.end()) { m->mothurOut("[ERROR]: this command doesn't create a " + type + " output file.\n"); }
69 if (type == "fasta") { outputFileName = "filter.fasta"; }
70 else if (type == "filter") { outputFileName = "filter"; }
71 else { m->mothurOut("[ERROR]: No definition for type " + type + " output file tag.\n"); m->control_pressed = true; }
73 return outputFileName;
76 m->errorOut(e, "FilterSeqsCommand", "getOutputFileNameTag");
80 //**********************************************************************************************************************
81 FilterSeqsCommand::FilterSeqsCommand(){
83 abort = true; calledHelp = true;
85 vector<string> tempOutNames;
86 outputTypes["fasta"] = tempOutNames;
87 outputTypes["filter"] = tempOutNames;
90 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
94 /**************************************************************************************/
95 FilterSeqsCommand::FilterSeqsCommand(string option) {
97 abort = false; calledHelp = false;
100 //allow user to run help
101 if(option == "help") { help(); abort = true; calledHelp = true; }
102 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
105 vector<string> myArray = setParameters();
107 OptionParser parser(option);
108 map<string,string> parameters = parser.getParameters();
110 ValidParameters validParameter("filter.seqs");
111 map<string,string>::iterator it;
113 //check to make sure all parameters are valid for command
114 for (it = parameters.begin(); it != parameters.end(); it++) {
115 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
118 //initialize outputTypes
119 vector<string> tempOutNames;
120 outputTypes["fasta"] = tempOutNames;
121 outputTypes["filter"] = tempOutNames;
123 //if the user changes the input directory command factory will send this info to us in the output parameter
124 string inputDir = validParameter.validFile(parameters, "inputdir", false);
125 if (inputDir == "not found"){ inputDir = ""; }
128 it = parameters.find("fasta");
129 //user has given a template file
130 if(it != parameters.end()){
131 path = m->hasPath(it->second);
132 //if the user has not given a path then, add inputdir. else leave path alone.
133 if (path == "") { parameters["fasta"] = inputDir + it->second; }
136 it = parameters.find("hard");
137 //user has given a template file
138 if(it != parameters.end()){
139 path = m->hasPath(it->second);
140 //if the user has not given a path then, add inputdir. else leave path alone.
141 if (path == "") { parameters["hard"] = inputDir + it->second; }
145 //check for required parameters
146 fasta = validParameter.validFile(parameters, "fasta", false);
147 if (fasta == "not found") {
148 fasta = m->getFastaFile();
150 fastafileNames.push_back(fasta);
151 m->mothurOut("Using " + fasta + " as input file for the fasta parameter."); m->mothurOutEndLine();
152 string simpleName = m->getSimpleName(fasta);
153 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
155 else { m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
158 m->splitAtDash(fasta, fastafileNames);
160 //go through files and make sure they are good, if not, then disregard them
161 for (int i = 0; i < fastafileNames.size(); i++) {
164 if (fastafileNames[i] == "current") {
165 fastafileNames[i] = m->getFastaFile();
166 if (fastafileNames[i] != "") { m->mothurOut("Using " + fastafileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
168 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
169 //erase from file list
170 fastafileNames.erase(fastafileNames.begin()+i);
176 if (inputDir != "") {
177 string path = m->hasPath(fastafileNames[i]);
178 //if the user has not given a path then, add inputdir. else leave path alone.
179 if (path == "") { fastafileNames[i] = inputDir + fastafileNames[i]; }
183 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
185 //if you can't open it, try default location
186 if (ableToOpen == 1) {
187 if (m->getDefaultPath() != "") { //default path is set
188 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
189 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
191 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
193 fastafileNames[i] = tryPath;
197 //if you can't open it, try default location
198 if (ableToOpen == 1) {
199 if (m->getOutputDir() != "") { //default path is set
200 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
201 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
203 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
205 fastafileNames[i] = tryPath;
211 if (ableToOpen == 1) {
212 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
213 //erase from file list
214 fastafileNames.erase(fastafileNames.begin()+i);
217 string simpleName = m->getSimpleName(fastafileNames[i]);
218 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
219 m->setFastaFile(fastafileNames[i]);
225 //make sure there is at least one valid file left
226 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
230 //if the user changes the output directory command factory will send this info to us in the output parameter
231 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
233 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
236 //check for optional parameter and set defaults
237 // ...at some point should added some additional type checking...
240 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
241 else if (hard == "not open") { hard = ""; abort = true; }
243 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
246 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
247 else { soft = (float)atoi(temp.c_str()) / 100.0; }
249 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
250 m->setProcessors(temp);
251 m->mothurConvert(temp, processors);
253 vertical = validParameter.validFile(parameters, "vertical", false);
254 if (vertical == "not found") {
255 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
256 else { vertical = "F"; }
263 catch(exception& e) {
264 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
268 /**************************************************************************************/
270 int FilterSeqsCommand::execute() {
273 if (abort == true) { if (calledHelp) { return 0; } return 2; }
276 m->openInputFile(fastafileNames[0], inFASTA);
278 Sequence testSeq(inFASTA);
279 alignmentLength = testSeq.getAlignLength();
282 ////////////create filter/////////////////
283 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
285 filter = createFilter();
287 m->mothurOutEndLine(); m->mothurOutEndLine();
289 if (m->control_pressed) { outputTypes.clear(); return 0; }
293 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
295 if (pid == 0) { //only one process should output the filter
300 //prevent giantic file name
302 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge." + getOutputFileNameTag("filter"); }
303 else { filterFile = outputDir + filterFileName + "." + getOutputFileNameTag("filter"); }
305 m->openOutputFile(filterFile, outFilter);
306 outFilter << filter << endl;
308 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
314 ////////////run filter/////////////////
316 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
320 m->mothurOutEndLine(); m->mothurOutEndLine();
322 int filteredLength = 0;
323 for(int i=0;i<alignmentLength;i++){
324 if(filter[i] == '1'){ filteredLength++; }
327 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); } return 0; }
330 m->mothurOutEndLine();
331 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
332 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
333 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
334 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
336 //set fasta file as new current fastafile
338 itTypes = outputTypes.find("fasta");
339 if (itTypes != outputTypes.end()) {
340 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
343 m->mothurOutEndLine();
344 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
345 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
346 m->mothurOutEndLine();
351 catch(exception& e) {
352 m->errorOut(e, "FilterSeqsCommand", "execute");
356 /**************************************************************************************/
357 int FilterSeqsCommand::filterSequences() {
362 for (int s = 0; s < fastafileNames.size(); s++) {
364 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
366 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + getOutputFileNameTag("fasta");
368 int pid, numSeqsPerProcessor, num;
370 vector<unsigned long long>MPIPos;
373 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
374 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
378 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
379 int inMode=MPI_MODE_RDONLY;
381 char outFilename[1024];
382 strcpy(outFilename, filteredFasta.c_str());
384 char inFileName[1024];
385 strcpy(inFileName, fastafileNames[s].c_str());
387 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
388 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
390 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
392 if (pid == 0) { //you are the root process
394 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
397 //send file positions to all processes
398 for(int i = 1; i < processors; i++) {
399 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
400 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
403 //figure out how many sequences you have to do
404 numSeqsPerProcessor = num / processors;
405 int startIndex = pid * numSeqsPerProcessor;
406 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
410 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
412 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
415 for(int i = 1; i < processors; i++) {
417 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
420 }else { //you are a child process
421 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
422 MPIPos.resize(num+1);
424 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
426 //figure out how many sequences you have to align
427 numSeqsPerProcessor = num / processors;
428 int startIndex = pid * numSeqsPerProcessor;
429 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
433 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
435 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
440 //tell parent you are done.
441 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
444 MPI_File_close(&outMPI);
445 MPI_File_close(&inMPI);
446 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
450 vector<unsigned long long> positions;
451 if (savedPositions.size() != 0) { positions = savedPositions[s]; }
453 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
454 positions = m->divideFile(fastafileNames[s], processors);
457 int numFastaSeqs = 0;
458 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
459 if (positions.size() < processors) { processors = positions.size(); }
463 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
464 //vector<unsigned long long> positions = m->divideFile(fastafileNames[s], processors);
466 for (int i = 0; i < (positions.size()-1); i++) {
467 lines.push_back(new linePair(positions[i], positions[(i+1)]));
471 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
472 numSeqs += numFastaSeqs;
474 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
475 numSeqs += numFastaSeqs;
478 if (m->control_pressed) { return 1; }
481 lines.push_back(new linePair(0, 1000));
482 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
483 numSeqs += numFastaSeqs;
485 int numFastaSeqs = positions.size()-1;
486 //positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
488 //figure out how many sequences you have to process
489 int numSeqsPerProcessor = numFastaSeqs / processors;
490 for (int i = 0; i < processors; i++) {
491 int startIndex = i * numSeqsPerProcessor;
492 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
493 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
496 numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
497 numSeqs += numFastaSeqs;
500 if (m->control_pressed) { return 1; }
503 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
508 catch(exception& e) {
509 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
514 /**************************************************************************************/
515 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long long>& MPIPos) {
517 string outputString = "";
521 for(int i=0;i<num;i++){
523 if (m->control_pressed) { return 0; }
526 int length = MPIPos[start+i+1] - MPIPos[start+i];
527 char* buf4 = new char[length];
528 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
530 string tempBuf = buf4;
531 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
532 istringstream iss (tempBuf,istringstream::in);
535 Sequence seq(iss); m->gobble(iss);
537 if (seq.getName() != "") {
538 string align = seq.getAligned();
539 string filterSeq = "";
541 for(int j=0;j<alignmentLength;j++){
542 if(filter[j] == '1'){
543 filterSeq += align[j];
548 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
550 if(count % 10 == 0){ //output to file
551 //send results to parent
552 int length = outputString.length();
553 char* buf = new char[length];
554 memcpy(buf, outputString.c_str(), length);
556 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
563 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
566 if(outputString != ""){ //output to file
567 //send results to parent
568 int length = outputString.length();
569 char* buf = new char[length];
570 memcpy(buf, outputString.c_str(), length);
572 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
577 if((num) % 100 != 0){ cout << (num) << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
581 catch(exception& e) {
582 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
587 /**************************************************************************************/
588 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
591 m->openOutputFile(outputFilename, out);
594 m->openInputFile(inputFilename, in);
596 in.seekg(filePos->start);
603 if (m->control_pressed) { in.close(); out.close(); return 0; }
605 Sequence seq(in); m->gobble(in);
606 if (seq.getName() != "") {
607 string align = seq.getAligned();
608 string filterSeq = "";
610 for(int j=0;j<alignmentLength;j++){
611 if(filter[j] == '1'){
612 filterSeq += align[j];
616 out << '>' << seq.getName() << endl << filterSeq << endl;
620 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
621 unsigned long long pos = in.tellg();
622 if ((pos == -1) || (pos >= filePos->end)) { break; }
624 if (in.eof()) { break; }
628 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
631 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
639 catch(exception& e) {
640 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
644 /**************************************************************************************************/
646 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename, string filteredFastaName) {
653 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
656 //loop through and create all the processes you want
657 while (process != processors) {
661 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
664 string filteredFasta = filename + toString(getpid()) + ".temp";
665 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
667 //pass numSeqs to parent
669 string tempFile = filename + toString(getpid()) + ".num.temp";
670 m->openOutputFile(tempFile, out);
676 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
677 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
682 num = driverRunFilter(F, filteredFastaName, filename, lines[0]);
684 //force parent to wait until all the processes are done
685 for (int i=0;i<processIDS.size();i++) {
686 int temp = processIDS[i];
690 for (int i = 0; i < processIDS.size(); i++) {
692 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
693 m->openInputFile(tempFile, in);
694 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
695 in.close(); m->mothurRemove(tempFile);
697 m->appendFiles((filename + toString(processIDS[i]) + ".temp"), filteredFastaName);
698 m->mothurRemove((filename + toString(processIDS[i]) + ".temp"));
703 //////////////////////////////////////////////////////////////////////////////////////////////////////
704 //Windows version shared memory, so be careful when passing variables through the filterData struct.
705 //Above fork() will clone, so memory is separate, but that's not the case with windows,
706 //Taking advantage of shared memory to allow both threads to add info to F.
707 //////////////////////////////////////////////////////////////////////////////////////////////////////
709 vector<filterRunData*> pDataArray;
710 DWORD dwThreadIdArray[processors-1];
711 HANDLE hThreadArray[processors-1];
713 //Create processor worker threads.
714 for( int i=0; i<processors-1; i++){
716 string extension = "";
717 if (i != 0) { extension = toString(i) + ".temp"; }
719 filterRunData* tempFilter = new filterRunData(filter, filename, (filteredFastaName + extension), m, lines[i]->start, lines[i]->end, alignmentLength, i);
720 pDataArray.push_back(tempFilter);
721 processIDS.push_back(i);
723 hThreadArray[i] = CreateThread(NULL, 0, MyRunFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
726 num = driverRunFilter(F, (filteredFastaName + toString(processors-1) + ".temp"), filename, lines[processors-1]);
728 //Wait until all threads have terminated.
729 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
731 //Close all thread handles and free memory allocations.
732 for(int i=0; i < pDataArray.size(); i++){
733 num += pDataArray[i]->count;
734 CloseHandle(hThreadArray[i]);
735 delete pDataArray[i];
738 for (int i = 1; i < processors; i++) {
739 m->appendFiles((filteredFastaName + toString(i) + ".temp"), filteredFastaName);
740 m->mothurRemove((filteredFastaName + toString(i) + ".temp"));
747 catch(exception& e) {
748 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
752 /**************************************************************************************/
753 string FilterSeqsCommand::createFilter() {
755 string filterString = "";
758 if (soft != 0) { F.setSoft(soft); }
759 if (trump != '*') { F.setTrump(trump); }
761 F.setLength(alignmentLength);
763 if(trump != '*' || m->isTrue(vertical) || soft != 0){
767 if(hard.compare("") != 0) { F.doHard(hard); }
768 else { F.setFilter(string(alignmentLength, '1')); }
771 if(trump != '*' || m->isTrue(vertical) || soft != 0){
772 for (int s = 0; s < fastafileNames.size(); s++) {
774 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
777 int pid, numSeqsPerProcessor, num;
779 vector<unsigned long long> MPIPos;
783 MPI_Comm_size(MPI_COMM_WORLD, &processors);
784 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
786 //char* tempFileName = new char(fastafileNames[s].length());
787 //tempFileName = &(fastafileNames[s][0]);
789 char tempFileName[1024];
790 strcpy(tempFileName, fastafileNames[s].c_str());
792 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
794 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
796 if (pid == 0) { //you are the root process
797 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
800 //send file positions to all processes
801 for(int i = 1; i < processors; i++) {
802 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
803 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
806 //figure out how many sequences you have to do
807 numSeqsPerProcessor = num / processors;
808 int startIndex = pid * numSeqsPerProcessor;
809 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
813 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
815 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
817 }else { //i am the child process
818 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
819 MPIPos.resize(num+1);
821 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
823 //figure out how many sequences you have to align
824 numSeqsPerProcessor = num / processors;
825 int startIndex = pid * numSeqsPerProcessor;
826 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
830 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
832 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
835 MPI_File_close(&inMPI);
836 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
840 vector<unsigned long long> positions;
841 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
842 positions = m->divideFile(fastafileNames[s], processors);
843 for (int i = 0; i < (positions.size()-1); i++) {
844 lines.push_back(new linePair(positions[i], positions[(i+1)]));
848 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
849 numSeqs += numFastaSeqs;
851 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
852 numSeqs += numFastaSeqs;
856 lines.push_back(new linePair(0, 1000));
857 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
858 numSeqs += numFastaSeqs;
860 int numFastaSeqs = 0;
861 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
862 if (positions.size() < processors) { processors = positions.size(); }
864 //figure out how many sequences you have to process
865 int numSeqsPerProcessor = numFastaSeqs / processors;
866 for (int i = 0; i < processors; i++) {
867 int startIndex = i * numSeqsPerProcessor;
868 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
869 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
872 numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
873 numSeqs += numFastaSeqs;
876 //save the file positions so we can reuse them in the runFilter function
877 savedPositions[s] = positions;
879 if (m->control_pressed) { return filterString; }
888 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
891 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
893 if(trump != '*' || m->isTrue(vertical) || soft != 0){
895 if (pid == 0) { //only one process should output the filter
897 vector<int> temp; temp.resize(alignmentLength+1);
899 //get the frequencies from the child processes
900 for(int i = 1; i < processors; i++) {
902 for (int j = 0; j < 5; j++) {
904 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
905 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
907 if (receiveTag == Atag) { //you are recieveing the A frequencies
908 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
909 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
910 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
911 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
912 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
913 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
914 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
915 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
916 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
922 //send my fequency counts
924 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
926 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
928 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
930 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
931 F.gap.push_back(Gaptag);
932 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
937 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
939 if (pid == 0) { //only one process should output the filter
942 F.setNumSeqs(numSeqs);
943 if(m->isTrue(vertical) == 1) { F.doVertical(); }
944 if(soft != 0) { F.doSoft(); }
945 filterString = F.getFilter();
948 //send filter string to kids
949 //for(int i = 1; i < processors; i++) {
950 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
952 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
954 //recieve filterString
955 char* tempBuf = new char[alignmentLength];
956 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
957 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
959 filterString = tempBuf;
960 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
964 MPI_Barrier(MPI_COMM_WORLD);
969 catch(exception& e) {
970 m->errorOut(e, "FilterSeqsCommand", "createFilter");
974 /**************************************************************************************/
975 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
979 m->openInputFile(filename, in);
981 in.seekg(filePos->start);
988 if (m->control_pressed) { in.close(); return 1; }
990 Sequence seq(in); m->gobble(in);
991 if (seq.getName() != "") {
992 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
994 if(trump != '*') { F.doTrump(seq); }
995 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
1000 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1001 unsigned long long pos = in.tellg();
1002 if ((pos == -1) || (pos >= filePos->end)) { break; }
1004 if (in.eof()) { break; }
1008 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
1011 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
1016 catch(exception& e) {
1017 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
1022 /**************************************************************************************/
1023 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long long>& MPIPos) {
1028 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
1030 for(int i=0;i<num;i++){
1032 if (m->control_pressed) { return 0; }
1034 //read next sequence
1035 int length = MPIPos[start+i+1] - MPIPos[start+i];
1037 char* buf4 = new char[length];
1038 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
1040 string tempBuf = buf4;
1041 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
1042 istringstream iss (tempBuf,istringstream::in);
1047 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
1049 if(trump != '*'){ F.doTrump(seq); }
1050 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
1054 if((i+1) % 100 == 0){ cout << (i+1) << endl; m->mothurOutJustToLog(toString(i+1) + "\n"); }
1058 if((num) % 100 != 0){ cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n"); }
1062 catch(exception& e) {
1063 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
1068 /**************************************************************************************************/
1070 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
1076 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1078 //loop through and create all the processes you want
1079 while (process != processors) {
1083 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
1085 }else if (pid == 0){
1086 //reset child's filter counts to 0;
1087 F.a.clear(); F.a.resize(alignmentLength, 0);
1088 F.t.clear(); F.t.resize(alignmentLength, 0);
1089 F.g.clear(); F.g.resize(alignmentLength, 0);
1090 F.c.clear(); F.c.resize(alignmentLength, 0);
1091 F.gap.clear(); F.gap.resize(alignmentLength, 0);
1093 num = driverCreateFilter(F, filename, lines[process]);
1095 //write out filter counts to file
1096 filename += toString(getpid()) + "filterValues.temp";
1098 m->openOutputFile(filename, out);
1101 out << F.getFilter() << endl;
1102 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
1103 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
1104 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
1105 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
1106 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
1108 //cout << F.getFilter() << endl;
1113 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
1114 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1119 //parent do your part
1120 num = driverCreateFilter(F, filename, lines[0]);
1122 //force parent to wait until all the processes are done
1123 for (int i=0;i<(processors-1);i++) {
1124 int temp = processIDS[i];
1128 //parent reads in and combines Filter info
1129 for (int i = 0; i < processIDS.size(); i++) {
1130 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1132 m->openInputFile(tempFilename, in);
1135 string tempFilterString;
1137 in >> tempNum; m->gobble(in); num += tempNum;
1139 in >> tempFilterString;
1140 F.mergeFilter(tempFilterString);
1142 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1143 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1144 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1145 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1146 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1149 m->mothurRemove(tempFilename);
1155 //////////////////////////////////////////////////////////////////////////////////////////////////////
1156 //Windows version shared memory, so be careful when passing variables through the filterData struct.
1157 //Above fork() will clone, so memory is separate, but that's not the case with windows,
1158 //Taking advantage of shared memory to allow both threads to add info to F.
1159 //////////////////////////////////////////////////////////////////////////////////////////////////////
1161 vector<filterData*> pDataArray;
1162 DWORD dwThreadIdArray[processors];
1163 HANDLE hThreadArray[processors];
1165 //Create processor worker threads.
1166 for( int i=0; i<processors; i++ ){
1168 filterData* tempFilter = new filterData(filename, m, lines[i]->start, lines[i]->end, alignmentLength, trump, vertical, soft, hard, i);
1169 pDataArray.push_back(tempFilter);
1170 processIDS.push_back(i);
1172 hThreadArray[i] = CreateThread(NULL, 0, MyCreateFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
1175 //Wait until all threads have terminated.
1176 WaitForMultipleObjects(processors, hThreadArray, TRUE, INFINITE);
1178 //Close all thread handles and free memory allocations.
1179 for(int i=0; i < pDataArray.size(); i++){
1180 num += pDataArray[i]->count;
1181 F.mergeFilter(pDataArray[i]->F.getFilter());
1183 for (int k = 0; k < alignmentLength; k++) { F.a[k] += pDataArray[i]->F.a[k]; }
1184 for (int k = 0; k < alignmentLength; k++) { F.t[k] += pDataArray[i]->F.t[k]; }
1185 for (int k = 0; k < alignmentLength; k++) { F.g[k] += pDataArray[i]->F.g[k]; }
1186 for (int k = 0; k < alignmentLength; k++) { F.c[k] += pDataArray[i]->F.c[k]; }
1187 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += pDataArray[i]->F.gap[k]; }
1189 CloseHandle(hThreadArray[i]);
1190 delete pDataArray[i];
1197 catch(exception& e) {
1198 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1202 /**************************************************************************************/