5 * Created by westcott on 11/4/09.
6 * Copyright 2009 Schloss Lab. All rights reserved.
10 #include "alignmentdb.h"
12 #include "suffixdb.hpp"
13 #include "blastdb.hpp"
14 #include "referencedb.h"
17 /**************************************************************************************************/
18 AlignmentDB::AlignmentDB(string fastaFileName, string s, int kmerSize, float gapOpen, float gapExtend, float match, float misMatch){ // This assumes that the template database is in fasta format, may
19 try { // need to alter this in the future?
20 m = MothurOut::getInstance();
23 bool needToGenerate = true;
24 ReferenceDB* rdb = ReferenceDB::getInstance();
26 if (fastaFileName == "saved") {
27 int start = time(NULL);
28 m->mothurOutEndLine(); m->mothurOut("Using sequences from " + rdb->getSavedReference() + " that are saved in memory."); m->mothurOutEndLine();
30 for (int i = 0; i < rdb->referenceSeqs.size(); i++) {
31 templateSequences.push_back(rdb->referenceSeqs[i]);
33 if (rdb->referenceSeqs[i].getUnaligned().length() >= longest) { longest = (rdb->referenceSeqs[i].getUnaligned().length()+1); }
35 fastaFileName = rdb->getSavedReference();
37 numSeqs = templateSequences.size();
38 m->mothurOut("It took " + toString(time(NULL) - start) + " to load " + toString(rdb->referenceSeqs.size()) + " sequences.");m->mothurOutEndLine();
41 int start = time(NULL);
42 m->mothurOutEndLine();
43 m->mothurOut("Reading in the " + fastaFileName + " template sequences...\t"); cout.flush();
47 vector<unsigned long int> positions;
51 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
52 MPI_Comm_size(MPI_COMM_WORLD, &processors);
55 char inFileName[1024];
56 strcpy(inFileName, fastaFileName.c_str());
58 MPI_File_open(MPI_COMM_WORLD, inFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
61 positions = m->setFilePosFasta(fastaFileName, numSeqs); //fills MPIPos, returns numSeqs
63 //send file positions to all processes
64 for(int i = 1; i < processors; i++) {
65 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
66 MPI_Send(&positions[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
69 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
70 positions.resize(numSeqs+1);
71 MPI_Recv(&positions[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
75 for(int i=0;i<numSeqs;i++){
77 if (m->control_pressed) { templateSequences.clear(); break; }
80 int length = positions[i+1] - positions[i];
81 char* buf4 = new char[length];
83 MPI_File_read_at(inMPI, positions[i], buf4, length, MPI_CHAR, &status);
85 string tempBuf = buf4;
86 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
89 istringstream iss (tempBuf,istringstream::in);
92 if (temp.getName() != "") {
93 templateSequences.push_back(temp);
95 if (rdb->save) { rdb->referenceSeqs.push_back(temp); }
98 if (temp.getUnaligned().length() >= longest) { longest = temp.getUnaligned().length()+1; }
102 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
104 MPI_File_close(&inMPI);
108 m->openInputFile(fastaFileName, fastaFile);
110 while (!fastaFile.eof()) {
111 Sequence temp(fastaFile); m->gobble(fastaFile);
113 if (m->control_pressed) { templateSequences.clear(); break; }
115 if (temp.getName() != "") {
116 templateSequences.push_back(temp);
118 if (rdb->save) { rdb->referenceSeqs.push_back(temp); }
121 if (temp.getUnaligned().length() >= longest) { longest = (temp.getUnaligned().length()+1); }
127 numSeqs = templateSequences.size();
128 //all of this is elsewhere already!
130 m->mothurOut("DONE.");
131 m->mothurOutEndLine(); cout.flush();
132 m->mothurOut("It took " + toString(time(NULL) - start) + " to read " + toString(rdb->referenceSeqs.size()) + " sequences."); m->mothurOutEndLine();
137 //in case you delete the seqs and then ask for them
138 emptySequence = Sequence();
139 emptySequence.setName("no_match");
140 emptySequence.setUnaligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
141 emptySequence.setAligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
145 if(method == "kmer") {
146 search = new KmerDB(fastaFileName, kmerSize);
150 kmerDBName = fastaFileName.substr(0,fastaFileName.find_last_of(".")+1) + char('0'+ kmerSize) + "mer";
151 ifstream kmerFileTest(kmerDBName.c_str());
154 bool GoodFile = m->checkReleaseVersion(kmerFileTest, m->getVersion());
155 if (GoodFile) { needToGenerate = false; }
159 else if(method == "suffix") { search = new SuffixDB(numSeqs); }
160 else if(method == "blast") { search = new BlastDB(fastaFileName.substr(0,fastaFileName.find_last_of(".")+1), gapOpen, gapExtend, match, misMatch); }
162 m->mothurOut(method + " is not a valid search option. I will run the command using kmer, ksize=8.");
163 m->mothurOutEndLine();
164 search = new KmerDB(fastaFileName, 8);
167 if (!(m->control_pressed)) {
168 if (needToGenerate) {
169 //add sequences to search
170 for (int i = 0; i < templateSequences.size(); i++) {
171 search->addSequence(templateSequences[i]);
173 if (m->control_pressed) { templateSequences.clear(); break; }
176 if (m->control_pressed) { templateSequences.clear(); }
178 search->generateDB();
180 }else if ((method == "kmer") && (!needToGenerate)) {
181 ifstream kmerFileTest(kmerDBName.c_str());
182 search->readKmerDB(kmerFileTest);
185 search->setNumSeqs(numSeqs);
189 catch(exception& e) {
190 m->errorOut(e, "AlignmentDB", "AlignmentDB");
194 /**************************************************************************************************/
195 AlignmentDB::AlignmentDB(string s){
197 m = MothurOut::getInstance();
200 if(method == "suffix") { search = new SuffixDB(); }
201 else if(method == "blast") { search = new BlastDB(); }
202 else { search = new KmerDB(); }
205 //in case you delete the seqs and then ask for them
206 emptySequence = Sequence();
207 emptySequence.setName("no_match");
208 emptySequence.setUnaligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
209 emptySequence.setAligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
212 catch(exception& e) {
213 m->errorOut(e, "AlignmentDB", "AlignmentDB");
217 /**************************************************************************************************/
218 AlignmentDB::~AlignmentDB() { delete search; }
219 /**************************************************************************************************/
220 Sequence AlignmentDB::findClosestSequence(Sequence* seq) {
223 vector<int> spot = search->findClosestSequences(seq, 1);
225 if (spot.size() != 0) { return templateSequences[spot[0]]; }
226 else { return emptySequence; }
229 catch(exception& e) {
230 m->errorOut(e, "AlignmentDB", "findClosestSequence");
234 /**************************************************************************************************/