Provided by: libbio-perl-run-perl_1.6.9-2_all bug

NAME

       Bio::Tools::Run::Alignment::Clustalw - Object for the calculation of a multiple sequence
       alignment from a set of unaligned sequences or alignments using the Clustalw program

SYNOPSIS

         #  Build a clustalw alignment factory
         @params = ('ktuple' => 2, 'matrix' => 'BLOSUM');
         $factory = Bio::Tools::Run::Alignment::Clustalw->new(@params);

         #  Pass the factory a list of sequences to be aligned.
         $inputfilename = 't/data/cysprot.fa';
         $aln = $factory->align($inputfilename); # $aln is a SimpleAlign object.
         # or
         $seq_array_ref = \@seq_array;
         # where @seq_array is an array of Bio::Seq objects
         $aln = $factory->align($seq_array_ref);

         # Or one can pass the factory a pair of (sub)alignments
         #to be aligned against each other, e.g.:
         $aln = $factory->profile_align($aln1,$aln2);
         # where $aln1 and $aln2 are Bio::SimpleAlign objects.

         # Or one can pass the factory an alignment and one or more unaligned
         # sequences to be added to the alignment. For example:
         $aln = $factory->profile_align($aln1,$seq); # $seq is a Bio::Seq object.

         # Get a tree of the sequences
         $tree = $factory->tree(\@seq_array);

         # Get both an alignment and a tree
         ($aln, $tree) = $factory->run(\@seq_array);

         # Do a footprinting analysis on the supplied sequences, getting back the
         # most conserved sub-alignments
         my @results = $factory->footprint(\@seq_array);
         foreach my $result (@results) {
           print $result->consensus_string, "\n";
         }

         # There are various additional options and input formats available.
         # See the DESCRIPTION section that follows for additional details.

