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","fasta-filter",false,true, 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, true); parameters.push_back(ptrump);
20 CommandParameter psoft("soft", "Number", "", "0", "", "", "","",false,false); parameters.push_back(psoft);
21 CommandParameter pvertical("vertical", "Boolean", "", "T", "", "", "","",false,false, true); parameters.push_back(pvertical);
22 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "","",false,false, true); 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::getOutputPattern(string type) {
64 if (type == "fasta") { pattern = "[filename],filter.fasta"; }
65 else if (type == "filter") { pattern = "[filename],filter"; }
66 else { m->mothurOut("[ERROR]: No definition for type " + type + " output pattern.\n"); m->control_pressed = true; }
71 m->errorOut(e, "FilterSeqsCommand", "getOutputPattern");
75 //**********************************************************************************************************************
76 FilterSeqsCommand::FilterSeqsCommand(){
78 abort = true; calledHelp = true;
80 vector<string> tempOutNames;
81 outputTypes["fasta"] = tempOutNames;
82 outputTypes["filter"] = tempOutNames;
85 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
89 /**************************************************************************************/
90 FilterSeqsCommand::FilterSeqsCommand(string option) {
92 abort = false; calledHelp = false;
95 //allow user to run help
96 if(option == "help") { help(); abort = true; calledHelp = true; }
97 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
100 vector<string> myArray = setParameters();
102 OptionParser parser(option);
103 map<string,string> parameters = parser.getParameters();
105 ValidParameters validParameter("filter.seqs");
106 map<string,string>::iterator it;
108 //check to make sure all parameters are valid for command
109 for (it = parameters.begin(); it != parameters.end(); it++) {
110 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
113 //initialize outputTypes
114 vector<string> tempOutNames;
115 outputTypes["fasta"] = tempOutNames;
116 outputTypes["filter"] = tempOutNames;
118 //if the user changes the input directory command factory will send this info to us in the output parameter
119 string inputDir = validParameter.validFile(parameters, "inputdir", false);
120 if (inputDir == "not found"){ inputDir = ""; }
123 it = parameters.find("fasta");
124 //user has given a template file
125 if(it != parameters.end()){
126 path = m->hasPath(it->second);
127 //if the user has not given a path then, add inputdir. else leave path alone.
128 if (path == "") { parameters["fasta"] = inputDir + it->second; }
131 it = parameters.find("hard");
132 //user has given a template file
133 if(it != parameters.end()){
134 path = m->hasPath(it->second);
135 //if the user has not given a path then, add inputdir. else leave path alone.
136 if (path == "") { parameters["hard"] = inputDir + it->second; }
140 //check for required parameters
141 fasta = validParameter.validFile(parameters, "fasta", false);
142 if (fasta == "not found") {
143 fasta = m->getFastaFile();
145 fastafileNames.push_back(fasta);
146 m->mothurOut("Using " + fasta + " as input file for the fasta parameter."); m->mothurOutEndLine();
147 string simpleName = m->getSimpleName(fasta);
148 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
150 else { m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
153 m->splitAtDash(fasta, fastafileNames);
155 //go through files and make sure they are good, if not, then disregard them
156 for (int i = 0; i < fastafileNames.size(); i++) {
159 if (fastafileNames[i] == "current") {
160 fastafileNames[i] = m->getFastaFile();
161 if (fastafileNames[i] != "") { m->mothurOut("Using " + fastafileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
163 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
164 //erase from file list
165 fastafileNames.erase(fastafileNames.begin()+i);
171 if (inputDir != "") {
172 string path = m->hasPath(fastafileNames[i]);
173 //if the user has not given a path then, add inputdir. else leave path alone.
174 if (path == "") { fastafileNames[i] = inputDir + fastafileNames[i]; }
178 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
180 //if you can't open it, try default location
181 if (ableToOpen == 1) {
182 if (m->getDefaultPath() != "") { //default path is set
183 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
184 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
186 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
188 fastafileNames[i] = tryPath;
192 //if you can't open it, try default location
193 if (ableToOpen == 1) {
194 if (m->getOutputDir() != "") { //default path is set
195 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
196 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
198 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
200 fastafileNames[i] = tryPath;
206 if (ableToOpen == 1) {
207 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
208 //erase from file list
209 fastafileNames.erase(fastafileNames.begin()+i);
212 string simpleName = m->getSimpleName(fastafileNames[i]);
213 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
214 m->setFastaFile(fastafileNames[i]);
220 //make sure there is at least one valid file left
221 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
225 //if the user changes the output directory command factory will send this info to us in the output parameter
226 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
228 outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it
231 //check for optional parameter and set defaults
232 // ...at some point should added some additional type checking...
235 hard = validParameter.validFile(parameters, "hard", true); if (hard == "not found") { hard = ""; }
236 else if (hard == "not open") { hard = ""; abort = true; }
238 temp = validParameter.validFile(parameters, "trump", false); if (temp == "not found") { temp = "*"; }
241 temp = validParameter.validFile(parameters, "soft", false); if (temp == "not found") { soft = 0; }
242 else { soft = (float)atoi(temp.c_str()) / 100.0; }
244 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
245 m->setProcessors(temp);
246 m->mothurConvert(temp, processors);
248 vertical = validParameter.validFile(parameters, "vertical", false);
249 if (vertical == "not found") {
250 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
251 else { vertical = "F"; }
258 catch(exception& e) {
259 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
263 /**************************************************************************************/
265 int FilterSeqsCommand::execute() {
268 if (abort == true) { if (calledHelp) { return 0; } return 2; }
271 m->openInputFile(fastafileNames[0], inFASTA);
273 Sequence testSeq(inFASTA);
274 alignmentLength = testSeq.getAlignLength();
277 ////////////create filter/////////////////
278 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
280 filter = createFilter();
282 m->mothurOutEndLine(); m->mothurOutEndLine();
284 if (m->control_pressed) { outputTypes.clear(); return 0; }
288 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
290 if (pid == 0) { //only one process should output the filter
295 //prevent giantic file name
296 map<string, string> variables;
297 variables["[filename]"] = outputDir + filterFileName + ".";
298 if (fastafileNames.size() > 3) { variables["[filename]"] = outputDir + "merge."; }
299 string filterFile = getOutputFileName("filter", variables);
301 m->openOutputFile(filterFile, outFilter);
302 outFilter << filter << endl;
304 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
310 ////////////run filter/////////////////
312 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
316 m->mothurOutEndLine(); m->mothurOutEndLine();
318 int filteredLength = 0;
319 for(int i=0;i<alignmentLength;i++){
320 if(filter[i] == '1'){ filteredLength++; }
323 if (m->control_pressed) { outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); } return 0; }
326 m->mothurOutEndLine();
327 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
328 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
329 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
330 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
332 //set fasta file as new current fastafile
334 itTypes = outputTypes.find("fasta");
335 if (itTypes != outputTypes.end()) {
336 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
339 m->mothurOutEndLine();
340 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
341 for(int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
342 m->mothurOutEndLine();
347 catch(exception& e) {
348 m->errorOut(e, "FilterSeqsCommand", "execute");
352 /**************************************************************************************/
353 int FilterSeqsCommand::filterSequences() {
358 for (int s = 0; s < fastafileNames.size(); s++) {
360 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
362 map<string, string> variables;
363 variables["[filename]"] = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s]));
364 string filteredFasta = getOutputFileName("fasta", variables);
366 int pid, numSeqsPerProcessor, num;
368 vector<unsigned long long>MPIPos;
371 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
372 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
376 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
377 int inMode=MPI_MODE_RDONLY;
379 char outFilename[1024];
380 strcpy(outFilename, filteredFasta.c_str());
382 char inFileName[1024];
383 strcpy(inFileName, fastafileNames[s].c_str());
385 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
386 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
388 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
390 if (pid == 0) { //you are the root process
392 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
395 //send file positions to all processes
396 for(int i = 1; i < processors; i++) {
397 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
398 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
401 //figure out how many sequences you have to do
402 numSeqsPerProcessor = num / processors;
403 int startIndex = pid * numSeqsPerProcessor;
404 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
408 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
410 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
413 for(int i = 1; i < processors; i++) {
415 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);
418 }else { //you are a child process
419 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
420 MPIPos.resize(num+1);
422 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
424 //figure out how many sequences you have to align
425 numSeqsPerProcessor = num / processors;
426 int startIndex = pid * numSeqsPerProcessor;
427 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
431 driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
433 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
438 //tell parent you are done.
439 MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
442 MPI_File_close(&outMPI);
443 MPI_File_close(&inMPI);
444 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
448 vector<unsigned long long> positions;
449 if (savedPositions.size() != 0) { positions = savedPositions[s]; }
451 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
452 positions = m->divideFile(fastafileNames[s], processors);
455 int numFastaSeqs = 0;
456 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
457 if (positions.size() < processors) { processors = positions.size(); }
461 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
462 //vector<unsigned long long> positions = m->divideFile(fastafileNames[s], processors);
464 for (int i = 0; i < (positions.size()-1); i++) {
465 lines.push_back(new linePair(positions[i], positions[(i+1)]));
469 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
470 numSeqs += numFastaSeqs;
472 int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
473 numSeqs += numFastaSeqs;
476 if (m->control_pressed) { return 1; }
479 lines.push_back(new linePair(0, 1000));
480 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
481 numSeqs += numFastaSeqs;
483 int numFastaSeqs = positions.size()-1;
484 //positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
486 //figure out how many sequences you have to process
487 int numSeqsPerProcessor = numFastaSeqs / processors;
488 for (int i = 0; i < processors; i++) {
489 int startIndex = i * numSeqsPerProcessor;
490 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
491 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
494 numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta);
495 numSeqs += numFastaSeqs;
498 if (m->control_pressed) { return 1; }
501 outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
506 catch(exception& e) {
507 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
512 /**************************************************************************************/
513 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long long>& MPIPos) {
515 string outputString = "";
519 for(int i=0;i<num;i++){
521 if (m->control_pressed) { return 0; }
524 int length = MPIPos[start+i+1] - MPIPos[start+i];
525 char* buf4 = new char[length];
526 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
528 string tempBuf = buf4;
529 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
530 istringstream iss (tempBuf,istringstream::in);
533 Sequence seq(iss); m->gobble(iss);
535 if (seq.getName() != "") {
536 string align = seq.getAligned();
537 string filterSeq = "";
539 for(int j=0;j<alignmentLength;j++){
540 if(filter[j] == '1'){
541 filterSeq += align[j];
546 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
548 if(count % 10 == 0){ //output to file
549 //send results to parent
550 int length = outputString.length();
551 char* buf = new char[length];
552 memcpy(buf, outputString.c_str(), length);
554 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
561 if((i+1) % 100 == 0){ cout << (i+1) << endl; }
564 if(outputString != ""){ //output to file
565 //send results to parent
566 int length = outputString.length();
567 char* buf = new char[length];
568 memcpy(buf, outputString.c_str(), length);
570 MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
575 if((num) % 100 != 0){ cout << (num) << endl; }
579 catch(exception& e) {
580 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
585 /**************************************************************************************/
586 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {
589 m->openOutputFile(outputFilename, out);
592 m->openInputFile(inputFilename, in);
594 in.seekg(filePos->start);
601 if (m->control_pressed) { in.close(); out.close(); return 0; }
603 Sequence seq(in); m->gobble(in);
604 if (seq.getName() != "") {
605 string align = seq.getAligned();
606 string filterSeq = "";
608 for(int j=0;j<alignmentLength;j++){
609 if(filter[j] == '1'){
610 filterSeq += align[j];
614 out << '>' << seq.getName() << endl << filterSeq << endl;
618 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
619 unsigned long long pos = in.tellg();
620 if ((pos == -1) || (pos >= filePos->end)) { break; }
622 if (in.eof()) { break; }
626 if((count) % 100 == 0){ m->mothurOutJustToScreen(toString(count)+"\n"); }
629 if((count) % 100 != 0){ m->mothurOutJustToScreen(toString(count)+"\n"); }
637 catch(exception& e) {
638 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
642 /**************************************************************************************************/
644 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename, string filteredFastaName) {
651 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
654 //loop through and create all the processes you want
655 while (process != processors) {
659 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
662 string filteredFasta = filename + toString(getpid()) + ".temp";
663 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
665 //pass numSeqs to parent
667 string tempFile = filename + toString(getpid()) + ".num.temp";
668 m->openOutputFile(tempFile, out);
674 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
675 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
680 num = driverRunFilter(F, filteredFastaName, filename, lines[0]);
682 //force parent to wait until all the processes are done
683 for (int i=0;i<processIDS.size();i++) {
684 int temp = processIDS[i];
688 for (int i = 0; i < processIDS.size(); i++) {
690 string tempFile = filename + toString(processIDS[i]) + ".num.temp";
691 m->openInputFile(tempFile, in);
692 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
693 in.close(); m->mothurRemove(tempFile);
695 m->appendFiles((filename + toString(processIDS[i]) + ".temp"), filteredFastaName);
696 m->mothurRemove((filename + toString(processIDS[i]) + ".temp"));
701 //////////////////////////////////////////////////////////////////////////////////////////////////////
702 //Windows version shared memory, so be careful when passing variables through the filterData struct.
703 //Above fork() will clone, so memory is separate, but that's not the case with windows,
704 //Taking advantage of shared memory to allow both threads to add info to F.
705 //////////////////////////////////////////////////////////////////////////////////////////////////////
707 vector<filterRunData*> pDataArray;
708 DWORD dwThreadIdArray[processors-1];
709 HANDLE hThreadArray[processors-1];
711 //Create processor worker threads.
712 for( int i=0; i<processors-1; i++){
714 string extension = "";
715 if (i != 0) { extension = toString(i) + ".temp"; }
717 filterRunData* tempFilter = new filterRunData(filter, filename, (filteredFastaName + extension), m, lines[i]->start, lines[i]->end, alignmentLength, i);
718 pDataArray.push_back(tempFilter);
719 processIDS.push_back(i);
721 hThreadArray[i] = CreateThread(NULL, 0, MyRunFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
724 num = driverRunFilter(F, (filteredFastaName + toString(processors-1) + ".temp"), filename, lines[processors-1]);
726 //Wait until all threads have terminated.
727 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
729 //Close all thread handles and free memory allocations.
730 for(int i=0; i < pDataArray.size(); i++){
731 num += pDataArray[i]->count;
732 if (pDataArray[i]->count != pDataArray[i]->end) {
733 m->mothurOut("[ERROR]: process " + toString(i) + " only processed " + toString(pDataArray[i]->count) + " of " + toString(pDataArray[i]->end) + " sequences assigned to it, quitting. \n"); m->control_pressed = true;
735 CloseHandle(hThreadArray[i]);
736 delete pDataArray[i];
739 for (int i = 1; i < processors; i++) {
740 m->appendFiles((filteredFastaName + toString(i) + ".temp"), filteredFastaName);
741 m->mothurRemove((filteredFastaName + toString(i) + ".temp"));
748 catch(exception& e) {
749 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
753 /**************************************************************************************/
754 string FilterSeqsCommand::createFilter() {
756 string filterString = "";
759 if (soft != 0) { F.setSoft(soft); }
760 if (trump != '*') { F.setTrump(trump); }
762 F.setLength(alignmentLength);
764 if(trump != '*' || m->isTrue(vertical) || soft != 0){
768 if(hard.compare("") != 0) { F.doHard(hard); }
769 else { F.setFilter(string(alignmentLength, '1')); }
772 if(trump != '*' || m->isTrue(vertical) || soft != 0){
773 for (int s = 0; s < fastafileNames.size(); s++) {
775 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
778 int pid, numSeqsPerProcessor, num;
780 vector<unsigned long long> MPIPos;
784 MPI_Comm_size(MPI_COMM_WORLD, &processors);
785 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
787 //char* tempFileName = new char(fastafileNames[s].length());
788 //tempFileName = &(fastafileNames[s][0]);
790 char tempFileName[1024];
791 strcpy(tempFileName, fastafileNames[s].c_str());
793 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
795 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
797 if (pid == 0) { //you are the root process
798 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
801 //send file positions to all processes
802 for(int i = 1; i < processors; i++) {
803 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
804 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
807 //figure out how many sequences you have to do
808 numSeqsPerProcessor = num / processors;
809 int startIndex = pid * numSeqsPerProcessor;
810 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
814 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
816 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
818 }else { //i am the child process
819 MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
820 MPIPos.resize(num+1);
822 MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
824 //figure out how many sequences you have to align
825 numSeqsPerProcessor = num / processors;
826 int startIndex = pid * numSeqsPerProcessor;
827 if(pid == (processors - 1)){ numSeqsPerProcessor = num - pid * numSeqsPerProcessor; }
831 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
833 if (m->control_pressed) { MPI_File_close(&inMPI); return 0; }
836 MPI_File_close(&inMPI);
837 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
841 vector<unsigned long long> positions;
842 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
843 positions = m->divideFile(fastafileNames[s], processors);
844 for (int i = 0; i < (positions.size()-1); i++) {
845 lines.push_back(new linePair(positions[i], positions[(i+1)]));
849 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
850 numSeqs += numFastaSeqs;
852 int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
853 numSeqs += numFastaSeqs;
857 lines.push_back(new linePair(0, 1000));
858 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
859 numSeqs += numFastaSeqs;
861 int numFastaSeqs = 0;
862 positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs);
863 if (positions.size() < processors) { processors = positions.size(); }
865 //figure out how many sequences you have to process
866 int numSeqsPerProcessor = numFastaSeqs / processors;
867 for (int i = 0; i < processors; i++) {
868 int startIndex = i * numSeqsPerProcessor;
869 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
870 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
873 numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]);
874 numSeqs += numFastaSeqs;
877 //save the file positions so we can reuse them in the runFilter function
878 savedPositions[s] = positions;
880 if (m->control_pressed) { return filterString; }
889 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
892 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
894 if(trump != '*' || m->isTrue(vertical) || soft != 0){
896 if (pid == 0) { //only one process should output the filter
898 vector<int> temp; temp.resize(alignmentLength+1);
900 //get the frequencies from the child processes
901 for(int i = 1; i < processors; i++) {
903 for (int j = 0; j < 5; j++) {
905 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
906 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what letter count this is for
908 if (receiveTag == Atag) { //you are recieveing the A frequencies
909 for (int k = 0; k < alignmentLength; k++) { F.a[k] += temp[k]; }
910 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
911 for (int k = 0; k < alignmentLength; k++) { F.t[k] += temp[k]; }
912 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
913 for (int k = 0; k < alignmentLength; k++) { F.c[k] += temp[k]; }
914 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
915 for (int k = 0; k < alignmentLength; k++) { F.g[k] += temp[k]; }
916 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
917 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += temp[k]; }
923 //send my fequency counts
925 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
927 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
929 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
931 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
932 F.gap.push_back(Gaptag);
933 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
938 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
940 if (pid == 0) { //only one process should output the filter
943 F.setNumSeqs(numSeqs);
944 if(m->isTrue(vertical) == 1) { F.doVertical(); }
945 if(soft != 0) { F.doSoft(); }
946 filterString = F.getFilter();
949 //send filter string to kids
950 //for(int i = 1; i < processors; i++) {
951 // MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
953 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
955 //recieve filterString
956 char* tempBuf = new char[alignmentLength];
957 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
958 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
960 filterString = tempBuf;
961 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength); }
965 MPI_Barrier(MPI_COMM_WORLD);
970 catch(exception& e) {
971 m->errorOut(e, "FilterSeqsCommand", "createFilter");
975 /**************************************************************************************/
976 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {
980 m->openInputFile(filename, in);
982 in.seekg(filePos->start);
989 if (m->control_pressed) { in.close(); return 1; }
991 Sequence seq(in); m->gobble(in);
992 if (seq.getName() != "") {
993 if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
995 if(trump != '*') { F.doTrump(seq); }
996 if(m->isTrue(vertical) || soft != 0) { F.getFreqs(seq); }
1001 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1002 unsigned long long pos = in.tellg();
1003 if ((pos == -1) || (pos >= filePos->end)) { break; }
1005 if (in.eof()) { break; }
1009 if((count) % 100 == 0){ m->mothurOutJustToScreen(toString(count)+"\n"); }
1012 if((count) % 100 != 0){ m->mothurOutJustToScreen(toString(count)+"\n"); }
1017 catch(exception& e) {
1018 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
1023 /**************************************************************************************/
1024 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long long>& MPIPos) {
1029 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
1031 for(int i=0;i<num;i++){
1033 if (m->control_pressed) { return 0; }
1035 //read next sequence
1036 int length = MPIPos[start+i+1] - MPIPos[start+i];
1038 char* buf4 = new char[length];
1039 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
1041 string tempBuf = buf4;
1042 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
1043 istringstream iss (tempBuf,istringstream::in);
1048 if (seq.getAligned().length() != alignmentLength) { cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1); }
1050 if(trump != '*'){ F.doTrump(seq); }
1051 if(m->isTrue(vertical) || soft != 0){ F.getFreqs(seq); }
1055 if((i+1) % 100 == 0){ cout << (i+1) << endl; }
1059 if((num) % 100 != 0){ cout << num << endl; }
1063 catch(exception& e) {
1064 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
1069 /**************************************************************************************************/
1071 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
1077 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1079 //loop through and create all the processes you want
1080 while (process != processors) {
1084 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
1086 }else if (pid == 0){
1087 //reset child's filter counts to 0;
1088 F.a.clear(); F.a.resize(alignmentLength, 0);
1089 F.t.clear(); F.t.resize(alignmentLength, 0);
1090 F.g.clear(); F.g.resize(alignmentLength, 0);
1091 F.c.clear(); F.c.resize(alignmentLength, 0);
1092 F.gap.clear(); F.gap.resize(alignmentLength, 0);
1094 num = driverCreateFilter(F, filename, lines[process]);
1096 //write out filter counts to file
1097 filename += toString(getpid()) + "filterValues.temp";
1099 m->openOutputFile(filename, out);
1102 out << F.getFilter() << endl;
1103 for (int k = 0; k < alignmentLength; k++) { out << F.a[k] << '\t'; } out << endl;
1104 for (int k = 0; k < alignmentLength; k++) { out << F.t[k] << '\t'; } out << endl;
1105 for (int k = 0; k < alignmentLength; k++) { out << F.g[k] << '\t'; } out << endl;
1106 for (int k = 0; k < alignmentLength; k++) { out << F.c[k] << '\t'; } out << endl;
1107 for (int k = 0; k < alignmentLength; k++) { out << F.gap[k] << '\t'; } out << endl;
1109 //cout << F.getFilter() << endl;
1114 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
1115 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1120 //parent do your part
1121 num = driverCreateFilter(F, filename, lines[0]);
1123 //force parent to wait until all the processes are done
1124 for (int i=0;i<(processors-1);i++) {
1125 int temp = processIDS[i];
1129 //parent reads in and combines Filter info
1130 for (int i = 0; i < processIDS.size(); i++) {
1131 string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1133 m->openInputFile(tempFilename, in);
1136 string tempFilterString;
1138 in >> tempNum; m->gobble(in); num += tempNum;
1140 in >> tempFilterString;
1141 F.mergeFilter(tempFilterString);
1143 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.a[k] += temp; } m->gobble(in);
1144 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.t[k] += temp; } m->gobble(in);
1145 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.g[k] += temp; } m->gobble(in);
1146 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.c[k] += temp; } m->gobble(in);
1147 for (int k = 0; k < alignmentLength; k++) { in >> temp; F.gap[k] += temp; } m->gobble(in);
1150 m->mothurRemove(tempFilename);
1156 //////////////////////////////////////////////////////////////////////////////////////////////////////
1157 //Windows version shared memory, so be careful when passing variables through the filterData struct.
1158 //Above fork() will clone, so memory is separate, but that's not the case with windows,
1159 //Taking advantage of shared memory to allow both threads to add info to F.
1160 //////////////////////////////////////////////////////////////////////////////////////////////////////
1162 vector<filterData*> pDataArray;
1163 DWORD dwThreadIdArray[processors];
1164 HANDLE hThreadArray[processors];
1166 //Create processor worker threads.
1167 for( int i=0; i<processors; i++ ){
1169 filterData* tempFilter = new filterData(filename, m, lines[i]->start, lines[i]->end, alignmentLength, trump, vertical, soft, hard, i);
1170 pDataArray.push_back(tempFilter);
1171 processIDS.push_back(i);
1173 hThreadArray[i] = CreateThread(NULL, 0, MyCreateFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
1176 //Wait until all threads have terminated.
1177 WaitForMultipleObjects(processors, hThreadArray, TRUE, INFINITE);
1179 //Close all thread handles and free memory allocations.
1180 for(int i=0; i < pDataArray.size(); i++){
1181 num += pDataArray[i]->count;
1182 if (pDataArray[i]->count != pDataArray[i]->end) {
1183 m->mothurOut("[ERROR]: process " + toString(i) + " only processed " + toString(pDataArray[i]->count) + " of " + toString(pDataArray[i]->end) + " sequences assigned to it, quitting. \n"); m->control_pressed = true;
1185 F.mergeFilter(pDataArray[i]->F.getFilter());
1187 for (int k = 0; k < alignmentLength; k++) { F.a[k] += pDataArray[i]->F.a[k]; }
1188 for (int k = 0; k < alignmentLength; k++) { F.t[k] += pDataArray[i]->F.t[k]; }
1189 for (int k = 0; k < alignmentLength; k++) { F.g[k] += pDataArray[i]->F.g[k]; }
1190 for (int k = 0; k < alignmentLength; k++) { F.c[k] += pDataArray[i]->F.c[k]; }
1191 for (int k = 0; k < alignmentLength; k++) { F.gap[k] += pDataArray[i]->F.gap[k]; }
1193 CloseHandle(hThreadArray[i]);
1194 delete pDataArray[i];
1201 catch(exception& e) {
1202 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1206 /**************************************************************************************/