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 //**********************************************************************************************************************
14 vector<string> TrimSeqsCommand::getValidParameters(){
16 string Array[] = {"fasta", "flip", "oligos", "maxambig", "maxhomop", "group","minlength", "maxlength", "qfile",
17 "qthreshold", "qwindowaverage", "qstepsize", "qwindowsize", "qaverage", "rollaverage",
18 "keepfirst", "removelast",
19 "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
20 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
24 m->errorOut(e, "TrimSeqsCommand", "getValidParameters");
28 //**********************************************************************************************************************
29 TrimSeqsCommand::TrimSeqsCommand(){
32 //initialize outputTypes
33 vector<string> tempOutNames;
34 outputTypes["fasta"] = tempOutNames;
35 outputTypes["qual"] = tempOutNames;
36 outputTypes["group"] = tempOutNames;
39 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
43 //**********************************************************************************************************************
44 vector<string> TrimSeqsCommand::getRequiredParameters(){
46 string Array[] = {"fasta"};
47 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
51 m->errorOut(e, "TrimSeqsCommand", "getRequiredParameters");
55 //**********************************************************************************************************************
56 vector<string> TrimSeqsCommand::getRequiredFiles(){
58 vector<string> myArray;
62 m->errorOut(e, "TrimSeqsCommand", "getRequiredFiles");
66 //***************************************************************************************************************
68 TrimSeqsCommand::TrimSeqsCommand(string option) {
74 //allow user to run help
75 if(option == "help") { help(); abort = true; }
78 //valid paramters for this command
79 string AlignArray[] = { "fasta", "flip", "oligos", "maxambig", "maxhomop", "group","minlength", "maxlength", "qfile",
80 "qthreshold", "qwindowaverage", "qstepsize", "qwindowsize", "qaverage", "rollaverage",
81 "keepfirst", "removelast",
82 "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
84 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
86 OptionParser parser(option);
87 map<string,string> parameters = parser.getParameters();
89 ValidParameters validParameter;
90 map<string,string>::iterator it;
92 //check to make sure all parameters are valid for command
93 for (it = parameters.begin(); it != parameters.end(); it++) {
94 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
97 //initialize outputTypes
98 vector<string> tempOutNames;
99 outputTypes["fasta"] = tempOutNames;
100 outputTypes["qual"] = tempOutNames;
101 outputTypes["group"] = tempOutNames;
103 //if the user changes the input directory command factory will send this info to us in the output parameter
104 string inputDir = validParameter.validFile(parameters, "inputdir", false);
105 if (inputDir == "not found"){ inputDir = ""; }
108 it = parameters.find("fasta");
109 //user has given a template file
110 if(it != parameters.end()){
111 path = m->hasPath(it->second);
112 //if the user has not given a path then, add inputdir. else leave path alone.
113 if (path == "") { parameters["fasta"] = inputDir + it->second; }
116 it = parameters.find("oligos");
117 //user has given a template file
118 if(it != parameters.end()){
119 path = m->hasPath(it->second);
120 //if the user has not given a path then, add inputdir. else leave path alone.
121 if (path == "") { parameters["oligos"] = inputDir + it->second; }
124 it = parameters.find("qfile");
125 //user has given a template file
126 if(it != parameters.end()){
127 path = m->hasPath(it->second);
128 //if the user has not given a path then, add inputdir. else leave path alone.
129 if (path == "") { parameters["qfile"] = inputDir + it->second; }
132 it = parameters.find("group");
133 //user has given a template file
134 if(it != parameters.end()){
135 path = m->hasPath(it->second);
136 //if the user has not given a path then, add inputdir. else leave path alone.
137 if (path == "") { parameters["group"] = inputDir + it->second; }
142 //check for required parameters
143 fastaFile = validParameter.validFile(parameters, "fasta", true);
144 if (fastaFile == "not found") { m->mothurOut("fasta is a required parameter for the trim.seqs command."); m->mothurOutEndLine(); abort = true; }
145 else if (fastaFile == "not open") { abort = true; }
147 //if the user changes the output directory command factory will send this info to us in the output parameter
148 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
150 outputDir += m->hasPath(fastaFile); //if user entered a file with a path then preserve it
154 //check for optional parameter and set defaults
155 // ...at some point should added some additional type checking...
157 temp = validParameter.validFile(parameters, "flip", false);
158 if (temp == "not found"){ flip = 0; }
159 else if(m->isTrue(temp)) { flip = 1; }
161 temp = validParameter.validFile(parameters, "oligos", true);
162 if (temp == "not found"){ oligoFile = ""; }
163 else if(temp == "not open"){ abort = true; }
164 else { oligoFile = temp; }
166 temp = validParameter.validFile(parameters, "group", true);
167 if (temp == "not found"){ groupfile = ""; }
168 else if(temp == "not open"){ abort = true; }
169 else { groupfile = temp; }
171 temp = validParameter.validFile(parameters, "maxambig", false); if (temp == "not found") { temp = "-1"; }
172 convert(temp, maxAmbig);
174 temp = validParameter.validFile(parameters, "maxhomop", false); if (temp == "not found") { temp = "0"; }
175 convert(temp, maxHomoP);
177 temp = validParameter.validFile(parameters, "minlength", false); if (temp == "not found") { temp = "0"; }
178 convert(temp, minLength);
180 temp = validParameter.validFile(parameters, "maxlength", false); if (temp == "not found") { temp = "0"; }
181 convert(temp, maxLength);
183 temp = validParameter.validFile(parameters, "bdiffs", false); if (temp == "not found") { temp = "0"; }
184 convert(temp, bdiffs);
186 temp = validParameter.validFile(parameters, "pdiffs", false); if (temp == "not found") { temp = "0"; }
187 convert(temp, pdiffs);
189 temp = validParameter.validFile(parameters, "tdiffs", false); if (temp == "not found") { int tempTotal = pdiffs + bdiffs; temp = toString(tempTotal); }
190 convert(temp, tdiffs);
192 if(tdiffs == 0){ tdiffs = bdiffs + pdiffs; }
194 temp = validParameter.validFile(parameters, "qfile", true);
195 if (temp == "not found") { qFileName = ""; }
196 else if(temp == "not open") { abort = true; }
197 else { qFileName = temp; }
199 temp = validParameter.validFile(parameters, "qthreshold", false); if (temp == "not found") { temp = "0"; }
200 convert(temp, qThreshold);
202 temp = validParameter.validFile(parameters, "qtrim", false); if (temp == "not found") { temp = "F"; }
203 qtrim = m->isTrue(temp);
205 temp = validParameter.validFile(parameters, "rollaverage", false); if (temp == "not found") { temp = "0"; }
206 convert(temp, qRollAverage);
208 temp = validParameter.validFile(parameters, "qwindowaverage", false);if (temp == "not found") { temp = "0"; }
209 convert(temp, qWindowAverage);
211 temp = validParameter.validFile(parameters, "qwindowsize", false); if (temp == "not found") { temp = "50"; }
212 convert(temp, qWindowSize);
214 temp = validParameter.validFile(parameters, "qstepsize", false); if (temp == "not found") { temp = "1"; }
215 convert(temp, qWindowStep);
217 temp = validParameter.validFile(parameters, "qaverage", false); if (temp == "not found") { temp = "0"; }
218 convert(temp, qAverage);
220 temp = validParameter.validFile(parameters, "keepfirst", false); if (temp == "not found") { temp = "0"; }
221 convert(temp, keepFirst);
223 temp = validParameter.validFile(parameters, "removelast", false); if (temp == "not found") { temp = "0"; }
224 convert(temp, removeLast);
226 temp = validParameter.validFile(parameters, "allfiles", false); if (temp == "not found") { temp = "F"; }
227 allFiles = m->isTrue(temp);
229 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found") { temp = "1"; }
230 convert(temp, processors);
232 if ((oligoFile != "") && (groupfile != "")) {
233 m->mothurOut("You given both a oligos file and a groupfile, only one is allowed."); m->mothurOutEndLine(); abort = true;
237 if(allFiles && (oligoFile == "") && (groupfile == "")){
238 m->mothurOut("You selected allfiles, but didn't enter an oligos or group file. Ignoring the allfiles request."); m->mothurOutEndLine();
240 if((qAverage != 0 && qThreshold != 0) && qFileName == ""){
241 m->mothurOut("You didn't provide a quality file name, quality criteria will be ignored."); m->mothurOutEndLine();
245 if(!flip && oligoFile=="" && !maxLength && !minLength && (maxAmbig==-1) && !maxHomoP && qFileName == ""){
246 m->mothurOut("You didn't set any options... quiting command."); m->mothurOutEndLine();
252 catch(exception& e) {
253 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
257 //**********************************************************************************************************************
259 void TrimSeqsCommand::help(){
261 m->mothurOut("The trim.seqs command reads a fastaFile and creates .....\n");
262 m->mothurOut("The trim.seqs command parameters are fasta, flip, oligos, group, maxambig, maxhomop, minlength, maxlength, qfile, qthreshold, qaverage, diffs, qtrim and allfiles.\n");
263 m->mothurOut("The fasta parameter is required.\n");
264 m->mothurOut("The group parameter allows you to enter a group file for your fasta file.\n");
265 m->mothurOut("The flip parameter will output the reverse compliment of your trimmed sequence. The default is false.\n");
266 m->mothurOut("The oligos parameter .... The default is "".\n");
267 m->mothurOut("The maxambig parameter .... The default is -1.\n");
268 m->mothurOut("The maxhomop parameter .... The default is 0.\n");
269 m->mothurOut("The minlength parameter .... The default is 0.\n");
270 m->mothurOut("The maxlength parameter .... The default is 0.\n");
271 m->mothurOut("The tdiffs parameter is used to specify the total number of differences allowed in the sequence. The default is pdiffs + bdiffs.\n");
272 m->mothurOut("The bdiffs parameter is used to specify the number of differences allowed in the barcode. The default is 0.\n");
273 m->mothurOut("The pdiffs parameter is used to specify the number of differences allowed in the primer. The default is 0.\n");
274 m->mothurOut("The qfile parameter .....\n");
275 m->mothurOut("The qthreshold parameter .... The default is 0.\n");
276 m->mothurOut("The qaverage parameter .... The default is 0.\n");
277 m->mothurOut("The allfiles parameter .... The default is F.\n");
278 m->mothurOut("The qtrim parameter .... The default is F.\n");
279 m->mothurOut("The trim.seqs command should be in the following format: \n");
280 m->mothurOut("trim.seqs(fasta=yourFastaFile, flip=yourFlip, oligos=yourOligos, maxambig=yourMaxambig, \n");
281 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength) \n");
282 m->mothurOut("Example trim.seqs(fasta=abrecovery.fasta, flip=..., oligos=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
283 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n");
284 m->mothurOut("For more details please check out the wiki http://www.mothur.org/wiki/Trim.seqs .\n\n");
287 catch(exception& e) {
288 m->errorOut(e, "TrimSeqsCommand", "help");
294 //***************************************************************************************************************
296 TrimSeqsCommand::~TrimSeqsCommand(){ /* do nothing */ }
298 //***************************************************************************************************************
300 int TrimSeqsCommand::execute(){
303 if (abort == true) { return 0; }
305 numFPrimers = 0; //this needs to be initialized
307 vector<string> fastaFileNames;
308 vector<string> qualFileNames;
310 string trimSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.fasta";
311 outputNames.push_back(trimSeqFile); outputTypes["fasta"].push_back(trimSeqFile);
312 string scrapSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.fasta";
313 outputNames.push_back(scrapSeqFile); outputTypes["fasta"].push_back(scrapSeqFile);
314 string trimQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.qual";
315 string scrapQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.qual";
316 if (qFileName != "") { outputNames.push_back(trimQualFile); outputNames.push_back(scrapQualFile); outputTypes["qual"].push_back(trimQualFile); outputTypes["qual"].push_back(scrapQualFile); }
317 string groupFile = "";
318 if (groupfile == "") { groupFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "groups"; }
320 groupFile = outputDir + m->getRootName(m->getSimpleName(groupfile)) + "trim.groups";
321 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
322 groupMap = new GroupMap(groupfile);
326 for (int i = 0; i < groupMap->namesOfGroups.size(); i++) {
327 groupToIndex[groupMap->namesOfGroups[i]] = i;
328 groupVector.push_back(groupMap->namesOfGroups[i]);
329 fastaFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupMap->namesOfGroups[i] + ".fasta"));
331 //we append later, so we want to clear file
333 m->openOutputFile(fastaFileNames[i], outRemove);
336 qualFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupMap->namesOfGroups[i] + ".qual"));
338 m->openOutputFile(qualFileNames[i], outRemove2);
343 comboStarts = fastaFileNames.size()-1;
347 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
348 getOligos(fastaFileNames, qualFileNames);
351 vector<unsigned long int> fastaFilePos;
352 vector<unsigned long int> qFilePos;
354 setLines(fastaFile, qFileName, fastaFilePos, qFilePos);
356 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
357 lines.push_back(new linePair(fastaFilePos[i], fastaFilePos[(i+1)]));
358 if (qFileName != "") { qLines.push_back(new linePair(qFilePos[i], qFilePos[(i+1)])); }
360 if(qFileName == "") { qLines = lines; } //files with duds
362 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
364 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
366 createProcessesCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames);
369 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
372 if (m->control_pressed) { return 0; }
375 for(int i=0;i<fastaFileNames.size();i++){
376 if (m->isBlank(fastaFileNames[i])) { blanks.insert(fastaFileNames[i]); }
377 else if (filesToRemove.count(fastaFileNames[i]) > 0) { remove(fastaFileNames[i].c_str()); }
381 m->openInputFile(fastaFileNames[i], inFASTA);
383 string outGroupFilename = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[i])) + "groups";
385 //if the fastafile is on the blanks list then the groups file should be as well
386 if (blanks.count(fastaFileNames[i]) != 0) { blanks.insert(outGroupFilename); }
388 m->openOutputFile(outGroupFilename, outGroups);
389 outputNames.push_back(outGroupFilename); outputTypes["group"].push_back(outGroupFilename);
391 string thisGroup = "";
392 if (i > comboStarts) {
393 map<string, int>::iterator itCombo;
394 for(itCombo=combos.begin();itCombo!=combos.end(); itCombo++){
395 if(itCombo->second == i){ thisGroup = itCombo->first; combos.erase(itCombo); break; }
397 }else{ thisGroup = groupVector[i]; }
399 while(!inFASTA.eof()){
400 if(inFASTA.get() == '>'){
402 outGroups << seqName << '\t' << thisGroup << endl;
404 while (!inFASTA.eof()) { char c = inFASTA.get(); if (c == 10 || c == 13){ break; } }
411 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
415 for(int i=0;i<qualFileNames.size();i++){
416 if (m->isBlank(qualFileNames[i])) { blanks.insert(qualFileNames[i]); }
417 else if (filesToRemove.count(qualFileNames[i]) > 0) { remove(qualFileNames[i].c_str()); }
421 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
423 if (m->control_pressed) {
424 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
428 m->mothurOutEndLine();
429 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
430 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
431 m->mothurOutEndLine();
436 catch(exception& e) {
437 m->errorOut(e, "TrimSeqsCommand", "execute");
442 /**************************************************************************************/
444 int TrimSeqsCommand::driverCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string trimQFile, string scrapQFile, string groupFile, vector<string> fastaNames, vector<string> qualNames, linePair* line, linePair* qline) {
449 int able = m->openOutputFile(trimFile, outFASTA);
452 m->openOutputFile(scrapFile, scrapFASTA);
457 m->openOutputFile(trimQFile, outQual);
458 m->openOutputFile(scrapQFile, scrapQual);
463 if (oligoFile != "") {
464 m->openOutputFile(groupFile, outGroups);
468 m->openInputFile(filename, inFASTA);
469 inFASTA.seekg(line->start);
472 if(qFileName != "") { m->openInputFile(qFileName, qFile); qFile.seekg(qline->start); }
475 for (int i = 0; i < fastaNames.size(); i++) { //clears old file
477 m->openOutputFile(fastaNames[i], temp);
480 for (int i = 0; i < qualNames.size(); i++) { //clears old file
482 m->openOutputFile(qualNames[i], temp);
492 if (m->control_pressed) {
493 inFASTA.close(); outFASTA.close(); scrapFASTA.close();
494 if (oligoFile != "") { outGroups.close(); }
499 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
507 Sequence currSeq(inFASTA); m->gobble(inFASTA);
509 QualityScores currQual;
511 currQual = QualityScores(qFile, currSeq.getNumBases()); m->gobble(qFile);
514 string origSeq = currSeq.getUnaligned();
516 int groupBar, groupPrime;
517 string trashCode = "";
518 int currentSeqsDiffs = 0;
520 if(barcodes.size() != 0){
521 success = stripBarcode(currSeq, currQual, groupBar);
522 if(success > bdiffs) { trashCode += 'b'; }
523 else{ currentSeqsDiffs += success; }
526 if(numFPrimers != 0){
527 success = stripForward(currSeq, currQual, groupPrime);
528 if(success > pdiffs) { trashCode += 'f'; }
529 else{ currentSeqsDiffs += success; }
532 if (currentSeqsDiffs > tdiffs) { trashCode += 't'; }
534 if(numRPrimers != 0){
535 success = stripReverse(currSeq, currQual);
536 if(!success) { trashCode += 'r'; }
540 success = keepFirstTrim(currSeq, currQual);
544 success = removeLastTrim(currSeq, currQual);
545 if(!success) { trashCode += 'l'; }
551 if(qThreshold != 0) { success = currQual.stripQualThreshold(currSeq, qThreshold); }
552 else if(qAverage != 0) { success = currQual.cullQualAverage(currSeq, qAverage); }
553 else if(qRollAverage != 0) { success = currQual.stripQualRollingAverage(currSeq, qRollAverage); }
554 else if(qWindowAverage != 0){ success = currQual.stripQualWindowAverage(currSeq, qWindowStep, qWindowSize, qWindowAverage); }
555 else { success = 1; }
557 // if (qtrim == 1 && (origSeq.length() != currSeq.getUnaligned().length())) {
558 // success = 0; //if you don't want to trim and the sequence does not meet quality requirements, move to scrap
561 if(!success) { trashCode += 'q'; }
564 if(minLength > 0 || maxLength > 0){
565 success = cullLength(currSeq);
566 if(!success) { trashCode += 'l'; }
569 success = cullHomoP(currSeq);
570 if(!success) { trashCode += 'h'; }
573 success = cullAmbigs(currSeq);
574 if(!success) { trashCode += 'n'; }
577 if(flip){ // should go last
578 currSeq.reverseComplement();
579 currQual.flipQScores();
582 if(trashCode.length() == 0){
583 currSeq.setAligned(currSeq.getUnaligned());
584 currSeq.printSequence(outFASTA);
585 currQual.printQScores(outQual);
587 if(barcodes.size() != 0){
588 string thisGroup = groupVector[groupBar];
589 int indexToFastaFile = groupBar;
590 if (primers.size() != 0){
591 //does this primer have a group
592 if (groupVector[groupPrime] != "") {
593 thisGroup += "." + groupVector[groupPrime];
594 indexToFastaFile = combos[thisGroup];
597 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
600 m->openOutputFileAppend(fastaNames[indexToFastaFile], outTemp);
601 //currSeq.printSequence(*fastaFileNames[indexToFastaFile]);
602 currSeq.printSequence(outTemp);
606 //currQual.printQScores(*qualFileNames[indexToFastaFile]);
608 m->openOutputFileAppend(qualNames[indexToFastaFile], outTemp2);
609 currQual.printQScores(outTemp2);
615 if (groupfile != "") {
616 string thisGroup = groupMap->getGroup(currSeq.getName());
618 if (thisGroup != "not found") {
619 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
622 m->openOutputFileAppend(fastaNames[groupToIndex[thisGroup]], outTemp);
623 currSeq.printSequence(outTemp);
627 m->openOutputFileAppend(qualNames[groupToIndex[thisGroup]], outTemp2);
628 currQual.printQScores(outTemp2);
633 m->mothurOut(currSeq.getName() + " is not in your groupfile, adding to group XXX."); m->mothurOutEndLine();
634 outGroups << currSeq.getName() << '\t' << "XXX" << endl;
636 m->mothurOut("[ERROR]: " + currSeq.getName() + " will not be added to any .group.fasta or .group.qual file."); m->mothurOutEndLine();
642 currSeq.setName(currSeq.getName() + '|' + trashCode);
643 currSeq.setUnaligned(origSeq);
644 currSeq.setAligned(origSeq);
645 currSeq.printSequence(scrapFASTA);
646 currQual.printQScores(scrapQual);
651 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
652 unsigned long int pos = inFASTA.tellg();
653 if ((pos == -1) || (pos >= line->end)) { break; }
655 if (inFASTA.eof()) { break; }
659 if((count) % 1000 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
663 if((count) % 1000 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
669 if (oligoFile != "") { outGroups.close(); }
670 if(qFileName != "") { qFile.close(); scrapQual.close(); outQual.close(); }
674 catch(exception& e) {
675 m->errorOut(e, "TrimSeqsCommand", "driverCreateTrim");
680 /**************************************************************************************************/
682 int TrimSeqsCommand::createProcessesCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string trimQFile, string scrapQFile, string groupFile, vector<string> fastaNames, vector<string> qualNames) {
684 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
689 //loop through and create all the processes you want
690 while (process != processors) {
694 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
697 for (int i = 0; i < fastaNames.size(); i++) {
698 fastaNames[i] = (fastaNames[i] + toString(getpid()) + ".temp");
699 //clear old file if it exists
701 m->openOutputFile(fastaNames[i], temp);
704 qualNames[i] = (qualNames[i] + toString(getpid()) + ".temp");
705 //clear old file if it exists
707 m->openOutputFile(qualNames[i], temp2);
712 driverCreateTrim(filename, qFileName, (trimFile + toString(getpid()) + ".temp"), (scrapFile + toString(getpid()) + ".temp"), (trimQFile + toString(getpid()) + ".temp"), (scrapQFile + toString(getpid()) + ".temp"), (groupFile + toString(getpid()) + ".temp"), fastaNames, qualNames, lines[process], qLines[process]);
715 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
716 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
722 for (int i = 0; i < fastaNames.size(); i++) {
723 //clear old file if it exists
725 m->openOutputFile(fastaNames[i], temp);
728 //clear old file if it exists
730 m->openOutputFile(qualNames[i], temp2);
735 driverCreateTrim(filename, qFileName, trimFile, scrapFile, trimQFile, scrapQFile, groupFile, fastaNames, qualNames, lines[0], qLines[0]);
738 //force parent to wait until all the processes are done
739 for (int i=0;i<processIDS.size();i++) {
740 int temp = processIDS[i];
745 for(int i=0;i<processIDS.size();i++){
746 m->mothurOut("Appending files from process " + processIDS[i]); m->mothurOutEndLine();
748 m->appendFiles((trimFile + toString(processIDS[i]) + ".temp"), trimFile);
749 remove((trimFile + toString(processIDS[i]) + ".temp").c_str());
750 m->appendFiles((scrapFile + toString(processIDS[i]) + ".temp"), scrapFile);
751 remove((scrapFile + toString(processIDS[i]) + ".temp").c_str());
753 m->mothurOut("Done with fasta files"); m->mothurOutEndLine();
756 m->appendFiles((trimQFile + toString(processIDS[i]) + ".temp"), trimQFile);
757 remove((trimQFile + toString(processIDS[i]) + ".temp").c_str());
758 m->appendFiles((scrapQFile + toString(processIDS[i]) + ".temp"), scrapQFile);
759 remove((scrapQFile + toString(processIDS[i]) + ".temp").c_str());
761 m->mothurOut("Done with quality files"); m->mothurOutEndLine();
764 m->appendFiles((groupFile + toString(processIDS[i]) + ".temp"), groupFile);
765 remove((groupFile + toString(processIDS[i]) + ".temp").c_str());
767 m->mothurOut("Done with group file"); m->mothurOutEndLine();
769 for (int j = 0; j < fastaNames.size(); j++) {
770 m->appendFiles((fastaNames[j] + toString(processIDS[i]) + ".temp"), fastaNames[j]);
771 remove((fastaNames[j] + toString(processIDS[i]) + ".temp").c_str());
775 for (int j = 0; j < qualNames.size(); j++) {
776 m->appendFiles((qualNames[j] + toString(processIDS[i]) + ".temp"), qualNames[j]);
777 remove((qualNames[j] + toString(processIDS[i]) + ".temp").c_str());
781 if (allFiles) { m->mothurOut("Done with allfiles"); m->mothurOutEndLine(); }
787 catch(exception& e) {
788 m->errorOut(e, "TrimSeqsCommand", "createProcessesCreateTrim");
793 /**************************************************************************************************/
795 int TrimSeqsCommand::setLines(string filename, string qfilename, vector<unsigned long int>& fastaFilePos, vector<unsigned long int>& qfileFilePos) {
798 //set file positions for fasta file
799 fastaFilePos = m->divideFile(filename, processors);
801 if (qfilename == "") { return processors; }
803 //get name of first sequence in each chunk
804 map<string, int> firstSeqNames;
805 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
807 m->openInputFile(filename, in);
808 in.seekg(fastaFilePos[i]);
811 firstSeqNames[temp.getName()] = i;
816 //seach for filePos of each first name in the qfile and save in qfileFilePos
818 m->openInputFile(qfilename, inQual);
821 while(!inQual.eof()){
822 input = m->getline(inQual);
824 if (input.length() != 0) {
825 if(input[0] == '>'){ //this is a sequence name line
826 istringstream nameStream(input);
828 string sname = ""; nameStream >> sname;
829 sname = sname.substr(1);
831 map<string, int>::iterator it = firstSeqNames.find(sname);
833 if(it != firstSeqNames.end()) { //this is the start of a new chunk
834 unsigned long int pos = inQual.tellg();
835 qfileFilePos.push_back(pos - input.length() - 1);
836 firstSeqNames.erase(it);
841 if (firstSeqNames.size() == 0) { break; }
846 if (firstSeqNames.size() != 0) {
847 for (map<string, int>::iterator it = firstSeqNames.begin(); it != firstSeqNames.end(); it++) {
848 m->mothurOut(it->first + " is in your fasta file and not in your quality file, not using quality file."); m->mothurOutEndLine();
854 //get last file position of qfile
856 unsigned long int size;
858 //get num bytes in file
859 pFile = fopen (qfilename.c_str(),"rb");
860 if (pFile==NULL) perror ("Error opening file");
862 fseek (pFile, 0, SEEK_END);
867 qfileFilePos.push_back(size);
871 catch(exception& e) {
872 m->errorOut(e, "TrimSeqsCommand", "setLines");
876 //***************************************************************************************************************
878 void TrimSeqsCommand::getOligos(vector<string>& outFASTAVec, vector<string>& outQualVec){
881 m->openInputFile(oligoFile, inOligos);
885 string type, oligo, group;
887 //int indexPrimer = 0;
889 while(!inOligos.eof()){
890 inOligos >> type; m->gobble(inOligos);
893 while (!inOligos.eof()) { char c = inOligos.get(); if (c == 10 || c == 13){ break; } } // get rest of line if there's any crap there
896 //make type case insensitive
897 for(int i=0;i<type.length();i++){ type[i] = toupper(type[i]); }
901 for(int i=0;i<oligo.length();i++){
902 oligo[i] = toupper(oligo[i]);
903 if(oligo[i] == 'U') { oligo[i] = 'T'; }
906 if(type == "FORWARD"){
909 // get rest of line in case there is a primer name
910 while (!inOligos.eof()) {
911 char c = inOligos.get();
912 if (c == 10 || c == 13){ break; }
913 else if (c == 32 || c == 9){;} //space or tab
917 //check for repeat barcodes
918 map<string, int>::iterator itPrime = primers.find(oligo);
919 if (itPrime != primers.end()) { m->mothurOut("primer " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
921 primers[oligo]=index; index++;
922 groupVector.push_back(group);
925 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
927 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
929 if (group == "") { //if there is not a group for this primer, then this file will not get written to, but we add it to keep the indexes correct
930 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
932 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
935 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
936 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
938 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
939 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
945 else if(type == "REVERSE"){
946 Sequence oligoRC("reverse", oligo);
947 oligoRC.reverseComplement();
948 revPrimer.push_back(oligoRC.getUnaligned());
950 else if(type == "BARCODE"){
953 //check for repeat barcodes
954 map<string, int>::iterator itBar = barcodes.find(oligo);
955 if (itBar != barcodes.end()) { m->mothurOut("barcode " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
957 barcodes[oligo]=index; index++;
958 groupVector.push_back(group);
961 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
962 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
963 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
965 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
966 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
967 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
971 }else{ m->mothurOut(type + " is not recognized as a valid type. Choices are forward, reverse, and barcode. Ignoring " + oligo + "."); m->mothurOutEndLine(); }
978 //add in potential combos
980 comboStarts = outFASTAVec.size()-1;
981 for (map<string, int>::iterator itBar = barcodes.begin(); itBar != barcodes.end(); itBar++) {
982 for (map<string, int>::iterator itPrime = primers.begin(); itPrime != primers.end(); itPrime++) {
983 if (groupVector[itPrime->second] != "") { //there is a group for this primer
984 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
985 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
986 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
987 combos[(groupVector[itBar->second] + "." + groupVector[itPrime->second])] = outFASTAVec.size()-1;
990 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
991 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
992 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
999 numFPrimers = primers.size();
1000 numRPrimers = revPrimer.size();
1003 catch(exception& e) {
1004 m->errorOut(e, "TrimSeqsCommand", "getOligos");
1008 //***************************************************************************************************************
1010 int TrimSeqsCommand::stripBarcode(Sequence& seq, QualityScores& qual, int& group){
1013 string rawSequence = seq.getUnaligned();
1014 int success = bdiffs + 1; //guilty until proven innocent
1016 //can you find the barcode
1017 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
1018 string oligo = it->first;
1019 if(rawSequence.length() < oligo.length()){ //let's just assume that the barcodes are the same length
1020 success = bdiffs + 10; //if the sequence is shorter than the barcode then bail out
1024 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1026 seq.setUnaligned(rawSequence.substr(oligo.length()));
1028 if(qual.getName() != ""){
1029 qual.trimQScores(oligo.length(), -1);
1037 //if you found the barcode or if you don't want to allow for diffs
1039 if ((bdiffs == 0) || (success == 0)) { return success; }
1041 else { //try aligning and see if you can find it
1046 Alignment* alignment;
1047 if (barcodes.size() > 0) {
1048 map<string,int>::iterator it=barcodes.begin();
1050 for(it;it!=barcodes.end();it++){
1051 if(it->first.length() > maxLength){
1052 maxLength = it->first.length();
1055 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+bdiffs+1));
1057 }else{ alignment = NULL; }
1059 //can you find the barcode
1065 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
1066 string oligo = it->first;
1067 // int length = oligo.length();
1069 if(rawSequence.length() < maxLength){ //let's just assume that the barcodes are the same length
1070 success = bdiffs + 10;
1074 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1075 alignment->align(oligo, rawSequence.substr(0,oligo.length()+bdiffs));
1076 oligo = alignment->getSeqAAln();
1077 string temp = alignment->getSeqBAln();
1079 int alnLength = oligo.length();
1081 for(int i=oligo.length()-1;i>=0;i--){
1082 if(oligo[i] != '-'){ alnLength = i+1; break; }
1084 oligo = oligo.substr(0,alnLength);
1085 temp = temp.substr(0,alnLength);
1088 int numDiff = countDiffs(oligo, temp);
1090 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1092 if(numDiff < minDiff){
1095 minGroup = it->second;
1097 for(int i=0;i<alnLength;i++){
1103 else if(numDiff == minDiff){
1109 if(minDiff > bdiffs) { success = minDiff; } //no good matches
1110 else if(minCount > 1) { success = bdiffs + 100; } //can't tell the difference between multiple barcodes
1111 else{ //use the best match
1113 seq.setUnaligned(rawSequence.substr(minPos));
1115 if(qual.getName() != ""){
1116 qual.trimQScores(minPos, -1);
1121 if (alignment != NULL) { delete alignment; }
1124 // cout << success << endl;
1129 catch(exception& e) {
1130 m->errorOut(e, "TrimSeqsCommand", "stripBarcode");
1136 //***************************************************************************************************************
1138 int TrimSeqsCommand::stripForward(Sequence& seq, QualityScores& qual, int& group){
1140 string rawSequence = seq.getUnaligned();
1141 int success = pdiffs + 1; //guilty until proven innocent
1143 //can you find the primer
1144 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1145 string oligo = it->first;
1146 if(rawSequence.length() < oligo.length()){ //let's just assume that the primers are the same length
1147 success = pdiffs + 10; //if the sequence is shorter than the barcode then bail out
1151 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1153 seq.setUnaligned(rawSequence.substr(oligo.length()));
1154 if(qual.getName() != ""){
1155 qual.trimQScores(oligo.length(), -1);
1162 //if you found the barcode or if you don't want to allow for diffs
1164 if ((pdiffs == 0) || (success == 0)) { return success; }
1166 else { //try aligning and see if you can find it
1171 Alignment* alignment;
1172 if (primers.size() > 0) {
1173 map<string,int>::iterator it=primers.begin();
1175 for(it;it!=primers.end();it++){
1176 if(it->first.length() > maxLength){
1177 maxLength = it->first.length();
1180 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+pdiffs+1));
1182 }else{ alignment = NULL; }
1184 //can you find the barcode
1190 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1191 string oligo = it->first;
1192 // int length = oligo.length();
1194 if(rawSequence.length() < maxLength){
1195 success = pdiffs + 100;
1199 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1200 alignment->align(oligo, rawSequence.substr(0,oligo.length()+pdiffs));
1201 oligo = alignment->getSeqAAln();
1202 string temp = alignment->getSeqBAln();
1204 int alnLength = oligo.length();
1206 for(int i=oligo.length()-1;i>=0;i--){
1207 if(oligo[i] != '-'){ alnLength = i+1; break; }
1209 oligo = oligo.substr(0,alnLength);
1210 temp = temp.substr(0,alnLength);
1213 int numDiff = countDiffs(oligo, temp);
1215 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1217 if(numDiff < minDiff){
1220 minGroup = it->second;
1222 for(int i=0;i<alnLength;i++){
1228 else if(numDiff == minDiff){
1234 if(minDiff > pdiffs) { success = minDiff; } //no good matches
1235 else if(minCount > 1) { success = pdiffs + 10; } //can't tell the difference between multiple primers
1236 else{ //use the best match
1238 seq.setUnaligned(rawSequence.substr(minPos));
1239 if(qual.getName() != ""){
1240 qual.trimQScores(minPos, -1);
1245 if (alignment != NULL) { delete alignment; }
1252 catch(exception& e) {
1253 m->errorOut(e, "TrimSeqsCommand", "stripForward");
1258 //***************************************************************************************************************
1260 bool TrimSeqsCommand::stripReverse(Sequence& seq, QualityScores& qual){
1262 string rawSequence = seq.getUnaligned();
1263 bool success = 0; //guilty until proven innocent
1265 for(int i=0;i<numRPrimers;i++){
1266 string oligo = revPrimer[i];
1268 if(rawSequence.length() < oligo.length()){
1273 if(compareDNASeq(oligo, rawSequence.substr(rawSequence.length()-oligo.length(),oligo.length()))){
1274 seq.setUnaligned(rawSequence.substr(0,rawSequence.length()-oligo.length()));
1275 if(qual.getName() != ""){
1276 qual.trimQScores(-1, rawSequence.length()-oligo.length());
1285 catch(exception& e) {
1286 m->errorOut(e, "TrimSeqsCommand", "stripReverse");
1291 //***************************************************************************************************************
1293 bool TrimSeqsCommand::keepFirstTrim(Sequence& sequence, QualityScores& qscores){
1296 if(qscores.getName() != ""){
1297 qscores.trimQScores(-1, keepFirst);
1299 sequence.trim(keepFirst);
1302 catch(exception& e) {
1303 m->errorOut(e, "keepFirstTrim", "countDiffs");
1309 //***************************************************************************************************************
1311 bool TrimSeqsCommand::removeLastTrim(Sequence& sequence, QualityScores& qscores){
1315 int length = sequence.getNumBases() - removeLast;
1318 if(qscores.getName() != ""){
1319 qscores.trimQScores(-1, length);
1321 sequence.trim(length);
1330 catch(exception& e) {
1331 m->errorOut(e, "removeLastTrim", "countDiffs");
1337 //***************************************************************************************************************
1339 bool TrimSeqsCommand::cullLength(Sequence& seq){
1342 int length = seq.getNumBases();
1343 bool success = 0; //guilty until proven innocent
1345 if(length >= minLength && maxLength == 0) { success = 1; }
1346 else if(length >= minLength && length <= maxLength) { success = 1; }
1347 else { success = 0; }
1352 catch(exception& e) {
1353 m->errorOut(e, "TrimSeqsCommand", "cullLength");
1359 //***************************************************************************************************************
1361 bool TrimSeqsCommand::cullHomoP(Sequence& seq){
1363 int longHomoP = seq.getLongHomoPolymer();
1364 bool success = 0; //guilty until proven innocent
1366 if(longHomoP <= maxHomoP){ success = 1; }
1367 else { success = 0; }
1371 catch(exception& e) {
1372 m->errorOut(e, "TrimSeqsCommand", "cullHomoP");
1378 //***************************************************************************************************************
1380 bool TrimSeqsCommand::cullAmbigs(Sequence& seq){
1382 int numNs = seq.getAmbigBases();
1383 bool success = 0; //guilty until proven innocent
1385 if(numNs <= maxAmbig) { success = 1; }
1386 else { success = 0; }
1390 catch(exception& e) {
1391 m->errorOut(e, "TrimSeqsCommand", "cullAmbigs");
1397 //***************************************************************************************************************
1399 bool TrimSeqsCommand::compareDNASeq(string oligo, string seq){
1402 int length = oligo.length();
1404 for(int i=0;i<length;i++){
1406 if(oligo[i] != seq[i]){
1407 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C') { success = 0; }
1408 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { success = 0; }
1409 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { success = 0; }
1410 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { success = 0; }
1411 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { success = 0; }
1412 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1413 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { success = 0; }
1414 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1415 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1416 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1417 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { success = 0; }
1418 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1420 if(success == 0) { break; }
1429 catch(exception& e) {
1430 m->errorOut(e, "TrimSeqsCommand", "compareDNASeq");
1435 //***************************************************************************************************************
1437 int TrimSeqsCommand::countDiffs(string oligo, string seq){
1440 int length = oligo.length();
1443 for(int i=0;i<length;i++){
1445 if(oligo[i] != seq[i]){
1446 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C' || oligo[i] == '-' || oligo[i] == '.') { countDiffs++; }
1447 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { countDiffs++; }
1448 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { countDiffs++; }
1449 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { countDiffs++; }
1450 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { countDiffs++; }
1451 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1452 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { countDiffs++; }
1453 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1454 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1455 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1456 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { countDiffs++; }
1457 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1464 catch(exception& e) {
1465 m->errorOut(e, "TrimSeqsCommand", "countDiffs");
1471 //***************************************************************************************************************