pynanalogue Python API Reference
Note: This file is auto-generated.
Functions
peek
peek(bam_path, treat_as_url=False)
Peeks at a BAM file to extract contig information and detected modifications.
This function reads the BAM header and examines up to 100 records to determine the contigs present in the file and any DNA/RNA modifications detected.
Args
-
bam_path (str): Path to the BAM file.
-
treat_as_url (optional, bool): If True, treat
bam_pathas a URL, default False.
Returns
A dictionary with two keys:
contigs: A dictionary mapping contig names to their lengthsmodifications: A list of modifications, where each modification is a list of three strings: [base, strand,mod_code]. For example: [[“T”, “+”, “T”], [“G”, “-”, “7200”]]
Example
import pynanalogue
result = pynanalogue.peek("example.bam")
print(result["contigs"]) # {"chr1": 248956422, "chr2": 242193529, ...}
print(result["modifications"]) # [["C", "+", "m"], ["A", "+", "28871"]]
Errors
If the BAM file cannot be read or parsed
polars_bam_mods
polars_bam_mods(
bam_path,
treat_as_url=False,
min_seq_len=0,
min_align_len=0,
read_ids=Ellipsis,
threads=2,
include_zero_len=False,
read_filter='',
sample_fraction=1.0,
mapq_filter=0,
exclude_mapq_unavail=False,
region='',
full_region=False,
tag='',
mod_strand='',
min_mod_qual=0,
reject_mod_qual_non_inclusive=Ellipsis,
trim_read_ends_mod=0,
base_qual_filter_mod=0,
mod_region=''
)
Converts modification data in mod BAM files into a Polars Dataframe.
Columns are read_id, seq_len, alignment_type, align_start, align_end,
contig, contig_id, base, is_strand_plus, mod_code, position, ref_position,
and mod_quality.
Sets various options through builder functions before running
the function and capturing the output.
Parses records, collects them, and runs nanalogue_core::read_utils::curr_reads_to_dataframe.
Args
-
bam_path (str): Path to the BAM file. Must be associated with a BAM index.
-
treat_as_url (optional, bool): If True, treat
bam_pathas a URL, default False. -
min_seq_len (optional, int): Only retain sequences above this length, default 0.
-
min_align_len (optional, int): Only retain sequences with an alignment length above this value. Defaults to unused.
-
read_ids (optional, set of str): Only retrieve these read ids, defaults to unused.
-
threads (optional, int): Number of threads used in some aspects of program execution, defaults to 2.
-
include_zero_len (optional, bool): Include sequences of zero length. WARNING: our program may crash if you do this. Defaults to False. Helps to check if sequences of zero length exist in our BAM file.
-
read_filter (optional, str): Comma-separated sequence of one to many of the following strings: primary_forward, primary_reverse, secondary_forward, secondary_reverse, supplementary_forward, supplementary_reverse, unmapped. If specified, only reads with a mapping belonging to this set are retained. Defaults to no filter.
-
sample_fraction (optional, float): Set to between 0 and 1 to subsample BAM file. WARNING: seeds are not set, so you may get a new set of reads every time. WARNING: we sample every read with the given probability, so the total number of reads fluctuates according to standard counting statistics.
-
mapq_filter (optional, int): Exclude reads with mapping quality below this number. defaults to unused.
-
exclude_mapq_unavail (optional, bool): Exclude reads where mapping quality is unavailable. defaults to false.
-
region (optional, str): Only include reads with at least one mapped base from this region. Use the format “contig”, “contig:start-”, or “contig:start-end”. These are 0-based, half-open intervals. Defaults to read entire BAM file. Can be used in combination with
mod_region. -
full_region (optional, bool): Only include reads if they pass through the region above in full. Defaults to false.
-
tag (optional, str): If set, only this type of modification is processed. Input is a string type, for example a single-letter code “m”, a number as a string “76792” etc. Defaults to processing all modifications.
-
mod_strand (optional, str): Set this to
bcorbc_compto retrieve information about mods only from the basecalled strand or only from its complement. Some sequencing technologies likePacBioorONT duplexrecord mod information both on a strand and its complement. It may be useful in some scenarios to separate this information. Defaults to not filter. -
min_mod_qual (optional, int): Set to a number 0-255. Reject modification calls whose probability is below this value (0, 255 correspond to a probability of 0 and 1 respectively). Defaults to 0.
-
reject_mod_qual_non_inclusive (optional, (int, int)): Reject modification calls whose probability is such that int_low < prob < int_high. Set both to a number between 0-255 and such that the first entry is <= the second (if they are equal, no filtering is performed). Defaults to no filtering. Also see comments under
min_mod_qual. -
trim_read_ends_mod (optional, int): Reject modification information within so many bp of either end of the read. Defaults to 0.
-
base_qual_filter_mod (optional, int): Reject modification information on any base whose basecalling quality is below this number. Defaults to 0.
-
mod_region (optional, str): Genomic region in the format “contig”, “contig:start-” or “contig:start-end”. Reject any modification information outside this region. These are half-open, 0-based intervals. Can be used in combination with
region.
Returns
A Polars dataframe.
Example output
If the polars dataframe were converted to a plain-text table, it might look like the following:
read_id seq_len alignment_type align_start align_end contig contig_id base is_strand_plus mod_code position ref_position mod_quality
5d10eb9a-aae1-4db8-8ec6-7ebb34d32575 8 primary_forward 9 17 dummyI 0 T true T 0 9 4
5d10eb9a-aae1-4db8-8ec6-7ebb34d32575 8 primary_forward 9 17 dummyI 0 T true T 3 12 7
5d10eb9a-aae1-4db8-8ec6-7ebb34d32575 8 primary_forward 9 17 dummyI 0 T true T 4 13 9
5d10eb9a-aae1-4db8-8ec6-7ebb34d32575 8 primary_forward 9 17 dummyI 0 T true T 7 16 6
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 primary_forward 23 71 dummyIII 2 T true T 3 26 221
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 primary_forward 23 71 dummyIII 2 T true T 8 31 242
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 primary_forward 23 71 dummyIII 2 T true T 27 50 3
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 primary_forward 23 71 dummyIII 2 T true T 39 62 47
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 primary_forward 23 71 dummyIII 2 T true T 47 70 239
fffffff1-10d2-49cb-8ca3-e8d48979001b 33 primary_reverse 3 36 dummyII 1 T true T 12 15 3
fffffff1-10d2-49cb-8ca3-e8d48979001b 33 primary_reverse 3 36 dummyII 1 T true T 13 16 3
fffffff1-10d2-49cb-8ca3-e8d48979001b 33 primary_reverse 3 36 dummyII 1 T true T 16 19 4
fffffff1-10d2-49cb-8ca3-e8d48979001b 33 primary_reverse 3 36 dummyII 1 T true T 19 22 3
fffffff1-10d2-49cb-8ca3-e8d48979001b 33 primary_reverse 3 36 dummyII 1 T true T 20 23 182
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 28 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 29 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 30 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 32 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 43 -1 77
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped G false 7200 44 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped T true T 3 -1 221
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped T true T 8 -1 242
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped T true T 27 -1 0
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped T true T 39 -1 47
a4f36092-b4d5-47a9-813e-c22c3b477a0c 48 unmapped T true T 47 -1 239
Errors
If building of option-related structs fails, if BAM input
cannot be obtained, if preparing records fails, or running the
nanalogue_core commands fail
read_info
read_info(
bam_path,
treat_as_url=False,
min_seq_len=0,
min_align_len=0,
read_ids=Ellipsis,
threads=2,
include_zero_len=False,
read_filter='',
sample_fraction=1.0,
mapq_filter=0,
exclude_mapq_unavail=False,
region='',
full_region=False,
tag='',
mod_strand='',
min_mod_qual=0,
reject_mod_qual_non_inclusive=Ellipsis,
trim_read_ends_mod=0,
base_qual_filter_mod=0,
mod_region=''
)
Produces bytes which can be decoded to JSON format, with one JSON record per BAM record with information per record such as alignment length, sequence length, read id, mod counts etc.
Sets various options through builder functions before running
the function and capturing the output as a stream of bytes
which can be decoded to JSON (see the Example output section below).
Runs nanalogue_core::read_info::run.
This function can be used to count how many reads are above a threshold length or modification count or are primary mappings. The function can also be used to analyze relationships such as alignment length vs basecalled length. The arguments can be used to filter the BAM data (e.g. passing through a specific region etc.).
Args
-
bam_path (str): Path to the BAM file. Must be associated with a BAM index.
-
treat_as_url (optional, bool): If True, treat
bam_pathas a URL, default False. -
min_seq_len (optional, int): Only retain sequences above this length, default 0.
-
min_align_len (optional, int): Only retain sequences with an alignment length above this value. Defaults to unused.
-
read_ids (optional, set of str): Only retrieve these read ids, defaults to unused.
-
threads (optional, int): Number of threads used in some aspects of program execution, defaults to 2.
-
include_zero_len (optional, bool): Include sequences of zero length. WARNING: our program may crash if you do this. Defaults to False. Helps to check if sequences of zero length exist in our BAM file.
-
read_filter (optional, str): Comma-separated sequence of one to many of the following strings: primary_forward, primary_reverse, secondary_forward, secondary_reverse, supplementary_forward, supplementary_reverse, unmapped. If specified, only reads with a mapping belonging to this set are retained. Defaults to no filter.
-
sample_fraction (optional, float): Set to between 0 and 1 to subsample BAM file. WARNING: seeds are not set, so you may get a new set of reads every time. WARNING: we sample every read with the given probability, so the total number of reads fluctuates according to standard counting statistics.
-
mapq_filter (optional, int): Exclude reads with mapping quality below this number. defaults to unused.
-
exclude_mapq_unavail (optional, bool): Exclude reads where mapping quality is unavailable. defaults to false.
-
region (optional, str): Only include reads with at least one mapped base from this region. Use the format “contig”, “contig:start-”, or “contig:start-end”. These are 0-based, half-open intervals. Defaults to read entire BAM file. Can be used in combination with
mod_region. -
full_region (optional, bool): Only include reads if they pass through the region above in full. Defaults to false.
-
tag (optional, str): If set, only this type of modification is processed. Input is a string type, for example a single-letter code “m”, a number as a string “76792” etc. Defaults to processing all modifications.
-
mod_strand (optional, str): Set this to
bcorbc_compto retrieve information about mods only from the basecalled strand or only from its complement. Some sequencing technologies likePacBioorONT duplexrecord mod information both on a strand and its complement. It may be useful in some scenarios to separate this information. Defaults to not filter. -
min_mod_qual (optional, int): Set to a number 0-255. Reject modification calls whose probability is below this value (0, 255 correspond to a probability of 0 and 1 respectively). Defaults to 0.
-
reject_mod_qual_non_inclusive (optional, (int, int)): Reject modification calls whose probability is such that int_low < prob < int_high. Set both to a number between 0-255 and such that the first entry is <= the second (if they are equal, no filtering is performed). Defaults to no filtering. Also see comments under
min_mod_qual. -
trim_read_ends_mod (optional, int): Reject modification information within so many bp of either end of the read. Defaults to 0.
-
base_qual_filter_mod (optional, int): Reject modification information on any base whose basecalling quality is below this number. Defaults to 0.
-
mod_region (optional, str): Genomic region in the format “contig”, “contig:start-” or “contig:start-end”. Reject any modification information outside this region. These are half-open, 0-based intervals. Can be used in combination with
region.
Returns
A stream of bytes that can be decoded to JSON (See the snippet from Example output).
Example output
You’ve to decode the output of the function using something like:
import json
# Assume the function output is in 'result_bytes'
decoded_output = json.loads(result_bytes)
A record from the decoded output might look like
[
{
"read_id": "cd623d4a-510d-4c6c-9d88-10eb475ac59d",
"sequence_length": 2104,
"contig": "contig_0",
"reference_start": 7369,
"reference_end": 9473,
"alignment_length": 2104,
"alignment_type": "primary_reverse",
"mod_count": "C-m:263;N+N:2104;(probabilities >= 0.5020, PHRED base qual >= 0)"
}
]
When mods are not available, you will see NA in the mod_count field.
Errors
If building of option-related structs fails, if BAM input
cannot be obtained, if preparing records fails, or running the
nanalogue_core::read_info::run function fails
seq_table
seq_table(
bam_path,
region,
treat_as_url=False,
min_seq_len=0,
min_align_len=0,
read_ids=Ellipsis,
threads=2,
include_zero_len=False,
read_filter='',
sample_fraction=1.0,
mapq_filter=0,
exclude_mapq_unavail=False,
tag='',
mod_strand='',
min_mod_qual=0,
reject_mod_qual_non_inclusive=Ellipsis,
trim_read_ends_mod=0,
base_qual_filter_mod=0
)
Returns a Polars DataFrame with read sequence and quality information for a genomic region.
Represents insertions as lowercase, deletions as periods, and modifications as ‘Z’.
Basecalling qualities represented as a string of numbers separated by periods, with value set
to ‘255’ at a deleted base.
Calls nanalogue_core::reads_table::run_df to extract read sequences and qualities
from a specified genomic region. The output DataFrame contains columns for read ID,
sequence (with modifications shown as ‘Z’), and base qualities.
Args
-
bam_path (str): Path to the BAM file. Must be associated with a BAM index.
-
region (str): Genomic region from which to extract sequences. Required. Use the format “contig”, “contig:start-”, or “contig:start-end”. These are 0-based, half-open intervals.
-
treat_as_url (optional, bool): If True, treat
bam_pathas a URL, default False. -
min_seq_len (optional, int): Only retain sequences above this length, default 0.
-
min_align_len (optional, int): Only retain sequences with an alignment length above this value. Defaults to unused.
-
read_ids (optional, set of str): Only retrieve these read ids, defaults to unused.
-
threads (optional, int): Number of threads used in some aspects of program execution, defaults to 2.
-
include_zero_len (optional, bool): Include sequences of zero length. WARNING: our program may crash if you do this. Defaults to False.
-
read_filter (optional, str): Comma-separated sequence of one to many of the following strings: primary_forward, primary_reverse, secondary_forward, secondary_reverse, supplementary_forward, supplementary_reverse, unmapped. If specified, only reads with a mapping belonging to this set are retained. Defaults to no filter.
-
sample_fraction (optional, float): Set to between 0 and 1 to subsample BAM file. WARNING: seeds are not set, so you may get a new set of reads every time.
-
mapq_filter (optional, int): Exclude reads with mapping quality below this number. defaults to unused.
-
exclude_mapq_unavail (optional, bool): Exclude reads where mapping quality is unavailable. defaults to false.
-
tag (optional, str): If set, only this type of modification is processed. Input is a string type, for example a single-letter code “m”, a number as a string “76792” etc. Defaults to processing all modifications.
-
mod_strand (optional, str): Set this to
bcorbc_compto retrieve information about mods only from the basecalled strand or only from its complement. Defaults to not filter. -
min_mod_qual (optional, int): Set to a number 0-255. Reject modification calls whose probability is below this value. Defaults to 0.
-
reject_mod_qual_non_inclusive (optional, (int, int)): Reject modification calls whose probability is such that int_low < prob < int_high. Defaults to no filtering.
-
trim_read_ends_mod (optional, int): Reject modification information within so many bp of either end of the read. Defaults to 0.
-
base_qual_filter_mod (optional, int): Reject modification information on any base whose basecalling quality is below this number. Defaults to 0.
Returns
A Polars dataframe with columns: read_id, sequence, qualities.
Sequence column conventions:
- Insertions are shown in lowercase
- Deletions are shown as periods (
.) - Modified bases on the reference are shown as
Z - Modified bases in an insertion are shown as
z
Qualities column: a period-separated string of integer quality scores, with one score per base in the sequence.
Example output
If the polars dataframe were converted to a plain-text table, it might look like:
read_id sequence qualities
5d10eb9a-aae1-4db8-8ec6-7ebb34d32575 ACGTZac..TZgt 30.35.40.38.42.20.22.255.255.41.25.28.30
a4f36092-b4d5-47a9-813e-c22c3b477a0c TGCAZz.ATGCA 28.33.39.41.40.18.255.36.42.38.35.30
Errors
If the region parameter is empty, if building of option-related structs fails,
if BAM input cannot be obtained, if preparing records fails, or if running the
nanalogue_core::reads_table::run_df function fails
simulate_mod_bam
simulate_mod_bam(json_config, bam_path, fasta_path)
Simulates a BAM file with or without modifications based on a JSON configuration. Creates both a BAM file and a corresponding FASTA reference file.
This function takes a JSON configuration string that specifies how to generate synthetic sequencing data, including contig specifications, read parameters, and optional modification patterns. The output includes a sorted, indexed BAM file and a FASTA reference file.
Args
-
json_config (str): JSON string containing the simulation configuration. Must conform to the
SimulationConfigschema. The JSON should specify -
(some inputs may be optional): -
contigs: Configuration for generating reference sequences -number: Number of contigs to generate -len_range: Range of contig lengths [min, max] -repeated_seq: Sequence pattern to repeat for contig generation. If not specified, random DNA sequences are generated. -reads: Array of read group specifications, each containing: -number: Number of reads to generate -mapq_range: Range of mapping quality values [min, max] -base_qual_range: Range of base quality values [min, max] -len_range: Range of read lengths as fraction of contig [min, max] -barcode: Optional barcode sequence -mods: Optional array of modification specifications -
bam_path (str): Output path for the BAM file. If the file exists, it will be overwritten. A corresponding
.baiindex file will be created automatically. -
fasta_path (str): Output path for the FASTA reference file. If the file exists, it will be overwritten.
Returns
Returns None on success. The function creates two files on disk:
- A sorted, indexed BAM file at
bam_path(with accompanying.baiindex) - A FASTA reference file at
fasta_path
Example
import pynanalogue
json_config = '''
{
"contigs": {
"number": 2,
"len_range": [100, 200],
"repeated_seq": "ACGTACGT"
},
"reads": [
{
"number": 10,
"mapq_range": [10, 30],
"base_qual_range": [20, 40],
"len_range": [0.1, 0.9],
"barcode": "ACGTAA",
"mods": [{
"base": "C",
"is_strand_plus": true,
"mod_code": "m",
"win": [5, 3],
"mod_range": [[0.3, 0.7], [0.1, 0.5]]
}]
}
]
}
'''
pynanalogue.simulate_mod_bam(
json_config=json_config,
bam_path="output.bam",
fasta_path="output.fasta"
)
Errors
Returns a Python exception if:
- The JSON configuration is invalid or cannot be parsed
- The JSON structure doesn’t match the expected
SimulationConfigschema - File I/O operations fail (e.g., permission issues, disk full)
- BAM or FASTA generation fails due to invalid configuration parameters
window_reads
window_reads(
bam_path,
win,
step,
win_op='density',
treat_as_url=False,
min_seq_len=0,
min_align_len=0,
read_ids=Ellipsis,
threads=2,
include_zero_len=False,
read_filter='',
sample_fraction=1.0,
mapq_filter=0,
exclude_mapq_unavail=False,
region='',
full_region=False,
tag='',
mod_strand='',
min_mod_qual=0,
reject_mod_qual_non_inclusive=Ellipsis,
trim_read_ends_mod=0,
base_qual_filter_mod=0,
mod_region=''
)
Window modification density across single molecules and return a Polars DataFrame.
With the gradient option, the gradient in mod density within each window is reported.
The output is a BED format, with windowed densities per window per read.
Details: runs nanalogue_core::window_reads::run_df.
Sets various options through builder functions before running
the function and capturing the output (see Example Output).
Args
-
bam_path (str): Path to the BAM file. Must be associated with a BAM index.
-
win (int): Size of window in number of bases whose mod is being queried. i.e. let’s say a read contains cytosine mods and win is set to 10, then each window is chosen so that there are 10 cytosines in it. If a read has multiple mods, then multiple windows are set up such that each window has the specified number of bases of that type in it.
-
step (int): Length by which the window is slid in the same units as win above.
-
win_op (optional, str): Type of windowing operation to use, allows “density” and “grad_density” i.e. measure modification density or the gradient of it within each window. Default is “density”.
-
treat_as_url (optional, bool): If True, treat
bam_pathas a URL, default False. -
min_seq_len (optional, int): Only retain sequences above this length, default 0.
-
min_align_len (optional, int): Only retain sequences with an alignment length above this value. Defaults to unused.
-
read_ids (optional, set of str): Only retrieve these read ids, defaults to unused.
-
threads (optional, int): Number of threads used in some aspects of program execution, defaults to 2.
-
include_zero_len (optional, bool): Include sequences of zero length. WARNING: our program may crash if you do this. Defaults to False. Helps to check if sequences of zero length exist in our BAM file.
-
read_filter (optional, str): Comma-separated sequence of one to many of the following strings: primary_forward, primary_reverse, secondary_forward, secondary_reverse, supplementary_forward, supplementary_reverse, unmapped. If specified, only reads with a mapping belonging to this set are retained. Defaults to no filter.
-
sample_fraction (optional, float): Set to between 0 and 1 to subsample BAM file. WARNING: seeds are not set, so you may get a new set of reads every time. WARNING: we sample every read with the given probability, so the total number of reads fluctuates according to standard counting statistics.
-
mapq_filter (optional, int): Exclude reads with mapping quality below this number. defaults to unused.
-
exclude_mapq_unavail (optional, bool): Exclude reads where mapping quality is unavailable. defaults to false.
-
region (optional, str): Only include reads with at least one mapped base from this region. Use the format “contig”, “contig:start-”, or “contig:start-end”. These are 0-based, half-open intervals. Defaults to read entire BAM file. Can be used in combination with
mod_region. -
full_region (optional, bool): Only include reads if they pass through the region above in full. Defaults to false.
-
tag (optional, str): If set, only this type of modification is processed. Input is a string type, for example a single-letter code “m”, a number as a string “76792” etc. Defaults to processing all modifications.
-
mod_strand (optional, str): Set this to
bcorbc_compto retrieve information about mods only from the basecalled strand or only from its complement. Some sequencing technologies likePacBioorONT duplexrecord mod information both on a strand and its complement. It may be useful in some scenarios to separate this information. Defaults to not filter. -
min_mod_qual (optional, int): Set to a number 0-255. Reject modification calls whose probability is below this value (0, 255 correspond to a probability of 0 and 1 respectively). Defaults to 0.
-
reject_mod_qual_non_inclusive (optional, (int, int)): Reject modification calls whose probability is such that int_low < prob < int_high. Set both to a number between 0-255 and such that the first entry is <= the second (if they are equal, no filtering is performed). Defaults to no filtering. Also see comments under
min_mod_qual. -
trim_read_ends_mod (optional, int): Reject modification information within so many bp of either end of the read. Defaults to 0.
-
base_qual_filter_mod (optional, int): Reject modification information on any base whose basecalling quality is below this number. Defaults to 0.
-
mod_region (optional, str): Genomic region in the format “contig”, “contig:start-” or “contig:start-end”. Reject any modification information outside this region. These are half-open, 0-based intervals. Can be used in combination with
region.
Returns
A Polars dataframe (See the snippet from Example output).
Example output
If the dataframe were converted to a TSV format, it might look like the following. The basecalling qualities are all 255 i.e. unknown because this is a BAM file where basecalling qualities haven’t been recorded.
#contig ref_win_start ref_win_end read_id win_val strand base mod_strand mod_type win_start win_end basecall_qual
dummyIII 26 32 a4f36092-b4d5-47a9-813e-c22c3b477a0c 1 + T + T 3 9 255
dummyIII 31 51 a4f36092-b4d5-47a9-813e-c22c3b477a0c 0.5 + T + T 8 28 255
dummyIII 62 71 a4f36092-b4d5-47a9-813e-c22c3b477a0c 0.5 + T + T 39 48 255
dummyII 22 24 fffffff1-10d2-49cb-8ca3-e8d48979001b 0.5 - T + T 19 21 255
. -1 -1 a4f36092-b4d5-47a9-813e-c22c3b477a0c 1 . T + T 3 9 255
. -1 -1 a4f36092-b4d5-47a9-813e-c22c3b477a0c 0.5 . T + T 8 28 255
. -1 -1 a4f36092-b4d5-47a9-813e-c22c3b477a0c 0.5 . T + T 39 48 255
Under the gradient option, win_val reports the gradient in modification density within each window.
#contig ref_win_start ref_win_end read_id win_val strand base mod_strand mod_type win_start win_end basecall_qual
dummyIII 40 50 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.054545455 + N + N 17 27 255
dummyIII 41 51 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.096969694 + N + N 18 28 255
dummyIII 42 52 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.12727273 + N + N 19 29 255
dummyIII 43 53 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.14545454 + N + N 20 30 255
dummyIII 44 54 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.15151516 + N + N 21 31 255
dummyIII 45 55 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.14545454 + N + N 22 32 255
dummyIII 46 56 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.12727273 + N + N 23 33 255
dummyIII 47 57 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.096969694 + N + N 24 34 255
dummyIII 48 58 c4f36092-b4d5-47a9-813e-c22c3b477a0c 0.054545455 + N + N 25 35 255
Errors
If building of option-related structs fails, if BAM input
cannot be obtained, if preparing records fails, or running the
nanalogue_core::window_reads::run_df function fails