File: BAMFile.java

package info (click to toggle)
fastqc 0.11.9%2Bdfsg-6
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 1,388 kB
  • sloc: java: 6,900; perl: 285; xml: 56; sh: 33; makefile: 7
file content (224 lines) | stat: -rw-r--r-- 5,917 bytes parent folder | download | duplicates (3)
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
/**
 * Copyright Copyright 2010-12 Simon Andrews
 *
 *    This file is part of FastQC.
 *
 *    FastQC is free software; you can redistribute it and/or modify
 *    it under the terms of the GNU General Public License as published by
 *    the Free Software Foundation; either version 3 of the License, or
 *    (at your option) any later version.
 *
 *    FastQC is distributed in the hope that it will be useful,
 *    but WITHOUT ANY WARRANTY; without even the implied warranty of
 *    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *    GNU General Public License for more details.
 *
 *    You should have received a copy of the GNU General Public License
 *    along with FastQC; if not, write to the Free Software
 *    Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
 */
package uk.ac.babraham.FastQC.Sequence;

import java.io.File;
import java.io.FileInputStream;
import java.io.IOException;
import java.util.Iterator;
import java.util.List;

import htsjdk.samtools.CigarElement;
import htsjdk.samtools.CigarOperator;
import htsjdk.samtools.SamInputResource;
import htsjdk.samtools.SamReader;
import htsjdk.samtools.SamReaderFactory;
import htsjdk.samtools.SAMFormatException;
import htsjdk.samtools.SAMRecord;
import htsjdk.samtools.ValidationStringency;

public class BAMFile implements SequenceFile {

	private File file;
	private boolean onlyMapped;
	private long fileSize = 0;
	private long recordSize = 0;
	
	// We keep the file stream around just so we can see how far through
	// the file we've got.  We don't read from this directly, but it's the
	// only way to access the file pointer.
	private FileInputStream fis;

	private SamReader br;
	private String name;
	private Sequence nextSequence = null;
	Iterator<SAMRecord> it;
	
	
	protected BAMFile (File file, boolean onlyMapped) throws SequenceFormatException, IOException {
		this.file = file;
		fileSize = file.length();
		name = file.getName();
		this.onlyMapped = onlyMapped;

		fis = new FileInputStream(file);

		SamInputResource sir = SamInputResource.of(fis);
		SamReaderFactory srf = SamReaderFactory.makeDefault().validationStringency(ValidationStringency.SILENT);

		br = srf.open(sir);
		
		it = br.iterator();
		readNext();
	}
	
	public String name () {
		return name;
	}
		
	public int getPercentComplete() {
		if (!hasNext()) return 100;
		
		try {
			int percent = (int) (((double)fis.getChannel().position()/ fileSize)*100);
			return percent;
		} 
		catch (IOException e) {
			e.printStackTrace();
		}
		return 0;
	}

	public boolean isColorspace () {
		return false;
	}
		
	public boolean hasNext() {
		return nextSequence != null;
	}

	public Sequence next () throws SequenceFormatException {
		Sequence returnSeq = nextSequence;
		readNext();
		return returnSeq;
	}
	
	private void readNext() throws SequenceFormatException {
		
		SAMRecord record;
		
		while (true) {
			
			if (!it.hasNext()) {
				nextSequence = null;
				try {
					br.close();
					fis.close();
				}
				catch (IOException ioe) {
					ioe.printStackTrace();
				}
				return;
			}
		
			try {
				record = it.next();
			}
			catch (SAMFormatException sfe) {
				throw new SequenceFormatException(sfe.getMessage());
			}
		
			// We skip over entries with no mapping if that's what the user asked for
			if (onlyMapped && record.getReadUnmappedFlag()) {
				continue;
			}
			else {
				break;
			}
		}
		
		// This is a very rough calculation of the record size so we can approximately track progress
		// through the file.
		if (recordSize == 0) {
			recordSize = (record.getReadLength()*2)+150;
			if (br.type() == SamReader.Type.BAM_TYPE || br.type() == SamReader.Type.CRAM_TYPE) {
				recordSize /= 4;
			}
		}
				

		String sequence = record.getReadString();
		String qualities = record.getBaseQualityString();
		
		
		// TODO: TEST THIS!!!
		// If we're only working with mapped data then we need to exclude any regions which have been either
		// hard or soft clipped by our aligner.
		if (onlyMapped) {
			List<CigarElement> elements = record.getCigar().getCigarElements();

			// We need to clip the 3' end first otherwise the numbers at the 5' end won't be right.
			if (elements.get(elements.size()-1).getOperator().equals(CigarOperator.S)) {
				int value = elements.get(elements.size()-1).getLength();
				sequence = sequence.substring(0,sequence.length()-value);
				qualities = qualities.substring(0,qualities.length()-value);
				
			}

			
			if (elements.get(0).getOperator().equals(CigarOperator.S)) {
				int value = elements.get(0).getLength();
				sequence = sequence.substring(value);
				qualities = qualities.substring(value);
			}
			
			
		}
		

		// BAM/SAM files always show sequence relative to the top strand of
		// the mapped reference so if this sequence maps to the reverse strand
		// we need to reverse complement the sequence and reverse the qualities
		// to get the original orientation of the read.
		if (record.getReadNegativeStrandFlag()) {
			sequence = reverseComplement(sequence);
			qualities = reverse(qualities);
		}

		nextSequence = new Sequence(this, sequence, qualities, record.getReadName());

	}

	
	private String reverseComplement (String sequence) {
		
		char [] letters = reverse(sequence).toUpperCase().toCharArray();
		char [] rc = new char[letters.length];
		
		for (int i=0;i<letters.length;i++) {
			switch(letters[i]) {
			case 'G': rc[i] = 'C';break;
			case 'A': rc[i] = 'T';break;
			case 'T': rc[i] = 'A';break;
			case 'C': rc[i] = 'G';break;
			default: rc[i] = letters[i];
			}
		}
	
		return new String(rc);

	}
	
	private String reverse (String sequence) {
		char [] starting = sequence.toCharArray();
		char [] reversed = new char[starting.length];
		
		for (int i=0;i<starting.length;i++) {
			reversed[reversed.length-(1+i)] = starting[i];
		}
		
		return new String(reversed);
	}

	public File getFile() {
		return file;
	}
	
}