DESCRIPTION

       Note: this DESCRIPTION only documents the Bioperl interface to Clustalw.  Clustalw,
       itself, is a large & complex program - for more information regarding clustalw, please see
       the clustalw documentation which accompanies the clustalw distribution. Clustalw is
       available from (among others) ftp://ftp.ebi.ac.uk/pub/software/. Clustalw.pm has only been
       tested using version 1.8 of clustalw.  Compatibility with earlier versions of the clustalw
       program is currently unknown. Before running Clustalw successfully it will be necessary:
       to install clustalw on your system, and to ensure that users have execute privileges for
       the clustalw program.

   Helping the module find your executable
       You will need to enable Clustalw to find the clustalw program. This can be done in (at
       least) three ways:

        1. Make sure the clustalw executable is in your path so that
           which clustalw
           returns a clustalw executable on your system.

        2. Define an environmental variable CLUSTALDIR which is a
           directory which contains the 'clustalw' application:
           In bash:

           export CLUSTALDIR=/home/username/clustalw1.8

           In csh/tcsh:

           setenv CLUSTALDIR /home/username/clustalw1.8

        3. Include a definition of an environmental variable CLUSTALDIR in
           every script that will use this Clustalw wrapper module, e.g.:

           BEGIN { $ENV{CLUSTALDIR} = '/home/username/clustalw1.8/' }
           use Bio::Tools::Run::Alignment::Clustalw;

       If you are running an application on a webserver make sure the webserver environment has
       the proper PATH set or use the options 2 or 3 to set the variables.

   How it works
       Bio::Tools::Run::Alignment::Clustalw is an object for performing a multiple sequence
       alignment from a set of unaligned sequences and/or sub-alignments by means of the clustalw
       program.

       Initially, a clustalw "factory object" is created. Optionally, the factory may be passed
       most of the parameters or switches of the clustalw program, e.g.:

               @params = ('ktuple' => 2, 'matrix' => 'BLOSUM');
               $factory = Bio::Tools::Run::Alignment::Clustalw->new(@params);

       Any parameters not explicitly set will remain as the defaults of the clustalw program.
       Additional parameters and switches (not available in clustalw) may also be set.
       Currently, the only such parameter is "quiet", which when set to a non-zero value,
       suppresses clustalw terminal output. Not all clustalw parameters are supported at this
       stage.

       By default, Clustalw output is returned solely in a the form of a Bio::SimpleAlign object
       which can then be printed and/or saved in multiple formats using the AlignIO.pm module.
       Optionally the raw clustalw output file can be saved if the calling script specifies an
       output file (with the clustalw parameter OUTFILE).  Currently only the GCG-MSF output file
       formats is supported.

       Not all parameters and features have been implemented yet in Perl format.

       Alignment parameters can be changed and/or examined at any time after the factory has been
       created.  The program checks that any parameter/switch being set/read is valid.  However,
       currently no additional checks are included to check that parameters are of the proper
       type (eg string or numeric) or that their values are within the proper range.  As an
       example, to change the value of the clustalw parameter ktuple to 3 and subsequently to
       check its value one would write:

               $ktuple = 3;
               $factory->ktuple($ktuple);
               $get_ktuple = $factory->ktuple();

       Once the factory has been created and the appropriate parameters set, one can call the
       method align() to align a set of unaligned sequences, or call profile_align() to add one
       or more sequences or a second alignment to an initial alignment.

       Input to align() may consist of a set of unaligned sequences in the form of the name of
       file containing the sequences. For example,

         $inputfilename = 't/data/cysprot.fa';
         $aln = $factory-E<gt>align($inputfilename);

       Alternately one can create an array of Bio::Seq objects somehow

               $str = Bio::SeqIO->new(-file=> 't/data/cysprot.fa', -format => 'Fasta');
               @seq_array =();
               while ( my $seq = $str->next_seq() ) {push (@seq_array, $seq) ;}

       and pass the factory a reference to that array

               $seq_array_ref = \@seq_array;
               $aln = $factory->align($seq_array_ref);

       In either case, align() returns a reference to a SimpleAlign object which can then used
       (see Bio::SimpleAlign).

       Once an initial alignment exists, one can pass the factory additional sequence(s) to be
       added (ie aligned) to the original alignment.  The alignment can be passed as either an
       alignment file or a Bio:SimpleAlign object.  The unaligned sequence(s) can be passed as a
       filename or as an array of BioPerl sequence objects or as a single BioPerl Seq object.
       For example (to add a single sequence to an alignment),

               $str = Bio::AlignIO->new(-file=> 't/data/cysprot1a.msf');
               $aln = $str->next_aln();
               $str1 = Bio::SeqIO->new(-file=> 't/data/cysprot1b.fa');
               $seq = $str1->next_seq();
               $aln = $factory->profile_align($aln,$seq);

       In either case, profile_align() returns a reference to a SimpleAlign object containing a
       new SimpleAlign object of the alignment with the additional sequence(s) added in.

       Finally one can pass the factory a pair of (sub)alignments to be aligned against each
       other.  The alignments can be passed in the form of either a pair of alignment files or a
       pair of Bio:SimpleAlign objects. For example,

               $profile1 = 't/data/cysprot1a.msf';
               $profile2 = 't/data/cysprot1b.msf';
               $aln = $factory->profile_align($profile1,$profile2);

       or

               $str1 = Bio::AlignIO->new(-file=> 't/data/cysprot1a.msf');
               $aln1 = $str1->next_aln();
               $str2 = Bio::AlignIO->new(-file=> 't/data/cysprot1b.msf');
               $aln2 = $str2->next_aln();
               $aln = $factory->profile_align($aln1,$aln2);

       In either case, profile_align() returns a reference to a SimpleAlign object containing an
       (super)alignment of the two input alignments.

       For more examples of syntax and use of Clustalw, the user is encouraged to look at the
       script Clustalw.t in the t/ directory.

       Note: Clustalw is still under development. Various features of the clustalw program have
       not yet been implemented.  If you would like that a specific clustalw feature be added to
       this perl contact bioperl-l@bioperl.org.

       These can be specified as parameters when instantiating a new Clustalw object, or through
       get/set methods of the same name (lowercase).

