File: Molphy.pm

package info (click to toggle)
bioperl 1.7.2-3
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 49,564 kB
  • sloc: perl: 170,474; xml: 22,869; lisp: 2,034; sh: 1,990; makefile: 22
file content (308 lines) | stat: -rw-r--r-- 8,424 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
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
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
#
# BioPerl module for Bio::Tools::Phylo::Molphy
#
# Please direct questions and support issues to <bioperl-l@bioperl.org> 
#
# Cared for by Jason Stajich <jason-at-bioperl.org>
#
# Copyright Jason Stajich
#
# You may distribute this module under the same terms as perl itself

# POD documentation - main docs before the code

=head1 NAME

Bio::Tools::Phylo::Molphy - parser for Molphy output

=head1 SYNOPSIS

  use Bio::Tools::Phylo::Molphy;
  my $parser = Bio::Tools::Phylo::Molphy->new(-file => 'output.protml');
  while( my $r = $parser->next_result ) {
    # r is a Bio::Tools::Phylo::Molphy::Result object

    # print the model name
    print $r->model, "\n";

    # get the substitution matrix
    # this is a hash of 3letter aa codes -> 3letter aa codes representing
    # substitution rate
    my $smat = $r->substitution_matrix;
    print "Arg -> Gln substitution rate is %d\n", 
          $smat->{'Arg'}->{'Gln'}, "\n";

    # get the transition probablity matrix
    # this is a hash of 3letter aa codes -> 3letter aa codes representing
    # transition probabilty
    my $tmat = $r->transition_probability_matrix;
    print "Arg -> Gln transition probablity is %.2f\n", 
          $tmat->{'Arg'}->{'Gln'}, "\n";

    # get the frequency for each of the residues
    my $rfreqs = $r->residue_frequencies;

    foreach my $residue ( keys %{$rfreqs} ) {
       printf "residue %s  expected freq: %.2f observed freq: %.2f\n",
              $residue,$rfreqs->{$residue}->[0], $rfreqs->{$residue}->[1];     
    }

    my @trees;
    while( my $t = $r->next_tree ) {
        push @trees, $t;
    }

    print "search space is ", $r->search_space, "\n",
          "1st tree score is ", $trees[0]->score, "\n";

    # writing to STDOUT, use -file => '>filename' to specify a file
    my $out = Bio::TreeIO->new(-format => "newick");
    $out->write_tree($trees[0]); # writing only the 1st tree
  }

=head1 DESCRIPTION

A parser for Molphy output (protml,dnaml)

=head1 FEEDBACK

=head2 Mailing Lists

User feedback is an integral part of the evolution of this and other
Bioperl modules. Send your comments and suggestions preferably to
the Bioperl mailing list.  Your participation is much appreciated.

  bioperl-l@bioperl.org                  - General discussion
  http://bioperl.org/wiki/Mailing_lists  - About the mailing lists

=head2 Support 

Please direct usage questions or support issues to the mailing list:

I<bioperl-l@bioperl.org>

rather than to the module maintainer directly. Many experienced and 
reponsive experts will be able look at the problem and quickly 
address it. Please include a thorough description of the problem 
with code and data examples if at all possible.

=head2 Reporting Bugs

Report bugs to the Bioperl bug tracking system to help us keep track
of the bugs and their resolution. Bug reports can be submitted via the
web:

  https://github.com/bioperl/bioperl-live/issues

=head1 AUTHOR - Jason Stajich

Email jason-at-bioperl.org

=head1 APPENDIX

The rest of the documentation details each of the object methods.
Internal methods are usually preceded with a _

=cut


# Let the code begin...


package Bio::Tools::Phylo::Molphy;
use strict;

use Bio::Tools::Phylo::Molphy::Result;
use Bio::TreeIO;
use IO::String;

use base qw(Bio::Root::Root Bio::Root::IO);

=head2 new

 Title   : new
 Usage   : my $obj = Bio::Tools::Phylo::Molphy->new();
 Function: Builds a new Bio::Tools::Phylo::Molphy object 
 Returns : Bio::Tools::Phylo::Molphy
 Args    : -fh/-file => $val, # for initing input, see Bio::Root::IO


=cut

sub new {
  my($class,@args) = @_;

  my $self = $class->SUPER::new(@args);
  $self->_initialize_io(@args);

  return $self;
}

=head2 next_result

 Title   : next_result
 Usage   : my $r = $molphy->next_result
 Function: Get the next result set from parser data
 Returns : Bio::Tools::Phylo::Molphy::Result object
 Args    : none


=cut

