Introduction

Pybwa is a python module that makes it easy to align sequence data. It is a lightweight wrapper of bwa.

This page provides task-oriented examples followed by the full API reference.

Examples

Align short reads (bwa aln)

Use BwaAln for aligning short reads (<100bp, e.g. Illumina). The align() method returns one AlignedSegment per input read.

from pybwa import BwaAln

aln = BwaAln(prefix="/path/to/genome.fasta")
for rec in aln.align(queries=["GATTACA"]):
    print(rec)

gives:

read.1       0       chr1    1       37      7M      *       0       0       GATTACA *       XT:A:U  NM:i:0  X0:i:1  X1:i:0  XM:i:0  XO:i:0  XG:i:0  MD:Z:7

Align reads (bwa mem)

Use BwaMem for aligning longer reads (>70bp, Illumina 100bp+, PacBio, ONT). The align() method returns a list of AlignedSegment per input read, since a single query may produce primary, supplementary, and secondary alignments.

from pybwa import BwaMem

mem = BwaMem(prefix="/path/to/genome.fasta")
for recs in mem.align(queries=["CTCAAGGTTGTTGCAAGGGGGTCTATGTGAACAAA"]):
    for rec in recs:
        print(rec)

gives:

chr1 0       chr1    1       60      35M     *       0       0       CTCAAGGTTGTTGCAAGGGGGTCTATGTGAACAAA     *       NM:i:0  MD:Z:35 AS:i:35 XS:i:0

Customize alignment options

Both aligners accept an options object to control alignment behavior.

For BwaAln, construct a BwaAlnOptions and pass it to align():

from pybwa import BwaAln, BwaAlnOptions

aln = BwaAln(prefix="/path/to/genome.fasta")
opt = BwaAlnOptions()
opt.max_mismatches = 5
recs = aln.align(queries=["GATTACA"], opt=opt)

For BwaMem, construct a BwaMemOptions and pass it to align():

from pybwa import BwaMem, BwaMemOptions

mem = BwaMem(prefix="/path/to/genome.fasta")
opt = BwaMemOptions()
opt.min_seed_len = 32
recs = mem.align(queries=["CTCAAGGTTGTTGCAAGGGGGTCTATGTGAACAAA"], opt=opt)

Note: the finalize() method will both apply the presets as specified by the mode option, as well as scale various other options (-TdBOELU) based on the match_score. The presets and scaling will only be applied to other options that have not been modified from their defaults. After calling the finalize() method, the options are immutable, unless copy=True is passed to finalize() method, in which case a copy of the options are returned by the method. Regardless, the finalize() method does not need to be called before the pybwa.BwaMem.align() is invoked, as the latter will do so (making a local copy).

Reuse an index across multiple alignments

Load a BwaIndex once and pass it to multiple aligners to avoid re-loading the index into memory for each aligner:

from pybwa import BwaIndex, BwaAln, BwaMem

index = BwaIndex(prefix="/path/to/genome.fasta")
aln = BwaAln(index=index)
mem = BwaMem(index=index)

This is also useful when aligning multiple batches of reads; construct the aligner once and call align() or align() repeatedly.

Parse alternative hits from the XA tag

BWA reports alternative alignment locations in the XA SAM tag. Use to_xa_hits() to parse these into AuxHit objects:

from pybwa import BwaAln, to_xa_hits

aln = BwaAln(prefix="/path/to/genome.fasta")
rec = aln.align(queries=["GATTACA"])[0]
hits = to_xa_hits(rec)
for hit in hits:
    print(f"{hit.refname}:{hit.start}-{hit.end} edits={hit.edits}")

You can also parse an XA tag value directly (without the leading XA:Z: prefix):

hits = to_xa_hits("chr4,-97592047,24M,3;chr8,-32368749,24M,3;")

Control BWA verbosity

By default BWA outputs informational and warning messages to stderr. Use set_bwa_verbosity() to suppress or increase output:

from pybwa import BwaVerbosity, set_bwa_verbosity

set_bwa_verbosity(BwaVerbosity.QUIET)  # suppress all output

Use long-read presets (PacBio/ONT)

Use BwaMemMode to apply presets for PacBio or Oxford Nanopore reads. This adjusts multiple options (penalties, seed length, etc.) to be appropriate for the given read type:

from pybwa import BwaMem, BwaMemOptions, BwaMemMode