PARAMETER FOR ALIGNMENT COMPUTATION

   KTUPLE
        Title       : KTUPLE
        Description : (optional) set the word size to be used in the alignment
                      This is the size of exactly matching fragment that is used.
                      INCREASE for speed (max= 2 for proteins; 4 for DNA),
                      DECREASE for sensitivity.
                      For longer sequences (e.g. >1000 residues) you may
                      need to increase the default

   TOPDIAGS
        Title       : TOPDIAGS
        Description : (optional) number of best diagonals to use
                      The number of k-tuple matches on each diagonal
                      (in an imaginary dot-matrix plot) is calculated.
                      Only the best ones (with most matches) are used in
                      the alignment.  This parameter specifies how many.
                      Decrease for speed; increase for sensitivity.

   WINDOW
        Title       : WINDOW
        Description : (optional) window size
                      This is the number of diagonals around each of the 'best'
                      diagonals that will be used.  Decrease for speed;
                      increase for sensitivity.

   PAIRGAP
        Title       : PAIRGAP
        Description : (optional) gap penalty for pairwise alignments
                      This is a penalty for each gap in the fast alignments.
                      It has little affect on the speed or sensitivity except
                      for extreme values.

   FIXEDGAP
        Title       : FIXEDGAP
        Description : (optional) fixed length gap penalty

   FLOATGAP
        Title       : FLOATGAP
        Description : (optional) variable length gap penalty

   MATRIX
        Title       : MATRIX
        Default     : PAM100 for DNA - PAM250 for protein alignment
        Description : (optional) substitution matrix used in the multiple
                      alignments. Depends on the version of clustalw as to
                      what default matrix will be used

                      PROTEIN WEIGHT MATRIX leads to a new menu where you are
                      offered a choice of weight matrices. The default for
                      proteins in version 1.8 is the PAM series derived by
                      Gonnet and colleagues. Note, a series is used! The
                      actual matrix that is used depends on how similar the
                      sequences to be aligned at this alignment step
                      are. Different matrices work differently at each
                      evolutionary distance.

                      DNA WEIGHT MATRIX leads to a new menu where a single
                      matrix (not a series) can be selected. The default is
                      the matrix used by BESTFIT for comparison of nucleic
                      acid sequences.

   TYPE
        Title       : TYPE
        Description : (optional) sequence type: protein or DNA. This allows
                           you to explicitly overide the programs attempt at
                           guessing the type of the sequence.  It is only useful
                           if you are using sequences with a VERY strange
                           composition.

   OUTPUT
        Title       : OUTPUT
        Description : (optional) clustalw supports GCG or PHYLIP or PIR or
                       Clustal format.  See the Bio::AlignIO modules for
                       which formats are supported by bioperl.

   OUTFILE
        Title       : OUTFILE
        Description : (optional) Name of clustalw output file. If not set
                           module will erase output file.  In any case alignment will
                           be returned in the form of SimpleAlign objects

   TRANSMIT
        Title       : TRANSMIT
        Description : (optional) transitions not weighted.  The default is to
                           weight transitions as more favourable than other
                           mismatches in DNA alignments.  This switch makes all
                           nucleotide mismatches equally weighted.

FEEDBACK

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

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

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

       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.

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

         http://redmine.open-bio.org/projects/bioperl/

AUTHOR - Peter Schattner

       Email schattner@alum.mit.edu

CONTRIBUTORS

       Jason Stajich jason-AT-bioperl_DOT_org Sendu Bala    bix@sendu.me.uk

