File: alignmentdb.cpp

package info (click to toggle)
mothur 1.33.3%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 11,248 kB
  • ctags: 12,231
  • sloc: cpp: 152,046; fortran: 665; makefile: 74; sh: 34
file content (256 lines) | stat: -rw-r--r-- 8,529 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
/*
 *  alignmentdb.cpp
 *  Mothur
 *
 *  Created by westcott on 11/4/09.
 *  Copyright 2009 Schloss Lab. All rights reserved.
 *
 */

#include "alignmentdb.h"
#include "kmerdb.hpp"
#include "suffixdb.hpp"
#include "blastdb.hpp"
#include "referencedb.h"

/**************************************************************************************************/
AlignmentDB::AlignmentDB(string fastaFileName, string s, int kmerSize, float gapOpen, float gapExtend, float match, float misMatch, int tid){		//	This assumes that the template database is in fasta format, may 
	try {											//	need to alter this in the future?
		m = MothurOut::getInstance();
		longest = 0;
		method = s;
		bool needToGenerate = true;
		ReferenceDB* rdb = ReferenceDB::getInstance();
		bool silent = false;
		threadID = tid;
		
		if (fastaFileName == "saved-silent") {
			fastaFileName = "saved"; silent = true;
		}
		
		if (fastaFileName == "saved") {
			int start = time(NULL);
			
			if (!silent) { m->mothurOutEndLine();  m->mothurOut("Using sequences from " + rdb->getSavedReference() + " that are saved in memory.");	m->mothurOutEndLine(); }

			for (int i = 0; i < rdb->referenceSeqs.size(); i++) {
				templateSequences.push_back(rdb->referenceSeqs[i]);
				//save longest base
				if (rdb->referenceSeqs[i].getUnaligned().length() >= longest)  { longest = (rdb->referenceSeqs[i].getUnaligned().length()+1); }
			}
			fastaFileName = rdb->getSavedReference();
			
			numSeqs = templateSequences.size();
			if (!silent) { m->mothurOut("It took " + toString(time(NULL) - start) + " to load " + toString(rdb->referenceSeqs.size()) + " sequences.");m->mothurOutEndLine();  }
            
		}else {
			int start = time(NULL);
			m->mothurOutEndLine();
			m->mothurOut("Reading in the " + fastaFileName + " template sequences...\t");	cout.flush();
			//bool aligned = false;
            int tempLength = 0;
            
			#ifdef USE_MPI	
				int pid, processors;
				vector<unsigned long long> positions;
			
				MPI_Status status; 
				MPI_File inMPI;
				MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
				MPI_Comm_size(MPI_COMM_WORLD, &processors);
				int tag = 2001;
		
				char inFileName[1024];
				strcpy(inFileName, fastaFileName.c_str());
		
				MPI_File_open(MPI_COMM_WORLD, inFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
				
				if (pid == 0) {
					positions = m->setFilePosFasta(fastaFileName, numSeqs); //fills MPIPos, returns numSeqs

					//send file positions to all processes
					for(int i = 1; i < processors; i++) { 
						MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
						MPI_Send(&positions[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
					}
				}else{
					MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
					positions.resize(numSeqs+1);
					MPI_Recv(&positions[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
				}
			
				//read file 
				for(int i=0;i<numSeqs;i++){
					
					if (m->control_pressed) {  templateSequences.clear(); break;  }
					
					//read next sequence
					int length = positions[i+1] - positions[i];
					char* buf4 = new char[length];
				
					MPI_File_read_at(inMPI, positions[i], buf4, length, MPI_CHAR, &status);
			
					string tempBuf = buf4;
					if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
					delete buf4;

					istringstream iss (tempBuf,istringstream::in);
			
					Sequence temp(iss);  
					if (temp.getName() != "") {
						templateSequences.push_back(temp);
						
						if (rdb->save) { rdb->referenceSeqs.push_back(temp); }
						
						//save longest base
						if (temp.getUnaligned().length() >= longest)  { longest = temp.getUnaligned().length()+1; }
                        if (tempLength != 0) {
                            if (tempLength != temp.getAligned().length()) { m->mothurOut("[ERROR]: template is not aligned, aborting.\n"); m->control_pressed=true; }
                        }else { tempLength = temp.getAligned().length(); }
					}
				}
				
				MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
				
				MPI_File_close(&inMPI);
			
		#else
			ifstream fastaFile;
			m->openInputFile(fastaFileName, fastaFile);

			while (!fastaFile.eof()) {
				Sequence temp(fastaFile);  m->gobble(fastaFile);
				
				if (m->control_pressed) {  templateSequences.clear(); break;  }
				
				if (temp.getName() != "") {
					templateSequences.push_back(temp);
					
					if (rdb->save) { rdb->referenceSeqs.push_back(temp); }
					
					//save longest base
					if (temp.getUnaligned().length() >= longest)  { longest = (temp.getUnaligned().length()+1); }
                    
                    if (tempLength != 0) {
                        if (tempLength != temp.getAligned().length()) { m->mothurOut("[ERROR]: template is not aligned, aborting.\n"); m->control_pressed=true; }
                    }else { tempLength = temp.getAligned().length(); }
				}
			}
			fastaFile.close();
		#endif
		
			numSeqs = templateSequences.size();
			//all of this is elsewhere already!
			
			m->mothurOut("DONE.");
			m->mothurOutEndLine();	cout.flush();
			m->mothurOut("It took " + toString(time(NULL) - start) + " to read  " + toString(templateSequences.size()) + " sequences."); m->mothurOutEndLine();  

		}
		
		
		//in case you delete the seqs and then ask for them
		emptySequence = Sequence();
		emptySequence.setName("no_match");
		emptySequence.setUnaligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
		emptySequence.setAligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
		
		
		string kmerDBName;
		if(method == "kmer")			{	
			search = new KmerDB(fastaFileName, kmerSize);			
			
			#ifdef USE_MPI
			#else
				kmerDBName = fastaFileName.substr(0,fastaFileName.find_last_of(".")+1) + char('0'+ kmerSize) + "mer";
				ifstream kmerFileTest(kmerDBName.c_str());
				
				if(kmerFileTest){	
					bool GoodFile = m->checkReleaseVersion(kmerFileTest, m->getVersion());
					if (GoodFile) {  needToGenerate = false;	}
				}
			#endif
		}
		else if(method == "suffix")		{	search = new SuffixDB(numSeqs);								}
		else if(method == "blast")		{	search = new BlastDB(fastaFileName.substr(0,fastaFileName.find_last_of(".")+1), gapOpen, gapExtend, match, misMatch, "", threadID);	}
		else {
			method = "kmer";
			m->mothurOut(method + " is not a valid search option. I will run the command using kmer, ksize=8.");
			m->mothurOutEndLine();
			search = new KmerDB(fastaFileName, 8);
		}
		
		if (!(m->control_pressed)) {
			if (needToGenerate) {
				//add sequences to search 
				for (int i = 0; i < templateSequences.size(); i++) {
					search->addSequence(templateSequences[i]);
					
					if (m->control_pressed) {  templateSequences.clear(); break;  }
				}
				
				if (m->control_pressed) {  templateSequences.clear();  }
				
				search->generateDB();
				
			}else if ((method == "kmer") && (!needToGenerate)) {
				ifstream kmerFileTest(kmerDBName.c_str());
				search->readKmerDB(kmerFileTest);	
			}
		
			search->setNumSeqs(numSeqs);
		}
		
	}
	catch(exception& e) {
		m->errorOut(e, "AlignmentDB", "AlignmentDB");
		exit(1);
	}
}
/**************************************************************************************************/
AlignmentDB::AlignmentDB(string s){		 
	try {											
		m = MothurOut::getInstance();
		method = s;
		
		if(method == "suffix")		{	search = new SuffixDB();	}
		else if(method == "blast")	{	search = new BlastDB("", 0);		}
		else						{	search = new KmerDB();		}

				
		//in case you delete the seqs and then ask for them
		emptySequence = Sequence();
		emptySequence.setName("no_match");
		emptySequence.setUnaligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
		emptySequence.setAligned("XXXXXXXXXXXXXXXXXXXXXXXXXXXXX");
		
	}
	catch(exception& e) {
		m->errorOut(e, "AlignmentDB", "AlignmentDB");
		exit(1);
	}
}
/**************************************************************************************************/
AlignmentDB::~AlignmentDB() {  delete search;	}
/**************************************************************************************************/
Sequence AlignmentDB::findClosestSequence(Sequence* seq) {
	try{
	
		vector<int> spot = search->findClosestSequences(seq, 1);
	
		if (spot.size() != 0)	{		return templateSequences[spot[0]];	}
		else					{		return emptySequence;				}
		
	}
	catch(exception& e) {
		m->errorOut(e, "AlignmentDB", "findClosestSequence");
		exit(1);
	}
}
/**************************************************************************************************/