mem = BwaMem(prefix="/path/to/genome.fasta")
opt = BwaMemOptions(mode=BwaMemMode.PACBIO)
recs = mem.align(queries=["ACGT" * 100], opt=opt)

Available presets:

  • PACBIO – PacBio CLR reads

  • ONT2D – Oxford Nanopore 2D reads

  • INTRACTG – intra-species contig alignment

API versus Command-line Differences

The reported alignments from pybwa may differ from those reported by the bwa command line. This is true when bwa is run with a different number of threads (see bwa aln -t and bwa mem -t), or when using a different number of chunks (see bwa mem -K).

Finally, the following additions have been made to bwa aln/samse:

  1. The standard SAM tag HN is added. This is useful if we find too many hits (max_hits) and therefore no hits are reported in the XA tag, we can still know how many were found.

  2. The with_md option will add the standard SAM tag MD to the XA tag, otherwise . will be used. This provides additional information on the quality of alternative alignments.

API

Bwa Aln

class pybwa.BwaAlnOptions(int max_mismatches: int = -1, int max_gap_opens: int = 1, int max_gap_extensions: int = 6, int min_indel_to_end_distance: int = 5, int max_occurrences_for_extending_long_deletion: int = 10, int seed_length: int = 32, int max_mismatches_in_seed: int = 2, int mismatch_penalty: int = 3, int gap_open_penalty: int = 11, int gap_extension_penalty: int = 4, int stop_at_max_best_hits: int = 30, int max_hits: int = 3, bool log_scaled_gap_penalty: bool = False, bool find_all_hits: bool = False, bool with_md: bool = False, int max_entries: int = 2000000, int threads: int = 1)

The container for options for pybwa.BwaAln.

Parameters:
  • max_mismatches – Maximum number of mismatches allowed. When set to -1, the maximum edit distance is calculated automatically using a fraction of missing alignments (fnr=0.04). -n <int>

  • max_gap_opens – Maximum number of gap opens allowed in an alignment. -o <int>

  • max_gap_extensions – Maximum number of gap extensions allowed. When set to the default (6), gap extensions count against the max mismatches limit. -e <int>

  • min_indel_to_end_distance – Indels within this many bases of the read ends are not allowed. -i <int>

  • max_occurrences_for_extending_long_deletion – Maximum occurrences for extending a long deletion. -d <int>

  • seed_length – Length of the seed used in the first alignment phase. Only the first seed_length bases of each read are used for seeding. -l <int>

  • max_mismatches_in_seed – Maximum number of mismatches allowed in the seed region. -k <int>

  • mismatch_penalty – Mismatch penalty in the alignment scoring. A higher penalty makes the aligner more strict about mismatches. -M <int>

  • gap_open_penalty – Penalty for opening a gap in the alignment. -O <int>

  • gap_extension_penalty – Penalty for each base of gap extension. -E <int>

  • stop_at_max_best_hits – Stop searching when there are more than this many best hits. -R <int>

  • max_hits – Maximum number of alternative alignments to report in the XA tag. If a read has more hits than this, no alternative hits are reported. bwa samse -n <int>

  • log_scaled_gap_penalty – Use logarithmic-scaled gap penalty for long deletions. -L

  • find_all_hits – Find all hits with at most max_mismatches differences, without an upper limit on the number of hits. -N

  • with_md – Output the MD string for each alignment in the XA tag, otherwise use ".". bwa samse -d

  • max_entries – Maximum number of entries in the occurrence queue during alignment; limits memory usage for highly repetitive sequences. -m

  • threads – The number of threads to use for alignment. -t

find_all_hits

bool

bwa aln -N

Type:

BwaAlnOptions.find_all_hits

gap_extension_penalty

int

bwa aln -E <int>

Type:

BwaAlnOptions.gap_extension_penalty

gap_open_penalty

int

bwa aln -O <int>

Type:

BwaAlnOptions.gap_open_penalty

log_scaled_gap_penalty

bool

bwa aln -L

Type:

BwaAlnOptions.log_scaled_gap_penalty

max_entries

int

bwa aln -m

Type:

BwaAlnOptions.max_entries

max_gap_extensions

int

bwa aln -e <int>

Type:

BwaAlnOptions.max_gap_extensions

max_gap_opens

int

bwa aln -o <int>

Type:

BwaAlnOptions.max_gap_opens

max_hits

int

bwa samse -n <int>

Type:

BwaAlnOptions.max_hits

max_mismatches

int

bwa aln -n <int>

Type:

BwaAlnOptions.max_mismatches

