]> git.donarmstrong.com Git - samtools.git/blob - misc/sam2vcf.pl
Added possibility to select indels only and fixed a bug in reporting homozygous indels.
[samtools.git] / misc / sam2vcf.pl
1 #!/usr/bin/perl -w
2
3 # VCF specs: http://www.1000genomes.org/wiki/doku.php?id=1000_genomes:analysis:vcf3.3
4
5 # Contact: pd3@sanger
6 # Version: 2009-12-08
7
8 use strict;
9 use warnings;
10 use Carp;
11
12 my $opts = parse_params();
13 do_pileup_to_vcf($opts);
14
15 exit;
16
17 #---------------
18
19 sub error
20 {
21     my (@msg) = @_;
22     if ( scalar @msg ) { croak(@msg); }
23     die
24         "Usage: sam2vcf.pl [OPTIONS] < in.pileup > out.vcf\n",
25         "Options:\n",
26         "   -i, --indels-only                Ignore SNPs.\n",
27         "   -r, --refseq <file.fa>           The reference sequence, required when indels are present.\n",
28         "   -s, --snps-only                  Ignore indels.\n",
29         "   -t, --column-title <string>      The column title.\n",
30         "   -h, -?, --help                   This help message.\n",
31         "\n";
32 }
33
34
35 sub parse_params
36 {
37     my %opts = ();
38
39     $opts{fh_in}  = *STDIN;
40     $opts{fh_out} = *STDOUT;
41
42     while (my $arg=shift(@ARGV))
43     {
44         if ( $arg eq '-r' || $arg eq '--refseq' ) { $opts{refseq}=shift(@ARGV); next; }
45         if ( $arg eq '-t' || $arg eq '--column-title' ) { $opts{title}=shift(@ARGV); next; }
46         if ( $arg eq '-s' || $arg eq '--snps-only' ) { $opts{snps_only}=1; next; }
47         if ( $arg eq '-i' || $arg eq '--indels-only' ) { $opts{indels_only}=1; next; }
48         if ( $arg eq '-?' || $arg eq '-h' || $arg eq '--help' ) { error(); }
49
50         error("Unknown parameter \"$arg\". Run -h for help.\n");
51     }
52     return \%opts;
53 }
54
55 sub iupac_to_gtype
56 {
57     my ($ref,$base) = @_;
58     my %iupac = (
59             'K' => ['G','T'],
60             'M' => ['A','C'],
61             'S' => ['C','G'],
62             'R' => ['A','G'],
63             'W' => ['A','T'],
64             'Y' => ['C','T'],
65             );
66     if ( !exists($iupac{$base}) ) 
67     { 
68         if ( $base ne 'A' && $base ne 'C' && $base ne 'G' && $base ne 'T' ) { error("FIXME: what is this [$base]?\n"); }
69         if ( $ref eq $base ) { return ('.','0/0'); }
70         return ($base,'1/1');
71     }
72     my $gt = $iupac{$base};
73     if ( $$gt[0] eq $ref  ) { return ($$gt[1],'0/1'); }
74     elsif ( $$gt[1] eq $ref ) { return ($$gt[0],'0/1'); }
75     return ("$$gt[0],$$gt[1]",'1/2');
76 }
77
78
79 sub parse_indel
80 {
81     my ($cons) = @_;
82     if ( $cons=~/^-/ ) 
83     { 
84         my $len = length($');
85         return "D$len"; 
86     }
87     elsif ( $cons=~/^\+/ ) { return "I$'"; }
88     elsif ( $cons eq '*' ) { return undef; }
89     error("FIXME: could not parse [$cons]\n");
90 }
91
92
93 # An example of the pileup format:
94 #   1       3000011 C       C       32      0       98      1       ^~,     A
95 #   1       3002155 *       +T/+T   53      119     52      5       +T      *       4       1       0
96 #   1       3003094 *       -TT/-TT 31      164     60      11      -TT     *       5       6       0
97 #   1       3073986 *       */-AAAAAAAAAAAAAA       3       3       45      9       *       -AAAAAAAAAAAAAA 7       2       0
98 #
99 sub do_pileup_to_vcf
100 {
101     my ($opts) = @_;
102
103     my $fh_in  = $$opts{fh_in};
104     my $fh_out = $$opts{fh_out};
105     my ($prev_chr,$prev_pos,$prev_ref);
106     my $refseq;
107     my $ignore_indels = $$opts{snps_only} ? 1 : 0;
108     my $ignore_snps   = $$opts{indels_only} ? 1 : 0;
109     my $title = exists($$opts{title}) ? $$opts{title} : 'data';
110
111     print $fh_out 
112         qq[##fileformat=VCFv3.3\n],
113         qq[##INFO=DP,1,Integer,"Total Depth"\n],
114         qq[##FORMAT=GT,1,String,"Genotype"\n],
115         qq[##FORMAT=GQ,1,Integer,"Genotype Quality"\n],
116         qq[##FORMAT=DP,1,Integer,"Read Depth"\n],
117         qq[#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO\tFORMAT\t$title\n]
118         ;
119
120     while (my $line=<$fh_in>)
121     {
122         chomp($line);
123         my (@items) = split(/\t/,$line);
124         if ( scalar @items<8 ) 
125         { 
126             error("\nToo few columns, does not look like output of 'samtools pileup -c': $line\n"); 
127         }
128         my ($chr,$pos,$ref,$cons,$cons_qual,$snp_qual,$rms_qual,$depth,$a1,$a2) = @items;
129
130         my ($alt,$gt);
131         if ( $ref eq '*' )
132         {
133             # An indel is involved.
134             if ( $ignore_indels ) { next; }
135
136             if (!defined $prev_chr || $chr ne $prev_chr || $pos ne $prev_pos) 
137             {
138                 if ( !$$opts{refseq} ) { error("Cannot do indels without the reference.\n"); }
139                 if ( !$refseq ) { $refseq = Fasta->new(file=>$$opts{refseq}); }
140                 $ref = $refseq->get_base($chr,$pos);
141             }
142             else { $ref = $prev_ref; }
143
144             # One of the alleles can be a reference and it can come in arbitrary order. In some
145             #   cases */* can be encountered. In such a case, look in the additional columns.
146             my ($al1,$al2) = split(m{/},$cons);
147             if ( $al1 eq $al2 && $al1 eq '*' ) { $al1=$a1; $al2=$a2; }
148             my $alt1 = parse_indel($al1);
149             my $alt2 = parse_indel($al2);
150             if ( !$alt1 && !$alt2 ) { error("FIXME: could not parse indel:\n", $line); }
151             if ( !$alt1 ) 
152             { 
153                 $alt=$alt2; 
154                 $gt='0/1'; 
155             }
156             elsif ( !$alt2 ) 
157             { 
158                 $alt=$alt1; 
159                 $gt='0/1'; 
160             }
161             elsif ( $alt1 eq $alt2 )
162             { 
163                 $alt="$alt1"; 
164                 $gt='1/1'; 
165             }
166             else
167             { 
168                 $alt="$alt1,$alt2"; 
169                 $gt='1/2'; 
170             }
171         }
172         else
173         {
174             if ( $ignore_snps ) { next; }
175
176             # SNP
177             ($alt,$gt) = iupac_to_gtype($ref,$cons);
178         }
179
180         print $fh_out "$chr\t$pos\t.\t$ref\t$alt\t$snp_qual\t0\tDP=$depth\tGT:GQ:DP\t$gt:$cons_qual:$depth\n";
181
182         $prev_ref = $ref;
183         $prev_pos = $pos;
184         $prev_chr = $chr;
185     }
186 }
187
188
189 #------------- Fasta --------------------
190 #
191 # Uses samtools to get a requested base from a fasta file. For efficiency, preloads
192 #   a chunk to memory. The size of the cached sequence can be controlled by the 'size'
193 #   parameter.
194 #
195 package Fasta;
196
197 use strict;
198 use warnings;
199 use Carp;
200
201 sub Fasta::new
202 {
203     my ($class,@args) = @_;
204     my $self = {@args};
205     bless $self, ref($class) || $class;
206     if ( !$$self{file} ) { $self->throw(qq[Missing the parameter "file"\n]); }
207     $$self{chr}  = undef;
208     $$self{from} = undef;
209     $$self{to}   = undef;
210     if ( !$$self{size} ) { $$self{size}=10_000_000; }
211     bless $self, ref($class) || $class;
212     return $self;
213 }
214
215 sub read_chunk
216 {
217     my ($self,$chr,$pos) = @_;
218     my $to = $pos + $$self{size};
219     my $cmd = "samtools faidx $$self{file} $chr:$pos-$to";
220     my @out = `$cmd`;
221     if ( $? ) { $self->throw("$cmd: $!"); }
222     my $line = shift(@out);
223     if ( !($line=~/^>$chr:(\d+)-(\d+)/) ) { $self->throw("Could not parse: $line"); }
224     $$self{chr}  = $chr;
225     $$self{from} = $1;
226     $$self{to}   = $2;
227     my $chunk = '';
228     while ($line=shift(@out))
229     {
230         chomp($line);
231         $chunk .= $line;
232     }
233     $$self{chunk} = $chunk;
234     return;
235 }
236
237 sub get_base
238 {
239     my ($self,$chr,$pos) = @_;
240     if ( !$$self{chr} || $chr ne $$self{chr} || $pos<$$self{from} || $pos>$$self{to} )
241     {
242         $self->read_chunk($chr,$pos);
243     }
244     my $idx = $pos - $$self{from};
245     return substr($$self{chunk},$idx,1);
246 }
247
248 sub throw
249 {
250     my ($self,@msg) = @_;
251     croak(@msg);
252 }