APPENDIX

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

   program_name
        Title   : program_name
        Usage   : $factory>program_name()
        Function: holds the program name
        Returns:  string
        Args    : None

   program_dir
        Title   : program_dir
        Usage   : $factory->program_dir(@params)
        Function: returns the program directory, obtained from ENV variable.
        Returns:  string
        Args    :

   version
        Title   : version
        Usage   : exit if $prog->version() < 1.8
        Function: Determine the version number of the program
        Example :
        Returns : float or undef
        Args    : none

   run
        Title   : run
        Usage   : ($aln, $tree) = $factory->run($inputfilename);
                  ($aln, $tree) = $factory->run($seq_array_ref);
        Function: Perform a multiple sequence alignment, generating a tree at the same
                  time. (Like align() and tree() combined.)
        Returns : A SimpleAlign object containing the sequence alignment and a
                  Bio::Tree::Tree object with the tree relating the sequences.
        Args    : Name of a file containing a set of unaligned fasta sequences
                  or else an array of references to Bio::Seq objects.

   align
        Title   : align
        Usage   : $inputfilename = 't/data/cysprot.fa';
                  $aln = $factory->align($inputfilename);
                  or
                  $seq_array_ref = \@seq_array; # @seq_array is array of Seq objs
                  $aln = $factory->align($seq_array_ref);
        Function: Perform a multiple sequence alignment
        Returns : Reference to a SimpleAlign object containing the
                  sequence alignment.
        Args    : Name of a file containing a set of unaligned fasta sequences
                  or else an array of references to Bio::Seq objects.

        Throws an exception if argument is not either a string (eg a
        filename) or a reference to an array of Bio::Seq objects.  If
        argument is string, throws exception if file corresponding to string
        name can not be found. If argument is Bio::Seq array, throws
        exception if less than two sequence objects are in array.

   profile_align
        Title   : profile_align
        Usage   : $aln = $factory->profile_align(@simple_aligns);
                  or
                  $aln = $factory->profile_align(@subalignment_filenames);
        Function: Perform an alignment of 2 (sub)alignments
        Returns : Reference to a SimpleAlign object containing the (super)alignment.
        Args    : Names of 2 files containing the subalignments
                  or references to 2 Bio::SimpleAlign objects.

       Throws an exception if arguments are not either strings (eg filenames) or references to
       SimpleAlign objects.

   add_sequences
        Title   : add_sequences
        Usage   :
        Function: Align and add sequences into an alignment
        Example :
        Returns : Reference to a SimpleAlign object containing the (super)alignment.
        Args    : Names of 2 files, the first one containing an alignment and the second one containing sequences to be added
                or references to 2 Bio::SimpleAlign objects.

       Throws an exception if arguments are not either strings (eg filenames) or references to
       SimpleAlign objects.

   tree
        Title   : tree
        Usage   : @params = ('bootstrap' => 1000,
                                   'tossgaps'  => 1,
                                   'kimura'    => 1,
                                   'seed'      => 121,
                                   'bootlabels'=> 'nodes',
                                   'quiet'     => 1);
                  $factory = Bio::Tools::Run::Alignment::Clustalw->new(@params);
                  $tree_obj = $factory->tree($aln_obj);
                  or
                  $tree_obj = $factory->tree($treefilename);
        Function: Retrieve a tree corresponding to the input
        Returns : Bio::TreeIO object
        Args    : Bio::SimpleAlign or filename of a tree

   footprint
        Title   : footprint
        Usage   : @alns = $factory->footprint($treefilename, $window_size, $diff);
                  @alns = $factory->footprint($seqs_array_ref);
        Function: Aligns all the supplied sequences and slices out of the alignment
                  those regions along a sliding window who's tree length differs
                  significantly from the total average tree length.
        Returns : list of Bio::SimpleAlign objects
        Args    : first argument as per run(), optional second argument to specify
                  the size of the sliding window (default 5 bp) and optional third
                  argument to specify the % difference from the total tree length
                  needed for a window to be considered a footprint (default 33%).

   _run
        Title   : _run
        Usage   : Internal function, not to be called directly
        Function: makes actual system call to clustalw program
        Returns : nothing; clustalw output is written to a
                  temporary file
        Args    : Name of a file containing a set of unaligned fasta sequences
                  and hash of parameters to be passed to clustalw

   _setinput()
        Title   : _setinput
        Usage   : Internal function, not to be called directly
        Function: Create input file for clustalw program
        Returns : name of file containing clustalw data input
        Args    : Seq or Align object reference or input file name

   _setparams()
        Title   : _setparams
        Usage   : Internal function, not to be called directly
        Function: Create parameter inputs for clustalw program
        Returns : parameter string to be passed to clustalw
                  during align or profile_align
        Args    : name of calling object