max_mismatches_in_seed

int

bwa aln -k <int>

Type:

BwaAlnOptions.max_mismatches_in_seed

max_occurrences_for_extending_long_deletion

int

bwa aln -d <int>

Type:

BwaAlnOptions.max_occurrences_for_extending_long_deletion

min_indel_to_end_distance

int

bwa aln -i <int>

Type:

BwaAlnOptions.min_indel_to_end_distance

mismatch_penalty

int

bwa aln -M <int>

Type:

BwaAlnOptions.mismatch_penalty

seed_length

int

bwa aln -l <int>

Type:

BwaAlnOptions.seed_length

stop_at_max_best_hits

int

bwa aln -R <int>

Type:

BwaAlnOptions.stop_at_max_best_hits

threads

int

bwa aln -t

Type:

BwaAlnOptions.threads

validate() None

Validate the options and raise ValueError if any are invalid.

Raises:

ValueError – If any option values are invalid

Return type:

None

with_md

bool

bwa samse -d

Output the MD to each alignment in the XA tag, otherwise use ".".

Type:

BwaAlnOptions.with_md

class pybwa.BwaAln(prefix: str | Path | None = None, BwaIndex index: BwaIndex | None = None)

Aligner for short reads using the bwa aln algorithm.

This wraps bwa’s short-read alignment algorithm, which is based on backward search with the Burrows-Wheeler Transform. Best suited for Illumina reads shorter than ~100bp.

This class can be used as a context manager for automatic resource cleanup:

Example

>>> with BwaAln(prefix="genome.fa") as aligner:
...     results = aligner.align(queries)
align(queries: Sequence[FastxRecord | str] | None = None, opt: BwaAlnOptions | None = None) List[AlignedSegment]

Align one or more queries with bwa aln.

Parameters:
  • queries (Sequence[FastxRecord | str] | None, default: None) – the queries to align, as either strings or FastxRecord objects

  • opt (BwaAlnOptions | None, default: None) – the alignment options, or None to use the default options

Return type:

List[AlignedSegment]

Returns:

one alignment (AlignedSegment) per query

reinitialize_seed(self)

Re-initializes the random seed.

class pybwa.AuxHit(refname: str, start: int, negative: bool, cigar: Cigar, edits: int, md: str | None = None, rest: str | None = None) None

Represents a single hit or alignment of a sequence to a location in the genome.

Typically, this is parsed from the XA SAM tag.

refname

the reference name of the hit.

start

the start position of the hit (1-based inclusive).

negative

whether the hit is on the negative strand.

cigar

the cigar returned by BWA.

edits

the number of edits (mismatches only, no indels) between the read and the reference. Calculated as the NM value in the XA SAM tag minus the number of indel bases from the cigar.

md

if present, the MD value that’s appended to the entry in the XA SAM tag.

rest

if present, any string that’s appended to the entry in the XA SAM tag (after the MD).

property end: int

The end position of the hit (1-based inclusive).

property mismatches: int

The number of mismatches for the hit (not including indels).

pybwa.to_xa_hits(rec: AlignedSegment | str | bytes, max_hits: int | None = None)

Parses the value of the XA SAM tag, returning a list of AuxHits.

Parameters:

rec (AlignedSegment | str | bytes) – one of an AlignedSegment, string, or bytes. If a string or bytes are provided, the _value_ of the XA SAM tag must be provided (i.e. not including the leading XA:Z:).

Returns:

An empty list if a provided AlignedSegment has no “XA” tag, otherwise a list of AuxHits parsed from the XA SAM tag.

:raises A ValueError if the XA SAM tag could not be parsed.:

Bwa Mem

class pybwa.BwaMemOptions(int min_seed_len: int = 19, mode: BwaMemMode = None, int band_width: int = 100, int match_score: int = 1, int mismatch_penalty: int = 4, int minimum_score: int = 30, int unpaired_penalty: int = 17, bool skip_pairing: bool = False, bool output_all_for_fragments: bool = False, bool interleaved_paired_end: bool = False, bool short_split_as_secondary: bool = False, bool skip_mate_rescue: bool = False, bool soft_clip_supplementary: bool = False, bool with_xr_tag: bool = False, bool query_coord_as_primary: bool = False, bool keep_mapq_for_supplementary: bool = False, bool with_xb_tag: bool = False, int max_occurrences: int = 500, int off_diagonal_x_dropoff: int = 100, bool ignore_alternate_contigs: bool = False, double internal_seed_split_factor: float = 1.5, double drop_chain_fraction: float = 0.50, int max_mate_rescue_rounds: int = 50, int min_seeded_bases_in_chain: int = 0, int seed_occurrence_in_3rd_round: int = 20, xa_max_hits: int | tuple[int, int] = (5, 200), double xa_drop_ratio: float = 0.80, gap_open_penalty: int | tuple[int, int] = 6, gap_extension_penalty: int | tuple[int, int] = 1, clipping_penalty: int | tuple[int, int] = 5, int threads: int = 1, int chunk_size: int = 10000000)

