Skip to content
Open
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
7 changes: 7 additions & 0 deletions doc/howto/write_an_alphabet/index.md
Original file line number Diff line number Diff line change
Expand Up @@ -188,3 +188,10 @@ from '1', 't' and 'T' for value 1 as well as from '0', 'f' and 'F' for value 0.
\note
You should really make your alphabet types [no-throw-default-constructible](\ref std::is_nothrow_default_constructible)
if you can!

## Adaptation of seqan2::Dna

In case you want to use `seqan2::Dna` in your SeqAn3 application, you can adapt it to the SeqAn3 alphabet
concepts by specialising the `seqan3::custom::alphabet` class template.

\include test/unit/alphabet/nucleotide/custom_alphabet_seqan2_dna.hpp
1 change: 1 addition & 0 deletions test/unit/alphabet/nucleotide/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -2,6 +2,7 @@
# SPDX-FileCopyrightText: 2016-2025 Knut Reinert & MPI für molekulare Genetik
# SPDX-License-Identifier: CC0-1.0

seqan3_test (custom_alphabet_seqan2_dna_test.cpp)
seqan3_test (dna4_test.cpp)
seqan3_test (dna5_test.cpp)
seqan3_test (dna15_test.cpp)
Expand Down
72 changes: 72 additions & 0 deletions test/unit/alphabet/nucleotide/custom_alphabet_seqan2_dna.hpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,72 @@
// SPDX-FileCopyrightText: 2006-2025 Knut Reinert & Freie Universität Berlin
// SPDX-FileCopyrightText: 2016-2025 Knut Reinert & MPI für molekulare Genetik
// SPDX-License-Identifier: CC0-1.0

#pragma once

#include <seqan3/alphabet/concept.hpp>

#include <seqan/basic.h>
#include <seqan/modifier/modifier_functors.h>

// Ranks and letters of `seqan2::Dna` and `seqan3::dna4` are the same:
// 'A' = 0, 'C' = 1, 'G' = 2, 'T' = 3.
template <>
struct seqan3::custom::alphabet<seqan2::Dna>
{
using alphabet_t = seqan2::Dna;

static constexpr size_t alphabet_size = 4u;

static uint8_t to_rank(alphabet_t const a) noexcept
{
return seqan2::ordValue(a);
}

static alphabet_t & assign_rank_to(uint8_t const r, alphabet_t & a) noexcept
{
seqan2::assign(a, r);
return a;
}

static char to_char(alphabet_t const a) noexcept
{
char c;
seqan2::assign(c, a);
return c;
}

static alphabet_t & assign_char_to(char const c, alphabet_t & a) noexcept
{
seqan2::assign(a, c);
return a;
}

static alphabet_t complement(alphabet_t const a) noexcept
{
static seqan2::FunctorComplement<alphabet_t> func;
return func(a);
}
};

// Disable `seqan2::begin` overload for `seqan::Dna`.
// Otherwise, `seqan2::Dna` would be considered a range.
// `seqan2::Dna` is indeed not a range, and assuming otherwise breaks many views.
namespace seqan2
{

void begin(seqan2::Dna) = delete;

} // namespace seqan2

// This is only needed if you want to use `seqan3::debug_stream` with `seqan2::Dna`.
namespace seqan3
{

inline debug_stream_type<char> & operator<<(debug_stream_type<char> & stream, seqan2::Dna const & data)
{
stream << seqan2::convert<char>(data);
return stream;
}

} // namespace seqan3
28 changes: 28 additions & 0 deletions test/unit/alphabet/nucleotide/custom_alphabet_seqan2_dna_test.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,28 @@
// SPDX-FileCopyrightText: 2006-2025 Knut Reinert & Freie Universität Berlin
// SPDX-FileCopyrightText: 2016-2025 Knut Reinert & MPI für molekulare Genetik
// SPDX-License-Identifier: BSD-3-Clause

#include <gtest/gtest.h>

#ifdef SEQAN3_HAS_SEQAN2

// # include "../alphabet_constexpr_test_template.hpp"
# include "../alphabet_test_template.hpp"
// # include "../semi_alphabet_constexpr_test_template.hpp"
# include "../semi_alphabet_test_template.hpp"
# include "custom_alphabet_seqan2_dna.hpp"
# include "nucleotide_test_template.hpp"

template <>
struct nucleotide<seqan2::Dna> : public ::testing::Test
{
static constexpr bool skip_trivial_thirdparty = true;
};

INSTANTIATE_TYPED_TEST_SUITE_P(seqan2_dna, alphabet, seqan2::Dna, );
INSTANTIATE_TYPED_TEST_SUITE_P(seqan2_dna, semi_alphabet_test, seqan2::Dna, );
// INSTANTIATE_TYPED_TEST_SUITE_P(seqan2_dna, alphabet_constexpr, seqan2::Dna, );
// INSTANTIATE_TYPED_TEST_SUITE_P(seqan2_dna, semi_alphabet_constexpr, seqan2::Dna, );
INSTANTIATE_TYPED_TEST_SUITE_P(seqan2_dna, nucleotide, seqan2::Dna, );

#endif
19 changes: 13 additions & 6 deletions test/unit/alphabet/nucleotide/nucleotide_test_template.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -10,13 +10,20 @@
#include <seqan3/alphabet/nucleotide/concept.hpp>

template <typename t>
using nucleotide = ::testing::Test;
struct nucleotide : public ::testing::Test
{
// Do not check the `trivial` concept, because `t` is a third party class and not trivial.
static constexpr bool skip_trivial_thirdparty = false;
};

TYPED_TEST_SUITE_P(nucleotide);

TYPED_TEST_P(nucleotide, concept_check)
{
EXPECT_TRUE(seqan3::trivial<TypeParam>);
if constexpr (!TestFixture::skip_trivial_thirdparty)
{
EXPECT_TRUE(seqan3::trivial<TypeParam>);
}

EXPECT_TRUE(seqan3::nucleotide_alphabet<TypeParam>);
EXPECT_TRUE(seqan3::nucleotide_alphabet<TypeParam &>);
Expand All @@ -26,10 +33,10 @@ TYPED_TEST_P(nucleotide, concept_check)

TYPED_TEST_P(nucleotide, complement)
{
EXPECT_EQ(seqan3::complement(TypeParam{}.assign_char('A')), TypeParam{}.assign_char('T'));
EXPECT_EQ(seqan3::complement(TypeParam{}.assign_char('C')), TypeParam{}.assign_char('G'));
EXPECT_EQ(seqan3::complement(TypeParam{}.assign_char('G')), TypeParam{}.assign_char('C'));
EXPECT_EQ(seqan3::complement(TypeParam{}.assign_char('T')), TypeParam{}.assign_char('A'));
EXPECT_EQ(seqan3::complement(seqan3::assign_char_to('A', TypeParam{})), seqan3::assign_char_to('T', TypeParam{}));
EXPECT_EQ(seqan3::complement(seqan3::assign_char_to('C', TypeParam{})), seqan3::assign_char_to('G', TypeParam{}));
EXPECT_EQ(seqan3::complement(seqan3::assign_char_to('G', TypeParam{})), seqan3::assign_char_to('C', TypeParam{}));
EXPECT_EQ(seqan3::complement(seqan3::assign_char_to('T', TypeParam{})), seqan3::assign_char_to('A', TypeParam{}));

using vsize_t = std::decay_t<decltype(seqan3::alphabet_size<TypeParam>)>;

Expand Down
Loading