sub next_result{
   my ($self) = @_;

   # A little statemachine for the parser here
   my ($state,$transition_ct,
       @transition_matrix, %transition_mat, @resloc,) = ( 0,0);
   my ( %subst_matrix, @treelines, @treedata, %frequencies);
   my ( $treenum,$possible_trees, $model);
   my ($trans_type,$trans_amount);
   my $parsed = 0;
   while( defined ( $_ = $self->_readline()) ) {
       $parsed = 1;
       if( /^Relative Substitution Rate Matrix/ ) {
	   if( %subst_matrix ) { 
	       $self->_pushback($_);
	       last;
	   }
	   $state = 0;
	   my ( @tempdata);
	   @resloc = ();
	   while( defined ($_ = $self->_readline) ) {
	       last if (/^\s+$/);
	       # remove leading/trailing spaces
	       s/^\s+//;
	       s/\s+$//;
	       my @data = split;
	       my $i = 0;
	       for my $l ( @data ) {
		   if( $l =~ /\D+/ ) { 
		       push @resloc, $l;
		   }
		   $i++;
	       }
	       push @tempdata, \@data;
	   }
	   my $i = 0;
	   for my $row ( @tempdata ) {
	       my $j = 0;
	       for my $col ( @$row ) {
		   if( $i == $j ) {
		       # empty string for diagonals
		       $subst_matrix{$resloc[$i]}->{$resloc[$j]} = '';
		   } else {
		       $subst_matrix{$resloc[$i]}->{$resloc[$j]} = $col;
		   }
		   $j++;
	       }
	       $i++;
	   }
       } elsif( /^Transition Probability Matrix/ ) {	   
	   if( /(1\.0e(5|7))\)\s+(\S+)/ ) {
	       $state = 1;
	       my $newtrans_type = "$3-$1";
	       $trans_amount = $1;
	       if( defined $trans_type ) {
		   # finish processing the transition_matrix
		   my $i =0;
		   foreach my $row ( @transition_matrix ) {
		       my $j = 0;
		       foreach my $col ( @$row ) {
			   $transition_mat{$trans_type}->{$resloc[$i]}->{$resloc[$j]} = $col;
			   $j++;
		       }
		       $i++;
		   }
	       }
	       $trans_type = $newtrans_type;
	       $transition_ct = 0;
	       @transition_matrix = ();
	   }
       } elsif ( /Acid Frequencies/ ) {
	   $state = 0;
	   $self->_readline(); # skip the next line
	   while( defined( $_ = $self->_readline) ) {
	       unless( /^\s+/) {
		   $self->_pushback($_);
		   last;
	       }
	       s/^\s+//;
	       s/\s+$//;
	       my ($index,$res,$model,$data) = split;
	       $frequencies{$res} = [ $model,$data];
	   }
       } elsif( /^(\d+)\s*\/\s*(\d+)\s+(.+)\s+model/ ) {
	   my @save = ($1,$2,$3);	   
	   # finish processing the transition_matrix
	   my $i =0;
	   foreach my $row ( @transition_matrix ) {
	       my $j = 0;
	       foreach my $col ( @$row ) {
		   $transition_mat{$trans_type}->{$resloc[$i]}->{$resloc[$j]} = $col;
		   $j++;
	       }
	       $i++;
	   }	   
	   if( defined $treenum ) { 	       
	       $self->_pushback($_);
	       last;
	   }
	   
	   $state = 2;	   
	   ($treenum,$possible_trees, $model) = @save;
	   $model =~ s/\s+/ /g;
       } elsif( $state == 1 ) {
	   next if( /^\s+$/ || /^\s+Ala/);
	   s/^\s+//;
	   s/\s+$//;
	   if( $trans_type eq '1PAM-1.0e7' ) {
	       # because the matrix is split up into 2-10 column sets 
	       push @{$transition_matrix[$transition_ct++]}, split ;	   
	       $transition_ct = 0 if $transition_ct % 20 == 0;
	   } elsif( $trans_type eq '1PAM-1.0e5' ) {
	       # because the matrix is split up into 2-10 column sets 
	       my ($res,@row) = split;
	       next if $transition_ct >= 20; # skip last 
	       push @{$transition_matrix[$transition_ct++]}, @row;	   	       
	   }
       } elsif( $state == 2 ) {
	   if( s/^(\d+)\s+(\-?\d+(\.\d+)?)\s+// ) {
	       push @treedata, [ $1,$2];
	   }
	   # save this for the end so that we can 
	   # be efficient and only open one tree parser
	   push @treelines, $_;
       }
   }
   # waiting till the end to do this, is it better
   my @trees;
   if( @treelines ) {
       my $strdat = IO::String->new(join('',@treelines));
       my $treeio = Bio::TreeIO->new(-fh => $strdat,
				    -format => 'newick');
       while( my $tree = $treeio->next_tree ) {
	   if( @treedata ) {
	       my $dat = shift @treedata;
	       # set the associated information
	       $tree->id($dat->[0]);
	       $tree->score($dat->[1]);
	   }
	   push @trees, $tree;
       }
   }
   return unless( $parsed );
   my $result = Bio::Tools::Phylo::Molphy::Result->new
       (-trees => \@trees,
	-substitution_matrix => \%subst_matrix,
	-frequencies         => \%frequencies,
	-model               => $model,
	-search_space        => $possible_trees,
	);
   while( my ($type,$mat) = each %transition_mat ) {
       $result->transition_probability_matrix( $type,$mat);
   }
   $result;
}

1;