]> git.donarmstrong.com Git - mothur.git/blob - trimseqscommand.cpp
cbf5e13603ced805d185196fad9129eacb573d22
[mothur.git] / trimseqscommand.cpp
1 /*
2  *  trimseqscommand.cpp
3  *  Mothur
4  *
5  *  Created by Pat Schloss on 6/6/09.
6  *  Copyright 2009 Patrick D. Schloss. All rights reserved.
7  *
8  */
9
10 #include "trimseqscommand.h"
11 #include "needlemanoverlap.hpp"
12
13 //***************************************************************************************************************
14
15 TrimSeqsCommand::TrimSeqsCommand(string option)  {
16         try {
17                 
18                 abort = false;
19                 comboStarts = 0;
20                 
21                 //allow user to run help
22                 if(option == "help") { help(); abort = true; }
23                 
24                 else {
25                         //valid paramters for this command
26                         string AlignArray[] =  {"fasta", "flip", "oligos", "maxambig", "maxhomop", "minlength", "maxlength", "qfile", 
27                                                                         "qthreshold", "qaverage", "allfiles", "qtrim","tdiffs", "pdiffs", "bdiffs", "processors", "outputdir","inputdir"};
28                         
29                         vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
30                         
31                         OptionParser parser(option);
32                         map<string,string> parameters = parser.getParameters();
33                         
34                         ValidParameters validParameter;
35                         map<string,string>::iterator it;
36                         
37                         //check to make sure all parameters are valid for command
38                         for (it = parameters.begin(); it != parameters.end(); it++) { 
39                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
40                         }
41                         
42                         //if the user changes the input directory command factory will send this info to us in the output parameter 
43                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
44                         if (inputDir == "not found"){   inputDir = "";          }
45                         else {
46                                 string path;
47                                 it = parameters.find("fasta");
48                                 //user has given a template file
49                                 if(it != parameters.end()){ 
50                                         path = hasPath(it->second);
51                                         //if the user has not given a path then, add inputdir. else leave path alone.
52                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
53                                 }
54                                 
55                                 it = parameters.find("oligos");
56                                 //user has given a template file
57                                 if(it != parameters.end()){ 
58                                         path = hasPath(it->second);
59                                         //if the user has not given a path then, add inputdir. else leave path alone.
60                                         if (path == "") {       parameters["oligos"] = inputDir + it->second;           }
61                                 }
62                                 
63                                 it = parameters.find("qfile");
64                                 //user has given a template file
65                                 if(it != parameters.end()){ 
66                                         path = hasPath(it->second);
67                                         //if the user has not given a path then, add inputdir. else leave path alone.
68                                         if (path == "") {       parameters["qfile"] = inputDir + it->second;            }
69                                 }
70                         }
71
72                         
73                         //check for required parameters
74                         fastaFile = validParameter.validFile(parameters, "fasta", true);
75                         if (fastaFile == "not found") { m->mothurOut("fasta is a required parameter for the screen.seqs command."); m->mothurOutEndLine(); abort = true; }
76                         else if (fastaFile == "not open") { abort = true; }     
77                         
78                         //if the user changes the output directory command factory will send this info to us in the output parameter 
79                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
80                                 outputDir = ""; 
81                                 outputDir += hasPath(fastaFile); //if user entered a file with a path then preserve it  
82                         }
83                 
84                         //check for optional parameter and set defaults
85                         // ...at some point should added some additional type checking...
86                         string temp;
87                         temp = validParameter.validFile(parameters, "flip", false);
88                         if (temp == "not found"){       flip = 0;       }
89                         else if(isTrue(temp))   {       flip = 1;       }
90                 
91                         temp = validParameter.validFile(parameters, "oligos", true);
92                         if (temp == "not found"){       oligoFile = "";         }
93                         else if(temp == "not open"){    abort = true;   } 
94                         else                                    {       oligoFile = temp;               }
95                         
96                         temp = validParameter.validFile(parameters, "maxambig", false);         if (temp == "not found") { temp = "-1"; }
97                         convert(temp, maxAmbig);  
98
99                         temp = validParameter.validFile(parameters, "maxhomop", false);         if (temp == "not found") { temp = "0"; }
100                         convert(temp, maxHomoP);  
101
102                         temp = validParameter.validFile(parameters, "minlength", false);        if (temp == "not found") { temp = "0"; }
103                         convert(temp, minLength); 
104                         
105                         temp = validParameter.validFile(parameters, "maxlength", false);        if (temp == "not found") { temp = "0"; }
106                         convert(temp, maxLength);
107                         
108                         
109                         temp = validParameter.validFile(parameters, "bdiffs", false);           if (temp == "not found") { temp = "0"; }
110                         convert(temp, bdiffs);
111                         
112                         temp = validParameter.validFile(parameters, "pdiffs", false);           if (temp == "not found") { temp = "0"; }
113                         convert(temp, pdiffs);
114                         
115                         temp = validParameter.validFile(parameters, "tdiffs", false);           if (temp == "not found") { int tempTotal = pdiffs + bdiffs;  temp = toString(tempTotal); }
116                         convert(temp, tdiffs);
117                         
118                         if(tdiffs == 0){        tdiffs = bdiffs + pdiffs;       }
119                         
120                         temp = validParameter.validFile(parameters, "qfile", true);     
121                         if (temp == "not found")        {       qFileName = "";         }
122                         else if(temp == "not open")     {       abort = true;           }
123                         else                                            {       qFileName = temp;       }
124                         
125                         temp = validParameter.validFile(parameters, "qthreshold", false);       if (temp == "not found") { temp = "0"; }
126                         convert(temp, qThreshold);
127                         
128                         temp = validParameter.validFile(parameters, "qtrim", false);    if (temp == "not found") { temp = "F"; }
129                         qtrim = isTrue(temp);
130
131                         temp = validParameter.validFile(parameters, "qaverage", false);         if (temp == "not found") { temp = "0"; }
132                         convert(temp, qAverage);
133                         
134                         temp = validParameter.validFile(parameters, "allfiles", false);         if (temp == "not found") { temp = "F"; }
135                         allFiles = isTrue(temp);
136                         
137                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
138                         convert(temp, processors); 
139                         
140                         if(allFiles && oligoFile == ""){
141                                 m->mothurOut("You selected allfiles, but didn't enter an oligos file.  Ignoring the allfiles request."); m->mothurOutEndLine();
142                         }
143                         if((qAverage != 0 && qThreshold != 0) && qFileName == ""){
144                                 m->mothurOut("You didn't provide a quality file name, quality criteria will be ignored."); m->mothurOutEndLine();
145                                 qAverage=0;
146                                 qThreshold=0;
147                         }
148                         if(!flip && oligoFile=="" && !maxLength && !minLength && (maxAmbig==-1) && !maxHomoP && qFileName == ""){               
149                                 m->mothurOut("You didn't set any options... quiting command."); m->mothurOutEndLine();
150                                 abort = true;
151                         }
152                 }
153
154         }
155         catch(exception& e) {
156                 m->errorOut(e, "TrimSeqsCommand", "TrimSeqsCommand");
157                 exit(1);
158         }
159 }
160 //**********************************************************************************************************************
161
162 void TrimSeqsCommand::help(){
163         try {
164                 m->mothurOut("The trim.seqs command reads a fastaFile and creates .....\n");
165                 m->mothurOut("The trim.seqs command parameters are fasta, flip, oligos, maxambig, maxhomop, minlength, maxlength, qfile, qthreshold, qaverage, diffs, qtrim and allfiles.\n");
166                 m->mothurOut("The fasta parameter is required.\n");
167                 m->mothurOut("The flip parameter will output the reverse compliment of your trimmed sequence. The default is false.\n");
168                 m->mothurOut("The oligos parameter .... The default is "".\n");
169                 m->mothurOut("The maxambig parameter .... The default is -1.\n");
170                 m->mothurOut("The maxhomop parameter .... The default is 0.\n");
171                 m->mothurOut("The minlength parameter .... The default is 0.\n");
172                 m->mothurOut("The maxlength parameter .... The default is 0.\n");
173                 m->mothurOut("The tdiffs parameter is used to specify the total number of differences allowed in the sequence. The default is pdiffs + bdiffs.\n");
174                 m->mothurOut("The bdiffs parameter is used to specify the number of differences allowed in the barcode. The default is 0.\n");
175                 m->mothurOut("The pdiffs parameter is used to specify the number of differences allowed in the primer. The default is 0.\n");
176                 m->mothurOut("The qfile parameter .....\n");
177                 m->mothurOut("The qthreshold parameter .... The default is 0.\n");
178                 m->mothurOut("The qaverage parameter .... The default is 0.\n");
179                 m->mothurOut("The allfiles parameter .... The default is F.\n");
180                 m->mothurOut("The qtrim parameter .... The default is F.\n");
181                 m->mothurOut("The trim.seqs command should be in the following format: \n");
182                 m->mothurOut("trim.seqs(fasta=yourFastaFile, flip=yourFlip, oligos=yourOligos, maxambig=yourMaxambig,  \n");
183                 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength)  \n");    
184                 m->mothurOut("Example trim.seqs(fasta=abrecovery.fasta, flip=..., oligos=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
185                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n");
186                 m->mothurOut("For more details please check out the wiki http://www.mothur.org/wiki/Trim.seqs .\n\n");
187
188         }
189         catch(exception& e) {
190                 m->errorOut(e, "TrimSeqsCommand", "help");
191                 exit(1);
192         }
193 }
194
195
196 //***************************************************************************************************************
197
198 TrimSeqsCommand::~TrimSeqsCommand(){    /*      do nothing      */      }
199
200 //***************************************************************************************************************
201
202 int TrimSeqsCommand::execute(){
203         try{
204         
205                 if (abort == true) { return 0; }
206                 
207                 numFPrimers = 0;  //this needs to be initialized
208                 numRPrimers = 0;
209                 
210                 string trimSeqFile = outputDir + getRootName(getSimpleName(fastaFile)) + "trim.fasta";
211                 outputNames.push_back(trimSeqFile);
212                 string scrapSeqFile = outputDir + getRootName(getSimpleName(fastaFile)) + "scrap.fasta";
213                 outputNames.push_back(scrapSeqFile);
214                 string groupFile = outputDir + getRootName(getSimpleName(fastaFile)) + "groups";
215                 
216                 vector<string> fastaFileNames;
217                 if(oligoFile != ""){
218                         outputNames.push_back(groupFile);
219                         getOligos(fastaFileNames);
220                 }
221                 
222                 if(qFileName != "")     {       setLines(qFileName, qLines);    }
223
224                 
225                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
226                                 if(processors == 1){
227                                         ifstream inFASTA;
228                                         int numSeqs;
229                                         openInputFile(fastaFile, inFASTA);
230                                         getNumSeqs(inFASTA, numSeqs);
231                                         inFASTA.close();
232                                         
233                                         lines.push_back(new linePair(0, numSeqs));
234                                         
235                                         driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames, lines[0], lines[0]);
236                                         
237                                         for (int j = 0; j < fastaFileNames.size(); j++) {
238                                                 rename((fastaFileNames[j] + toString(getpid()) + ".temp").c_str(), fastaFileNames[j].c_str());
239                                         }
240
241                                 }else{
242                                         setLines(fastaFile, lines);     
243                                         if(qFileName == "")     {       qLines = lines; }       
244                                                                 
245                                         createProcessesCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames); 
246                                         
247                                         rename((trimSeqFile + toString(processIDS[0]) + ".temp").c_str(), trimSeqFile.c_str());
248                                         rename((scrapSeqFile + toString(processIDS[0]) + ".temp").c_str(), scrapSeqFile.c_str());
249                                         rename((groupFile + toString(processIDS[0]) + ".temp").c_str(), groupFile.c_str());
250                                         for (int j = 0; j < fastaFileNames.size(); j++) {
251                                                 rename((fastaFileNames[j] + toString(processIDS[0]) + ".temp").c_str(), fastaFileNames[j].c_str());
252                                         }
253                                         //append files
254                                         for(int i=1;i<processors;i++){
255                                                 appendFiles((trimSeqFile + toString(processIDS[i]) + ".temp"), trimSeqFile);
256                                                 remove((trimSeqFile + toString(processIDS[i]) + ".temp").c_str());
257                                                 appendFiles((scrapSeqFile + toString(processIDS[i]) + ".temp"), scrapSeqFile);
258                                                 remove((scrapSeqFile + toString(processIDS[i]) + ".temp").c_str());
259                                                 appendFiles((groupFile + toString(processIDS[i]) + ".temp"), groupFile);
260                                                 remove((groupFile + toString(processIDS[i]) + ".temp").c_str());
261                                                 for (int j = 0; j < fastaFileNames.size(); j++) {
262                                                         appendFiles((fastaFileNames[j] + toString(processIDS[i]) + ".temp"), fastaFileNames[j]);
263                                                         remove((fastaFileNames[j] + toString(processIDS[i]) + ".temp").c_str());
264                                                 }
265                                         }
266                                 }
267                                 
268                                 if (m->control_pressed) {  return 0; }
269                 #else
270                                 ifstream inFASTA;
271                                 int numSeqs;
272                                 openInputFile(fastaFile, inFASTA);
273                                 getNumSeqs(inFASTA, numSeqs);
274                                 inFASTA.close();
275                                 
276                                 lines.push_back(new linePair(0, numSeqs));
277                                 
278                                 driverCreateTrim(fastaFile, qFileName, trimSeqFile, scrapSeqFile, groupFile, fastaFileNames, lines[0], lines[0]);
279                                 
280                                 if (m->control_pressed) {  return 0; }
281                 #endif
282                                                 
283                                                                                 
284                 for(int i=0;i<fastaFileNames.size();i++){
285                         if (isBlank(fastaFileNames[i])) { remove(fastaFileNames[i].c_str()); }
286                         else if (filesToRemove.count(fastaFileNames[i]) > 0) { remove(fastaFileNames[i].c_str()); }
287                         else {
288                                 ifstream inFASTA;
289                                 string seqName;
290                                 //openInputFile(getRootName(fastaFile) +  groupVector[i] + ".fasta", inFASTA);
291                                 openInputFile(fastaFileNames[i], inFASTA);
292                                 ofstream outGroups;
293                                 string outGroupFilename = outputDir + getRootName(getSimpleName(fastaFileNames[i])) + "groups";
294                                 openOutputFile(outGroupFilename, outGroups);
295                                 //openOutputFile(outputDir + getRootName(getSimpleName(fastaFile)) + groupVector[i] + ".groups", outGroups);
296                                 outputNames.push_back(outGroupFilename);
297                                 
298                                 string thisGroup = "";
299                                 if (i > comboStarts) {
300                                         map<string, int>::iterator itCombo;
301                                         for(itCombo=combos.begin();itCombo!=combos.end(); itCombo++){
302                                                 if(itCombo->second == i){       thisGroup = itCombo->first;     combos.erase(itCombo);  break;  }
303                                         }
304                                 }else{ thisGroup = groupVector[i]; }
305                                 
306                                 while(!inFASTA.eof()){
307                                         if(inFASTA.get() == '>'){
308                                                 inFASTA >> seqName;
309                                                 outGroups << seqName << '\t' << thisGroup << endl;
310                                         }
311                                         while (!inFASTA.eof())  {       char c = inFASTA.get(); if (c == 10 || c == 13){        break;  }       }
312                                 }
313                                 outGroups.close();
314                                 inFASTA.close();
315                         }
316                 }
317                 
318                 if (m->control_pressed) { 
319                         for (int i = 0; i < outputNames.size(); i++) {  remove(outputNames[i].c_str()); }
320                         return 0;
321                 }
322
323                 m->mothurOutEndLine();
324                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
325                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }
326                 m->mothurOutEndLine();
327                 
328                 return 0;       
329                         
330         }
331         catch(exception& e) {
332                 m->errorOut(e, "TrimSeqsCommand", "execute");
333                 exit(1);
334         }
335 }
336                 
337 /**************************************************************************************/
338 int TrimSeqsCommand::driverCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string groupFile, vector<string> fastaNames, linePair* line, linePair* qline) {     
339         try {
340                 
341                 ofstream outFASTA;
342                 int able = openOutputFile(trimFile, outFASTA);
343                 
344                 ofstream scrapFASTA;
345                 openOutputFile(scrapFile, scrapFASTA);
346                 
347                 ofstream outGroups;
348                 vector<ofstream*> fastaFileNames;
349                 if (oligoFile != "") {          
350                         openOutputFile(groupFile, outGroups);   
351                         for (int i = 0; i < fastaNames.size(); i++) {
352                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
353                                 fastaFileNames.push_back(new ofstream((fastaNames[i] + toString(getpid()) + ".temp").c_str(), ios::ate)); 
354                         #else
355                                 fastaFileNames.push_back(new ofstream((fastaNames[i] + toString(i) + ".temp").c_str(), ios::ate));                      
356                         #endif
357                         }
358                 }
359                 
360                 ifstream inFASTA;
361                 openInputFile(filename, inFASTA);
362                 
363                 ifstream qFile;
364                 if(qFileName != "")     {       openInputFile(qFileName, qFile);        }
365                 
366                 qFile.seekg(qline->start);
367                 inFASTA.seekg(line->start);
368                 
369                 for(int i=0;i<line->num;i++){
370                                 
371                         if (m->control_pressed) { 
372                                 inFASTA.close(); outFASTA.close(); scrapFASTA.close(); if (oligoFile != "") {    outGroups.close();   } if(qFileName != "")     {       qFile.close();  }
373                                 for(int i=0;i<fastaFileNames.size();i++){  fastaFileNames[i]->close(); delete fastaFileNames[i];  }     
374                                 for (int i = 0; i < outputNames.size(); i++) {  remove(outputNames[i].c_str()); }
375                                 return 0;
376                         }
377                         
378                         int success = 1;
379                         
380                         Sequence currSeq(inFASTA);
381
382                         string origSeq = currSeq.getUnaligned();
383                         if (origSeq != "") {
384                                 int groupBar, groupPrime;
385                                 string trashCode = "";
386                                 int currentSeqsDiffs = 0;
387                                 
388                                 if(qFileName != ""){
389                                         if(qThreshold != 0)             {       success = stripQualThreshold(currSeq, qFile);   }
390                                         else if(qAverage != 0)  {       success = cullQualAverage(currSeq, qFile);              }
391                                         
392                                         if (qtrim == 1 && (origSeq.length() != currSeq.getUnaligned().length())) { 
393                                                 success = 0; //if you don't want to trim and the sequence does not meet quality requirements, move to scrap
394                                         }
395
396                                         if(!success)                    {       trashCode += 'q';                                                               }
397                                 }
398                         
399                                 if(barcodes.size() != 0){
400                                         success = stripBarcode(currSeq, groupBar);
401                                         if(success > bdiffs){   trashCode += 'b';       }
402                                         else{ currentSeqsDiffs += success;  }
403                                 }
404
405                                 if(numFPrimers != 0){
406                                         success = stripForward(currSeq, groupPrime);
407                                         if(success > pdiffs){   trashCode += 'f';       }
408                                         else{ currentSeqsDiffs += success;  }
409                                 }
410                                 
411                                 if (currentSeqsDiffs > tdiffs) { trashCode += 't';   }
412
413                                 if(numRPrimers != 0){
414                                         success = stripReverse(currSeq);
415                                         if(!success){   trashCode += 'r';       }
416                                 }
417                 
418                                 if(minLength > 0 || maxLength > 0){
419                                         success = cullLength(currSeq);
420                                         if(!success){   trashCode += 'l'; }
421                                 }
422                                 if(maxHomoP > 0){
423                                         success = cullHomoP(currSeq);
424                                         if(!success){   trashCode += 'h';       }
425                                 }
426                                 if(maxAmbig != -1){
427                                         success = cullAmbigs(currSeq);
428                                         if(!success){   trashCode += 'n';       }
429                                 }
430                                 
431                                 if(flip){       currSeq.reverseComplement();    }               // should go last                       
432                                 
433                                 if(trashCode.length() == 0){
434                                         currSeq.setAligned(currSeq.getUnaligned());
435                                         currSeq.printSequence(outFASTA);
436                                         if(barcodes.size() != 0){
437                                                 string thisGroup = groupVector[groupBar];
438                                                 int indexToFastaFile = groupBar;
439                                                 if (primers.size() != 0){
440                                                         //does this primer have a group
441                                                         if (groupVector[groupPrime] != "") {  
442                                                                 thisGroup += "." + groupVector[groupPrime]; 
443                                                                 indexToFastaFile = combos[thisGroup];
444                                                         }
445                                                 }
446                                                 outGroups << currSeq.getName() << '\t' << thisGroup << endl;
447                                                 if(allFiles){
448                                                         currSeq.printSequence(*fastaFileNames[indexToFastaFile]);                                       
449                                                 }
450                                         }
451                                 }
452                                 else{
453                                         currSeq.setName(currSeq.getName() + '|' + trashCode);
454                                         currSeq.setUnaligned(origSeq);
455                                         currSeq.setAligned(origSeq);
456                                         currSeq.printSequence(scrapFASTA);
457                                 }
458                         }
459                         gobble(inFASTA);
460                 }
461                 
462                 inFASTA.close();
463                 outFASTA.close();
464                 scrapFASTA.close();
465                 if (oligoFile != "") {   outGroups.close();   }
466                 if(qFileName != "")     {       qFile.close();  }
467                 
468                 for(int i=0;i<fastaFileNames.size();i++){
469                         fastaFileNames[i]->close();
470                         delete fastaFileNames[i];
471                 }               
472                 
473                 return 0;
474         }
475         catch(exception& e) {
476                 m->errorOut(e, "TrimSeqsCommand", "driverCreateTrim");
477                 exit(1);
478         }
479 }
480 /**************************************************************************************************/
481 int TrimSeqsCommand::createProcessesCreateTrim(string filename, string qFileName, string trimFile, string scrapFile, string groupFile, vector<string> fastaNames) {
482         try {
483 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
484                 int process = 0;
485                 int exitCommand = 1;
486                 processIDS.clear();
487                 
488                 //loop through and create all the processes you want
489                 while (process != processors) {
490                         int pid = fork();
491                         
492                         if (pid > 0) {
493                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
494                                 process++;
495                         }else if (pid == 0){
496                                 driverCreateTrim(filename, qFileName, (trimFile + toString(getpid()) + ".temp"), (scrapFile + toString(getpid()) + ".temp"), (groupFile + toString(getpid()) + ".temp"), fastaNames, lines[process], qLines[process]);
497                                 exit(0);
498                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
499                 }
500                 
501                 //force parent to wait until all the processes are done
502                 for (int i=0;i<processors;i++) { 
503                         int temp = processIDS[i];
504                         wait(&temp);
505                 }
506                 
507                 return exitCommand;
508 #endif          
509         }
510         catch(exception& e) {
511                 m->errorOut(e, "TrimSeqsCommand", "createProcessesCreateTrim");
512                 exit(1);
513         }
514 }
515 /**************************************************************************************************/
516
517 int TrimSeqsCommand::setLines(string filename, vector<linePair*>& lines) {
518         try {
519                 
520                 lines.clear();
521                 
522                 vector<long int> positions;
523                 
524                 ifstream inFASTA;
525                 openInputFile(filename, inFASTA);
526                         
527                 string input;
528                 while(!inFASTA.eof()){  
529                         input = getline(inFASTA);
530
531                         if (input.length() != 0) {
532                                 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);     }
533                         }
534                 }
535                 inFASTA.close();
536                 
537                 int numFastaSeqs = positions.size();
538         
539                 FILE * pFile;
540                 long size;
541                 
542                 //get num bytes in file
543                 pFile = fopen (filename.c_str(),"rb");
544                 if (pFile==NULL) perror ("Error opening file");
545                 else{
546                         fseek (pFile, 0, SEEK_END);
547                         size=ftell (pFile);
548                         fclose (pFile);
549                 }
550                 
551                 int numSeqsPerProcessor = numFastaSeqs / processors;
552                 
553                 for (int i = 0; i < processors; i++) {
554
555                         long int startPos = positions[ i * numSeqsPerProcessor ];
556                         if(i == processors - 1){
557                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
558                         }else{  
559                                 long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
560                         }
561                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
562                 }
563                 
564                 return numFastaSeqs;
565         }
566         catch(exception& e) {
567                 m->errorOut(e, "TrimSeqsCommand", "setLines");
568                 exit(1);
569         }
570 }
571 //***************************************************************************************************************
572
573 void TrimSeqsCommand::getOligos(vector<string>& outFASTAVec){ //vector<ofstream*>& outFASTAVec
574         try {
575                 ifstream inOligos;
576                 openInputFile(oligoFile, inOligos);
577                 
578                 ofstream test;
579                 
580                 string type, oligo, group;
581                 int index=0;
582                 //int indexPrimer = 0;
583                 
584                 while(!inOligos.eof()){
585                         inOligos >> type;
586                                         
587                         if(type[0] == '#'){
588                                 while (!inOligos.eof()) {       char c = inOligos.get(); if (c == 10 || c == 13){       break;  }       } // get rest of line if there's any crap there
589                         }
590                         else{
591                                 //make type case insensitive
592                                 for(int i=0;i<type.length();i++){       type[i] = toupper(type[i]);  }
593                                 
594                                 inOligos >> oligo;
595                                 
596                                 for(int i=0;i<oligo.length();i++){
597                                         oligo[i] = toupper(oligo[i]);
598                                         if(oligo[i] == 'U')     {       oligo[i] = 'T'; }
599                                 }
600                                 
601                                 if(type == "FORWARD"){
602                                         group = "";
603                                         
604                                         // get rest of line in case there is a primer name
605                                         while (!inOligos.eof()) {       
606                                                 char c = inOligos.get(); 
607                                                 if (c == 10 || c == 13){        break;  }
608                                                 else if (c == 32 || c == 9){;} //space or tab
609                                                 else {  group += c;  }
610                                         } 
611                                         
612                                         //check for repeat barcodes
613                                         map<string, int>::iterator itPrime = primers.find(oligo);
614                                         if (itPrime != primers.end()) { m->mothurOut("primer " + oligo + " is in your oligos file already."); m->mothurOutEndLine();  }
615                                         
616                                         primers[oligo]=index; index++;
617                                         groupVector.push_back(group);
618                                         
619                                         if(allFiles){
620                                                 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
621                                                 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
622                                                         filesToRemove.insert((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
623                                                 }else {
624                                                         outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
625                                                 }
626                                         }
627
628                                 }
629                                 else if(type == "REVERSE"){
630                                         Sequence oligoRC("reverse", oligo);
631                                         oligoRC.reverseComplement();
632                                         revPrimer.push_back(oligoRC.getUnaligned());
633                                 }
634                                 else if(type == "BARCODE"){
635                                         inOligos >> group;
636                                         
637                                         //check for repeat barcodes
638                                         map<string, int>::iterator itBar = barcodes.find(oligo);
639                                         if (itBar != barcodes.end()) { m->mothurOut("barcode " + oligo + " is in your oligos file already."); m->mothurOutEndLine();  }
640                                         
641                                         barcodes[oligo]=index; index++;
642                                         groupVector.push_back(group);
643                                         
644                                         if(allFiles){
645                                                 outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
646                                                 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(index) + "." + group + ".fasta"));
647                                         }
648                                 }else{  m->mothurOut(type + " is not recognized as a valid type. Choices are forward, reverse, and barcode. Ignoring " + oligo + "."); m->mothurOutEndLine();  }
649                         }
650                         gobble(inOligos);
651                 }
652                 
653                 inOligos.close();
654                 
655                 //add in potential combos
656                 if(allFiles){
657                         comboStarts = outFASTAVec.size()-1;
658                         for (map<string, int>::iterator itBar = barcodes.begin(); itBar != barcodes.end(); itBar++) {
659                                 for (map<string, int>::iterator itPrime = primers.begin(); itPrime != primers.end(); itPrime++) {
660                                         if (groupVector[itPrime->second] != "") { //there is a group for this primer
661                                                 outputNames.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(itBar->second) + "." + groupVector[itBar->second] + "." + toString(itPrime->second) + "." + groupVector[itPrime->second] + ".fasta"));
662                                                 outFASTAVec.push_back((outputDir + getRootName(getSimpleName(fastaFile)) + toString(itBar->second) + "." + groupVector[itBar->second] + "." + toString(itPrime->second) + "." + groupVector[itPrime->second] + ".fasta"));
663                                                 combos[(groupVector[itBar->second] + "." + groupVector[itPrime->second])] = outFASTAVec.size()-1;
664                                         }
665                                 }
666                         }
667                 }
668                 
669                 numFPrimers = primers.size();
670                 numRPrimers = revPrimer.size();
671                 
672         }
673         catch(exception& e) {
674                 m->errorOut(e, "TrimSeqsCommand", "getOligos");
675                 exit(1);
676         }
677 }
678 //***************************************************************************************************************
679
680 int TrimSeqsCommand::stripBarcode(Sequence& seq, int& group){
681         try {
682                 
683                 string rawSequence = seq.getUnaligned();
684                 int success = bdiffs + 1;       //guilty until proven innocent
685                 
686                 //can you find the barcode
687                 for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
688                         string oligo = it->first;
689                         if(rawSequence.length() < oligo.length()){      //let's just assume that the barcodes are the same length
690                                 success = bdiffs + 10;                                  //if the sequence is shorter than the barcode then bail out
691                                 break;  
692                         }
693                         
694                         if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
695                                 group = it->second;
696                                 seq.setUnaligned(rawSequence.substr(oligo.length()));
697                                 success = 0;
698                                 break;
699                         }
700                 }
701                 
702                 //if you found the barcode or if you don't want to allow for diffs
703 //              cout << success;
704                 if ((bdiffs == 0) || (success == 0)) { return success;  }
705                 
706                 else { //try aligning and see if you can find it
707 //                      cout << endl;
708
709                         int maxLength = 0;
710
711                         Alignment* alignment;
712                         if (barcodes.size() > 0) {
713                                 map<string,int>::iterator it=barcodes.begin();
714
715                                 for(it;it!=barcodes.end();it++){
716                                         if(it->first.length() > maxLength){
717                                                 maxLength = it->first.length();
718                                         }
719                                 }
720                                 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+bdiffs+1));  
721
722                         }else{ alignment = NULL; } 
723                         
724                         //can you find the barcode
725                         int minDiff = 1e6;
726                         int minCount = 1;
727                         int minGroup = -1;
728                         int minPos = 0;
729                         
730                         for(map<string,int>::iterator it=barcodes.begin();it!=barcodes.end();it++){
731                                 string oligo = it->first;
732 //                              int length = oligo.length();
733                                 
734                                 if(rawSequence.length() < maxLength){   //let's just assume that the barcodes are the same length
735                                         success = bdiffs + 10;
736                                         break;
737                                 }
738                                 
739                                 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
740                                 alignment->align(oligo, rawSequence.substr(0,oligo.length()+bdiffs));
741                                 oligo = alignment->getSeqAAln();
742                                 string temp = alignment->getSeqBAln();
743                 
744                                 int alnLength = oligo.length();
745                                 
746                                 for(int i=oligo.length()-1;i>=0;i--){
747                                         if(oligo[i] != '-'){    alnLength = i+1;        break;  }
748                                 }
749                                 oligo = oligo.substr(0,alnLength);
750                                 temp = temp.substr(0,alnLength);
751                                 
752                                 int newStart=0;
753                                 int numDiff = countDiffs(oligo, temp);
754                                 
755 //                              cout << oligo << '\t' << temp << '\t' << numDiff << endl;                               
756                                 
757                                 if(numDiff < minDiff){
758                                         minDiff = numDiff;
759                                         minCount = 1;
760                                         minGroup = it->second;
761                                         minPos = 0;
762                                         for(int i=0;i<alnLength;i++){
763                                                 if(temp[i] != '-'){
764                                                         minPos++;
765                                                 }
766                                         }
767                                 }
768                                 else if(numDiff == minDiff){
769                                         minCount++;
770                                 }
771
772                         }
773
774                         if(minDiff > bdiffs)    {       success = minDiff;              }       //no good matches
775                         else if(minCount > 1)   {       success = bdiffs + 100; }       //can't tell the difference between multiple barcodes
776                         else{                                                                                                   //use the best match
777                                 group = minGroup;
778                                 seq.setUnaligned(rawSequence.substr(minPos));
779                                 success = minDiff;
780                         }
781                         
782                         if (alignment != NULL) {  delete alignment;  }
783                         
784                 }
785 //              cout << success << endl;
786                 
787                 return success;
788                 
789         }
790         catch(exception& e) {
791                 m->errorOut(e, "TrimSeqsCommand", "stripBarcode");
792                 exit(1);
793         }
794
795 }
796
797 //***************************************************************************************************************
798
799 int TrimSeqsCommand::stripForward(Sequence& seq, int& group){
800         try {
801                 string rawSequence = seq.getUnaligned();
802                 int success = pdiffs + 1;       //guilty until proven innocent
803                 
804                 //can you find the primer
805                 for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
806                         string oligo = it->first;
807                         if(rawSequence.length() < oligo.length()){      //let's just assume that the primers are the same length
808                                 success = pdiffs + 10;                                  //if the sequence is shorter than the barcode then bail out
809                                 break;  
810                         }
811                         
812                         if(compareDNASeq(oligo, rawSequence.substr(0,oligo.length()))){
813                                 group = it->second;
814                                 seq.setUnaligned(rawSequence.substr(oligo.length()));
815                                 success = 0;
816                                 break;
817                         }
818                 }
819
820                 //if you found the barcode or if you don't want to allow for diffs
821 //              cout << success;
822                 if ((pdiffs == 0) || (success == 0)) { return success;  }
823                 
824                 else { //try aligning and see if you can find it
825 //                      cout << endl;
826
827                         int maxLength = 0;
828
829                         Alignment* alignment;
830                         if (primers.size() > 0) {
831                                 map<string,int>::iterator it=primers.begin();
832
833                                 for(it;it!=primers.end();it++){
834                                         if(it->first.length() > maxLength){
835                                                 maxLength = it->first.length();
836                                         }
837                                 }
838                                 alignment = new NeedlemanOverlap(-1.0, 1.0, -1.0, (maxLength+pdiffs+1));  
839
840                         }else{ alignment = NULL; } 
841                         
842                         //can you find the barcode
843                         int minDiff = 1e6;
844                         int minCount = 1;
845                         int minGroup = -1;
846                         int minPos = 0;
847                         
848                         for(map<string,int>::iterator it=primers.begin();it!=primers.end();it++){
849                                 string oligo = it->first;
850 //                              int length = oligo.length();
851                                 
852                                 if(rawSequence.length() < maxLength){   
853                                         success = pdiffs + 100;
854                                         break;
855                                 }
856                                 
857                                 //use needleman to align first barcode.length()+numdiffs of sequence to each barcode
858                                 alignment->align(oligo, rawSequence.substr(0,oligo.length()+pdiffs));
859                                 oligo = alignment->getSeqAAln();
860                                 string temp = alignment->getSeqBAln();
861                 
862                                 int alnLength = oligo.length();
863                                 
864                                 for(int i=oligo.length()-1;i>=0;i--){
865                                         if(oligo[i] != '-'){    alnLength = i+1;        break;  }
866                                 }
867                                 oligo = oligo.substr(0,alnLength);
868                                 temp = temp.substr(0,alnLength);
869                                 
870                                 int newStart=0;
871                                 int numDiff = countDiffs(oligo, temp);
872                                 
873 //                              cout << oligo << '\t' << temp << '\t' << numDiff << endl;                               
874                                 
875                                 if(numDiff < minDiff){
876                                         minDiff = numDiff;
877                                         minCount = 1;
878                                         minGroup = it->second;
879                                         minPos = 0;
880                                         for(int i=0;i<alnLength;i++){
881                                                 if(temp[i] != '-'){
882                                                         minPos++;
883                                                 }
884                                         }
885                                 }
886                                 else if(numDiff == minDiff){
887                                         minCount++;
888                                 }
889
890                         }
891
892                         if(minDiff > pdiffs)    {       success = minDiff;              }       //no good matches
893                         else if(minCount > 1)   {       success = pdiffs + 10;  }       //can't tell the difference between multiple primers
894                         else{                                                                                                   //use the best match
895                                 group = minGroup;
896                                 seq.setUnaligned(rawSequence.substr(minPos));
897                                 success = minDiff;
898                         }
899                         
900                         if (alignment != NULL) {  delete alignment;  }
901                         
902                 }
903                 
904                 return success;
905
906         }
907         catch(exception& e) {
908                 m->errorOut(e, "TrimSeqsCommand", "stripForward");
909                 exit(1);
910         }
911 }
912
913 //***************************************************************************************************************
914
915 bool TrimSeqsCommand::stripReverse(Sequence& seq){
916         try {
917                 string rawSequence = seq.getUnaligned();
918                 bool success = 0;       //guilty until proven innocent
919                 
920                 for(int i=0;i<numRPrimers;i++){
921                         string oligo = revPrimer[i];
922                         
923                         if(rawSequence.length() < oligo.length()){
924                                 success = 0;
925                                 break;
926                         }
927                         
928                         if(compareDNASeq(oligo, rawSequence.substr(rawSequence.length()-oligo.length(),oligo.length()))){
929                                 seq.setUnaligned(rawSequence.substr(0,rawSequence.length()-oligo.length()));
930                                 success = 1;
931                                 break;
932                         }
933                 }       
934                 return success;
935                 
936         }
937         catch(exception& e) {
938                 m->errorOut(e, "TrimSeqsCommand", "stripReverse");
939                 exit(1);
940         }
941 }
942
943 //***************************************************************************************************************
944
945 bool TrimSeqsCommand::cullLength(Sequence& seq){
946         try {
947         
948                 int length = seq.getNumBases();
949                 bool success = 0;       //guilty until proven innocent
950                 
951                 if(length >= minLength && maxLength == 0)                       {       success = 1;    }
952                 else if(length >= minLength && length <= maxLength)     {       success = 1;    }
953                 else                                                                                            {       success = 0;    }
954                 
955                 return success;
956         
957         }
958         catch(exception& e) {
959                 m->errorOut(e, "TrimSeqsCommand", "cullLength");
960                 exit(1);
961         }
962         
963 }
964
965 //***************************************************************************************************************
966
967 bool TrimSeqsCommand::cullHomoP(Sequence& seq){
968         try {
969                 int longHomoP = seq.getLongHomoPolymer();
970                 bool success = 0;       //guilty until proven innocent
971                 
972                 if(longHomoP <= maxHomoP){      success = 1;    }
973                 else                                    {       success = 0;    }
974                 
975                 return success;
976         }
977         catch(exception& e) {
978                 m->errorOut(e, "TrimSeqsCommand", "cullHomoP");
979                 exit(1);
980         }
981         
982 }
983
984 //***************************************************************************************************************
985
986 bool TrimSeqsCommand::cullAmbigs(Sequence& seq){
987         try {
988                 int numNs = seq.getAmbigBases();
989                 bool success = 0;       //guilty until proven innocent
990                 
991                 if(numNs <= maxAmbig)   {       success = 1;    }
992                 else                                    {       success = 0;    }
993                 
994                 return success;
995         }
996         catch(exception& e) {
997                 m->errorOut(e, "TrimSeqsCommand", "cullAmbigs");
998                 exit(1);
999         }
1000         
1001 }
1002
1003 //***************************************************************************************************************
1004
1005 bool TrimSeqsCommand::compareDNASeq(string oligo, string seq){
1006         try {
1007                 bool success = 1;
1008                 int length = oligo.length();
1009                 
1010                 for(int i=0;i<length;i++){
1011                         
1012                         if(oligo[i] != seq[i]){
1013                                 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C')    {       success = 0;    }
1014                                 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N'))                                {       success = 0;    }
1015                                 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G'))                                    {       success = 0;    }
1016                                 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T'))                                    {       success = 0;    }
1017                                 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A'))                                    {       success = 0;    }
1018                                 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G'))                                    {       success = 0;    }
1019                                 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A'))                                    {       success = 0;    }
1020                                 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G'))                                    {       success = 0;    }
1021                                 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G'))   {       success = 0;    }
1022                                 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G'))   {       success = 0;    }
1023                                 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C'))   {       success = 0;    }
1024                                 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G'))   {       success = 0;    }                       
1025                                 
1026                                 if(success == 0)        {       break;   }
1027                         }
1028                         else{
1029                                 success = 1;
1030                         }
1031                 }
1032                 
1033                 return success;
1034         }
1035         catch(exception& e) {
1036                 m->errorOut(e, "TrimSeqsCommand", "compareDNASeq");
1037                 exit(1);
1038         }
1039
1040 }
1041 //***************************************************************************************************************
1042
1043 int TrimSeqsCommand::countDiffs(string oligo, string seq){
1044         try {
1045
1046                 int length = oligo.length();
1047                 int countDiffs = 0;
1048                 
1049                 for(int i=0;i<length;i++){
1050                                                                 
1051                         if(oligo[i] != seq[i]){
1052                                 if(oligo[i] == 'A' || oligo[i] == 'T' || oligo[i] == 'G' || oligo[i] == 'C' || oligo[i] == '-' || oligo[i] == '.')      {       countDiffs++;   }
1053                                 else if((oligo[i] == 'N' || oligo[i] == 'I') && (seq[i] == 'N'))                                {       countDiffs++;   }
1054                                 else if(oligo[i] == 'R' && (seq[i] != 'A' && seq[i] != 'G'))                                    {       countDiffs++;   }
1055                                 else if(oligo[i] == 'Y' && (seq[i] != 'C' && seq[i] != 'T'))                                    {       countDiffs++;   }
1056                                 else if(oligo[i] == 'M' && (seq[i] != 'C' && seq[i] != 'A'))                                    {       countDiffs++;   }
1057                                 else if(oligo[i] == 'K' && (seq[i] != 'T' && seq[i] != 'G'))                                    {       countDiffs++;   }
1058                                 else if(oligo[i] == 'W' && (seq[i] != 'T' && seq[i] != 'A'))                                    {       countDiffs++;   }
1059                                 else if(oligo[i] == 'S' && (seq[i] != 'C' && seq[i] != 'G'))                                    {       countDiffs++;   }
1060                                 else if(oligo[i] == 'B' && (seq[i] != 'C' && seq[i] != 'T' && seq[i] != 'G'))   {       countDiffs++;   }
1061                                 else if(oligo[i] == 'D' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'G'))   {       countDiffs++;   }
1062                                 else if(oligo[i] == 'H' && (seq[i] != 'A' && seq[i] != 'T' && seq[i] != 'C'))   {       countDiffs++;   }
1063                                 else if(oligo[i] == 'V' && (seq[i] != 'A' && seq[i] != 'C' && seq[i] != 'G'))   {       countDiffs++;   }                       
1064                                 
1065                         }
1066                 }
1067                 
1068                 return countDiffs;
1069         }
1070         catch(exception& e) {
1071                 m->errorOut(e, "TrimSeqsCommand", "countDiffs");
1072                 exit(1);
1073         }
1074
1075 }
1076 //***************************************************************************************************************
1077
1078 bool TrimSeqsCommand::stripQualThreshold(Sequence& seq, ifstream& qFile){
1079         try {
1080 //              string rawSequence = seq.getUnaligned();
1081 //              int seqLength;  // = rawSequence.length();
1082 //              string name, temp, temp2;
1083 //              
1084 //              qFile >> name;
1085 //              
1086 //              //get rest of line
1087 //              temp = "";
1088 //              while (!qFile.eof())    {       
1089 //                      char c = qFile.get(); 
1090 //                      if (c == 10 || c == 13){        break;  }       
1091 //                      else { temp += c; }
1092 //              } 
1093 //      
1094 //              int pos = temp.find("length");
1095 //              if (pos == temp.npos) { m->mothurOut("Cannot find length in qfile for " + seq.getName()); m->mothurOutEndLine();  seqLength = 0;  }
1096 //              else {
1097 //                      string tempLength = temp.substr(pos);
1098 //                      istringstream iss (tempLength,istringstream::in);
1099 //                      iss >> temp;
1100 //              }
1101 //              
1102 //              splitAtEquals(temp2, temp); //separates length=242, temp=length, temp2=242
1103 //              convert(temp, seqLength); //converts string to int
1104 //      
1105 //              if (name.length() != 0) {  if(name.substr(1) != seq.getName())  {       m->mothurOut("sequence name mismatch btwn fasta and qual file"); m->mothurOutEndLine(); }  } 
1106                 
1107                 string rawSequence = seq.getUnaligned();
1108                 int seqLength = seq.getNumBases();
1109                 bool success = 0;       //guilty until proven innocent
1110                 string name;
1111                 
1112                 qFile >> name;
1113                 if (name[0] == '>') {  if(name.substr(1) != seq.getName())      {       m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); }  }
1114                 
1115                 while (!qFile.eof())    {       char c = qFile.get(); if (c == 10 || c == 13){  break;  }       }
1116                 
1117                 int score;
1118                 int end = seqLength;
1119                 
1120                 for(int i=0;i<seqLength;i++){
1121                         qFile >> score;
1122                         
1123                         if(score < qThreshold){
1124                                 end = i;
1125                                 break;
1126                         }
1127                 }
1128                 for(int i=end+1;i<seqLength;i++){
1129                         qFile >> score;
1130                 }
1131                 
1132                 seq.setUnaligned(rawSequence.substr(0,end));
1133                 
1134                 return 1;
1135         }
1136         catch(exception& e) {
1137                 m->errorOut(e, "TrimSeqsCommand", "stripQualThreshold");
1138                 exit(1);
1139         }
1140 }
1141
1142 //***************************************************************************************************************
1143
1144 bool TrimSeqsCommand::cullQualAverage(Sequence& seq, ifstream& qFile){
1145         try {
1146                 string rawSequence = seq.getUnaligned();
1147                 int seqLength = seq.getNumBases();
1148                 bool success = 0;       //guilty until proven innocent
1149                 string name;
1150                 
1151                 qFile >> name;
1152                 if (name[0] == '>') {  if(name.substr(1) != seq.getName())      {       m->mothurOut("sequence name mismatch btwn fasta: " + seq.getName() + " and qual file: " + name); m->mothurOutEndLine(); }  }
1153                 
1154                 while (!qFile.eof())    {       char c = qFile.get(); if (c == 10 || c == 13){  break;  }       }
1155                 
1156                 float score;    
1157                 float average = 0;
1158                 
1159                 for(int i=0;i<seqLength;i++){
1160                         qFile >> score;
1161                         average += score;
1162                 }
1163                 average /= seqLength;
1164
1165                 if(average >= qAverage) {       success = 1;    }
1166                 else                                    {       success = 0;    }
1167                 
1168                 return success;
1169         }
1170         catch(exception& e) {
1171                 m->errorOut(e, "TrimSeqsCommand", "cullQualAverage");
1172                 exit(1);
1173         }
1174 }
1175
1176 //***************************************************************************************************************