Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix get_end VCF parser - ENSINT-2005 for main #172

Merged
merged 1 commit into from
Nov 15, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
35 changes: 13 additions & 22 deletions modules/Bio/EnsEMBL/IO/Parser/BaseVCF4.pm
Original file line number Diff line number Diff line change
Expand Up @@ -279,48 +279,39 @@ sub get_start {
# Sequence end

=head2 get_raw_end
Description : Return the end position of the feature
Description : Return the end position of the feature. Synonym to `get_end`
Returntype : Integer
=cut

sub get_raw_end {
my $self = shift;
my $info = $self->get_info;
my $end;
if (defined($info->{END})) {
$end = $info->{END};
}
elsif(defined($info->{SVLEN})) {
my $svlen = (split(',',$info->{SVLEN}))[0];
$end = $self->get_raw_start + abs($svlen);
}
else {
$end = $self->get_raw_start + length($self->get_raw_reference) - 1;
}
return $end;
return $self->get_end;
}


=head2 get_end
Description : Return the adjusted end position of the feature
Description : Return the end position of the feature
Returntype : Integer
=cut

sub get_end {
my $self = shift;
my $info = $self->get_info;
my $end;
if (defined($info->{END})) {
$end = $info->{END};
}
elsif(defined($info->{SVLEN})) {
my $alternatives = join(",", @{$self->get_alternatives});
if(defined($info->{SVLEN})) {
return unless $self->get_start;
my $svlen = (split(',',$info->{SVLEN}))[0];
$end = $self->get_start + abs($svlen)-1;
$end = $self->get_start + abs($svlen) - 1;
}
elsif (defined($info->{END})) {
$end = $info->{END};
}
elsif (($self->get_raw_info && $self->get_raw_info =~ /SVTYPE/) || ($alternatives && $alternatives =~ /\<|\[|\]|\>/)) {
$end = $self->get_start;
}
else {
return unless $self->get_start;
$end = $self->get_start + length($self->get_raw_reference) - 1;
$end = $self->get_raw_start + length($self->get_raw_reference) - 1;
}
return $end;
}
Expand Down
1 change: 1 addition & 0 deletions modules/t/input/data.vcf
Original file line number Diff line number Diff line change
Expand Up @@ -22,5 +22,6 @@
20 17330 . T A 3 q10 NS=3;DP=11;AF=0.017 GT:GQ:DP:HQ 0|0:49:3:58,50 0|1:3:5:65,3 0/0:41:3
20 1110696 rs6040355 A G,T 67 PASS NS=2;DP=10;AF=0.333,0.667;AA=T;DB GT:GQ:DP:HQ 1|2:21:6:23,27 2|1:2:0:18,2 2/2:35:4
20 1230237 . T . 47 PASS NS=3;DP=13;AA=T GT:GQ:DP:HQ 0|0:54:7:56,60 0|0:48:4:51,51 0/0:61:2
20 186302389 . TTA T 6 PASS . . . . .
20 1234567 microsat1 GTC G,GTCT 50 PASS NS=3;DP=9;AA=G GT:GQ:DP 0/1:35:4 0/2:17:2 1/1:40:3
20 1234567 INS0 C <ctg1> 6 PASS . . . . .
47 changes: 36 additions & 11 deletions modules/t/vcf.t
Original file line number Diff line number Diff line change
Expand Up @@ -35,7 +35,7 @@ is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
my $index = 0;
$test_sample = [$inds[$index], $test_row[9]];
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
Expand All @@ -49,7 +49,7 @@ is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
$index = 1;
$test_sample = [$inds[$index], $test_row[10]];
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
Expand All @@ -63,7 +63,7 @@ is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
$index = 2;
$test_sample = [$inds[$index], $test_row[11]];
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
Expand All @@ -77,30 +77,41 @@ is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
$index = 0;
$test_sample = [$inds[$index], $test_row[9]];
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
is_deeply($test_sample, $sample_info->[$index], 'Sample data');


note "Record 5";
ok ($parser->next(), "Loading fifth record");
ok ($parser->next(), "Loading the fifth record");
@test_row = (qw(20 186302389 . TTA T 6 PASS . . . . .));
is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
ok($parser->get_seqname eq '20', 'get_seqname');
ok($parser->get_start == 186302390, 'get_start');
ok($parser->get_end == 186302391, 'get_end');


note "Record 6";
ok ($parser->next(), "Loading the sixth record");
@test_row = (qw(20 1234567 microsat1 GTC), 'G,GTCT', 50, 'PASS', qw(NS=3;DP=9;AA=G GT:GQ:DP 0/1:35:4 0/2:17:2 1/1:40:3));
is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
$index = 0;
$test_sample = [$inds[$index], $test_row[9]];
$index = 1;
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
is_deeply($test_sample, $sample_info->[$index], 'Sample data');

print "\n> Testing the getters (only for the last record):\n";
print "> Testing the getters (only for the 6th record):\n";
ok($parser->get_seqname eq '20', 'get_seqname');
ok($parser->get_start == 1234568, 'get_start');
ok($parser->get_end == 1234570, 'get_end');
ok($parser->get_end == 1234569, 'get_end');
ok($parser->get_IDs->[0] eq 'microsat1', 'get_IDs');
ok($parser->get_reference eq 'GTC', 'get_reference');
ok($parser->get_alternatives->[0] eq 'G', 'get_alternatives');
Expand All @@ -116,7 +127,21 @@ ok($parser->get_individuals_genotypes($inds[$index])->{$inds[$index]} eq 'GTC/G'
ok($parser->get_samples_info($inds[$index])->{$inds[$index]}->{'GT'} eq '0/1', 'get_samples_info');
ok($parser->get_samples_genotypes($inds[$index])->{$inds[$index]} eq 'GTC/G', 'get_samples_genotypes');

note "> Testing metadata getters";

note "Record 7";
ok ($parser->next(), "Loading the seventh record");
@test_row = (qw(20 1234567 INS0 C <ctg1> 6 PASS . . . . .));
is_deeply($parser->{'record'},\@test_row,"Test basic parsing of a row");
note "> Testing each column of the row";
do_the_tests(\@test_row);
$index = 2;
$test_sample = [$inds[$index], $test_row[9+$index]];
$ind_info = $parser->get_raw_individuals_info($inds[$index]);
is_deeply($test_sample, $ind_info->[$index], 'Individual data (DEPRECATED)');
$sample_info = $parser->get_raw_samples_info($inds[$index]);
is_deeply($test_sample, $sample_info->[$index], 'Sample data');

note "\n> Testing metadata getters";
ok($parser->get_metadata_key_list eq 'FILTER, FORMAT, INFO, contig, fileDate, fileformat, header, phasing, reference, source', 'getMetadataKeyList');
ok($parser->get_metadata_by_pragma('fileDate') eq '20090805', 'getMetadataByPragma');
ok($parser->get_vcf_version eq 'VCFv4.2', 'getVCFversion');
Expand All @@ -125,7 +150,7 @@ ok($parser->get_metadata_field('contig', '20', 'length') eq 62435964, 'getMetaFi
ok(!defined($parser->get_metadata_field('contig', '20', 'URL')), 'getMetaField - non-existing field');
ok($parser->get_metadata_field('contig', 'ctg1', 'URL') =~ /ftp/, 'Multiple contigs');

note "> Testing format validation";
note "\n> Testing format validation";
$parser->reset();
$parser->shift_block;
ok ($parser->validate(), "Validating vcf format");
Expand Down