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 vector<string> TrimSeqsCommand::getValidParameters(){
17 string Array[] = {"fasta", "flip", "oligos", "maxambig", "maxhomop", "group","minlength", "maxlength", "qfile",
18 "qthreshold", "qwindowaverage", "qstepsize", "qwindowsize", "qaverage", "rollaverage",
19 "keepfirst", "removelast",
20 "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
21 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
25 m->errorOut(e, "TrimSeqsCommand", "getValidParameters");
30 //**********************************************************************************************************************
32 TrimSeqsCommand::TrimSeqsCommand(){
35 //initialize outputTypes
36 vector<string> tempOutNames;
37 outputTypes["fasta"] = tempOutNames;
38 outputTypes["qual"] = tempOutNames;
39 outputTypes["group"] = tempOutNames;
42 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
47 //**********************************************************************************************************************
49 vector<string> TrimSeqsCommand::getRequiredParameters(){
51 string Array[] = {"fasta"};
52 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
56 m->errorOut(e, "TrimSeqsCommand", "getRequiredParameters");
61 //**********************************************************************************************************************
63 vector<string> TrimSeqsCommand::getRequiredFiles(){
65 vector<string> myArray;
69 m->errorOut(e, "TrimSeqsCommand", "getRequiredFiles");
74 //***************************************************************************************************************
76 TrimSeqsCommand::TrimSeqsCommand(string option) {
82 //allow user to run help
83 if(option == "help") { help(); abort = true; }
86 //valid paramters for this command
87 string AlignArray[] = { "fasta", "flip", "oligos", "maxambig", "maxhomop", "group","minlength", "maxlength", "qfile",
88 "qthreshold", "qwindowaverage", "qstepsize", "qwindowsize", "qaverage", "rollaverage",
89 "keepfirst", "removelast",
90 "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
92 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
94 OptionParser parser(option);
95 map<string,string> parameters = parser.getParameters();
97 ValidParameters validParameter;
98 map<string,string>::iterator it;
100 //check to make sure all parameters are valid for command
101 for (it = parameters.begin(); it != parameters.end(); it++) {
102 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
105 //initialize outputTypes
106 vector<string> tempOutNames;
107 outputTypes["fasta"] = tempOutNames;
108 outputTypes["qual"] = tempOutNames;
109 outputTypes["group"] = tempOutNames;
111 //if the user changes the input directory command factory will send this info to us in the output parameter
112 string inputDir = validParameter.validFile(parameters, "inputdir", false);
113 if (inputDir == "not found"){ inputDir = ""; }
116 it = parameters.find("fasta");
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["fasta"] = inputDir + it->second; }
124 it = parameters.find("oligos");
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["oligos"] = inputDir + it->second; }
132 it = parameters.find("qfile");
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["qfile"] = inputDir + it->second; }
140 it = parameters.find("group");
141 //user has given a template file
142 if(it != parameters.end()){
143 path = m->hasPath(it->second);
144 //if the user has not given a path then, add inputdir. else leave path alone.
145 if (path == "") { parameters["group"] = inputDir + it->second; }
150 //check for required parameters
151 fastaFile = validParameter.validFile(parameters, "fasta", true);
152 if (fastaFile == "not found") { m->mothurOut("fasta is a required parameter for the trim.seqs command."); m->mothurOutEndLine(); abort = true; }
153 else if (fastaFile == "not open") { abort = true; }
155 //if the user changes the output directory command factory will send this info to us in the output parameter
156 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
158 outputDir += m->hasPath(fastaFile); //if user entered a file with a path then preserve it
162 //check for optional parameter and set defaults
163 // ...at some point should added some additional type checking...
165 temp = validParameter.validFile(parameters, "flip", false);
166 if (temp == "not found"){ flip = 0; }
167 else if(m->isTrue(temp)) { flip = 1; }
169 temp = validParameter.validFile(parameters, "oligos", true);
170 if (temp == "not found"){ oligoFile = ""; }
171 else if(temp == "not open"){ abort = true; }
172 else { oligoFile = temp; }
174 temp = validParameter.validFile(parameters, "group", true);
175 if (temp == "not found"){ groupfile = ""; }
176 else if(temp == "not open"){ abort = true; }
177 else { groupfile = temp; }
179 temp = validParameter.validFile(parameters, "maxambig", false); if (temp == "not found") { temp = "-1"; }
180 convert(temp, maxAmbig);
182 temp = validParameter.validFile(parameters, "maxhomop", false); if (temp == "not found") { temp = "0"; }
183 convert(temp, maxHomoP);
185 temp = validParameter.validFile(parameters, "minlength", false); if (temp == "not found") { temp = "0"; }
186 convert(temp, minLength);
188 temp = validParameter.validFile(parameters, "maxlength", false); if (temp == "not found") { temp = "0"; }
189 convert(temp, maxLength);
191 temp = validParameter.validFile(parameters, "bdiffs", false); if (temp == "not found") { temp = "0"; }
192 convert(temp, bdiffs);
194 temp = validParameter.validFile(parameters, "pdiffs", false); if (temp == "not found") { temp = "0"; }
195 convert(temp, pdiffs);
197 temp = validParameter.validFile(parameters, "tdiffs", false); if (temp == "not found") { int tempTotal = pdiffs + bdiffs; temp = toString(tempTotal); }
198 convert(temp, tdiffs);
200 if(tdiffs == 0){ tdiffs = bdiffs + pdiffs; }
202 temp = validParameter.validFile(parameters, "qfile", true);
203 if (temp == "not found") { qFileName = ""; }
204 else if(temp == "not open") { abort = true; }
205 else { qFileName = temp; }
207 temp = validParameter.validFile(parameters, "qthreshold", false); if (temp == "not found") { temp = "0"; }
208 convert(temp, qThreshold);
210 temp = validParameter.validFile(parameters, "qtrim", false); if (temp == "not found") { temp = "F"; }
211 qtrim = m->isTrue(temp);
213 temp = validParameter.validFile(parameters, "rollaverage", false); if (temp == "not found") { temp = "0"; }
214 convert(temp, qRollAverage);
216 temp = validParameter.validFile(parameters, "qwindowaverage", false);if (temp == "not found") { temp = "0"; }
217 convert(temp, qWindowAverage);
219 temp = validParameter.validFile(parameters, "qwindowsize", false); if (temp == "not found") { temp = "50"; }
220 convert(temp, qWindowSize);
222 temp = validParameter.validFile(parameters, "qstepsize", false); if (temp == "not found") { temp = "1"; }
223 convert(temp, qWindowStep);
225 temp = validParameter.validFile(parameters, "qaverage", false); if (temp == "not found") { temp = "0"; }
226 convert(temp, qAverage);
228 temp = validParameter.validFile(parameters, "keepfirst", false); if (temp == "not found") { temp = "0"; }
229 convert(temp, keepFirst);
231 temp = validParameter.validFile(parameters, "removelast", false); if (temp == "not found") { temp = "0"; }
232 convert(temp, removeLast);
234 temp = validParameter.validFile(parameters, "allfiles", false); if (temp == "not found") { temp = "F"; }
235 allFiles = m->isTrue(temp);
237 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found") { temp = "1"; }
238 convert(temp, processors);
240 if ((oligoFile != "") && (groupfile != "")) {
241 m->mothurOut("You given both a oligos file and a groupfile, only one is allowed."); m->mothurOutEndLine(); abort = true;
245 if(allFiles && (oligoFile == "") && (groupfile == "")){
246 m->mothurOut("You selected allfiles, but didn't enter an oligos or group file. Ignoring the allfiles request."); m->mothurOutEndLine();
248 if((qAverage != 0 && qThreshold != 0) && qFileName == ""){
249 m->mothurOut("You didn't provide a quality file name, quality criteria will be ignored."); m->mothurOutEndLine();
253 if(!flip && oligoFile=="" && !maxLength && !minLength && (maxAmbig==-1) && !maxHomoP && qFileName == ""){
254 m->mothurOut("You didn't set any options... quiting command."); m->mothurOutEndLine();
260 catch(exception& e) {
261 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
266 //**********************************************************************************************************************
268 void TrimSeqsCommand::help(){
270 m->mothurOut("The trim.seqs command reads a fastaFile and creates .....\n");
271 m->mothurOut("The trim.seqs command parameters are fasta, flip, oligos, group, maxambig, maxhomop, minlength, maxlength, qfile, qthreshold, qaverage, diffs, qtrim and allfiles.\n");
272 m->mothurOut("The fasta parameter is required.\n");
273 m->mothurOut("The group parameter allows you to enter a group file for your fasta file.\n");
274 m->mothurOut("The flip parameter will output the reverse compliment of your trimmed sequence. The default is false.\n");
275 m->mothurOut("The oligos parameter .... The default is "".\n");
276 m->mothurOut("The maxambig parameter .... The default is -1.\n");
277 m->mothurOut("The maxhomop parameter .... The default is 0.\n");
278 m->mothurOut("The minlength parameter .... The default is 0.\n");
279 m->mothurOut("The maxlength parameter .... The default is 0.\n");
280 m->mothurOut("The tdiffs parameter is used to specify the total number of differences allowed in the sequence. The default is pdiffs + bdiffs.\n");
281 m->mothurOut("The bdiffs parameter is used to specify the number of differences allowed in the barcode. The default is 0.\n");
282 m->mothurOut("The pdiffs parameter is used to specify the number of differences allowed in the primer. The default is 0.\n");
283 m->mothurOut("The qfile parameter .....\n");
284 m->mothurOut("The qthreshold parameter .... The default is 0.\n");
285 m->mothurOut("The qaverage parameter .... The default is 0.\n");
286 m->mothurOut("The allfiles parameter .... The default is F.\n");
287 m->mothurOut("The qtrim parameter .... The default is F.\n");
288 m->mothurOut("The trim.seqs command should be in the following format: \n");
289 m->mothurOut("trim.seqs(fasta=yourFastaFile, flip=yourFlip, oligos=yourOligos, maxambig=yourMaxambig, \n");
290 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength) \n");
291 m->mothurOut("Example trim.seqs(fasta=abrecovery.fasta, flip=..., oligos=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
292 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n");
293 m->mothurOut("For more details please check out the wiki http://www.mothur.org/wiki/Trim.seqs .\n\n");
296 catch(exception& e) {
297 m->errorOut(e, "TrimSeqsCommand", "help");
303 //***************************************************************************************************************
305 TrimSeqsCommand::~TrimSeqsCommand(){ /* do nothing */ }
307 //***************************************************************************************************************
309 int TrimSeqsCommand::execute(){
312 if (abort == true) { return 0; }
314 numFPrimers = 0; //this needs to be initialized
316 vector<string> fastaFileNames;
317 vector<string> qualFileNames;
319 string trimSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.fasta";
320 outputNames.push_back(trimSeqFile); outputTypes["fasta"].push_back(trimSeqFile);
321 string scrapSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.fasta";
322 outputNames.push_back(scrapSeqFile); outputTypes["fasta"].push_back(scrapSeqFile);
323 string trimQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.qual";
324 string scrapQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.qual";
325 if (qFileName != "") { outputNames.push_back(trimQualFile); outputNames.push_back(scrapQualFile); outputTypes["qual"].push_back(trimQualFile); outputTypes["qual"].push_back(scrapQualFile); }
326 string groupFile = "";
327 if (groupfile == "") { groupFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "groups"; }
329 groupFile = outputDir + m->getRootName(m->getSimpleName(groupfile)) + "trim.groups";
330 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
331 groupMap = new GroupMap(groupfile);
335 for (int i = 0; i < groupMap->namesOfGroups.size(); i++) {
336 groupToIndex[groupMap->namesOfGroups[i]] = i;
337 groupVector.push_back(groupMap->namesOfGroups[i]);
338 fastaFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupMap->namesOfGroups[i] + ".fasta"));
340 //we append later, so we want to clear file
342 m->openOutputFile(fastaFileNames[i], outRemove);
345 qualFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupMap->namesOfGroups[i] + ".qual"));
347 m->openOutputFile(qualFileNames[i], outRemove2);
352 comboStarts = fastaFileNames.size()-1;
356 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
357 getOligos(fastaFileNames, qualFileNames);
360 vector<unsigned long int> fastaFilePos;
361 vector<unsigned long int> qFilePos;
363 setLines(fastaFile, qFileName, fastaFilePos, qFilePos);
365 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
366 lines.push_back(new linePair(fastaFilePos[i], fastaFilePos[(i+1)]));
367 if (qFileName != "") { qLines.push_back(new linePair(qFilePos[i], qFilePos[(i+1)])); }
369 if(qFileName == "") { qLines = lines; } //files with duds
371 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
373 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
375 createProcessesCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames);
378 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
381 if (m->control_pressed) { return 0; }
384 for(int i=0;i<fastaFileNames.size();i++){
385 if (m->isBlank(fastaFileNames[i])) { blanks.insert(fastaFileNames[i]); }
386 else if (filesToRemove.count(fastaFileNames[i]) > 0) { remove(fastaFileNames[i].c_str()); }
390 m->openInputFile(fastaFileNames[i], inFASTA);
392 string outGroupFilename = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[i])) + "groups";
394 //if the fastafile is on the blanks list then the groups file should be as well
395 if (blanks.count(fastaFileNames[i]) != 0) { blanks.insert(outGroupFilename); }
397 m->openOutputFile(outGroupFilename, outGroups);
398 outputNames.push_back(outGroupFilename); outputTypes["group"].push_back(outGroupFilename);
400 string thisGroup = "";
401 if (i > comboStarts) {
402 map<string, int>::iterator itCombo;
403 for(itCombo=combos.begin();itCombo!=combos.end(); itCombo++){
404 if(itCombo->second == i){ thisGroup = itCombo->first; combos.erase(itCombo); break; }
406 }else{ thisGroup = groupVector[i]; }
408 while(!inFASTA.eof()){
409 if(inFASTA.get() == '>'){
411 outGroups << seqName << '\t' << thisGroup << endl;
413 while (!inFASTA.eof()) { char c = inFASTA.get(); if (c == 10 || c == 13){ break; } }
420 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
424 for(int i=0;i<qualFileNames.size();i++){
425 if (m->isBlank(qualFileNames[i])) { blanks.insert(qualFileNames[i]); }
426 else if (filesToRemove.count(qualFileNames[i]) > 0) { remove(qualFileNames[i].c_str()); }
430 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
432 if (m->control_pressed) {
433 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
437 m->mothurOutEndLine();
438 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
439 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
440 m->mothurOutEndLine();
445 catch(exception& e) {
446 m->errorOut(e, "TrimSeqsCommand", "execute");
451 /**************************************************************************************/
453 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) {
458 int able = m->openOutputFile(trimFile, outFASTA);
461 m->openOutputFile(scrapFile, scrapFASTA);
466 m->openOutputFile(trimQFile, outQual);
467 m->openOutputFile(scrapQFile, scrapQual);
472 if (oligoFile != "") {
473 m->openOutputFile(groupFile, outGroups);
477 m->openInputFile(filename, inFASTA);
478 inFASTA.seekg(line->start);
481 if(qFileName != "") { m->openInputFile(qFileName, qFile); qFile.seekg(qline->start); }
484 for (int i = 0; i < fastaNames.size(); i++) { //clears old file
486 m->openOutputFile(fastaNames[i], temp);
489 for (int i = 0; i < qualNames.size(); i++) { //clears old file
491 m->openOutputFile(qualNames[i], temp);
501 if (m->control_pressed) {
502 inFASTA.close(); outFASTA.close(); scrapFASTA.close();
503 if (oligoFile != "") { outGroups.close(); }
508 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
516 Sequence currSeq(inFASTA); m->gobble(inFASTA);
518 QualityScores currQual;
520 currQual = QualityScores(qFile, currSeq.getNumBases()); m->gobble(qFile);
523 string origSeq = currSeq.getUnaligned();
525 int groupBar, groupPrime;
526 string trashCode = "";
527 int currentSeqsDiffs = 0;
529 if(barcodes.size() != 0){
530 success = stripBarcode(currSeq, currQual, groupBar);
531 if(success > bdiffs) { trashCode += 'b'; }
532 else{ currentSeqsDiffs += success; }
535 if(numFPrimers != 0){
536 success = stripForward(currSeq, currQual, groupPrime);
537 if(success > pdiffs) { trashCode += 'f'; }
538 else{ currentSeqsDiffs += success; }
541 if (currentSeqsDiffs > tdiffs) { trashCode += 't'; }
543 if(numRPrimers != 0){
544 success = stripReverse(currSeq, currQual);
545 if(!success) { trashCode += 'r'; }
549 success = keepFirstTrim(currSeq, currQual);
553 success = removeLastTrim(currSeq, currQual);
554 if(!success) { trashCode += 'l'; }
560 if(qThreshold != 0) { success = currQual.stripQualThreshold(currSeq, qThreshold); }
561 else if(qAverage != 0) { success = currQual.cullQualAverage(currSeq, qAverage); }
562 else if(qRollAverage != 0) { success = currQual.stripQualRollingAverage(currSeq, qRollAverage); }
563 else if(qWindowAverage != 0){ success = currQual.stripQualWindowAverage(currSeq, qWindowStep, qWindowSize, qWindowAverage); }
564 else { success = 1; }
566 // if (qtrim == 1 && (origSeq.length() != currSeq.getUnaligned().length())) {
567 // success = 0; //if you don't want to trim and the sequence does not meet quality requirements, move to scrap
570 if(!success) { trashCode += 'q'; }
573 if(minLength > 0 || maxLength > 0){
574 success = cullLength(currSeq);
575 if(!success) { trashCode += 'l'; }
578 success = cullHomoP(currSeq);
579 if(!success) { trashCode += 'h'; }
582 success = cullAmbigs(currSeq);
583 if(!success) { trashCode += 'n'; }
586 if(flip){ // should go last
587 currSeq.reverseComplement();
588 currQual.flipQScores();
591 if(trashCode.length() == 0){
592 currSeq.setAligned(currSeq.getUnaligned());
593 currSeq.printSequence(outFASTA);
594 currQual.printQScores(outQual);
596 if(barcodes.size() != 0){
597 string thisGroup = groupVector[groupBar];
598 int indexToFastaFile = groupBar;
599 if (primers.size() != 0){
600 //does this primer have a group
601 if (groupVector[groupPrime] != "") {
602 thisGroup += "." + groupVector[groupPrime];
603 indexToFastaFile = combos[thisGroup];
606 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
609 m->openOutputFileAppend(fastaNames[indexToFastaFile], outTemp);
610 //currSeq.printSequence(*fastaFileNames[indexToFastaFile]);
611 currSeq.printSequence(outTemp);
615 //currQual.printQScores(*qualFileNames[indexToFastaFile]);
617 m->openOutputFileAppend(qualNames[indexToFastaFile], outTemp2);
618 currQual.printQScores(outTemp2);
624 if (groupfile != "") {
625 string thisGroup = groupMap->getGroup(currSeq.getName());
627 if (thisGroup != "not found") {
628 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
631 m->openOutputFileAppend(fastaNames[groupToIndex[thisGroup]], outTemp);
632 currSeq.printSequence(outTemp);
636 m->openOutputFileAppend(qualNames[groupToIndex[thisGroup]], outTemp2);
637 currQual.printQScores(outTemp2);
642 m->mothurOut(currSeq.getName() + " is not in your groupfile, adding to group XXX."); m->mothurOutEndLine();
643 outGroups << currSeq.getName() << '\t' << "XXX" << endl;
645 m->mothurOut("[ERROR]: " + currSeq.getName() + " will not be added to any .group.fasta or .group.qual file."); m->mothurOutEndLine();
651 currSeq.setName(currSeq.getName() + '|' + trashCode);
652 currSeq.setUnaligned(origSeq);
653 currSeq.setAligned(origSeq);
654 currSeq.printSequence(scrapFASTA);
655 currQual.printQScores(scrapQual);
660 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
661 unsigned long int pos = inFASTA.tellg();
662 if ((pos == -1) || (pos >= line->end)) { break; }
664 if (inFASTA.eof()) { break; }
668 if((count) % 1000 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
672 if((count) % 1000 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
678 if (oligoFile != "") { outGroups.close(); }
679 if(qFileName != "") { qFile.close(); scrapQual.close(); outQual.close(); }
683 catch(exception& e) {
684 m->errorOut(e, "TrimSeqsCommand", "driverCreateTrim");
689 /**************************************************************************************************/
691 int TrimSeqsCommand::createProcessesCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string trimQFile, string scrapQFile, string groupFile, vector<string> fastaNames, vector<string> qualNames) {
693 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
698 //loop through and create all the processes you want
699 while (process != processors) {
703 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
706 for (int i = 0; i < fastaNames.size(); i++) {
707 fastaNames[i] = (fastaNames[i] + toString(getpid()) + ".temp");
708 //clear old file if it exists
710 m->openOutputFile(fastaNames[i], temp);
713 qualNames[i] = (qualNames[i] + toString(getpid()) + ".temp");
714 //clear old file if it exists
716 m->openOutputFile(qualNames[i], temp2);
721 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]);
724 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
725 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
731 for (int i = 0; i < fastaNames.size(); i++) {
732 //clear old file if it exists
734 m->openOutputFile(fastaNames[i], temp);
737 //clear old file if it exists
739 m->openOutputFile(qualNames[i], temp2);
744 driverCreateTrim(filename, qFileName, trimFile, scrapFile, trimQFile, scrapQFile, groupFile, fastaNames, qualNames, lines[0], qLines[0]);
747 //force parent to wait until all the processes are done
748 for (int i=0;i<processIDS.size();i++) {
749 int temp = processIDS[i];
754 for(int i=0;i<processIDS.size();i++){
756 m->mothurOut("Appending files from process " + toString(processIDS[i])); m->mothurOutEndLine();
758 m->appendFiles((trimFile + toString(processIDS[i]) + ".temp"), trimFile);
759 remove((trimFile + toString(processIDS[i]) + ".temp").c_str());
760 m->appendFiles((scrapFile + toString(processIDS[i]) + ".temp"), scrapFile);
761 remove((scrapFile + toString(processIDS[i]) + ".temp").c_str());
763 m->mothurOut("Done with fasta files"); m->mothurOutEndLine();
766 m->appendFiles((trimQFile + toString(processIDS[i]) + ".temp"), trimQFile);
767 remove((trimQFile + toString(processIDS[i]) + ".temp").c_str());
768 m->appendFiles((scrapQFile + toString(processIDS[i]) + ".temp"), scrapQFile);
769 remove((scrapQFile + toString(processIDS[i]) + ".temp").c_str());
771 m->mothurOut("Done with quality files"); m->mothurOutEndLine();
774 m->appendFiles((groupFile + toString(processIDS[i]) + ".temp"), groupFile);
775 remove((groupFile + toString(processIDS[i]) + ".temp").c_str());
777 m->mothurOut("Done with group file"); m->mothurOutEndLine();
779 for (int j = 0; j < fastaNames.size(); j++) {
780 m->appendFiles((fastaNames[j] + toString(processIDS[i]) + ".temp"), fastaNames[j]);
781 remove((fastaNames[j] + toString(processIDS[i]) + ".temp").c_str());
785 for (int j = 0; j < qualNames.size(); j++) {
786 m->appendFiles((qualNames[j] + toString(processIDS[i]) + ".temp"), qualNames[j]);
787 remove((qualNames[j] + toString(processIDS[i]) + ".temp").c_str());
791 if (allFiles) { m->mothurOut("Done with allfiles"); m->mothurOutEndLine(); }
797 catch(exception& e) {
798 m->errorOut(e, "TrimSeqsCommand", "createProcessesCreateTrim");
803 /**************************************************************************************************/
805 int TrimSeqsCommand::setLines(string filename, string qfilename, vector<unsigned long int>& fastaFilePos, vector<unsigned long int>& qfileFilePos) {
808 //set file positions for fasta file
809 fastaFilePos = m->divideFile(filename, processors);
811 if (qfilename == "") { return processors; }
813 //get name of first sequence in each chunk
814 map<string, int> firstSeqNames;
815 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
817 m->openInputFile(filename, in);
818 in.seekg(fastaFilePos[i]);
821 firstSeqNames[temp.getName()] = i;
826 //seach for filePos of each first name in the qfile and save in qfileFilePos
828 m->openInputFile(qfilename, inQual);
831 while(!inQual.eof()){
832 input = m->getline(inQual);
834 if (input.length() != 0) {
835 if(input[0] == '>'){ //this is a sequence name line
836 istringstream nameStream(input);
838 string sname = ""; nameStream >> sname;
839 sname = sname.substr(1);
841 map<string, int>::iterator it = firstSeqNames.find(sname);
843 if(it != firstSeqNames.end()) { //this is the start of a new chunk
844 unsigned long int pos = inQual.tellg();
845 qfileFilePos.push_back(pos - input.length() - 1);
846 firstSeqNames.erase(it);
851 if (firstSeqNames.size() == 0) { break; }
856 if (firstSeqNames.size() != 0) {
857 for (map<string, int>::iterator it = firstSeqNames.begin(); it != firstSeqNames.end(); it++) {
858 m->mothurOut(it->first + " is in your fasta file and not in your quality file, not using quality file."); m->mothurOutEndLine();
864 //get last file position of qfile
866 unsigned long int size;
868 //get num bytes in file
869 pFile = fopen (qfilename.c_str(),"rb");
870 if (pFile==NULL) perror ("Error opening file");
872 fseek (pFile, 0, SEEK_END);
877 qfileFilePos.push_back(size);
881 catch(exception& e) {
882 m->errorOut(e, "TrimSeqsCommand", "setLines");
886 //***************************************************************************************************************
888 void TrimSeqsCommand::getOligos(vector<string>& outFASTAVec, vector<string>& outQualVec){
891 m->openInputFile(oligoFile, inOligos);
895 string type, oligo, group;
897 //int indexPrimer = 0;
899 while(!inOligos.eof()){
900 inOligos >> type; m->gobble(inOligos);
903 while (!inOligos.eof()) { char c = inOligos.get(); if (c == 10 || c == 13){ break; } } // get rest of line if there's any crap there
906 //make type case insensitive
907 for(int i=0;i<type.length();i++){ type[i] = toupper(type[i]); }
911 for(int i=0;i<oligo.length();i++){
912 oligo[i] = toupper(oligo[i]);
913 if(oligo[i] == 'U') { oligo[i] = 'T'; }
916 if(type == "FORWARD"){
919 // get rest of line in case there is a primer name
920 while (!inOligos.eof()) {
921 char c = inOligos.get();
922 if (c == 10 || c == 13){ break; }
923 else if (c == 32 || c == 9){;} //space or tab
927 //check for repeat barcodes
928 map<string, int>::iterator itPrime = primers.find(oligo);
929 if (itPrime != primers.end()) { m->mothurOut("primer " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
931 primers[oligo]=index; index++;
932 groupVector.push_back(group);
935 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
937 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
939 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
940 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
942 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
945 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
946 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
948 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
949 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
955 else if(type == "REVERSE"){
956 Sequence oligoRC("reverse", oligo);
957 oligoRC.reverseComplement();
958 revPrimer.push_back(oligoRC.getUnaligned());
960 else if(type == "BARCODE"){
963 //check for repeat barcodes
964 map<string, int>::iterator itBar = barcodes.find(oligo);
965 if (itBar != barcodes.end()) { m->mothurOut("barcode " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
967 barcodes[oligo]=index; index++;
968 groupVector.push_back(group);
971 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
972 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
973 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
975 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
976 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
977 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
981 }else{ m->mothurOut(type + " is not recognized as a valid type. Choices are forward, reverse, and barcode. Ignoring " + oligo + "."); m->mothurOutEndLine(); }
988 //add in potential combos
990 comboStarts = outFASTAVec.size()-1;
991 for (map<string, int>::iterator itBar = barcodes.begin(); itBar != barcodes.end(); itBar++) {
992 for (map<string, int>::iterator itPrime = primers.begin(); itPrime != primers.end(); itPrime++) {
993 if (groupVector[itPrime->second] != "") { //there is a group for this primer
994 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
995 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
996 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
997 combos[(groupVector[itBar->second] + "." + groupVector[itPrime->second])] = outFASTAVec.size()-1;
1000 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
1001 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
1002 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
1009 numFPrimers = primers.size();
1010 numRPrimers = revPrimer.size();
1013 catch(exception& e) {
1014 m->errorOut(e, "TrimSeqsCommand", "getOligos");
1018 //***************************************************************************************************************
1020 int TrimSeqsCommand::stripBarcode(Sequence& seq, QualityScores& qual, int& group){
1023 string rawSequence = seq.getUnaligned();
1024 int success = bdiffs + 1; //guilty until proven innocent
1026 //can you find the barcode
1027 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
1028 string oligo = it->first;
1029 if(rawSequence.length() < oligo.length()){ //let's just assume that the barcodes are the same length
1030 success = bdiffs + 10; //if the sequence is shorter than the barcode then bail out
1034 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1036 seq.setUnaligned(rawSequence.substr(oligo.length()));
1038 if(qual.getName() != ""){
1039 qual.trimQScores(oligo.length(), -1);
1047 //if you found the barcode or if you don't want to allow for diffs
1049 if ((bdiffs == 0) || (success == 0)) { return success; }
1051 else { //try aligning and see if you can find it
1056 Alignment* alignment;
1057 if (barcodes.size() > 0) {
1058 map<string,int>::iterator it=barcodes.begin();
1060 for(it;it!=barcodes.end();it++){
1061 if(it->first.length() > maxLength){
1062 maxLength = it->first.length();
1065 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+bdiffs+1));
1067 }else{ alignment = NULL; }
1069 //can you find the barcode
1075 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
1076 string oligo = it->first;
1077 // int length = oligo.length();
1079 if(rawSequence.length() < maxLength){ //let's just assume that the barcodes are the same length
1080 success = bdiffs + 10;
1084 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1085 alignment->align(oligo, rawSequence.substr(0,oligo.length()+bdiffs));
1086 oligo = alignment->getSeqAAln();
1087 string temp = alignment->getSeqBAln();
1089 int alnLength = oligo.length();
1091 for(int i=oligo.length()-1;i>=0;i--){
1092 if(oligo[i] != '-'){ alnLength = i+1; break; }
1094 oligo = oligo.substr(0,alnLength);
1095 temp = temp.substr(0,alnLength);
1098 int numDiff = countDiffs(oligo, temp);
1100 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1102 if(numDiff < minDiff){
1105 minGroup = it->second;
1107 for(int i=0;i<alnLength;i++){
1113 else if(numDiff == minDiff){
1119 if(minDiff > bdiffs) { success = minDiff; } //no good matches
1120 else if(minCount > 1) { success = bdiffs + 100; } //can't tell the difference between multiple barcodes
1121 else{ //use the best match
1123 seq.setUnaligned(rawSequence.substr(minPos));
1125 if(qual.getName() != ""){
1126 qual.trimQScores(minPos, -1);
1131 if (alignment != NULL) { delete alignment; }
1134 // cout << success << endl;
1139 catch(exception& e) {
1140 m->errorOut(e, "TrimSeqsCommand", "stripBarcode");
1146 //***************************************************************************************************************
1148 int TrimSeqsCommand::stripForward(Sequence& seq, QualityScores& qual, int& group){
1150 string rawSequence = seq.getUnaligned();
1151 int success = pdiffs + 1; //guilty until proven innocent
1153 //can you find the primer
1154 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1155 string oligo = it->first;
1156 if(rawSequence.length() < oligo.length()){ //let's just assume that the primers are the same length
1157 success = pdiffs + 10; //if the sequence is shorter than the barcode then bail out
1161 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1163 seq.setUnaligned(rawSequence.substr(oligo.length()));
1164 if(qual.getName() != ""){
1165 qual.trimQScores(oligo.length(), -1);
1172 //if you found the barcode or if you don't want to allow for diffs
1174 if ((pdiffs == 0) || (success == 0)) { return success; }
1176 else { //try aligning and see if you can find it
1181 Alignment* alignment;
1182 if (primers.size() > 0) {
1183 map<string,int>::iterator it=primers.begin();
1185 for(it;it!=primers.end();it++){
1186 if(it->first.length() > maxLength){
1187 maxLength = it->first.length();
1190 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+pdiffs+1));
1192 }else{ alignment = NULL; }
1194 //can you find the barcode
1200 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1201 string oligo = it->first;
1202 // int length = oligo.length();
1204 if(rawSequence.length() < maxLength){
1205 success = pdiffs + 100;
1209 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1210 alignment->align(oligo, rawSequence.substr(0,oligo.length()+pdiffs));
1211 oligo = alignment->getSeqAAln();
1212 string temp = alignment->getSeqBAln();
1214 int alnLength = oligo.length();
1216 for(int i=oligo.length()-1;i>=0;i--){
1217 if(oligo[i] != '-'){ alnLength = i+1; break; }
1219 oligo = oligo.substr(0,alnLength);
1220 temp = temp.substr(0,alnLength);
1223 int numDiff = countDiffs(oligo, temp);
1225 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1227 if(numDiff < minDiff){
1230 minGroup = it->second;
1232 for(int i=0;i<alnLength;i++){
1238 else if(numDiff == minDiff){
1244 if(minDiff > pdiffs) { success = minDiff; } //no good matches
1245 else if(minCount > 1) { success = pdiffs + 10; } //can't tell the difference between multiple primers
1246 else{ //use the best match
1248 seq.setUnaligned(rawSequence.substr(minPos));
1249 if(qual.getName() != ""){
1250 qual.trimQScores(minPos, -1);
1255 if (alignment != NULL) { delete alignment; }
1262 catch(exception& e) {
1263 m->errorOut(e, "TrimSeqsCommand", "stripForward");
1268 //***************************************************************************************************************
1270 bool TrimSeqsCommand::stripReverse(Sequence& seq, QualityScores& qual){
1272 string rawSequence = seq.getUnaligned();
1273 bool success = 0; //guilty until proven innocent
1275 for(int i=0;i<numRPrimers;i++){
1276 string oligo = revPrimer[i];
1278 if(rawSequence.length() < oligo.length()){
1283 if(compareDNASeq(oligo, rawSequence.substr(rawSequence.length()-oligo.length(),oligo.length()))){
1284 seq.setUnaligned(rawSequence.substr(0,rawSequence.length()-oligo.length()));
1285 if(qual.getName() != ""){
1286 qual.trimQScores(-1, rawSequence.length()-oligo.length());
1295 catch(exception& e) {
1296 m->errorOut(e, "TrimSeqsCommand", "stripReverse");
1301 //***************************************************************************************************************
1303 bool TrimSeqsCommand::keepFirstTrim(Sequence& sequence, QualityScores& qscores){
1306 if(qscores.getName() != ""){
1307 qscores.trimQScores(-1, keepFirst);
1309 sequence.trim(keepFirst);
1312 catch(exception& e) {
1313 m->errorOut(e, "keepFirstTrim", "countDiffs");
1319 //***************************************************************************************************************
1321 bool TrimSeqsCommand::removeLastTrim(Sequence& sequence, QualityScores& qscores){
1325 int length = sequence.getNumBases() - removeLast;
1328 if(qscores.getName() != ""){
1329 qscores.trimQScores(-1, length);
1331 sequence.trim(length);
1340 catch(exception& e) {
1341 m->errorOut(e, "removeLastTrim", "countDiffs");
1347 //***************************************************************************************************************
1349 bool TrimSeqsCommand::cullLength(Sequence& seq){
1352 int length = seq.getNumBases();
1353 bool success = 0; //guilty until proven innocent
1355 if(length >= minLength && maxLength == 0) { success = 1; }
1356 else if(length >= minLength && length <= maxLength) { success = 1; }
1357 else { success = 0; }
1362 catch(exception& e) {
1363 m->errorOut(e, "TrimSeqsCommand", "cullLength");
1369 //***************************************************************************************************************
1371 bool TrimSeqsCommand::cullHomoP(Sequence& seq){
1373 int longHomoP = seq.getLongHomoPolymer();
1374 bool success = 0; //guilty until proven innocent
1376 if(longHomoP <= maxHomoP){ success = 1; }
1377 else { success = 0; }
1381 catch(exception& e) {
1382 m->errorOut(e, "TrimSeqsCommand", "cullHomoP");
1388 //***************************************************************************************************************
1390 bool TrimSeqsCommand::cullAmbigs(Sequence& seq){
1392 int numNs = seq.getAmbigBases();
1393 bool success = 0; //guilty until proven innocent
1395 if(numNs <= maxAmbig) { success = 1; }
1396 else { success = 0; }
1400 catch(exception& e) {
1401 m->errorOut(e, "TrimSeqsCommand", "cullAmbigs");
1407 //***************************************************************************************************************
1409 bool TrimSeqsCommand::compareDNASeq(string oligo, string seq){
1412 int length = oligo.length();
1414 for(int i=0;i<length;i++){
1416 if(oligo[i] != seq[i]){
1417 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C') { success = 0; }
1418 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { success = 0; }
1419 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { success = 0; }
1420 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { success = 0; }
1421 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { success = 0; }
1422 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1423 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { success = 0; }
1424 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1425 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1426 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1427 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { success = 0; }
1428 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1430 if(success == 0) { break; }
1439 catch(exception& e) {
1440 m->errorOut(e, "TrimSeqsCommand", "compareDNASeq");
1446 //***************************************************************************************************************
1448 int TrimSeqsCommand::countDiffs(string oligo, string seq){
1451 int length = oligo.length();
1454 for(int i=0;i<length;i++){
1456 if(oligo[i] != seq[i]){
1457 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C' || oligo[i] == '-' || oligo[i] == '.') { countDiffs++; }
1458 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { countDiffs++; }
1459 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { countDiffs++; }
1460 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { countDiffs++; }
1461 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { countDiffs++; }
1462 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1463 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { countDiffs++; }
1464 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1465 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1466 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1467 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { countDiffs++; }
1468 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1475 catch(exception& e) {
1476 m->errorOut(e, "TrimSeqsCommand", "countDiffs");
1482 //***************************************************************************************************************