5 * Created by Pat Schloss on 6/6/09.
6 * Copyright 2009 Patrick D. Schloss. All rights reserved.
10 #include "trimseqscommand.h"
11 #include "needlemanoverlap.hpp"
13 //***************************************************************************************************************
15 TrimSeqsCommand::TrimSeqsCommand(string option) {
21 //allow user to run help
22 if(option == "help") { help(); abort = true; }
25 //valid paramters for this command
26 string AlignArray[] = {"fasta", "flip", "oligos", "maxambig", "maxhomop", "minlength", "maxlength", "qfile",
27 "qthreshold", "qaverage", "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
29 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
31 OptionParser parser(option);
32 map<string,string> parameters = parser.getParameters();
34 ValidParameters validParameter;
35 map<string,string>::iterator it;
37 //check to make sure all parameters are valid for command
38 for (it = parameters.begin(); it != parameters.end(); it++) {
39 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
42 //if the user changes the input directory command factory will send this info to us in the output parameter
43 string inputDir = validParameter.validFile(parameters, "inputdir", false);
44 if (inputDir == "not found"){ inputDir = ""; }
47 it = parameters.find("fasta");
48 //user has given a template file
49 if(it != parameters.end()){
50 path = hasPath(it->second);
51 //if the user has not given a path then, add inputdir. else leave path alone.
52 if (path == "") { parameters["fasta"] = inputDir + it->second; }
55 it = parameters.find("oligos");
56 //user has given a template file
57 if(it != parameters.end()){
58 path = hasPath(it->second);
59 //if the user has not given a path then, add inputdir. else leave path alone.
60 if (path == "") { parameters["oligos"] = inputDir + it->second; }
63 it = parameters.find("qfile");
64 //user has given a template file
65 if(it != parameters.end()){
66 path = hasPath(it->second);
67 //if the user has not given a path then, add inputdir. else leave path alone.
68 if (path == "") { parameters["qfile"] = inputDir + it->second; }
73 //check for required parameters
74 fastaFile = validParameter.validFile(parameters, "fasta", true);
75 if (fastaFile == "not found") { m->mothurOut("fasta is a required parameter for the screen.seqs command."); m->mothurOutEndLine(); abort = true; }
76 else if (fastaFile == "not open") { abort = true; }
78 //if the user changes the output directory command factory will send this info to us in the output parameter
79 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
81 outputDir += hasPath(fastaFile); //if user entered a file with a path then preserve it
84 //check for optional parameter and set defaults
85 // ...at some point should added some additional type checking...
87 temp = validParameter.validFile(parameters, "flip", false);
88 if (temp == "not found"){ flip = 0; }
89 else if(isTrue(temp)) { flip = 1; }
91 temp = validParameter.validFile(parameters, "oligos", true);
92 if (temp == "not found"){ oligoFile = ""; }
93 else if(temp == "not open"){ abort = true; }
94 else { oligoFile = temp; }
96 temp = validParameter.validFile(parameters, "maxambig", false); if (temp == "not found") { temp = "-1"; }
97 convert(temp, maxAmbig);
99 temp = validParameter.validFile(parameters, "maxhomop", false); if (temp == "not found") { temp = "0"; }
100 convert(temp, maxHomoP);
102 temp = validParameter.validFile(parameters, "minlength", false); if (temp == "not found") { temp = "0"; }
103 convert(temp, minLength);
105 temp = validParameter.validFile(parameters, "maxlength", false); if (temp == "not found") { temp = "0"; }
106 convert(temp, maxLength);
109 temp = validParameter.validFile(parameters, "bdiffs", false); if (temp == "not found") { temp = "0"; }
110 convert(temp, bdiffs);
112 temp = validParameter.validFile(parameters, "pdiffs", false); if (temp == "not found") { temp = "0"; }
113 convert(temp, pdiffs);
115 temp = validParameter.validFile(parameters, "tdiffs", false); if (temp == "not found") { int tempTotal = pdiffs + bdiffs; temp = toString(tempTotal); }
116 convert(temp, tdiffs);
118 if(tdiffs == 0){ tdiffs = bdiffs + pdiffs; }
120 temp = validParameter.validFile(parameters, "qfile", true);
121 if (temp == "not found") { qFileName = ""; }
122 else if(temp == "not open") { abort = true; }
123 else { qFileName = temp; }
125 temp = validParameter.validFile(parameters, "qthreshold", false); if (temp == "not found") { temp = "0"; }
126 convert(temp, qThreshold);
128 temp = validParameter.validFile(parameters, "qtrim", false); if (temp == "not found") { temp = "F"; }
129 qtrim = isTrue(temp);
131 temp = validParameter.validFile(parameters, "qaverage", false); if (temp == "not found") { temp = "0"; }
132 convert(temp, qAverage);
134 temp = validParameter.validFile(parameters, "allfiles", false); if (temp == "not found") { temp = "F"; }
135 allFiles = isTrue(temp);
137 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
138 convert(temp, processors);
140 if(allFiles && oligoFile == ""){
141 m->mothurOut("You selected allfiles, but didn't enter an oligos file. Ignoring the allfiles request."); m->mothurOutEndLine();
143 if((qAverage != 0 && qThreshold != 0) && qFileName == ""){
144 m->mothurOut("You didn't provide a quality file name, quality criteria will be ignored."); m->mothurOutEndLine();
148 if(!flip && oligoFile=="" && !maxLength && !minLength && (maxAmbig==-1) && !maxHomoP && qFileName == ""){
149 m->mothurOut("You didn't set any options... quiting command."); m->mothurOutEndLine();
155 catch(exception& e) {
156 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
160 //**********************************************************************************************************************
162 void TrimSeqsCommand::help(){
164 m->mothurOut("The trim.seqs command reads a fastaFile and creates .....\n");
165 m->mothurOut("The trim.seqs command parameters are fasta, flip, oligos, maxambig, maxhomop, minlength, maxlength, qfile, qthreshold, qaverage, diffs, qtrim and allfiles.\n");
166 m->mothurOut("The fasta parameter is required.\n");
167 m->mothurOut("The flip parameter will output the reverse compliment of your trimmed sequence. The default is false.\n");
168 m->mothurOut("The oligos parameter .... The default is "".\n");
169 m->mothurOut("The maxambig parameter .... The default is -1.\n");
170 m->mothurOut("The maxhomop parameter .... The default is 0.\n");
171 m->mothurOut("The minlength parameter .... The default is 0.\n");
172 m->mothurOut("The maxlength parameter .... The default is 0.\n");
173 m->mothurOut("The tdiffs parameter is used to specify the total number of differences allowed in the sequence. The default is pdiffs + bdiffs.\n");
174 m->mothurOut("The bdiffs parameter is used to specify the number of differences allowed in the barcode. The default is 0.\n");
175 m->mothurOut("The pdiffs parameter is used to specify the number of differences allowed in the primer. The default is 0.\n");
176 m->mothurOut("The qfile parameter .....\n");
177 m->mothurOut("The qthreshold parameter .... The default is 0.\n");
178 m->mothurOut("The qaverage parameter .... The default is 0.\n");
179 m->mothurOut("The allfiles parameter .... The default is F.\n");
180 m->mothurOut("The qtrim parameter .... The default is F.\n");
181 m->mothurOut("The trim.seqs command should be in the following format: \n");
182 m->mothurOut("trim.seqs(fasta=yourFastaFile, flip=yourFlip, oligos=yourOligos, maxambig=yourMaxambig, \n");
183 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength) \n");
184 m->mothurOut("Example trim.seqs(fasta=abrecovery.fasta, flip=..., oligos=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
185 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n");
186 m->mothurOut("For more details please check out the wiki http://www.mothur.org/wiki/Trim.seqs .\n\n");
189 catch(exception& e) {
190 m->errorOut(e, "TrimSeqsCommand", "help");
196 //***************************************************************************************************************
198 TrimSeqsCommand::~TrimSeqsCommand(){ /* do nothing */ }
200 //***************************************************************************************************************
202 int TrimSeqsCommand::execute(){
205 if (abort == true) { return 0; }
207 numFPrimers = 0; //this needs to be initialized
210 string trimSeqFile = outputDir + getRootName(getSimpleName(fastaFile)) + "trim.fasta";
211 outputNames.push_back(trimSeqFile);
212 string scrapSeqFile = outputDir + getRootName(getSimpleName(fastaFile)) + "scrap.fasta";
213 outputNames.push_back(scrapSeqFile);
214 string groupFile = outputDir + getRootName(getSimpleName(fastaFile)) + "groups";
216 vector<string> fastaFileNames;
218 outputNames.push_back(groupFile);
219 getOligos(fastaFileNames);
222 if(qFileName != "") { setLines(qFileName, qLines); }
225 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
229 openInputFile(fastaFile, inFASTA);
230 getNumSeqs(inFASTA, numSeqs);
233 lines.push_back(new linePair(0, numSeqs));
235 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames, lines[0], lines[0]);
237 for (int j = 0; j < fastaFileNames.size(); j++) {
238 rename((fastaFileNames[j] + toString(getpid()) + ".temp").c_str(), fastaFileNames[j].c_str());
242 setLines(fastaFile, lines);
243 if(qFileName == "") { qLines = lines; }
245 createProcessesCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames);
247 rename((trimSeqFile + toString(processIDS[0]) + ".temp").c_str(), trimSeqFile.c_str());
248 rename((scrapSeqFile + toString(processIDS[0]) + ".temp").c_str(), scrapSeqFile.c_str());
249 rename((groupFile + toString(processIDS[0]) + ".temp").c_str(), groupFile.c_str());
250 for (int j = 0; j < fastaFileNames.size(); j++) {
251 rename((fastaFileNames[j] + toString(processIDS[0]) + ".temp").c_str(), fastaFileNames[j].c_str());
254 for(int i=1;i<processors;i++){
255 appendFiles((trimSeqFile + toString(processIDS[i]) + ".temp"), trimSeqFile);
256 remove((trimSeqFile + toString(processIDS[i]) + ".temp").c_str());
257 appendFiles((scrapSeqFile + toString(processIDS[i]) + ".temp"), scrapSeqFile);
258 remove((scrapSeqFile + toString(processIDS[i]) + ".temp").c_str());
259 appendFiles((groupFile + toString(processIDS[i]) + ".temp"), groupFile);
260 remove((groupFile + toString(processIDS[i]) + ".temp").c_str());
261 for (int j = 0; j < fastaFileNames.size(); j++) {
262 appendFiles((fastaFileNames[j] + toString(processIDS[i]) + ".temp"), fastaFileNames[j]);
263 remove((fastaFileNames[j] + toString(processIDS[i]) + ".temp").c_str());
268 if (m->control_pressed) { return 0; }
272 openInputFile(fastaFile, inFASTA);
273 getNumSeqs(inFASTA, numSeqs);
276 lines.push_back(new linePair(0, numSeqs));
278 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames, lines[0], lines[0]);
280 if (m->control_pressed) { return 0; }
284 for(int i=0;i<fastaFileNames.size();i++){
285 if (isBlank(fastaFileNames[i])) { remove(fastaFileNames[i].c_str()); }
289 //openInputFile(getRootName(fastaFile) + groupVector[i] + ".fasta", inFASTA);
290 openInputFile(fastaFileNames[i], inFASTA);
292 string outGroupFilename = outputDir + getRootName(getSimpleName(fastaFileNames[i])) + "groups";
293 openOutputFile(outGroupFilename, outGroups);
294 //openOutputFile(outputDir + getRootName(getSimpleName(fastaFile)) + groupVector[i] + ".groups", outGroups);
295 outputNames.push_back(outGroupFilename);
297 string thisGroup = "";
298 if (i > comboStarts) {
299 map<string, int>::iterator itCombo;
300 for(itCombo=combos.begin();itCombo!=combos.end(); itCombo++){
301 if(itCombo->second == i){ thisGroup = itCombo->first; combos.erase(itCombo); break; }
303 }else{ thisGroup = groupVector[i]; }
305 while(!inFASTA.eof()){
306 if(inFASTA.get() == '>'){
308 outGroups << seqName << '\t' << thisGroup << endl;
310 while (!inFASTA.eof()) { char c = inFASTA.get(); if (c == 10 || c == 13){ break; } }
317 if (m->control_pressed) {
318 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
322 m->mothurOutEndLine();
323 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
324 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
325 m->mothurOutEndLine();
330 catch(exception& e) {
331 m->errorOut(e, "TrimSeqsCommand", "execute");
336 /**************************************************************************************/
337 int TrimSeqsCommand::driverCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string groupFile, vector<string> fastaNames, linePair* line, linePair* qline) {
341 int able = openOutputFile(trimFile, outFASTA);
344 openOutputFile(scrapFile, scrapFASTA);
347 vector<ofstream*> fastaFileNames;
348 if (oligoFile != "") {
349 openOutputFile(groupFile, outGroups);
350 for (int i = 0; i < fastaNames.size(); i++) {
351 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
352 fastaFileNames.push_back(new ofstream((fastaNames[i] + toString(getpid()) + ".temp").c_str(), ios::ate));
354 fastaFileNames.push_back(new ofstream((fastaNames[i] + toString(i) + ".temp").c_str(), ios::ate));
360 openInputFile(filename, inFASTA);
363 if(qFileName != "") { openInputFile(qFileName, qFile); }
365 qFile.seekg(qline->start);
366 inFASTA.seekg(line->start);
368 for(int i=0;i<line->num;i++){
370 if (m->control_pressed) {
371 inFASTA.close(); outFASTA.close(); scrapFASTA.close(); if (oligoFile != "") { outGroups.close(); } if(qFileName != "") { qFile.close(); }
372 for(int i=0;i<fastaFileNames.size();i++){ fastaFileNames[i]->close(); delete fastaFileNames[i]; }
373 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
379 Sequence currSeq(inFASTA);
381 string origSeq = currSeq.getUnaligned();
383 int groupBar, groupPrime;
384 string trashCode = "";
385 int currentSeqsDiffs = 0;
388 if(qThreshold != 0) { success = stripQualThreshold(currSeq, qFile); }
389 else if(qAverage != 0) { success = cullQualAverage(currSeq, qFile); }
391 if (qtrim == 1 && (origSeq.length() != currSeq.getUnaligned().length())) {
392 success = 0; //if you don't want to trim and the sequence does not meet quality requirements, move to scrap
395 if(!success) { trashCode += 'q'; }
398 if(barcodes.size() != 0){
399 success = stripBarcode(currSeq, groupBar);
400 if(success > bdiffs){ trashCode += 'b'; }
401 else{ currentSeqsDiffs += success; }
404 if(numFPrimers != 0){
405 success = stripForward(currSeq, groupPrime);
406 if(success > pdiffs){ trashCode += 'f'; }
407 else{ currentSeqsDiffs += success; }
410 if (currentSeqsDiffs > tdiffs) { trashCode += 't'; }
412 if(numRPrimers != 0){
413 success = stripReverse(currSeq);
414 if(!success){ trashCode += 'r'; }
417 if(minLength > 0 || maxLength > 0){
418 success = cullLength(currSeq);
419 if(!success){ trashCode += 'l'; }
422 success = cullHomoP(currSeq);
423 if(!success){ trashCode += 'h'; }
426 success = cullAmbigs(currSeq);
427 if(!success){ trashCode += 'n'; }
430 if(flip){ currSeq.reverseComplement(); } // should go last
432 if(trashCode.length() == 0){
433 currSeq.setAligned(currSeq.getUnaligned());
434 currSeq.printSequence(outFASTA);
435 if(barcodes.size() != 0){
436 string thisGroup = groupVector[groupBar];
437 int indexToFastaFile = groupBar;
438 if (primers.size() != 0){
439 //does this primer have a group
440 if (groupVector[groupPrime] != "") {
441 thisGroup += "." + groupVector[groupPrime];
442 indexToFastaFile = combos[thisGroup];
445 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
448 currSeq.printSequence(*fastaFileNames[indexToFastaFile]);
453 currSeq.setName(currSeq.getName() + '|' + trashCode);
454 currSeq.setUnaligned(origSeq);
455 currSeq.setAligned(origSeq);
456 currSeq.printSequence(scrapFASTA);
465 if (oligoFile != "") { outGroups.close(); }
466 if(qFileName != "") { qFile.close(); }
468 for(int i=0;i<fastaFileNames.size();i++){
469 fastaFileNames[i]->close();
470 delete fastaFileNames[i];
475 catch(exception& e) {
476 m->errorOut(e, "TrimSeqsCommand", "driverCreateTrim");
480 /**************************************************************************************************/
481 int TrimSeqsCommand::createProcessesCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string groupFile, vector<string> fastaNames) {
483 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
488 //loop through and create all the processes you want
489 while (process != processors) {
493 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
496 driverCreateTrim(filename, qFileName, (trimFile + toString(getpid()) + ".temp"), (scrapFile + toString(getpid()) + ".temp"), (groupFile + toString(getpid()) + ".temp"), fastaNames, lines[process], qLines[process]);
498 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
501 //force parent to wait until all the processes are done
502 for (int i=0;i<processors;i++) {
503 int temp = processIDS[i];
510 catch(exception& e) {
511 m->errorOut(e, "TrimSeqsCommand", "createProcessesCreateTrim");
515 /**************************************************************************************************/
517 int TrimSeqsCommand::setLines(string filename, vector<linePair*>& lines) {
522 vector<long int> positions;
525 openInputFile(filename, inFASTA);
528 while(!inFASTA.eof()){
529 input = getline(inFASTA);
531 if (input.length() != 0) {
532 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
537 int numFastaSeqs = positions.size();
542 //get num bytes in file
543 pFile = fopen (filename.c_str(),"rb");
544 if (pFile==NULL) perror ("Error opening file");
546 fseek (pFile, 0, SEEK_END);
551 int numSeqsPerProcessor = numFastaSeqs / processors;
553 for (int i = 0; i < processors; i++) {
555 long int startPos = positions[ i * numSeqsPerProcessor ];
556 if(i == processors - 1){
557 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
559 long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
561 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
566 catch(exception& e) {
567 m->errorOut(e, "TrimSeqsCommand", "setLines");
571 //***************************************************************************************************************
573 void TrimSeqsCommand::getOligos(vector<string>& outFASTAVec){ //vector<ofstream*>& outFASTAVec
576 openInputFile(oligoFile, inOligos);
580 string type, oligo, group;
582 //int indexPrimer = 0;
584 while(!inOligos.eof()){
588 while (!inOligos.eof()) { char c = inOligos.get(); if (c == 10 || c == 13){ break; } } // get rest of line if there's any crap there
593 for(int i=0;i<oligo.length();i++){
594 oligo[i] = toupper(oligo[i]);
595 if(oligo[i] == 'U') { oligo[i] = 'T'; }
598 if(type == "forward"){
601 // get rest of line in case there is a primer name
602 while (!inOligos.eof()) {
603 char c = inOligos.get();
604 if (c == 10 || c == 13){ break; }
605 else if (c == 32 || c == 9){;} //space or tab
609 //check for repeat barcodes
610 map<string, int>::iterator itPrime = primers.find(oligo);
611 if (itPrime != primers.end()) { m->mothurOut("primer " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
613 primers[oligo]=index; index++;
614 groupVector.push_back(group);
617 if (group != "") { //there is a group for this primer
618 outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
619 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
623 else if(type == "reverse"){
624 Sequence oligoRC("reverse", oligo);
625 oligoRC.reverseComplement();
626 revPrimer.push_back(oligoRC.getUnaligned());
628 else if(type == "barcode"){
631 //check for repeat barcodes
632 map<string, int>::iterator itBar = barcodes.find(oligo);
633 if (itBar != barcodes.end()) { m->mothurOut("barcode " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
635 barcodes[oligo]=index; index++;
636 groupVector.push_back(group);
639 outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
640 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
642 }else{ m->mothurOut(type + " is not recognized as a valid type. Choices are forward, reverse, and barcode. Ignoring " + oligo + "."); m->mothurOutEndLine(); }
649 //add in potential combos
651 comboStarts = outFASTAVec.size()-1;
652 for (map<string, int>::iterator itBar = barcodes.begin(); itBar != barcodes.end(); itBar++) {
653 for (map<string, int>::iterator itPrime = primers.begin(); itPrime != primers.end(); itPrime++) {
654 if (groupVector[itPrime->second] != "") { //there is a group for this primer
655 outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(itBar->second) + "." + groupVector[itBar->second] + "." + toString(itPrime->second) + "." + groupVector[itPrime->second] + ".fasta"));
656 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(itBar->second) + "." + groupVector[itBar->second] + "." + toString(itPrime->second) + "." + groupVector[itPrime->second] + ".fasta"));
657 combos[(groupVector[itBar->second] + "." + groupVector[itPrime->second])] = outFASTAVec.size()-1;
663 numFPrimers = primers.size();
664 numRPrimers = revPrimer.size();
667 catch(exception& e) {
668 m->errorOut(e, "TrimSeqsCommand", "getOligos");
672 //***************************************************************************************************************
674 int TrimSeqsCommand::stripBarcode(Sequence& seq, int& group){
677 string rawSequence = seq.getUnaligned();
678 int success = bdiffs + 1; //guilty until proven innocent
680 //can you find the barcode
681 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
682 string oligo = it->first;
683 if(rawSequence.length() < oligo.length()){ //let's just assume that the barcodes are the same length
684 success = bdiffs + 10; //if the sequence is shorter than the barcode then bail out
688 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
690 seq.setUnaligned(rawSequence.substr(oligo.length()));
696 //if you found the barcode or if you don't want to allow for diffs
698 if ((bdiffs == 0) || (success == 0)) { return success; }
700 else { //try aligning and see if you can find it
705 Alignment* alignment;
706 if (barcodes.size() > 0) {
707 map<string,int>::iterator it=barcodes.begin();
709 for(it;it!=barcodes.end();it++){
710 if(it->first.length() > maxLength){
711 maxLength = it->first.length();
714 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+bdiffs+1));
716 }else{ alignment = NULL; }
718 //can you find the barcode
724 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
725 string oligo = it->first;
726 // int length = oligo.length();
728 if(rawSequence.length() < maxLength){ //let's just assume that the barcodes are the same length
729 success = bdiffs + 10;
733 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
734 alignment->align(oligo, rawSequence.substr(0,oligo.length()+bdiffs));
735 oligo = alignment->getSeqAAln();
736 string temp = alignment->getSeqBAln();
738 int alnLength = oligo.length();
740 for(int i=oligo.length()-1;i>=0;i--){
741 if(oligo[i] != '-'){ alnLength = i+1; break; }
743 oligo = oligo.substr(0,alnLength);
744 temp = temp.substr(0,alnLength);
747 int numDiff = countDiffs(oligo, temp);
749 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
751 if(numDiff < minDiff){
754 minGroup = it->second;
756 for(int i=0;i<alnLength;i++){
762 else if(numDiff == minDiff){
768 if(minDiff > bdiffs) { success = minDiff; } //no good matches
769 else if(minCount > 1) { success = bdiffs + 100; } //can't tell the difference between multiple barcodes
770 else{ //use the best match
772 seq.setUnaligned(rawSequence.substr(minPos));
776 if (alignment != NULL) { delete alignment; }
779 // cout << success << endl;
784 catch(exception& e) {
785 m->errorOut(e, "TrimSeqsCommand", "stripBarcode");
791 //***************************************************************************************************************
793 int TrimSeqsCommand::stripForward(Sequence& seq, int& group){
795 string rawSequence = seq.getUnaligned();
796 int success = pdiffs + 1; //guilty until proven innocent
798 //can you find the primer
799 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
800 string oligo = it->first;
801 if(rawSequence.length() < oligo.length()){ //let's just assume that the primers are the same length
802 success = pdiffs + 10; //if the sequence is shorter than the barcode then bail out
806 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
808 seq.setUnaligned(rawSequence.substr(oligo.length()));
814 //if you found the barcode or if you don't want to allow for diffs
816 if ((pdiffs == 0) || (success == 0)) { return success; }
818 else { //try aligning and see if you can find it
823 Alignment* alignment;
824 if (primers.size() > 0) {
825 map<string,int>::iterator it=primers.begin();
827 for(it;it!=primers.end();it++){
828 if(it->first.length() > maxLength){
829 maxLength = it->first.length();
832 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+pdiffs+1));
834 }else{ alignment = NULL; }
836 //can you find the barcode
842 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
843 string oligo = it->first;
844 // int length = oligo.length();
846 if(rawSequence.length() < maxLength){
847 success = pdiffs + 100;
851 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
852 alignment->align(oligo, rawSequence.substr(0,oligo.length()+pdiffs));
853 oligo = alignment->getSeqAAln();
854 string temp = alignment->getSeqBAln();
856 int alnLength = oligo.length();
858 for(int i=oligo.length()-1;i>=0;i--){
859 if(oligo[i] != '-'){ alnLength = i+1; break; }
861 oligo = oligo.substr(0,alnLength);
862 temp = temp.substr(0,alnLength);
865 int numDiff = countDiffs(oligo, temp);
867 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
869 if(numDiff < minDiff){
872 minGroup = it->second;
874 for(int i=0;i<alnLength;i++){
880 else if(numDiff == minDiff){
886 if(minDiff > pdiffs) { success = minDiff; } //no good matches
887 else if(minCount > 1) { success = pdiffs + 10; } //can't tell the difference between multiple primers
888 else{ //use the best match
890 seq.setUnaligned(rawSequence.substr(minPos));
894 if (alignment != NULL) { delete alignment; }
901 catch(exception& e) {
902 m->errorOut(e, "TrimSeqsCommand", "stripForward");
907 //***************************************************************************************************************
909 bool TrimSeqsCommand::stripReverse(Sequence& seq){
911 string rawSequence = seq.getUnaligned();
912 bool success = 0; //guilty until proven innocent
914 for(int i=0;i<numRPrimers;i++){
915 string oligo = revPrimer[i];
917 if(rawSequence.length() < oligo.length()){
922 if(compareDNASeq(oligo, rawSequence.substr(rawSequence.length()-oligo.length(),oligo.length()))){
923 seq.setUnaligned(rawSequence.substr(0,rawSequence.length()-oligo.length()));
931 catch(exception& e) {
932 m->errorOut(e, "TrimSeqsCommand", "stripReverse");
937 //***************************************************************************************************************
939 bool TrimSeqsCommand::cullLength(Sequence& seq){
942 int length = seq.getNumBases();
943 bool success = 0; //guilty until proven innocent
945 if(length >= minLength && maxLength == 0) { success = 1; }
946 else if(length >= minLength && length <= maxLength) { success = 1; }
947 else { success = 0; }
952 catch(exception& e) {
953 m->errorOut(e, "TrimSeqsCommand", "cullLength");
959 //***************************************************************************************************************
961 bool TrimSeqsCommand::cullHomoP(Sequence& seq){
963 int longHomoP = seq.getLongHomoPolymer();
964 bool success = 0; //guilty until proven innocent
966 if(longHomoP <= maxHomoP){ success = 1; }
967 else { success = 0; }
971 catch(exception& e) {
972 m->errorOut(e, "TrimSeqsCommand", "cullHomoP");
978 //***************************************************************************************************************
980 bool TrimSeqsCommand::cullAmbigs(Sequence& seq){
982 int numNs = seq.getAmbigBases();
983 bool success = 0; //guilty until proven innocent
985 if(numNs <= maxAmbig) { success = 1; }
986 else { success = 0; }
990 catch(exception& e) {
991 m->errorOut(e, "TrimSeqsCommand", "cullAmbigs");
997 //***************************************************************************************************************
999 bool TrimSeqsCommand::compareDNASeq(string oligo, string seq){
1002 int length = oligo.length();
1004 for(int i=0;i<length;i++){
1006 if(oligo[i] != seq[i]){
1007 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C') { success = 0; }
1008 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { success = 0; }
1009 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { success = 0; }
1010 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { success = 0; }
1011 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { success = 0; }
1012 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1013 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { success = 0; }
1014 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1015 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1016 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1017 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { success = 0; }
1018 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1020 if(success == 0) { break; }
1029 catch(exception& e) {
1030 m->errorOut(e, "TrimSeqsCommand", "compareDNASeq");
1035 //***************************************************************************************************************
1037 int TrimSeqsCommand::countDiffs(string oligo, string seq){
1040 int length = oligo.length();
1043 for(int i=0;i<length;i++){
1045 if(oligo[i] != seq[i]){
1046 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C' || oligo[i] == '-' || oligo[i] == '.') { countDiffs++; }
1047 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { countDiffs++; }
1048 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { countDiffs++; }
1049 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { countDiffs++; }
1050 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { countDiffs++; }
1051 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1052 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { countDiffs++; }
1053 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1054 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1055 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1056 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { countDiffs++; }
1057 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1064 catch(exception& e) {
1065 m->errorOut(e, "TrimSeqsCommand", "countDiffs");
1070 //***************************************************************************************************************
1072 bool TrimSeqsCommand::stripQualThreshold(Sequence& seq, ifstream& qFile){
1074 // string rawSequence = seq.getUnaligned();
1075 // int seqLength; // = rawSequence.length();
1076 // string name, temp, temp2;
1080 // //get rest of line
1082 // while (!qFile.eof()) {
1083 // char c = qFile.get();
1084 // if (c == 10 || c == 13){ break; }
1085 // else { temp += c; }
1088 // int pos = temp.find("length");
1089 // if (pos == temp.npos) { m->mothurOut("Cannot find length in qfile for " + seq.getName()); m->mothurOutEndLine(); seqLength = 0; }
1091 // string tempLength = temp.substr(pos);
1092 // istringstream iss (tempLength,istringstream::in);
1096 // splitAtEquals(temp2, temp); //separates length=242, temp=length, temp2=242
1097 // convert(temp, seqLength); //converts string to int
1099 // if (name.length() != 0) { if(name.substr(1) != seq.getName()) { m->mothurOut("sequence name mismatch btwn fasta and qual file"); m->mothurOutEndLine(); } }
1101 string rawSequence = seq.getUnaligned();
1102 int seqLength = seq.getNumBases();
1103 bool success = 0; //guilty until proven innocent
1107 if (name[0] == '>') { if(name.substr(1) != seq.getName()) { m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); } }
1109 while (!qFile.eof()) { char c = qFile.get(); if (c == 10 || c == 13){ break; } }
1112 int end = seqLength;
1114 for(int i=0;i<seqLength;i++){
1117 if(score < qThreshold){
1122 for(int i=end+1;i<seqLength;i++){
1126 seq.setUnaligned(rawSequence.substr(0,end));
1130 catch(exception& e) {
1131 m->errorOut(e, "TrimSeqsCommand", "stripQualThreshold");
1136 //***************************************************************************************************************
1138 bool TrimSeqsCommand::cullQualAverage(Sequence& seq, ifstream& qFile){
1140 string rawSequence = seq.getUnaligned();
1141 int seqLength = seq.getNumBases();
1142 bool success = 0; //guilty until proven innocent
1146 if (name[0] == '>') { if(name.substr(1) != seq.getName()) { m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); } }
1148 while (!qFile.eof()) { char c = qFile.get(); if (c == 10 || c == 13){ break; } }
1153 for(int i=0;i<seqLength;i++){
1157 average /= seqLength;
1159 if(average >= qAverage) { success = 1; }
1160 else { success = 0; }
1164 catch(exception& e) {
1165 m->errorOut(e, "TrimSeqsCommand", "cullQualAverage");
1170 //***************************************************************************************************************