"""
Contact Binners
~~~~~~~~~~~~~~~
Binners are iterators that convert input data of various flavors into a
properly sorted, chunked stream of binned contacts.
"""
from __future__ import annotations
import itertools
import warnings
from bisect import bisect_left
from collections import Counter, OrderedDict
from collections.abc import Iterator
from functools import partial
from typing import Any, Callable
import h5py
import numpy as np
import pandas as pd
from pandas.api.types import is_integer_dtype
from .._logging import get_logger
from .._typing import MapFunctor
from ..util import (
GenomeSegmentation,
balanced_partition,
check_bins,
get_chromsizes,
partition,
rlencode,
)
logger = get_logger("cooler.create")
class BadInputError(ValueError):
pass
SANITIZE_PRESETS = {
"bg2": {
"decode_chroms": True,
"is_one_based": False,
"tril_action": "reflect",
"chrom_field": "chrom",
"anchor_field": "start",
"sided_fields": ("chrom", "start", "end"),
"suffixes": ("1", "2"),
"sort": True,
"validate": True,
},
"pairs": {
"decode_chroms": True,
"is_one_based": False,
"tril_action": "reflect",
"chrom_field": "chrom",
"anchor_field": "pos",
"sided_fields": ("chrom", "pos"),
"suffixes": ("1", "2"),
"sort": False,
"validate": True,
},
}
def _sanitize_records(
chunk,
gs,
decode_chroms,
is_one_based,
tril_action,
chrom_field,
anchor_field,
sided_fields,
suffixes,
sort,
validate,
):
# Get integer contig IDs
if decode_chroms:
# Unspecified chroms get assigned category = NaN and integer code = -1
chrom1_ids = np.array(
pd.Categorical(chunk["chrom1"], gs.contigs, ordered=True).codes
)
chrom2_ids = np.array(
pd.Categorical(chunk["chrom2"], gs.contigs, ordered=True).codes
)
else:
chrom1_ids = chunk["chrom1"].values
chrom2_ids = chunk["chrom2"].values
if validate:
for col, dt in [("chrom1", chrom1_ids.dtype), ("chrom2", chrom2_ids.dtype)]:
if not is_integer_dtype(dt):
raise BadInputError(
f"`{col}` column is non-integer. "
+ "If string, use `decode_chroms=True` to convert to enum"
)
# Drop records from non-requested chromosomes
to_drop = (chrom1_ids < 0) | (chrom2_ids < 0)
if np.any(to_drop):
mask = ~to_drop
chrom1_ids = chrom1_ids[mask]
chrom2_ids = chrom2_ids[mask]
chunk = chunk[mask].copy()
# Handle empty case
if not len(chunk):
chunk["bin1_id"] = []
chunk["bin2_id"] = []
return chunk
# Find positional anchor columns, convert to zero-based if needed
anchor1 = chunk[anchor_field + suffixes[0]].to_numpy(copy=True)
anchor2 = chunk[anchor_field + suffixes[1]].to_numpy(copy=True)
if is_one_based:
anchor1 -= 1
anchor2 -= 1
# Check types and bounds
if validate:
for dt in [anchor1.dtype, anchor2.dtype]:
if not is_integer_dtype(dt):
raise BadInputError("Found a non-integer anchor column")
is_neg = (anchor1 < 0) | (anchor2 < 0)
if np.any(is_neg):
err = chunk[is_neg]
raise BadInputError(
"Found an anchor position with negative value. Make sure your "
"coordinates are 1-based or use the --zero-based option "
"when loading. \n{}".format(err.head().to_csv(sep="\t"))
)
chromsizes1 = gs.chromsizes.values[chrom1_ids]
chromsizes2 = gs.chromsizes.values[chrom2_ids]
is_excess = (anchor1 > chromsizes1) | (anchor2 > chromsizes2)
if np.any(is_excess):
err = chunk[is_excess]
raise BadInputError(
"Found an anchor position exceeding chromosome length:\n{}".format(
err.head().to_csv(sep="\t")
)
)
# Handle lower triangle records
# Note: Swap assignment works as desired because boolean masks create copies
# See https://github.com/open2c/cooler/pull/229
if tril_action is not None:
is_tril = (chrom1_ids > chrom2_ids) | (
(chrom1_ids == chrom2_ids) & (anchor1 > anchor2)
)
if np.any(is_tril):
if tril_action == "reflect":
(chrom1_ids[is_tril], chrom2_ids[is_tril]) = (
chrom2_ids[is_tril],
chrom1_ids[is_tril],
)
anchor1[is_tril], anchor2[is_tril] = anchor2[is_tril], anchor1[is_tril]
for field in sided_fields:
(
chunk.loc[is_tril, field + suffixes[0]],
chunk.loc[is_tril, field + suffixes[1]],
) = (
chunk.loc[is_tril, field + suffixes[1]],
chunk.loc[is_tril, field + suffixes[0]],
)
elif tril_action == "drop":
mask = ~is_tril
chrom1_ids = chrom1_ids[mask]
chrom2_ids = chrom2_ids[mask]
anchor1 = anchor1[mask]
anchor2 = anchor2[mask]
chunk = chunk[mask].copy()
elif tril_action == "raise":
err = chunk[is_tril]
raise BadInputError(
"Found lower triangle pairs:\n{}".format(
err.head().to_csv(sep="\t")
)
)
else:
raise ValueError(f"Unknown tril_action value: '{tril_action}'")
# Assign bin IDs from bin table
chrom_binoffset = gs.chrom_binoffset
binsize = gs.binsize
if binsize is None:
chrom_abspos = gs.chrom_abspos
start_abspos = gs.start_abspos
bin1_ids = []
bin2_ids = []
for cid1, pos1, cid2, pos2 in zip(chrom1_ids, anchor1, chrom2_ids, anchor2):
lo = chrom_binoffset[cid1]
hi = chrom_binoffset[cid1 + 1]
bin1_ids.append(
lo
+ np.searchsorted(
start_abspos[lo:hi], chrom_abspos[cid1] + pos1, side="right"
)
- 1
)
lo = chrom_binoffset[cid2]
hi = chrom_binoffset[cid2 + 1]
bin2_ids.append(
lo
+ np.searchsorted(
start_abspos[lo:hi], chrom_abspos[cid2] + pos2, side="right"
)
- 1
)
chunk["bin1_id"] = bin1_ids
chunk["bin2_id"] = bin2_ids
else:
chunk["bin1_id"] = chrom_binoffset[chrom1_ids] + anchor1 // binsize
chunk["bin2_id"] = chrom_binoffset[chrom2_ids] + anchor2 // binsize
# Sort by bin IDs
if sort:
chunk = chunk.sort_values(["bin1_id", "bin2_id"])
# TODO: check for duplicate records and warn
return chunk
[docs]
def sanitize_records(
bins: pd.DataFrame,
schema: str | None = None,
**kwargs
) -> Callable[[pd.DataFrame], pd.DataFrame]:
"""
Builds a funtion to sanitize and assign bin IDs to a data frame of
paired genomic positions based on a provided genomic bin segmentation.
Parameters
----------
bins : DataFrame
Bin table to compare records against.
schema : str, optional
Use pre-defined parameters for a particular format. Any options can be
overriden via kwargs. If not provided, values for all the options below
must be given.
decode_chroms : bool
Convert string chromosome names to integer IDs based on the order given
in the bin table. Set to False if the chromosomes are already given as
an enumeration, starting at 0. Records with either chrom ID < 0 are
dropped.
is_one_based : bool
Whether the input anchor coordinates are one-based, rather than
zero-based. They will be converted to zero-based.
tril_action : 'reflect', 'drop', 'raise' or None
How to handle lower triangle ("tril") records.
If set to 'reflect', tril records will be flipped or "reflected"
to their mirror image: "sided" column pairs will have their values
swapped.
If set to 'drop', tril records will be discarded. This is useful if
your input data is symmetric, i.e. contains mirror duplicates of every
record.
If set to 'raise', an exception will be raised if any tril record is
encountered.
chrom_field : str
Base name of the two chromosome/scaffold/contig columns.
anchor_field : str
Base name of the positional anchor columns.
sided_fields : sequence of str
Base names of column pairs to swap values between when
mirror-reflecting records.
suffixes : pair of str
Suffixes used to identify pairs of sided columns. e.g.: ('1', '2'),
('_x', '_y'), etc.
sort : bool
Whether to sort the output dataframe by bin_id and bin2_id.
validate : bool
Whether to do type- and bounds-checking on the anchor position
columns. Raises BadInputError.
Returns
-------
callable :
Function of one argument that takes a raw dataframe and returns a
sanitized dataframe with bin IDs assigned.
"""
if schema is not None:
try:
options = SANITIZE_PRESETS[schema]
except KeyError:
raise ValueError(f"Unknown schema: '{schema}'") from None
else:
options = {}
options = options.copy()
options.update(**kwargs)
chromsizes = get_chromsizes(bins)
options["gs"] = GenomeSegmentation(chromsizes, bins)
return partial(_sanitize_records, **options)
def _sanitize_pixels(
chunk,
gs,
is_one_based=False,
tril_action="reflect",
bin1_field="bin1_id",
bin2_field="bin2_id",
sided_fields=(),
suffixes=("1", "2"),
sort=True,
):
if is_one_based:
chunk[bin1_field] -= 1
chunk[bin2_field] -= 1
# Note: Swap assignment syntax works as desired because boolean mask
# selection produces array copies rather than views.
# See https://github.com/open2c/cooler/pull/229.
if tril_action is not None:
is_tril = chunk[bin1_field] > chunk[bin2_field] # boolean mask
if np.any(is_tril):
if tril_action == "reflect":
(
chunk.loc[is_tril, bin1_field],
chunk.loc[is_tril, bin2_field],
) = (
chunk.loc[is_tril, bin2_field],
chunk.loc[is_tril, bin1_field],
)
for field in sided_fields:
(
chunk.loc[is_tril, field + suffixes[0]],
chunk.loc[is_tril, field + suffixes[1]],
) = (
chunk.loc[is_tril, field + suffixes[1]],
chunk.loc[is_tril, field + suffixes[0]],
)
elif tril_action == "drop":
chunk = chunk[~is_tril]
elif tril_action == "raise":
raise BadInputError("Found bin1_id greater than bin2_id")
else:
raise ValueError(f"Unknown tril_action value: '{tril_action}'")
return chunk.sort_values([bin1_field, bin2_field]) if sort else chunk
[docs]
def sanitize_pixels(
bins: pd.DataFrame,
**kwargs,
) -> Callable[[pd.DataFrame], pd.DataFrame]:
"""
Builds a function to sanitize an already-binned genomic data with
genomic bin assignments.
Parameters
----------
bins : DataFrame
Bin table to compare pixel records against.
is_one_based : bool, optional
Whether the input bin IDs are one-based, rather than zero-based.
They will be converted to zero-based.
tril_action : 'reflect', 'drop', 'raise' or None
How to handle lower triangle ("tril") pixels.
If set to 'reflect' [default], tril pixels will be flipped or
"reflected" to their mirror image: "sided" column pairs will have their
values swapped.
If set to 'drop', tril pixels will be discarded. This is useful if
your input data is duplexed, i.e. contains mirror duplicates of every
record.
If set to 'raise', an exception will be raised if any tril record is
encountered.
bin1_field : str
Name of the column representing ith (row) axis of the matrix.
Default is 'bin1_id'.
bin2_field : str
Name of the column representing jth (col) axis of the matrix.
Default is 'bin2_id'.
sided_fields : sequence of str
Base names of column pairs to swap values between when mirror-reflecting
pixels.
suffixes : pair of str
Suffixes used to identify pairs of sided columns. e.g.: ('1', '2'),
('_x', '_y'), etc.
sort : bool
Whether to sort the output dataframe by bin_id and bin2_id.
Returns
-------
callable :
Function of one argument that takes a raw dataframe and returns a
sanitized dataframe.
"""
chromsizes = get_chromsizes(bins)
kwargs["gs"] = GenomeSegmentation(chromsizes, bins)
return partial(_sanitize_pixels, **kwargs)
def _validate_pixels(chunk, n_bins, boundscheck, triucheck, dupcheck, ensure_sorted):
if boundscheck:
is_neg = (chunk["bin1_id"] < 0) | (chunk["bin2_id"] < 0)
if np.any(is_neg):
raise BadInputError("Found bin ID < 0")
is_excess = (chunk["bin1_id"] >= n_bins) | (chunk["bin2_id"] >= n_bins)
if np.any(is_excess):
raise BadInputError(
"Found a bin ID that exceeds the declared number of bins. "
"Check whether your bin table is correct."
)
if triucheck:
is_tril = chunk["bin1_id"] > chunk["bin2_id"]
if np.any(is_tril):
raise BadInputError("Found bin1_id greater than bin2_id")
if not isinstance(chunk, pd.DataFrame):
chunk = pd.DataFrame(chunk)
if dupcheck:
is_dup = chunk.duplicated(["bin1_id", "bin2_id"])
if is_dup.any():
err = chunk[is_dup]
raise BadInputError(
"Found duplicate pixels:\n{}".format(err.head().to_csv(sep="\t"))
)
if ensure_sorted:
chunk = chunk.sort_values(["bin1_id", "bin2_id"])
return chunk
def validate_pixels(
n_bins: pd.DataFrame,
boundscheck: bool,
triucheck: bool,
dupcheck: bool,
ensure_sorted: bool
) -> Callable[[pd.DataFrame], pd.DataFrame]:
return partial(
_validate_pixels,
n_bins=n_bins,
boundscheck=boundscheck,
triucheck=triucheck,
dupcheck=dupcheck,
ensure_sorted=ensure_sorted,
)
def aggregate_records(
sort: bool = True,
count: bool = True,
agg: dict[str, Any] | None = None,
rename: dict[str, str] | None = None
) -> Callable[[pd.DataFrame], pd.DataFrame]:
"""
Generates a function that aggregates bin-assigned records by pixel.
Parameters
----------
sort : bool, optional
Sort group keys. Get better performance by turning this off.
Note that this does not influence the order of observations within each
group.
count : bool, optional
Output the number of records per pixel. Default is True.
agg : dict, optional
Dict of column names -> functions or names.
rename : dict, optional
Dict to rename columns after aggregating.
Returns
-------
Function that takes a dataframe of records with bin IDs assigned, groups
them by pixel, counts them, and optionally aggregates other value columns.
Notes
-----
The GroupBy 'count' method ignores NaNs within groups, as opposed to 'size'.
"""
if agg is None:
agg = {}
if rename is None:
rename = {}
# We use one of the grouper columns to count the number of pairs per pixel.
# We always do count, even if 'count' isn't requested as output.
if count and "count" not in agg:
agg["bin1_id"] = "size"
rename["bin1_id"] = "count"
def _aggregate_records(chunk):
return (
chunk.groupby(["bin1_id", "bin2_id"], sort=sort)
.aggregate(agg)
.rename(columns=rename)
.reset_index()
)
return _aggregate_records
class ContactBinner:
"""
Base class for iterable contact binners.
"""
def __getstate__(self) -> dict:
d = self.__dict__.copy()
d.pop("_map", None)
return d
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
"""Iterator over chunks of binned contacts (i.e., nonzero pixels)
Chunks are expected to have the following format:
* dict of 1D arrays
* keys `bin1_id`, `bin2_id`, `count`
* arrays lexically sorted by `bin_id` then `bin2_id`
"""
raise NotImplementedError
class HDF5Aggregator(ContactBinner):
"""
Aggregate contacts from a hiclib-style HDF5 contacts file.
"""
def __init__(
self,
h5pairs: h5py.Group,
chromsizes: pd.Series,
bins: pd.DataFrame,
chunksize: int,
**kwargs
):
self.h5 = h5pairs
self.C1 = kwargs.pop("C1", "chrms1")
self.P1 = kwargs.pop("P1", "cuts1")
self.C2 = kwargs.pop("C2", "chrms2")
self.P2 = kwargs.pop("P2", "cuts2")
self.gs = GenomeSegmentation(chromsizes, bins)
self.chunksize = chunksize
self.partition = self._index_chroms()
def _index_chroms(self) -> dict[str, tuple[int, int]]:
# index extents of chromosomes on first axis of contact list
starts, lengths, values = rlencode(self.h5[self.C1], self.chunksize)
if len(set(values)) != len(values):
raise ValueError("Read pair coordinates are not sorted on the first axis")
return dict(zip(values, zip(starts, starts + lengths)))
def _load_chunk(self, lo, hi) -> dict[str, np.ndarray]:
data = OrderedDict(
[
("chrom_id1", self.h5[self.C1][lo:hi]),
("cut1", self.h5[self.P1][lo:hi]),
("chrom_id2", self.h5[self.C2][lo:hi]),
("cut2", self.h5[self.P2][lo:hi]),
]
)
return pd.DataFrame(data)
def aggregate(self, chrom: str) -> pd.DataFrame: # pragma: no cover
h5pairs = self.h5
C1, P1, C2, P2 = self.C1, self.P1, self.C2, self.P2
chunksize = self.chunksize
bins = self.gs.bins
binsize = self.gs.binsize
chrom_binoffset = self.gs.chrom_binoffset
chrom_abspos = self.gs.chrom_abspos
start_abspos = self.gs.start_abspos
cid = self.gs.idmap[chrom]
chrom_lo, chrom_hi = self.partition.get(cid, (-1, -1))
lo = chrom_lo
hi = lo
while hi < chrom_hi:
# update `hi` to make sure our selection doesn't split a bin1
lo, hi = hi, min(hi + chunksize, chrom_hi)
abspos = chrom_abspos[cid] + h5pairs[P1][hi - 1]
bin_id = int(np.searchsorted(start_abspos, abspos, side="right")) - 1
bin_end = bins["end"][bin_id]
hi = bisect_left(h5pairs[P1], bin_end, lo, chrom_hi)
if lo == hi:
hi = chrom_hi
logger.info(f"{lo} {hi}")
# load chunk and assign bin IDs to each read side
table = self._load_chunk(lo, hi)
abspos1 = chrom_abspos[h5pairs[C1][lo:hi]] + h5pairs[P1][lo:hi]
abspos2 = chrom_abspos[h5pairs[C2][lo:hi]] + h5pairs[P2][lo:hi]
if np.any(abspos1 > abspos2):
raise ValueError(
"Found a read pair that maps to the lower triangle of the "
"contact map (side1 > side2). Check that the provided "
"chromosome ordering and read pair file are consistent "
"such that all pairs map to the upper triangle with "
"respect to the given chromosome ordering."
)
if binsize is None:
table["bin1_id"] = (
np.searchsorted(start_abspos, abspos1, side="right") - 1
)
table["bin2_id"] = (
np.searchsorted(start_abspos, abspos2, side="right") - 1
)
else:
rel_bin1 = np.floor(table["cut1"] / binsize).astype(np.int64)
rel_bin2 = np.floor(table["cut2"] / binsize).astype(np.int64)
table["bin1_id"] = chrom_binoffset[table["chrom_id1"].values] + rel_bin1
table["bin2_id"] = chrom_binoffset[table["chrom_id2"].values] + rel_bin2
# reduce
gby = table.groupby(["bin1_id", "bin2_id"])
agg = (
gby["chrom_id1"]
.count()
.reset_index()
.rename(columns={"chrom_id1": "count"})
)
yield agg
def size(self) -> int:
return len(self.h5["chrms1"])
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
for chrom in self.gs.contigs:
for df in self.aggregate(chrom):
yield {k: v.values for k, v in df.items()}
class TabixAggregator(ContactBinner):
"""
Aggregate contacts from a sorted, BGZIP-compressed and tabix-indexed
tab-delimited text file.
"""
def __init__(
self,
filepath: str,
chromsizes: pd.Series,
bins: pd.DataFrame,
map: MapFunctor = map,
n_chunks: int = 1,
is_one_based: bool = False,
**kwargs,
):
try:
import pysam
except ImportError:
raise ImportError("pysam is required to read tabix files") from None
import pickle
import dill
dill.settings["protocol"] = pickle.HIGHEST_PROTOCOL
self._map = map
self.n_chunks = n_chunks
self.is_one_based = bool(is_one_based)
self.C2 = kwargs.pop("C2", 3)
self.P2 = kwargs.pop("P2", 4)
# all requested contigs will be placed in the output matrix
self.gs = GenomeSegmentation(chromsizes, bins)
# find available contigs in the contact list
self.filepath = filepath
self.n_records = None
with pysam.TabixFile(filepath, "r", encoding="ascii") as f:
try:
self.file_contigs = [c.decode("ascii") for c in f.contigs]
except AttributeError:
self.file_contigs = f.contigs
if not len(self.file_contigs):
raise RuntimeError("No reference sequences found.")
# warn about requested contigs not seen in the contact list
for chrom in self.gs.contigs:
if chrom not in self.file_contigs:
warnings.warn(
"Did not find contig " + f" '{chrom}' in contact list file.",
stacklevel=2,
)
warnings.warn(
"NOTE: When using the Tabix aggregator, make sure the order of "
"chromosomes in the provided chromsizes agrees with the chromosome "
"ordering of read ends in the contact list file.",
stacklevel=2,
)
def aggregate(
self, grange: tuple[str, int, int]
) -> pd.DataFrame | None: # pragma: no cover
chrom1, start, end = grange
import pysam
filepath = self.filepath
binsize = self.gs.binsize
idmap = self.gs.idmap
# chromsizes = self.gs.chromsizes
chrom_binoffset = self.gs.chrom_binoffset
chrom_abspos = self.gs.chrom_abspos
start_abspos = self.gs.start_abspos
decr = int(self.is_one_based)
C2 = self.C2
P2 = self.P2
logger.info(f"Binning {chrom1}:{start}-{end}|*")
these_bins = self.gs.fetch((chrom1, start, end))
rows = []
with pysam.TabixFile(filepath, "r", encoding="ascii") as f:
parser = pysam.asTuple()
accumulator = Counter()
for bin1_id, bin1 in these_bins.iterrows():
for line in f.fetch(chrom1, bin1.start, bin1.end, parser=parser):
chrom2 = line[C2]
pos2 = int(line[P2]) - decr
try:
cid2 = idmap[chrom2]
except KeyError:
# this chrom2 is not requested
continue
if binsize is None:
lo = chrom_binoffset[cid2]
hi = chrom_binoffset[cid2 + 1]
bin2_id = (
lo
+ np.searchsorted(
start_abspos[lo:hi],
chrom_abspos[cid2] + pos2,
side="right",
)
- 1
)
else:
bin2_id = chrom_binoffset[cid2] + (pos2 // binsize)
accumulator[bin2_id] += 1
if not accumulator:
continue
rows.append(
pd.DataFrame(
{
"bin1_id": bin1_id,
"bin2_id": list(accumulator.keys()),
"count": list(accumulator.values()),
},
columns=["bin1_id", "bin2_id", "count"],
).sort_values("bin2_id")
)
accumulator.clear()
logger.info(f"Finished {chrom1}:{start}-{end}|*")
return pd.concat(rows, axis=0) if len(rows) else None
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
granges = balanced_partition(self.gs, self.n_chunks, self.file_contigs)
for df in self._map(self.aggregate, granges):
if df is not None:
yield {k: v.values for k, v in df.items()}
class PairixAggregator(ContactBinner):
"""
Aggregate contacts from a sorted, BGZIP-compressed and pairix-indexed
tab-delimited text file.
"""
def __init__(
self,
filepath: str,
chromsizes: pd.Series,
bins: pd.DataFrame,
map: MapFunctor = map,
n_chunks: int = 1,
is_one_based: bool = False,
block_char: str = "|",
**kwargs,
):
try:
import pypairix
except ImportError:
raise ImportError(
"pypairix is required to read pairix-indexed files"
) from None
import pickle
import dill
dill.settings["protocol"] = pickle.HIGHEST_PROTOCOL
self._map = map
self.n_chunks = n_chunks
self.is_one_based = bool(is_one_based)
self.block_char = block_char
f = pypairix.open(filepath, "r")
self.C1 = f.get_chr1_col()
self.C2 = f.get_chr2_col()
self.P1 = f.get_startpos1_col()
self.P2 = f.get_startpos2_col()
blocknames = f.get_blocknames()
if block_char not in blocknames[0]:
raise ValueError(
f"The contig separator character `{block_char}` does not "
f"appear in the first block name `{blocknames[0]}`. Please "
"specify the correct character as `block_char`."
)
self.file_contigs = set(
itertools.chain.from_iterable([
blockname.split(block_char) for blockname in blocknames
])
)
if not len(self.file_contigs):
raise RuntimeError("No reference sequences found.")
for c1, c2 in itertools.combinations(self.file_contigs, 2):
if f.exists2(c1, c2) and f.exists2(c2, c1):
raise RuntimeError(
"Pairs are not triangular: found blocks "
+ f"'{c1}{block_char}{c2}'' and '{c2}{block_char}{c1}'"
)
# dumb heuristic to prevent excessively large chunks on one worker
if hasattr(f, "get_linecount"):
n_lines = f.get_linecount()
if n_lines < 0:
# correct int32 overflow bug
MAXINT32 = 2147483647
n_lines = MAXINT32 + MAXINT32 + n_lines
max_chunk = int(100e6)
n_chunks = n_lines // 2 // max_chunk
old_n = self.n_chunks
self.n_chunks = max(self.n_chunks, n_chunks)
if self.n_chunks > old_n:
logger.info(
f"Pairs file has {n_lines} lines. "
f"Increasing max-split to {self.n_chunks}."
)
# all requested contigs will be placed in the output matrix
self.gs = GenomeSegmentation(chromsizes, bins)
# find available contigs in the contact list
self.filepath = filepath
self.n_records = None
# warn about requested contigs not seen in the contact list
for chrom in self.gs.contigs:
if chrom not in self.file_contigs:
warnings.warn(
"Did not find contig " + f" '{chrom}' in contact list file.",
stacklevel=2,
)
def aggregate(
self, grange: tuple[str, int, int]
) -> pd.DataFrame | None: # pragma: no cover
chrom1, start, end = grange
import pypairix
filepath = self.filepath
binsize = self.gs.binsize
chromsizes = self.gs.chromsizes
chrom_binoffset = self.gs.chrom_binoffset
chrom_abspos = self.gs.chrom_abspos
start_abspos = self.gs.start_abspos
decr = int(self.is_one_based)
# C1 = self.C1
# C2 = self.C2
P1 = self.P1
P2 = self.P2
logger.info(f"Binning {chrom1}:{start}-{end}|*")
f = pypairix.open(filepath, "r")
these_bins = self.gs.fetch((chrom1, start, end))
remaining_chroms = self.gs.idmap[chrom1:]
# cid1 = self.gs.idmap[chrom1]
accumulator = Counter()
rows = []
for bin1_id, bin1 in these_bins.iterrows():
for chrom2, cid2 in remaining_chroms.items():
chrom2_size = chromsizes[chrom2]
if chrom1 != chrom2 and f.exists2(chrom2, chrom1): # flipped
iterator = f.query2D(
chrom2, 0, chrom2_size, chrom1, bin1.start, bin1.end
)
pos2_col = P1
else:
iterator = f.query2D(
chrom1, bin1.start, bin1.end, chrom2, 0, chrom2_size
)
pos2_col = P2
for line in iterator:
pos2 = int(line[pos2_col]) - decr
if binsize is None:
lo = chrom_binoffset[cid2]
hi = chrom_binoffset[cid2 + 1]
bin2_id = (
lo
+ np.searchsorted(
start_abspos[lo:hi],
chrom_abspos[cid2] + pos2,
side="right",
)
- 1
)
else:
bin2_id = chrom_binoffset[cid2] + (pos2 // binsize)
accumulator[bin2_id] += 1
if not accumulator:
continue
rows.append(
pd.DataFrame(
{
"bin1_id": bin1_id,
"bin2_id": list(accumulator.keys()),
"count": list(accumulator.values()),
},
columns=["bin1_id", "bin2_id", "count"],
).sort_values("bin2_id")
)
accumulator.clear()
logger.info(f"Finished {chrom1}:{start}-{end}|*")
return pd.concat(rows, axis=0) if len(rows) else None
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
granges = balanced_partition(self.gs, self.n_chunks, self.file_contigs)
for df in self._map(self.aggregate, granges):
if df is not None:
yield {k: v.values for k, v in df.items()}
class SparseBlockLoader(ContactBinner): # pragma: no cover
""" """
def __init__(self, chromsizes, bins, mapping, chunksize):
bins = check_bins(bins, chromsizes)
self.bins = bins
self.chromosomes = list(chromsizes.index)
self.chunksize = chunksize
n_chroms = len(chromsizes)
n_bins = len(bins)
chrom_ids = bins["chrom"].cat.codes
self.offsets = np.zeros(n_chroms + 1, dtype=int)
curr_val = 0
for start, _length, value in zip(*rlencode(chrom_ids)):
self.offsets[curr_val : value + 1] = start
curr_val = value + 1
self.offsets[curr_val:] = n_bins
self.mapping = mapping
def select_block(self, chrom1, chrom2):
try:
block = self.mapping[chrom1, chrom2]
except KeyError:
try:
block = self.mapping[chrom2, chrom1].T
except KeyError:
warnings.warn(
f"Block for {{{chrom1}, {chrom2}}} not found",
stacklevel=2,
)
raise
return block
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
# n_bins = len(self.bins)
import scipy.sparse
chromosomes = self.chromosomes
for cid1, chrom1 in enumerate(chromosomes):
offset = self.offsets[cid1]
chrom1_nbins = self.offsets[cid1 + 1] - offset
spans = partition(0, chrom1_nbins, self.chunksize)
for lo, hi in spans:
chunks = []
for chrom2 in chromosomes[cid1:]:
try:
block = self.select_block(chrom1, chrom2)
except KeyError:
continue
chunks.append(block.tocsr()[lo:hi, :])
X = scipy.sparse.hstack(chunks).tocsr().tocoo()
i, j, v = X.row, X.col, X.data
mask = (offset + i) <= (offset + j)
triu_i, triu_j, triu_v = i[mask], j[mask], v[mask]
yield {
"bin1_id": offset + triu_i,
"bin2_id": offset + triu_j,
"count": triu_v,
}
class ArrayLoader(ContactBinner):
"""
Load a dense genome-wide numpy array contact matrix.
Works with array-likes such as h5py.Dataset and memmapped arrays
See Also
--------
numpy.save, numpy.load (mmap_mode)
"""
def __init__(
self,
bins: pd.DataFrame,
array: np.ndarray | h5py.Dataset,
chunksize: int,
):
if len(bins) != array.shape[0]:
raise ValueError("Number of bins must equal the dimenion of the matrix")
self.array = array
self.chunksize = chunksize
def __iter__(self) -> Iterator[dict[str, np.ndarray]]:
n_bins = self.array.shape[0]
spans = partition(0, n_bins, self.chunksize)
# TRIU sparsify the matrix
for lo, hi in spans:
X = self.array[lo:hi, :]
i, j = np.nonzero(X)
mask = (lo + i) <= j
triu_i, triu_j = i[mask], j[mask]
yield {
"bin1_id": lo + triu_i,
"bin2_id": triu_j,
"count": X[triu_i, triu_j],
}
class ArrayBlockLoader(ContactBinner): # pragma: no cover
""" """
def __init__(self, chromsizes, bins, mapping, chunksize):
bins = check_bins(bins, chromsizes)
self.bins = bins
self.chromosomes = list(chromsizes.index)
self.chunksize = chunksize
n_chroms = len(chromsizes)
n_bins = len(bins)
chrom_ids = bins["chrom"].cat.codes
self.offsets = np.zeros(n_chroms + 1, dtype=int)
curr_val = 0
for start, _length, value in zip(*rlencode(chrom_ids)):
self.offsets[curr_val : value + 1] = start
curr_val = value + 1
self.offsets[curr_val:] = n_bins
self.mapping = mapping
def select_block(self, chrom1, chrom2):
try:
block = self.mapping[chrom1, chrom2]
except KeyError:
try:
block = self.mapping[chrom2, chrom1].T
except KeyError:
warnings.warn(
f"Block for {{{chrom1}, {chrom2}}} not found",
stacklevel=2,
)
raise
return block
def __iter__(self):
# n_bins = len(self.bins)
chromosomes = self.chromosomes
for cid1, chrom1 in enumerate(chromosomes):
offset = self.offsets[cid1]
chrom1_nbins = self.offsets[cid1 + 1] - offset
spans = partition(0, chrom1_nbins, self.chunksize)
for lo, hi in spans:
chunks = []
for chrom2 in chromosomes[cid1:]:
try:
block = self.select_block(chrom1, chrom2)
except KeyError:
continue
chunks.append(block[lo:hi, :])
X = np.concatenate(chunks, axis=1)
i, j = np.nonzero(X)
mask = (offset + i) <= (offset + j)
triu_i, triu_j = i[mask], j[mask]
yield {
"bin1_id": offset + triu_i,
"bin2_id": offset + triu_j,
"count": X[triu_i, triu_j],
}