The container for options for BwaMem.

Options may not be modified after finalize() has been called.

Parameters:
  • min_seed_len – Minimum seed length for the MEM algorithm. Shorter seeds increase sensitivity at the cost of speed. -k <int>

  • mode – Read type preset that overrides multiple options for specific sequencing platforms (PacBio, ONT, or intra-species contig alignment). -x <str>

  • band_width – Band width for banded Smith-Waterman alignment. Gaps longer than this are not found. -w <int>

  • match_score – Score for a matching base. Other scoring parameters are scaled by this if they are not explicitly set. -A <int>

  • mismatch_penalty – Penalty for a mismatching base. -B <int>

  • minimum_score – Minimum alignment score to report; alignments below this threshold are marked as unmapped. -T <int>

  • unpaired_penalty – Penalty for an unpaired read pair; used in paired-end mode only. -U <int>

  • skip_pairing – Skip read pairing; mate rescue is still performed unless also disabled. -P

  • output_all_for_fragments – Output all found alignments for single-end or unpaired reads, including those that are likely suboptimal. -a

  • interleaved_paired_end – Treat input queries as interleaved paired-end reads. -p

  • short_split_as_secondary – Mark shorter split hits as secondary rather than supplementary, for compatibility with Picard. -M

  • skip_mate_rescue – Skip mate rescue; faster but may reduce alignment quality for pairs where one end maps poorly. -S

  • soft_clip_supplementary – Use soft clipping for supplementary alignments instead of hard clipping. -Y

  • with_xr_tag – Add the XR tag with the full reference sequence name (including description) to each alignment. -V

  • query_coord_as_primary – For a multi-part alignment, use the alignment with the smallest query coordinate as the primary; also keeps mapping quality for supplementary alignments. -5

  • keep_mapq_for_supplementary – Don’t modify mapping quality of supplementary alignments. -q

  • with_xb_tag – Output the XB tag with base alignment quality (BAQ) scores. -u

  • max_occurrences – Discard a MEM if it has more than this many occurrences in the genome. -c <int>

  • off_diagonal_x_dropoff – Z-dropoff score for extending alignments: stop extension when the best score drops below the current by more than this value. -d <int>

  • ignore_alternate_contigs – Treat ALT contigs as part of the primary assembly (i.e., ignore the .alt file). -j

  • internal_seed_split_factor – Trigger re-seeding for a MEM longer than min_seed_len * internal_seed_split_factor. Larger values yield fewer seeds, faster alignment, but lower accuracy. -r <float>

  • drop_chain_fraction – Drop chains with a seed count below this fraction of the best overlapping chain. -D <float>

  • max_mate_rescue_rounds – Maximum rounds of mate-rescue Smith-Waterman per read. -m <int>

  • min_seeded_bases_in_chain – Discard a chain if the number of seeded bases is fewer than this value. -W <int>

  • seed_occurrence_in_3rd_round – Seed occurrence for the 3rd round of seeding; used when a read has very few seed hits. -y <int>

  • xa_max_hits – Maximum number of alternative hits to output in the XA tag, specified as a single int or a tuple of (primary, alt). -h <int<,int>>

  • xa_drop_ratio – Drop a chain in the XA tag if its score is below this fraction of the best chain’s score. -z <float>

  • gap_open_penalty – Gap open penalty, specified as a single int (same for deletions and insertions) or a tuple of (deletion, insertion). -O <int<,int>>

  • gap_extension_penalty – Gap extension penalty, specified as a single int or a tuple of (deletion, insertion). -E <int<,int>>

  • clipping_penalty – Penalty for 5’ and 3’ clipping, specified as a single int (same for both ends) or a tuple of (5’, 3’). -L <int<,int>>

  • threads – The number of threads to use for alignment. -t <int>

  • chunk_size – Number of bases processed in each batch; the actual batch size scales with the number of threads. -K <int>

