FU Logo
  • Startseite
  • Kontakt
  • Impressum
  • Home
  • Listenauswahl
  • Anleitungen

Re: [Seqan-dev] local alignment with blosum scoring matrix

<-- thread -->
<-- date -->
  • From: Young-ik Yang <maemukki@yahoo.com>
  • To: SeqAn Development <seqan-dev@lists.fu-berlin.de>
  • Date: Fri, 5 Aug 2011 16:01:05 -0700 (PDT)
  • Domainkey-signature: a=rsa-sha1; q=dns; c=nofws; s=s1024; d=yahoo.com; h=X-YMail-OSG:Received:X-Mailer:References:Message-ID:Date:From:Subject:To:In-Reply-To:MIME-Version:Content-Type:Content-Transfer-Encoding; b=dcbLY8N4CNMcTcTRhI85hqjJ51BCDpfdvE63kzf7zEhCIDck51RcbpixQBe5qGmmWHENPu2Bc5gSN8e9BCNlO4JB8f9lOLYBI0lXllUt1ECPPmN3PPDeNGBrFh7CuVJ1Qj9fVbfnQaFJnxDEY+nf8SnQDijFm7Uws/nfwZ6V0j8=;
  • Reply-to: SeqAn Development <seqan-dev@lists.fu-berlin.de>
  • Subject: Re: [Seqan-dev] local alignment with blosum scoring matrix

Dear Manuel,

Thank you for the reply.
After several tries and errors with looking up header file and example, I figure 
it out that SmithWaterman with scoring matrix seems to work only with Alignment 
Graph.

Is there any useful functions in Alignment Graph similar to 
clippedBeginPosition, countGaps, etc?
What I can find so far is converAlignment, which just converts alignment to 
string object. 

Do I need to write my own function to further process the string object to get 
the alignment start/end, gap op/ext, etc?
Or, is there any function that I directly utilize to get alignment summary?

Thank you,
Youngik.


#include <iostream>
#include <seqan/align.h>
#include <seqan/graph_msa.h>
using namespace seqan;

int main()
{
    typedef String<AminoAcid> TSequence;
    StringSet<TSequence> seq;
    appendValue(seq, "aphilologicaltheorem");
    appendValue(seq, "bizarreamphibology");
    Graph<Alignment<StringSet<TSequence, Dependent<> > > > aliG(seq);
    int score = localAlignment(aliG, Blosum62(-1, -11), SmithWaterman());
    std::cout << "score:" << score << "\n";
    ::std::cout << aliG << ::std::endl;
    ::std::cout << "Aligns Seq1[" << getFirstCoveredPosition(aliG,0) << ":" << 
getLastCoveredPosition(aliG,1) << "]";
    std::cout << stringSet(aliG)[0] << "\n";

    std::string mat;
    if (convertAlignment(aliG, mat)) {
      std::cout << mat << "\n";
    }
    return 0;
}








----- Original Message ----
From: Manuel Holtgrewe <manuel.holtgrewe@fu-berlin.de>
To: SeqAn Development <seqan-dev@lists.fu-berlin.de>
Sent: Fri, August 5, 2011 3:25:55 PM
Subject: Re: [Seqan-dev] local alignment with blosum scoring matrix

Dear Youngik,

without a working minimal example, we will not be able to help you with your 
code.

There is an example for working with Blosum matrices and global alignments in 
[1] as you probably already know. Maybe you can combine it with the local 
alignment example in [2] to achieve what you need?

Bests,
Manuel

[1] 
http://trac.mi.fu-berlin.de/seqan/wiki/Tutorial/Alignments#MultipleAlignments
[2] 
http://trac.mi.fu-berlin.de/seqan/browser/trunk/seqan/core/demos/alignment_local.cpp


Am 05.08.2011 um 17:30 schrieb Young-ik Yang:

> Would you let me know which local alignment algorithm support scoring matrix
> such as blosum62?
> If any, can you show me a function call or an example?
> 
> I am using seqan 1.3.
> I am attaching error messages from local alignment call at the end of the
> message.
> 
> Thanks,
> Youngik.
> 
> 
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h: In
> function 'typename seqan::Value<TRight, 0>::Type
> seqan::_alignSmithWaterman(TTrace&, const TStringSet&, const TScore&, typename
> seqan::Value<THolder, 0>::Type&, TIndexPair&, TForbidden&) [with TTrace =
> seqan::String<unsigned char, seqan::Alloc<void> >, TStringSet =
> seqan::Align<seqan::String<seqan::SimpleType<unsigned char, 
seqan::AminoAcid_>,
> seqan::Alloc<void> >, seqan::ArrayGaps>, TScore = seqan::Score<int,
> seqan::ScoreMatrix<seqan::SimpleType<unsigned char, seqan::AminoAcid_>,
> seqan::Blosum62_> >, TIndexPair = size_t [2], TForbidden = seqan::Nothing]':
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h:381:
> instantiated from 'typename seqan::Value<TLeft, 0>::Type
> seqan::_localAlignment(const TStringSet&, const TScore&, seqan::SmithWaterman)
> [with TStringSet = seqan::Align<seqan::String<seqan::SimpleType<unsigned char,
> seqan::AminoAcid_>, seqan::Alloc<void> >, seqan::ArrayGaps>, TScore =
> seqan::Score<int, seqan::ScoreMatrix<seqan::SimpleType<unsigned char,
> seqan::AminoAcid_>, seqan::Blosum62_> >]'
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_interface.h:378:
> instantiated from 'TScoreValue seqan::localAlignment(const TStringSet&, const
> seqan::Score<TScoreValue, TScoreSpec>&, TTag) [with TStringSet =
> seqan::Align<seqan::String<seqan::SimpleType<unsigned char, 
seqan::AminoAcid_>,
> seqan::Alloc<void> >, seqan::ArrayGaps>, TScoreValue = int, TSpec =
> seqan::ScoreMatrix<seqan::SimpleType<unsigned char, seqan::AminoAcid_>,
> seqan::Blosum62_>, TTag = seqan::Tag<seqan::SmithWaterman_>]'
> alignment_local.cpp:65:   instantiated from here
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h:233: 
>error:
> no match for 'operator[]' in 'str[0]'
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h:234: 
>error:
> no match for 'operator[]' in 'str[1]'
> /home/yyang/bin/seqan/seqan/score/score_base.h: In function 'TValue
> seqan::score(const seqan::Score<TValue, TSpec>&, TPos1, TPos2, const TSeq1&,
> const TSeq2&) [with TValue = int, TSpec =
> seqan::ScoreMatrix<seqan::SimpleType<unsigned char, seqan::AminoAcid_>,
> seqan::Blosum62_>, TPos1 = long unsigned int, TPos2 = long unsigned int, TSeq1 
>=
> seqan::SimpleType<unsigned char, seqan::AminoAcid_>, TSeq2 =
> seqan::SimpleType<unsigned char, seqan::AminoAcid_>]':
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h:291:
> instantiated from 'typename seqan::Value<TRight, 0>::Type
> seqan::_alignSmithWaterman(TTrace&, const TStringSet&, const TScore&, typename
> seqan::Value<THolder, 0>::Type&, TIndexPair&, TForbidden&) [with TTrace =
> seqan::String<unsigned char, seqan::Alloc<void> >, TStringSet =
> seqan::Align<seqan::String<seqan::SimpleType<unsigned char, 
seqan::AminoAcid_>,
> seqan::Alloc<void> >, seqan::ArrayGaps>, TScore = seqan::Score<int,
> seqan::ScoreMatrix<seqan::SimpleType<unsigned char, seqan::AminoAcid_>,
> seqan::Blosum62_> >, TIndexPair = size_t [2], TForbidden = seqan::Nothing]'
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_smith_waterman.h:381:
> instantiated from 'typename seqan::Value<TLeft, 0>::Type
> seqan::_localAlignment(const TStringSet&, const TScore&, seqan::SmithWaterman)
> [with TStringSet = seqan::Align<seqan::String<seqan::SimpleType<unsigned char,
> seqan::AminoAcid_>, seqan::Alloc<void> >, seqan::ArrayGaps>, TScore =
> seqan::Score<int, seqan::ScoreMatrix<seqan::SimpleType<unsigned char,
> seqan::AminoAcid_>, seqan::Blosum62_> >]'
> /home/yyang/bin/seqan/seqan/graph_align/graph_align_interface.h:378:
> instantiated from 'TScoreValue seqan::localAlignment(const TStringSet&, const
> seqan::Score<TScoreValue, TScoreSpec>&, TTag) [with TStringSet =
> seqan::Align<seqan::String<seqan::SimpleType<unsigned char, 
seqan::AminoAcid_>,
> seqan::Alloc<void> >, seqan::ArrayGaps>, TScoreValue = int, TSpec =
> seqan::ScoreMatrix<seqan::SimpleType<unsigned char, seqan::AminoAcid_>,
> seqan::Blosum62_>, TTag = seqan::Tag<seqan::SmithWaterman_>]'
> alignment_local.cpp:65:   instantiated from here
> /home/yyang/bin/seqan/seqan/score/score_base.h:206: error: no match for
> 'operator[]' in 'seq1[pos1]'
> /home/yyang/bin/seqan/seqan/score/score_base.h:206: error: no match for
> 'operator[]' in 'seq2[pos2]'
> 
> _______________________________________________
> seqan-dev mailing list
> seqan-dev@lists.fu-berlin.de
> https://lists.fu-berlin.de/listinfo/seqan-dev

--Manuel Holtgrewe            manuel.holtgrewe@fu-berlin.de
Freie Universität Berlin        http://www.inf.fu-berlin.de/
Institut für Informatik            Phone: +49 30 838 75246
Takustraße 9                Algorithmic Bioinformatics
14195 Berlin                Room 021


_______________________________________________
seqan-dev mailing list
seqan-dev@lists.fu-berlin.de
https://lists.fu-berlin.de/listinfo/seqan-dev




<-- thread -->
<-- date -->
  • References:
    • [Seqan-dev] local alignment with blosum scoring matrix
      • From: Young-ik Yang <maemukki@yahoo.com>
    • Re: [Seqan-dev] local alignment with blosum scoring matrix
      • From: Manuel Holtgrewe <manuel.holtgrewe@fu-berlin.de>
  • seqan-dev - August 2011 - Archives indexes sorted by:
    [ thread ] [ subject ] [ author ] [ date ]
  • Complete archive of the seqan-dev mailing list
  • More info on this list...

Hilfe

  • FAQ
  • Dienstbeschreibung
  • ZEDAT Beratung
  • postmaster@lists.fu-berlin.de

Service-Navigation

  • Startseite
  • Listenauswahl

Einrichtung Mailingliste

  • ZEDAT-Portal
  • Mailinglisten Portal