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", "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
18 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
22 m->errorOut(e, "TrimSeqsCommand", "getValidParameters");
26 //**********************************************************************************************************************
27 TrimSeqsCommand::TrimSeqsCommand(){
30 //initialize outputTypes
31 vector<string> tempOutNames;
32 outputTypes["fasta"] = tempOutNames;
33 outputTypes["qual"] = tempOutNames;
34 outputTypes["group"] = tempOutNames;
37 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
41 //**********************************************************************************************************************
42 vector<string> TrimSeqsCommand::getRequiredParameters(){
44 string Array[] = {"fasta"};
45 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
49 m->errorOut(e, "TrimSeqsCommand", "getRequiredParameters");
53 //**********************************************************************************************************************
54 vector<string> TrimSeqsCommand::getRequiredFiles(){
56 vector<string> myArray;
60 m->errorOut(e, "TrimSeqsCommand", "getRequiredFiles");
64 //***************************************************************************************************************
66 TrimSeqsCommand::TrimSeqsCommand(string option) {
72 //allow user to run help
73 if(option == "help") { help(); abort = true; }
76 //valid paramters for this command
77 string AlignArray[] = {"fasta", "flip", "group","oligos", "maxambig", "maxhomop", "minlength", "maxlength", "qfile",
78 "qthreshold", "qwindowaverage", "qstepsize", "qwindowsize", "qaverage", "rollaverage", "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
80 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
82 OptionParser parser(option);
83 map<string,string> parameters = parser.getParameters();
85 ValidParameters validParameter;
86 map<string,string>::iterator it;
88 //check to make sure all parameters are valid for command
89 for (it = parameters.begin(); it != parameters.end(); it++) {
90 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
93 //initialize outputTypes
94 vector<string> tempOutNames;
95 outputTypes["fasta"] = tempOutNames;
96 outputTypes["qual"] = tempOutNames;
97 outputTypes["group"] = tempOutNames;
99 //if the user changes the input directory command factory will send this info to us in the output parameter
100 string inputDir = validParameter.validFile(parameters, "inputdir", false);
101 if (inputDir == "not found"){ inputDir = ""; }
104 it = parameters.find("fasta");
105 //user has given a template file
106 if(it != parameters.end()){
107 path = m->hasPath(it->second);
108 //if the user has not given a path then, add inputdir. else leave path alone.
109 if (path == "") { parameters["fasta"] = inputDir + it->second; }
112 it = parameters.find("oligos");
113 //user has given a template file
114 if(it != parameters.end()){
115 path = m->hasPath(it->second);
116 //if the user has not given a path then, add inputdir. else leave path alone.
117 if (path == "") { parameters["oligos"] = inputDir + it->second; }
120 it = parameters.find("qfile");
121 //user has given a template file
122 if(it != parameters.end()){
123 path = m->hasPath(it->second);
124 //if the user has not given a path then, add inputdir. else leave path alone.
125 if (path == "") { parameters["qfile"] = inputDir + it->second; }
128 it = parameters.find("group");
129 //user has given a template file
130 if(it != parameters.end()){
131 path = m->hasPath(it->second);
132 //if the user has not given a path then, add inputdir. else leave path alone.
133 if (path == "") { parameters["group"] = inputDir + it->second; }
138 //check for required parameters
139 fastaFile = validParameter.validFile(parameters, "fasta", true);
140 if (fastaFile == "not found") { m->mothurOut("fasta is a required parameter for the trim.seqs command."); m->mothurOutEndLine(); abort = true; }
141 else if (fastaFile == "not open") { abort = true; }
143 //if the user changes the output directory command factory will send this info to us in the output parameter
144 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
146 outputDir += m->hasPath(fastaFile); //if user entered a file with a path then preserve it
150 //check for optional parameter and set defaults
151 // ...at some point should added some additional type checking...
153 temp = validParameter.validFile(parameters, "flip", false);
154 if (temp == "not found"){ flip = 0; }
155 else if(m->isTrue(temp)) { flip = 1; }
157 temp = validParameter.validFile(parameters, "oligos", true);
158 if (temp == "not found"){ oligoFile = ""; }
159 else if(temp == "not open"){ abort = true; }
160 else { oligoFile = temp; }
162 temp = validParameter.validFile(parameters, "group", true);
163 if (temp == "not found"){ groupfile = ""; }
164 else if(temp == "not open"){ abort = true; }
165 else { groupfile = temp; }
167 temp = validParameter.validFile(parameters, "maxambig", false); if (temp == "not found") { temp = "-1"; }
168 convert(temp, maxAmbig);
170 temp = validParameter.validFile(parameters, "maxhomop", false); if (temp == "not found") { temp = "0"; }
171 convert(temp, maxHomoP);
173 temp = validParameter.validFile(parameters, "minlength", false); if (temp == "not found") { temp = "0"; }
174 convert(temp, minLength);
176 temp = validParameter.validFile(parameters, "maxlength", false); if (temp == "not found") { temp = "0"; }
177 convert(temp, maxLength);
179 temp = validParameter.validFile(parameters, "bdiffs", false); if (temp == "not found") { temp = "0"; }
180 convert(temp, bdiffs);
182 temp = validParameter.validFile(parameters, "pdiffs", false); if (temp == "not found") { temp = "0"; }
183 convert(temp, pdiffs);
185 temp = validParameter.validFile(parameters, "tdiffs", false); if (temp == "not found") { int tempTotal = pdiffs + bdiffs; temp = toString(tempTotal); }
186 convert(temp, tdiffs);
188 if(tdiffs == 0){ tdiffs = bdiffs + pdiffs; }
190 temp = validParameter.validFile(parameters, "qfile", true);
191 if (temp == "not found") { qFileName = ""; }
192 else if(temp == "not open") { abort = true; }
193 else { qFileName = temp; }
195 temp = validParameter.validFile(parameters, "qthreshold", false); if (temp == "not found") { temp = "0"; }
196 convert(temp, qThreshold);
198 temp = validParameter.validFile(parameters, "qtrim", false); if (temp == "not found") { temp = "F"; }
199 qtrim = m->isTrue(temp);
201 temp = validParameter.validFile(parameters, "rollaverage", false); if (temp == "not found") { temp = "0"; }
202 convert(temp, qRollAverage);
204 temp = validParameter.validFile(parameters, "qwindowaverage", false);if (temp == "not found") { temp = "0"; }
205 convert(temp, qWindowAverage);
207 temp = validParameter.validFile(parameters, "qwindowsize", false); if (temp == "not found") { temp = "100"; }
208 convert(temp, qWindowSize);
210 temp = validParameter.validFile(parameters, "qstepsize", false); if (temp == "not found") { temp = "10"; }
211 convert(temp, qWindowStep);
213 temp = validParameter.validFile(parameters, "qaverage", false); if (temp == "not found") { temp = "0"; }
214 convert(temp, qAverage);
216 temp = validParameter.validFile(parameters, "allfiles", false); if (temp == "not found") { temp = "F"; }
217 allFiles = m->isTrue(temp);
219 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found") { temp = "1"; }
220 convert(temp, processors);
222 if ((oligoFile != "") && (groupfile != "")) {
223 m->mothurOut("You given both a oligos file and a groupfile, only one is allowed."); m->mothurOutEndLine(); abort = true;
227 if(allFiles && (oligoFile == "") && (groupfile == "")){
228 m->mothurOut("You selected allfiles, but didn't enter an oligos or group file. Ignoring the allfiles request."); m->mothurOutEndLine();
230 if((qAverage != 0 && qThreshold != 0) && qFileName == ""){
231 m->mothurOut("You didn't provide a quality file name, quality criteria will be ignored."); m->mothurOutEndLine();
235 if(!flip && oligoFile=="" && !maxLength && !minLength && (maxAmbig==-1) && !maxHomoP && qFileName == ""){
236 m->mothurOut("You didn't set any options... quiting command."); m->mothurOutEndLine();
242 catch(exception& e) {
243 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
247 //**********************************************************************************************************************
249 void TrimSeqsCommand::help(){
251 m->mothurOut("The trim.seqs command reads a fastaFile and creates .....\n");
252 m->mothurOut("The trim.seqs command parameters are fasta, flip, oligos, group, maxambig, maxhomop, minlength, maxlength, qfile, qthreshold, qaverage, diffs, qtrim and allfiles.\n");
253 m->mothurOut("The fasta parameter is required.\n");
254 m->mothurOut("The group parameter allows you to enter a group file for your fasta file.\n");
255 m->mothurOut("The flip parameter will output the reverse compliment of your trimmed sequence. The default is false.\n");
256 m->mothurOut("The oligos parameter .... The default is "".\n");
257 m->mothurOut("The maxambig parameter .... The default is -1.\n");
258 m->mothurOut("The maxhomop parameter .... The default is 0.\n");
259 m->mothurOut("The minlength parameter .... The default is 0.\n");
260 m->mothurOut("The maxlength parameter .... The default is 0.\n");
261 m->mothurOut("The tdiffs parameter is used to specify the total number of differences allowed in the sequence. The default is pdiffs + bdiffs.\n");
262 m->mothurOut("The bdiffs parameter is used to specify the number of differences allowed in the barcode. The default is 0.\n");
263 m->mothurOut("The pdiffs parameter is used to specify the number of differences allowed in the primer. The default is 0.\n");
264 m->mothurOut("The qfile parameter .....\n");
265 m->mothurOut("The qthreshold parameter .... The default is 0.\n");
266 m->mothurOut("The qaverage parameter .... The default is 0.\n");
267 m->mothurOut("The allfiles parameter .... The default is F.\n");
268 m->mothurOut("The qtrim parameter .... The default is F.\n");
269 m->mothurOut("The trim.seqs command should be in the following format: \n");
270 m->mothurOut("trim.seqs(fasta=yourFastaFile, flip=yourFlip, oligos=yourOligos, maxambig=yourMaxambig, \n");
271 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength) \n");
272 m->mothurOut("Example trim.seqs(fasta=abrecovery.fasta, flip=..., oligos=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
273 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n");
274 m->mothurOut("For more details please check out the wiki http://www.mothur.org/wiki/Trim.seqs .\n\n");
277 catch(exception& e) {
278 m->errorOut(e, "TrimSeqsCommand", "help");
284 //***************************************************************************************************************
286 TrimSeqsCommand::~TrimSeqsCommand(){ /* do nothing */ }
288 //***************************************************************************************************************
290 int TrimSeqsCommand::execute(){
293 if (abort == true) { return 0; }
295 numFPrimers = 0; //this needs to be initialized
297 vector<string> fastaFileNames;
298 vector<string> qualFileNames;
300 string trimSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.fasta";
301 outputNames.push_back(trimSeqFile); outputTypes["fasta"].push_back(trimSeqFile);
302 string scrapSeqFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.fasta";
303 outputNames.push_back(scrapSeqFile); outputTypes["fasta"].push_back(scrapSeqFile);
304 string trimQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "trim.qual";
305 string scrapQualFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "scrap.qual";
306 if (qFileName != "") { outputNames.push_back(trimQualFile); outputNames.push_back(scrapQualFile); outputTypes["qual"].push_back(trimQualFile); outputTypes["qual"].push_back(scrapQualFile); }
307 string groupFile = "";
308 if (groupfile == "") { groupFile = outputDir + m->getRootName(m->getSimpleName(fastaFile)) + "groups"; }
310 groupFile = outputDir + m->getRootName(m->getSimpleName(groupfile)) + "trim.groups";
311 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
312 groupMap = new GroupMap(groupfile);
316 for (int i = 0; i < groupMap->namesOfGroups.size(); i++) {
317 groupToIndex[groupMap->namesOfGroups[i]] = i;
318 groupVector.push_back(groupMap->namesOfGroups[i]);
319 fastaFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupMap->namesOfGroups[i] + ".fasta"));
321 //we append later, so we want to clear file
323 m->openOutputFile(fastaFileNames[i], outRemove);
326 qualFileNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupMap->namesOfGroups[i] + ".qual"));
328 m->openOutputFile(qualFileNames[i], outRemove2);
333 comboStarts = fastaFileNames.size()-1;
337 outputNames.push_back(groupFile); outputTypes["group"].push_back(groupFile);
338 getOligos(fastaFileNames, qualFileNames);
341 vector<unsigned long int> fastaFilePos;
342 vector<unsigned long int> qFilePos;
344 setLines(fastaFile, qFileName, fastaFilePos, qFilePos);
346 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
347 lines.push_back(new linePair(fastaFilePos[i], fastaFilePos[(i+1)]));
348 if (qFileName != "") { qLines.push_back(new linePair(qFilePos[i], qFilePos[(i+1)])); }
350 if(qFileName == "") { qLines = lines; } //files with duds
352 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
354 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
356 createProcessesCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames);
359 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, trimQualFile, scrapQualFile, groupFile, fastaFileNames, qualFileNames, lines[0], qLines[0]);
362 if (m->control_pressed) { return 0; }
365 for(int i=0;i<fastaFileNames.size();i++){
366 if (m->isBlank(fastaFileNames[i])) { blanks.insert(fastaFileNames[i]); }
367 else if (filesToRemove.count(fastaFileNames[i]) > 0) { remove(fastaFileNames[i].c_str()); }
371 m->openInputFile(fastaFileNames[i], inFASTA);
373 string outGroupFilename = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[i])) + "groups";
375 //if the fastafile is on the blanks list then the groups file should be as well
376 if (blanks.count(fastaFileNames[i]) != 0) { blanks.insert(outGroupFilename); }
378 m->openOutputFile(outGroupFilename, outGroups);
379 outputNames.push_back(outGroupFilename); outputTypes["group"].push_back(outGroupFilename);
381 string thisGroup = "";
382 if (i > comboStarts) {
383 map<string, int>::iterator itCombo;
384 for(itCombo=combos.begin();itCombo!=combos.end(); itCombo++){
385 if(itCombo->second == i){ thisGroup = itCombo->first; combos.erase(itCombo); break; }
387 }else{ thisGroup = groupVector[i]; }
389 while(!inFASTA.eof()){
390 if(inFASTA.get() == '>'){
392 outGroups << seqName << '\t' << thisGroup << endl;
394 while (!inFASTA.eof()) { char c = inFASTA.get(); if (c == 10 || c == 13){ break; } }
401 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
405 for(int i=0;i<qualFileNames.size();i++){
406 if (m->isBlank(qualFileNames[i])) { blanks.insert(qualFileNames[i]); }
407 else if (filesToRemove.count(qualFileNames[i]) > 0) { remove(qualFileNames[i].c_str()); }
411 for (set<string>::iterator itBlanks = blanks.begin(); itBlanks != blanks.end(); itBlanks++) { remove((*(itBlanks)).c_str()); }
413 if (m->control_pressed) {
414 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
418 m->mothurOutEndLine();
419 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
420 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
421 m->mothurOutEndLine();
426 catch(exception& e) {
427 m->errorOut(e, "TrimSeqsCommand", "execute");
432 /**************************************************************************************/
434 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) {
439 int able = m->openOutputFile(trimFile, outFASTA);
442 m->openOutputFile(scrapFile, scrapFASTA);
447 m->openOutputFile(trimQFile, outQual);
448 m->openOutputFile(scrapQFile, scrapQual);
453 if (oligoFile != "") {
454 m->openOutputFile(groupFile, outGroups);
458 m->openInputFile(filename, inFASTA);
459 inFASTA.seekg(line->start);
462 if(qFileName != "") { m->openInputFile(qFileName, qFile); qFile.seekg(qline->start); }
465 for (int i = 0; i < fastaNames.size(); i++) { //clears old file
467 m->openOutputFile(fastaNames[i], temp);
470 for (int i = 0; i < qualNames.size(); i++) { //clears old file
472 m->openOutputFile(qualNames[i], temp);
482 if (m->control_pressed) {
483 inFASTA.close(); outFASTA.close(); scrapFASTA.close();
484 if (oligoFile != "") { outGroups.close(); }
489 for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }
497 Sequence currSeq(inFASTA); m->gobble(inFASTA);
499 QualityScores currQual;
501 currQual = QualityScores(qFile, currSeq.getNumBases()); m->gobble(qFile);
504 string origSeq = currSeq.getUnaligned();
506 int groupBar, groupPrime;
507 string trashCode = "";
508 int currentSeqsDiffs = 0;
511 if(qThreshold != 0) { success = currQual.stripQualThreshold(currSeq, qThreshold); }
512 else if(qAverage != 0) { success = currQual.cullQualAverage(currSeq, qAverage); }
513 else if(qRollAverage != 0) { success = currQual.stripQualRollingAverage(currSeq, qRollAverage); }
514 else if(qWindowAverage != 0){ success = currQual.stripQualWindowAverage(currSeq, qWindowStep, qWindowSize, qWindowAverage); }
516 if (qtrim == 1 && (origSeq.length() != currSeq.getUnaligned().length())) {
517 success = 0; //if you don't want to trim and the sequence does not meet quality requirements, move to scrap
520 if(!success) { trashCode += 'q'; }
523 if(barcodes.size() != 0){
524 success = stripBarcode(currSeq, currQual, groupBar);
525 if(success > bdiffs) { trashCode += 'b'; }
526 else{ currentSeqsDiffs += success; }
529 if(numFPrimers != 0){
530 success = stripForward(currSeq, currQual, groupPrime);
531 if(success > pdiffs) { trashCode += 'f'; }
532 else{ currentSeqsDiffs += success; }
535 if (currentSeqsDiffs > tdiffs) { trashCode += 't'; }
537 if(numRPrimers != 0){
538 success = stripReverse(currSeq, currQual);
539 if(!success) { trashCode += 'r'; }
542 if(minLength > 0 || maxLength > 0){
543 success = cullLength(currSeq);
544 if(!success) { trashCode += 'l'; }
547 success = cullHomoP(currSeq);
548 if(!success) { trashCode += 'h'; }
551 success = cullAmbigs(currSeq);
552 if(!success) { trashCode += 'n'; }
555 if(flip){ currSeq.reverseComplement(); currQual.flipQScores(); } // should go last
557 if(trashCode.length() == 0){
558 currSeq.setAligned(currSeq.getUnaligned());
559 currSeq.printSequence(outFASTA);
560 currQual.printQScores(outQual);
562 if(barcodes.size() != 0){
563 string thisGroup = groupVector[groupBar];
564 int indexToFastaFile = groupBar;
565 if (primers.size() != 0){
566 //does this primer have a group
567 if (groupVector[groupPrime] != "") {
568 thisGroup += "." + groupVector[groupPrime];
569 indexToFastaFile = combos[thisGroup];
572 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
575 m->openOutputFileAppend(fastaNames[indexToFastaFile], outTemp);
576 //currSeq.printSequence(*fastaFileNames[indexToFastaFile]);
577 currSeq.printSequence(outTemp);
581 //currQual.printQScores(*qualFileNames[indexToFastaFile]);
583 m->openOutputFileAppend(qualNames[indexToFastaFile], outTemp2);
584 currQual.printQScores(outTemp2);
590 if (groupfile != "") {
591 string thisGroup = groupMap->getGroup(currSeq.getName());
593 if (thisGroup != "not found") {
594 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
597 m->openOutputFileAppend(fastaNames[groupToIndex[thisGroup]], outTemp);
598 currSeq.printSequence(outTemp);
602 m->openOutputFileAppend(qualNames[groupToIndex[thisGroup]], outTemp2);
603 currQual.printQScores(outTemp2);
608 m->mothurOut(currSeq.getName() + " is not in your groupfile, adding to group XXX."); m->mothurOutEndLine();
609 outGroups << currSeq.getName() << '\t' << "XXX" << endl;
611 m->mothurOut("[ERROR]: " + currSeq.getName() + " will not be added to any .group.fasta or .group.qual file."); m->mothurOutEndLine();
617 currSeq.setName(currSeq.getName() + '|' + trashCode);
618 currSeq.setUnaligned(origSeq);
619 currSeq.setAligned(origSeq);
620 currSeq.printSequence(scrapFASTA);
621 currQual.printQScores(scrapQual);
626 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
627 unsigned long int pos = inFASTA.tellg();
628 if ((pos == -1) || (pos >= line->end)) { break; }
630 if (inFASTA.eof()) { break; }
634 if((count) % 1000 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
638 if((count) % 1000 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
644 if (oligoFile != "") { outGroups.close(); }
645 if(qFileName != "") { qFile.close(); scrapQual.close(); outQual.close(); }
649 catch(exception& e) {
650 m->errorOut(e, "TrimSeqsCommand", "driverCreateTrim");
655 /**************************************************************************************************/
657 int TrimSeqsCommand::createProcessesCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string trimQFile, string scrapQFile, string groupFile, vector<string> fastaNames, vector<string> qualNames) {
659 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
664 //loop through and create all the processes you want
665 while (process != processors) {
669 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
672 for (int i = 0; i < fastaNames.size(); i++) {
673 fastaNames[i] = (fastaNames[i] + toString(getpid()) + ".temp");
674 //clear old file if it exists
676 m->openOutputFile(fastaNames[i], temp);
679 qualNames[i] = (qualNames[i] + toString(getpid()) + ".temp");
680 //clear old file if it exists
682 m->openOutputFile(qualNames[i], temp2);
687 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]);
690 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
691 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
697 for (int i = 0; i < fastaNames.size(); i++) {
698 //clear old file if it exists
700 m->openOutputFile(fastaNames[i], temp);
703 //clear old file if it exists
705 m->openOutputFile(qualNames[i], temp2);
710 driverCreateTrim(filename, qFileName, trimFile, scrapFile, trimQFile, scrapQFile, groupFile, fastaNames, qualNames, lines[0], qLines[0]);
713 //force parent to wait until all the processes are done
714 for (int i=0;i<processIDS.size();i++) {
715 int temp = processIDS[i];
720 for(int i=0;i<processIDS.size();i++){
721 m->mothurOut("Appending files from process " + processIDS[i]); m->mothurOutEndLine();
723 m->appendFiles((trimFile + toString(processIDS[i]) + ".temp"), trimFile);
724 remove((trimFile + toString(processIDS[i]) + ".temp").c_str());
725 m->appendFiles((scrapFile + toString(processIDS[i]) + ".temp"), scrapFile);
726 remove((scrapFile + toString(processIDS[i]) + ".temp").c_str());
728 m->mothurOut("Done with fasta files"); m->mothurOutEndLine();
731 m->appendFiles((trimQFile + toString(processIDS[i]) + ".temp"), trimQFile);
732 remove((trimQFile + toString(processIDS[i]) + ".temp").c_str());
733 m->appendFiles((scrapQFile + toString(processIDS[i]) + ".temp"), scrapQFile);
734 remove((scrapQFile + toString(processIDS[i]) + ".temp").c_str());
736 m->mothurOut("Done with quality files"); m->mothurOutEndLine();
739 m->appendFiles((groupFile + toString(processIDS[i]) + ".temp"), groupFile);
740 remove((groupFile + toString(processIDS[i]) + ".temp").c_str());
742 m->mothurOut("Done with group file"); m->mothurOutEndLine();
744 for (int j = 0; j < fastaNames.size(); j++) {
745 m->appendFiles((fastaNames[j] + toString(processIDS[i]) + ".temp"), fastaNames[j]);
746 remove((fastaNames[j] + toString(processIDS[i]) + ".temp").c_str());
750 for (int j = 0; j < qualNames.size(); j++) {
751 m->appendFiles((qualNames[j] + toString(processIDS[i]) + ".temp"), qualNames[j]);
752 remove((qualNames[j] + toString(processIDS[i]) + ".temp").c_str());
756 if (allFiles) { m->mothurOut("Done with allfiles"); m->mothurOutEndLine(); }
762 catch(exception& e) {
763 m->errorOut(e, "TrimSeqsCommand", "createProcessesCreateTrim");
768 /**************************************************************************************************/
770 int TrimSeqsCommand::setLines(string filename, string qfilename, vector<unsigned long int>& fastaFilePos, vector<unsigned long int>& qfileFilePos) {
773 //set file positions for fasta file
774 fastaFilePos = m->divideFile(filename, processors);
776 if (qfilename == "") { return processors; }
778 //get name of first sequence in each chunk
779 map<string, int> firstSeqNames;
780 for (int i = 0; i < (fastaFilePos.size()-1); i++) {
782 m->openInputFile(filename, in);
783 in.seekg(fastaFilePos[i]);
786 firstSeqNames[temp.getName()] = i;
791 //seach for filePos of each first name in the qfile and save in qfileFilePos
793 m->openInputFile(qfilename, inQual);
796 while(!inQual.eof()){
797 input = m->getline(inQual);
799 if (input.length() != 0) {
800 if(input[0] == '>'){ //this is a sequence name line
801 istringstream nameStream(input);
803 string sname = ""; nameStream >> sname;
804 sname = sname.substr(1);
806 map<string, int>::iterator it = firstSeqNames.find(sname);
808 if(it != firstSeqNames.end()) { //this is the start of a new chunk
809 unsigned long int pos = inQual.tellg();
810 qfileFilePos.push_back(pos - input.length() - 1);
811 firstSeqNames.erase(it);
816 if (firstSeqNames.size() == 0) { break; }
821 if (firstSeqNames.size() != 0) {
822 for (map<string, int>::iterator it = firstSeqNames.begin(); it != firstSeqNames.end(); it++) {
823 m->mothurOut(it->first + " is in your fasta file and not in your quality file, not using quality file."); m->mothurOutEndLine();
829 //get last file position of qfile
831 unsigned long int size;
833 //get num bytes in file
834 pFile = fopen (qfilename.c_str(),"rb");
835 if (pFile==NULL) perror ("Error opening file");
837 fseek (pFile, 0, SEEK_END);
842 qfileFilePos.push_back(size);
846 catch(exception& e) {
847 m->errorOut(e, "TrimSeqsCommand", "setLines");
851 //***************************************************************************************************************
853 void TrimSeqsCommand::getOligos(vector<string>& outFASTAVec, vector<string>& outQualVec){
856 m->openInputFile(oligoFile, inOligos);
860 string type, oligo, group;
862 //int indexPrimer = 0;
864 while(!inOligos.eof()){
865 inOligos >> type; m->gobble(inOligos);
868 while (!inOligos.eof()) { char c = inOligos.get(); if (c == 10 || c == 13){ break; } } // get rest of line if there's any crap there
871 //make type case insensitive
872 for(int i=0;i<type.length();i++){ type[i] = toupper(type[i]); }
876 for(int i=0;i<oligo.length();i++){
877 oligo[i] = toupper(oligo[i]);
878 if(oligo[i] == 'U') { oligo[i] = 'T'; }
881 if(type == "FORWARD"){
884 // get rest of line in case there is a primer name
885 while (!inOligos.eof()) {
886 char c = inOligos.get();
887 if (c == 10 || c == 13){ break; }
888 else if (c == 32 || c == 9){;} //space or tab
892 //check for repeat barcodes
893 map<string, int>::iterator itPrime = primers.find(oligo);
894 if (itPrime != primers.end()) { m->mothurOut("primer " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
896 primers[oligo]=index; index++;
897 groupVector.push_back(group);
900 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
902 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
904 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
905 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
907 filesToRemove.insert((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
910 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
911 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
913 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
914 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
920 else if(type == "REVERSE"){
921 Sequence oligoRC("reverse", oligo);
922 oligoRC.reverseComplement();
923 revPrimer.push_back(oligoRC.getUnaligned());
925 else if(type == "BARCODE"){
928 //check for repeat barcodes
929 map<string, int>::iterator itBar = barcodes.find(oligo);
930 if (itBar != barcodes.end()) { m->mothurOut("barcode " + oligo + " is in your oligos file already."); m->mothurOutEndLine(); }
932 barcodes[oligo]=index; index++;
933 groupVector.push_back(group);
936 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
937 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
938 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + group + ".fasta"));
940 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
941 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
942 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + group + ".qual"));
946 }else{ m->mothurOut(type + " is not recognized as a valid type. Choices are forward, reverse, and barcode. Ignoring " + oligo + "."); m->mothurOutEndLine(); }
953 //add in potential combos
955 comboStarts = outFASTAVec.size()-1;
956 for (map<string, int>::iterator itBar = barcodes.begin(); itBar != barcodes.end(); itBar++) {
957 for (map<string, int>::iterator itPrime = primers.begin(); itPrime != primers.end(); itPrime++) {
958 if (groupVector[itPrime->second] != "") { //there is a group for this primer
959 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
960 outputTypes["fasta"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
961 outFASTAVec.push_back((outputDir + m->getRootName(m->getSimpleName(fastaFile)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".fasta"));
962 combos[(groupVector[itBar->second] + "." + groupVector[itPrime->second])] = outFASTAVec.size()-1;
965 outQualVec.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
966 outputNames.push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
967 outputTypes["qual"].push_back((outputDir + m->getRootName(m->getSimpleName(qFileName)) + groupVector[itBar->second] + "." + groupVector[itPrime->second] + ".qual"));
974 numFPrimers = primers.size();
975 numRPrimers = revPrimer.size();
978 catch(exception& e) {
979 m->errorOut(e, "TrimSeqsCommand", "getOligos");
983 //***************************************************************************************************************
985 int TrimSeqsCommand::stripBarcode(Sequence& seq, QualityScores& qual, int& group){
988 string rawSequence = seq.getUnaligned();
989 int success = bdiffs + 1; //guilty until proven innocent
991 //can you find the barcode
992 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
993 string oligo = it->first;
994 if(rawSequence.length() < oligo.length()){ //let's just assume that the barcodes are the same length
995 success = bdiffs + 10; //if the sequence is shorter than the barcode then bail out
999 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1001 seq.setUnaligned(rawSequence.substr(oligo.length()));
1003 if(qual.getName() != ""){
1004 qual.trimQScores(oligo.length(), -1);
1012 //if you found the barcode or if you don't want to allow for diffs
1014 if ((bdiffs == 0) || (success == 0)) { return success; }
1016 else { //try aligning and see if you can find it
1021 Alignment* alignment;
1022 if (barcodes.size() > 0) {
1023 map<string,int>::iterator it=barcodes.begin();
1025 for(it;it!=barcodes.end();it++){
1026 if(it->first.length() > maxLength){
1027 maxLength = it->first.length();
1030 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+bdiffs+1));
1032 }else{ alignment = NULL; }
1034 //can you find the barcode
1040 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
1041 string oligo = it->first;
1042 // int length = oligo.length();
1044 if(rawSequence.length() < maxLength){ //let's just assume that the barcodes are the same length
1045 success = bdiffs + 10;
1049 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1050 alignment->align(oligo, rawSequence.substr(0,oligo.length()+bdiffs));
1051 oligo = alignment->getSeqAAln();
1052 string temp = alignment->getSeqBAln();
1054 int alnLength = oligo.length();
1056 for(int i=oligo.length()-1;i>=0;i--){
1057 if(oligo[i] != '-'){ alnLength = i+1; break; }
1059 oligo = oligo.substr(0,alnLength);
1060 temp = temp.substr(0,alnLength);
1063 int numDiff = countDiffs(oligo, temp);
1065 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1067 if(numDiff < minDiff){
1070 minGroup = it->second;
1072 for(int i=0;i<alnLength;i++){
1078 else if(numDiff == minDiff){
1084 if(minDiff > bdiffs) { success = minDiff; } //no good matches
1085 else if(minCount > 1) { success = bdiffs + 100; } //can't tell the difference between multiple barcodes
1086 else{ //use the best match
1088 seq.setUnaligned(rawSequence.substr(minPos));
1090 if(qual.getName() != ""){
1091 qual.trimQScores(minPos, -1);
1096 if (alignment != NULL) { delete alignment; }
1099 // cout << success << endl;
1104 catch(exception& e) {
1105 m->errorOut(e, "TrimSeqsCommand", "stripBarcode");
1111 //***************************************************************************************************************
1113 int TrimSeqsCommand::stripForward(Sequence& seq, QualityScores& qual, int& group){
1115 string rawSequence = seq.getUnaligned();
1116 int success = pdiffs + 1; //guilty until proven innocent
1118 //can you find the primer
1119 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1120 string oligo = it->first;
1121 if(rawSequence.length() < oligo.length()){ //let's just assume that the primers are the same length
1122 success = pdiffs + 10; //if the sequence is shorter than the barcode then bail out
1126 if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
1128 seq.setUnaligned(rawSequence.substr(oligo.length()));
1129 if(qual.getName() != ""){
1130 qual.trimQScores(oligo.length(), -1);
1138 //if you found the barcode or if you don't want to allow for diffs
1140 if ((pdiffs == 0) || (success == 0)) { return success; }
1142 else { //try aligning and see if you can find it
1147 Alignment* alignment;
1148 if (primers.size() > 0) {
1149 map<string,int>::iterator it=primers.begin();
1151 for(it;it!=primers.end();it++){
1152 if(it->first.length() > maxLength){
1153 maxLength = it->first.length();
1156 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+pdiffs+1));
1158 }else{ alignment = NULL; }
1160 //can you find the barcode
1166 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
1167 string oligo = it->first;
1168 // int length = oligo.length();
1170 if(rawSequence.length() < maxLength){
1171 success = pdiffs + 100;
1175 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
1176 alignment->align(oligo, rawSequence.substr(0,oligo.length()+pdiffs));
1177 oligo = alignment->getSeqAAln();
1178 string temp = alignment->getSeqBAln();
1180 int alnLength = oligo.length();
1182 for(int i=oligo.length()-1;i>=0;i--){
1183 if(oligo[i] != '-'){ alnLength = i+1; break; }
1185 oligo = oligo.substr(0,alnLength);
1186 temp = temp.substr(0,alnLength);
1189 int numDiff = countDiffs(oligo, temp);
1191 // cout << oligo << '\t' << temp << '\t' << numDiff << endl;
1193 if(numDiff < minDiff){
1196 minGroup = it->second;
1198 for(int i=0;i<alnLength;i++){
1204 else if(numDiff == minDiff){
1210 if(minDiff > pdiffs) { success = minDiff; } //no good matches
1211 else if(minCount > 1) { success = pdiffs + 10; } //can't tell the difference between multiple primers
1212 else{ //use the best match
1214 seq.setUnaligned(rawSequence.substr(minPos));
1215 if(qual.getName() != ""){
1216 qual.trimQScores(minPos, -1);
1221 if (alignment != NULL) { delete alignment; }
1228 catch(exception& e) {
1229 m->errorOut(e, "TrimSeqsCommand", "stripForward");
1234 //***************************************************************************************************************
1236 bool TrimSeqsCommand::stripReverse(Sequence& seq, QualityScores& qual){
1238 string rawSequence = seq.getUnaligned();
1239 bool success = 0; //guilty until proven innocent
1241 for(int i=0;i<numRPrimers;i++){
1242 string oligo = revPrimer[i];
1244 if(rawSequence.length() < oligo.length()){
1249 if(compareDNASeq(oligo, rawSequence.substr(rawSequence.length()-oligo.length(),oligo.length()))){
1250 seq.setUnaligned(rawSequence.substr(0,rawSequence.length()-oligo.length()));
1251 if(qual.getName() != ""){
1252 qual.trimQScores(-1, rawSequence.length()-oligo.length());
1261 catch(exception& e) {
1262 m->errorOut(e, "TrimSeqsCommand", "stripReverse");
1267 //***************************************************************************************************************
1269 bool TrimSeqsCommand::cullLength(Sequence& seq){
1272 int length = seq.getNumBases();
1273 bool success = 0; //guilty until proven innocent
1275 if(length >= minLength && maxLength == 0) { success = 1; }
1276 else if(length >= minLength && length <= maxLength) { success = 1; }
1277 else { success = 0; }
1282 catch(exception& e) {
1283 m->errorOut(e, "TrimSeqsCommand", "cullLength");
1289 //***************************************************************************************************************
1291 bool TrimSeqsCommand::cullHomoP(Sequence& seq){
1293 int longHomoP = seq.getLongHomoPolymer();
1294 bool success = 0; //guilty until proven innocent
1296 if(longHomoP <= maxHomoP){ success = 1; }
1297 else { success = 0; }
1301 catch(exception& e) {
1302 m->errorOut(e, "TrimSeqsCommand", "cullHomoP");
1308 //***************************************************************************************************************
1310 bool TrimSeqsCommand::cullAmbigs(Sequence& seq){
1312 int numNs = seq.getAmbigBases();
1313 bool success = 0; //guilty until proven innocent
1315 if(numNs <= maxAmbig) { success = 1; }
1316 else { success = 0; }
1320 catch(exception& e) {
1321 m->errorOut(e, "TrimSeqsCommand", "cullAmbigs");
1327 //***************************************************************************************************************
1329 bool TrimSeqsCommand::compareDNASeq(string oligo, string seq){
1332 int length = oligo.length();
1334 for(int i=0;i<length;i++){
1336 if(oligo[i] != seq[i]){
1337 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C') { success = 0; }
1338 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { success = 0; }
1339 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { success = 0; }
1340 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { success = 0; }
1341 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { success = 0; }
1342 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1343 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { success = 0; }
1344 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1345 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1346 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { success = 0; }
1347 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { success = 0; }
1348 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { success = 0; }
1350 if(success == 0) { break; }
1359 catch(exception& e) {
1360 m->errorOut(e, "TrimSeqsCommand", "compareDNASeq");
1365 //***************************************************************************************************************
1367 int TrimSeqsCommand::countDiffs(string oligo, string seq){
1370 int length = oligo.length();
1373 for(int i=0;i<length;i++){
1375 if(oligo[i] != seq[i]){
1376 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C' || oligo[i] == '-' || oligo[i] == '.') { countDiffs++; }
1377 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N')) { countDiffs++; }
1378 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G')) { countDiffs++; }
1379 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T')) { countDiffs++; }
1380 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A')) { countDiffs++; }
1381 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1382 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A')) { countDiffs++; }
1383 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1384 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1385 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G')) { countDiffs++; }
1386 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C')) { countDiffs++; }
1387 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G')) { countDiffs++; }
1394 catch(exception& e) {
1395 m->errorOut(e, "TrimSeqsCommand", "countDiffs");
1400 //***************************************************************************************************************
1402 //bool TrimSeqsCommand::stripQualThreshold(Sequence& seq, ifstream& qFile){
1405 // string rawSequence = seq.getUnaligned();
1406 // int seqLength = seq.getNumBases();
1407 // bool success = 0; //guilty until proven innocent
1411 // if (name[0] == '>') { if(name.substr(1) != seq.getName()) { m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); } }
1413 // while (!qFile.eof()) { char c = qFile.get(); if (c == 10 || c == 13){ break; } }
1416 // int end = seqLength;
1418 // for(int i=0;i<seqLength;i++){
1421 // if(score < qThreshold){
1426 // for(int i=end+1;i<seqLength;i++){
1430 // seq.setUnaligned(rawSequence.substr(0,end));
1434 // catch(exception& e) {
1435 // m->errorOut(e, "TrimSeqsCommand", "stripQualThreshold");
1440 //***************************************************************************************************************
1442 //bool TrimSeqsCommand::cullQualAverage(Sequence& seq, ifstream& qFile){
1444 // string rawSequence = seq.getUnaligned();
1445 // int seqLength = seq.getNumBases();
1446 // bool success = 0; //guilty until proven innocent
1450 // if (name[0] == '>') { if(name.substr(1) != seq.getName()) { m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); } }
1452 // while (!qFile.eof()) { char c = qFile.get(); if (c == 10 || c == 13){ break; } }
1455 // float average = 0;
1457 // for(int i=0;i<seqLength;i++){
1459 // average += score;
1461 // average /= seqLength;
1463 // if(average >= qAverage) { success = 1; }
1464 // else { success = 0; }
1468 // catch(exception& e) {
1469 // m->errorOut(e, "TrimSeqsCommand", "cullQualAverage");
1474 //***************************************************************************************************************