band_width

int

bwa mem -w <int>

Type:

BwaMemOptions.band_width

chunk_size

int

bwa mem -K <int>

Type:

BwaMemOptions.chunk_size

clipping_penalty

int | tuple[int, int]

bwa mem -L <int<,int>>

Type:

BwaMemOptions.clipping_penalty

drop_chain_fraction

float

bwa mem -D <float>

Type:

BwaMemOptions.drop_chain_fraction

finalize(copy: bool = False) BwaMemOptions

Performs final initialization of these options. The object returned may not be modified further.

If the mode is given, then the presets are applied to options that have not been explicitly set. Otherwise, if the match score has been set, the match score scales the (-TdBOELU) options if they have not been explicitly set.

Parameters:

copy (bool, default: False) – true to return a finalized copy of this object, false to finalize this object

Return type:

BwaMemOptions

finalized

bool

True if the options have been finalized with finalize().

Type:

BwaMemOptions.finalized

gap_extension_penalty

int | tuple[int, int]

bwa mem -E <int<,int>>

Type:

BwaMemOptions.gap_extension_penalty

gap_open_penalty

int | tuple[int, int]

bwa mem -O <int<,int>>

Type:

BwaMemOptions.gap_open_penalty

ignore_alternate_contigs

bool

bwa mem -j

Type:

BwaMemOptions.ignore_alternate_contigs

interleaved_paired_end

bool

bwa mem -p

Type:

BwaMemOptions.interleaved_paired_end

internal_seed_split_factor

float

bwa mem -r <float>

Type:

BwaMemOptions.internal_seed_split_factor

keep_mapq_for_supplementary

bool

bwa mem -q

Type:

BwaMemOptions.keep_mapq_for_supplementary

match_score

int

bwa mem -A <int>

Type:

BwaMemOptions.match_score

max_mate_rescue_rounds

int

bwa mem -m <int>

Type:

BwaMemOptions.max_mate_rescue_rounds

max_occurrences

int

bwa mem -c <int>

Type:

BwaMemOptions.max_occurrences

min_seed_len

int

bwa mem -k <int>

Type:

BwaMemOptions.min_seed_len

min_seeded_bases_in_chain

int

bwa mem -W <int>

Type:

BwaMemOptions.min_seeded_bases_in_chain

minimum_score

int

bwa mem -T <int>

Type:

BwaMemOptions.minimum_score

mismatch_penalty

int

bwa mem -B <int>

Type:

BwaMemOptions.mismatch_penalty

mode

BwaMemMode | None

bwa mem -x <str>

Type:

BwaMemOptions.mode

off_diagonal_x_dropoff

int

bwa mem -d <int>

Type:

BwaMemOptions.off_diagonal_x_dropoff

output_all_for_fragments

bool

bwa mem -a

Type:

BwaMemOptions.output_all_for_fragments

query_coord_as_primary

bool

bwa mem -5

Type:

BwaMemOptions.query_coord_as_primary

seed_occurrence_in_3rd_round

int

bwa mem -y <int>

Type:

BwaMemOptions.seed_occurrence_in_3rd_round

short_split_as_secondary

bool

bwa mem -M

Type:

BwaMemOptions.short_split_as_secondary

skip_mate_rescue

bool

bwa mem -S

Type:

BwaMemOptions.skip_mate_rescue

skip_pairing

bool

bwa mem -P

Type:

BwaMemOptions.skip_pairing

soft_clip_supplementary

bool

bwa mem -Y

Type:

BwaMemOptions.soft_clip_supplementary

threads

int

bwa mem -t <int>

Note

BWA uses internal multithreading for alignment. Each BwaMem instance maintains its own thread pool. Multiple BwaMem instances can be used concurrently from different Python threads safely.

Type:

BwaMemOptions.threads

unpaired_penalty

int

bwa mem -U <int>

Type:

BwaMemOptions.unpaired_penalty

validate() None

Validate the options and raise ValueError if any are invalid.

Raises:

ValueError – If any option values are invalid

Return type:

None

with_xb_tag

bool

bwa mem -u

Type:

BwaMemOptions.with_xb_tag

with_xr_tag

bool

bwa mem -V

Type:

BwaMemOptions.with_xr_tag

xa_drop_ratio

float

bwa mem -z <float>

Type:

BwaMemOptions.xa_drop_ratio

xa_max_hits

int | tuple[int, int]

