Cython bindings and Python interface to Opal, a SIMD-accelerated database search aligner.
Opal is a sequence aligner enabling fast sequence similarity search using either of the Smith-Waterman, semi-global or Needleman-Wunsch algorithms. It is used part of the SW#db method[1] to align a query sequence to multiple database sequences on CPU, using the multi-sequence vectorization method described in SWIPE[2]
PyOpal is a Python module that provides bindings to Opal using Cython. It implements a user-friendly, Pythonic interface to query a database of sequences and access the search results. It interacts with the Opal interface rather than with the CLI, which has the following advantages:
- no binary dependency: PyOpal is distributed as a Python package, so you can add it as a dependency to your project, and stop worrying about the Opal binary being present on the end-user machine.
- no intermediate files: Everything happens in memory, in a Python object you control, so you don't have to invoke the Opal CLI using a sub-process and temporary files.
- better portability: Opal uses SIMD to accelerate alignment scoring, but doesn't support dynamic dispatch, so it has to be compiled on the local machine to be able to use the full capabilities of the local CPU. PyOpal ships several versions of Opal instead, each compiled with different target features, and selects the best one for the local platform at runtime.
- wider platform support: The Opal code has been backported to work on SSE2 rather than SSE4.1, allowing PyOpal to run on older x86 CPUs (all x86 CPUs support it since 2003). In addition, Armv7 and Aarch64 CPUs are also supported if they implement NEON extensions. Finally, the C++ code of Opal has been modified to compile on Windows.
PyOpal is available for all modern versions (3.6+), optionally depending on
the lightweight Python package archspec
for runtime CPU feature detection.
It can be installed directly from PyPI, which hosts some pre-built x86-64 wheels for Linux, MacOS, and Windows, Aarch64 wheels for Linux and MacOS, as well as the code required to compile from source with Cython:
$ pip install pyopal
Otherwise, PyOpal is also available as a Bioconda package:
$ conda install -c bioconda pyopal
Check the install page of the documentation for other ways to install PyOpal on your machine.
All classes are imported in the main namespace pyopal
:
import pyopal
pyopal
can work with sequences passed as Python strings,
as well as with ASCII strings in bytes
objects:
query = "MAGFLKVVQLLAKYGSKAVQWAWANKGKILDWLNAGQAIDWVVSKIKQILGIK"
database = [
"MESILDLQELETSEEESALMAASTVSNNC",
"MKKAVIVENKGCATCSIGAACLVDGPIPDFEIAGATGLFGLWG",
"MAGFLKVVQILAKYGSKAVQWAWANKGKILDWINAGQAIDWVVEKIKQILGIK",
"MTQIKVPTALIASVHGEGQHLFEPMAARCTCTTIISSSSTF",
]
If you plan to reuse the database across several queries, you can store it in
a Database
,
which will keep sequences encoded according to
an Alphabet
:
database = pyopal.Database(database)
The top-level function pyopal.align
can be used to align a query
sequence against a database, using multithreading to process chunks
of the database in parallel:
for result in pyopal.align(query, database):
print(result.score, result.target_index, database[result.target_index])
See the API documentation for more examples, including how to use the internal API, and detailed reference of the parameters and result types.
Database
objects are thread safe through a
C++17 read/write lock
that prevents modification while the database is searched. In addition, the
Aligner.align
method is re-entrant and can be safely used to query the
same database in parallel with different queries across different threads:
import multiprocessing.pool
import pyopal
import Bio.SeqIO
queries = [
"MEQQIELDVLEISDLIAGAGENDDLAQVMAASCTTSSVSTSSSSSSS",
"MTQIKVPTALIASVHGEGQHLFEPMAARCTCTTIISSSSTF",
"MGAIAKLVAKFGWPIVKKYYKQIMQFIGEGWAINKIIDWIKKHI",
"MGPVVVFDCMTADFLNDDPNNAELSALEMEELESWGAWDGEATS",
]
database = pyopal.Database([
str(record.seq)
for record in Bio.SeqIO.parse("vendor/opal/test_data/db/uniprot_sprot12071.fasta", "fasta")
])
aligner = pyopal.Aligner()
with multiprocessing.pool.ThreadPool() as pool:
hits = dict(pool.map(lambda q: (q, aligner.align(q, database)), queries))
Found a bug ? Have an enhancement request ? Head over to the GitHub issue tracker if you need to report or ask something. If you are filing in on a bug, please include as much information as you can about the issue, and try to recreate the same bug in a simple, easily reproducible situation.
Contributions are more than welcome! See
CONTRIBUTING.md
for more details.
This project adheres to Semantic Versioning and provides a changelog in the Keep a Changelog format.
This library is provided under the MIT License.
Opal is developed by Martin Šošić and is distributed under the
terms of the MIT License as well. See vendor/opal/LICENSE
for more information.
This project is in no way not affiliated, sponsored, or otherwise endorsed by the Opal authors. It was developed by Martin Larralde during his PhD project at the European Molecular Biology Laboratory in the Zeller team.
- [1] Korpar Matija, Martin Šošić, Dino Blažeka, Mile Šikić. SW#db: ‘GPU-Accelerated Exact Sequence Similarity Database Search’. PLoS One. 2015 Dec 31;10(12):e0145857. doi:10.1371/journal.pone.0145857. PMID:26719890. PMC4699916.
- [2] Rognes Torbjørn. Faster Smith-Waterman database searches with inter-sequence SIMD parallelisation. BMC Bioinformatics. 2011 Jun 1;12:221. doi:10.1186/1471-2105-12-221. PMID:21631914.PMC3120707.