bwa mem -h <int<,int>>

Type:

BwaMemOptions.xa_max_hits

enum pybwa.BwaMemMode(value)

The read type preset for overriding multiple options (bwa mem -x).

Valid values are as follows:

PACBIO = <BwaMemMode.PACBIO: 1>
ONT2D = <BwaMemMode.ONT2D: 2>
INTRACTG = <BwaMemMode.INTRACTG: 3>
class pybwa.BwaMem(prefix: str | Path | None = None, BwaIndex index: BwaIndex | None = None)

Aligner for reads using the bwa mem algorithm.

This wraps bwa’s maximal exact match (MEM) alignment algorithm, which is the recommended algorithm for reads longer than ~70bp (e.g. Illumina 100bp+, PacBio, and Oxford Nanopore).

This class can be used as a context manager for automatic resource cleanup:

Example

>>> with BwaMem(prefix="genome.fa") as aligner:
...     results = aligner.align(queries)
align(queries: Sequence[FastxRecord | str], opt: BwaMemOptions | None = None) List[List[AlignedSegment]]

Align one or more queries with bwa mem.

Parameters:
  • queries (Sequence[FastxRecord | str]) – the queries to align, as either strings or FastxRecord objects

  • opt (BwaMemOptions | None, default: None) – the alignment options, or None to use the default options. If the options have not been finalized, a finalized copy will be made automatically.

Return type:

List[List[AlignedSegment]]

Returns:

a list of alignments (AlignedSegment) per query, since a single query may produce multiple alignments (primary, supplementary)

Bwa Index

class pybwa.BwaIndex(prefix: str | Path, bool bwt: bool = True, bool bns: bool = True, bool pac: bool = True)

Contains the index and nucleotide sequence for Bwa. Use bwa index on the command line to generate the bwa index.

Note: the accompanying sequence dictionary must exist (i.e. .dict file, which can be generated with samtools dict <fasta>).

Parameters:
  • prefix – the path prefix for the BWA index (typically a FASTA)

  • bwt – load the BWT (FM-index)

  • bns – load the BNS (reference sequence metadata)

  • pac – load the PAC (the actual 2-bit encoded reference sequences with ‘N’ converted to a random base)

header

object

Type:

header

classmethod index(fasta: str | Path, method: BwaIndexBuildMethod = BwaIndexBuildMethod.AUTO, prefix: str | Path | None = None, block_size: int = 10000000, out_64: bool = False) None

Indexes a given FASTA. Also builds the sequence dictionary (.dict).

Parameters:
  • fasta (str | Path) – the path to the FASTA to index

  • method (BwaIndexBuildMethod, default: <BwaIndexBuildMethod.AUTO: 1>) – the BWT construction algorithm (bwa index -a <str>)

  • prefix (str | Path | None, default: None) – the path prefix for the BWA index (typically a FASTA)

  • block_size (int, default: 10000000) – block size for the bwtsw algorithm (effective with -a bwtsw)

  • out_64 (bool, default: False) – index files named as <in.fasta>.64.* instead of <in.fasta>.*

Return type:

None

enum pybwa.BwaIndexBuildMethod(value)

The BWT construction algorithm (bwa index -a <str>)

Valid values are as follows:

AUTO = <BwaIndexBuildMethod.AUTO: 1>
RB2 = <BwaIndexBuildMethod.RB2: 2>
BWTSW = <BwaIndexBuildMethod.BWTSW: 3>
IS = <BwaIndexBuildMethod.IS: 4>

Verbosity

enum pybwa.BwaVerbosity(value)

The verbosity level for the BWA C-API.

Member Type:

int

Valid values are as follows:

QUIET = <BwaVerbosity.QUIET: 0>
ERROR = <BwaVerbosity.ERROR: 1>
WARNING = <BwaVerbosity.WARNING: 2>
INFO = <BwaVerbosity.INFO: 3>
DEBUG = <BwaVerbosity.DEBUG: 4>
pybwa.set_bwa_verbosity(level: BwaVerbosity) bool

Set the BWA C-API verbosity level.

By default BWA outputs informational and warning messages to stderr. Use this to suppress or increase the output level.

Warning

This function modifies a global C variable and is NOT thread-safe. If you are using pybwa in a multi-threaded environment, call this function before creating any threads, or ensure proper synchronization.

Parameters:

level (BwaVerbosity) – the desired verbosity level

Return type:

bool

Returns:

True if the verbosity level was changed, False otherwise