Skip to content

Commit

Permalink
Remove extraneous str casts (#159)
Browse files Browse the repository at this point in the history
* Remove extraneous str casts

* Edge case checking for empty tables

* imports and test was doing unnecessary things:

* kick tests

* kick tests

* the value of a local env
  • Loading branch information
wasade authored Nov 14, 2023
1 parent bd3cfac commit 8ffb515
Show file tree
Hide file tree
Showing 2 changed files with 48 additions and 9 deletions.
32 changes: 23 additions & 9 deletions unifrac/_methods.py
Original file line number Diff line number Diff line change
Expand Up @@ -50,6 +50,16 @@ def is_biom_v210(f, ids=None):
return True


def has_samples_biom_v210(f):
# assumes this is already checked to be biom v210
import h5py
with h5py.File(f, 'r') as fp:
if len(fp['sample/ids']) == 0:
return False

return True


def is_newick(f):
sniffer = skbio.io.format.newick.newick.sniffer_function
return sniffer(f)[0]
Expand All @@ -58,12 +68,16 @@ def is_newick(f):
def _validate(table, phylogeny, ids=None):
if not is_biom_v210(table, ids):
raise ValueError("Table does not appear to be a BIOM-Format v2.1")
if not has_samples_biom_v210(table):
raise ValueError("Table does not contain any samples")
if not is_newick(phylogeny):
raise ValueError("The phylogeny does not appear to be newick")


def _call_ssu(table, phylogeny, *args):
if isinstance(table, Table) and isinstance(phylogeny, (TreeNode, BP)):
if table.is_empty():
raise ValueError("Table does not contain any samples")
return qsu.ssu_inmem(table, phylogeny, *args)
elif isinstance(table, str) and isinstance(phylogeny, str):
ids = []
Expand Down Expand Up @@ -356,7 +370,7 @@ def weighted_normalized(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_normalized',
return _call_ssu(table, phylogeny, 'weighted_normalized',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -423,7 +437,7 @@ def weighted_normalized_fp64(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_normalized_fp64',
return _call_ssu(table, phylogeny, 'weighted_normalized_fp64',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -490,7 +504,7 @@ def weighted_normalized_fp32(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_normalized_fp32',
return _call_ssu(table, phylogeny, 'weighted_normalized_fp32',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -557,7 +571,7 @@ def weighted_unnormalized(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_unnormalized',
return _call_ssu(table, phylogeny, 'weighted_unnormalized',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -625,7 +639,7 @@ def weighted_unnormalized_fp64(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_unnormalized_fp64',
return _call_ssu(table, phylogeny, 'weighted_unnormalized_fp64',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -693,7 +707,7 @@ def weighted_unnormalized_fp32(table: Union[str, Table],
powerful beta diversity measure for comparing communities based on
phylogeny. BMC Bioinformatics 12:118 (2011).
"""
return _call_ssu(str(table), str(phylogeny), 'weighted_unnormalized_fp32',
return _call_ssu(table, phylogeny, 'weighted_unnormalized_fp32',
variance_adjusted, 1.0, bypass_tips, n_substeps)


Expand Down Expand Up @@ -779,7 +793,7 @@ def generalized(table: Union[str, Table],
return weighted_normalized(table, phylogeny, threads,
variance_adjusted, bypass_tips, n_substeps)
else:
return _call_ssu(str(table), str(phylogeny), 'generalized',
return _call_ssu(table, phylogeny, 'generalized',
variance_adjusted, alpha, bypass_tips, n_substeps)


Expand Down Expand Up @@ -866,7 +880,7 @@ def generalized_fp64(table: Union[str, Table],
variance_adjusted, bypass_tips,
n_substeps)
else:
return _call_ssu(str(table), str(phylogeny), 'generalized_fp64',
return _call_ssu(table, phylogeny, 'generalized_fp64',
variance_adjusted, alpha, bypass_tips, n_substeps)


Expand Down Expand Up @@ -953,7 +967,7 @@ def generalized_fp32(table: Union[str, Table],
variance_adjusted, bypass_tips,
n_substeps)
else:
return _call_ssu(str(table), str(phylogeny), 'generalized_fp32',
return _call_ssu(table, phylogeny, 'generalized_fp32',
variance_adjusted, alpha, bypass_tips, n_substeps)


Expand Down
25 changes: 25 additions & 0 deletions unifrac/tests/test_methods.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,11 +7,16 @@
# ----------------------------------------------------------------------------
import unittest
import pkg_resources
import os

import h5py
import biom
import skbio
import numpy as np
import numpy.testing as npt

from unifrac import meta
from unifrac._methods import _call_ssu, has_samples_biom_v210


class StateUnifracTests(unittest.TestCase):
Expand Down Expand Up @@ -107,6 +112,26 @@ def test_meta_unifrac_alpha_not_generalized(self):
meta((self.table1, ), (self.tree1, ), method='unweighted',
alpha=1, consolidation='skipping_missing_matrices')

def test_has_samples_biom_v210(self):
fp = self.get_data_path('crawford.biom')
self.assertTrue(has_samples_biom_v210(fp))

tmpfile = '/tmp/fake.biom'
empty = biom.Table([], [], [])
try:
with h5py.File(tmpfile, 'w') as fp:
empty.to_hdf5(fp, 'asd')
fp = self.get_data_path(tmpfile)
self.assertFalse(has_samples_biom_v210(tmpfile))
finally:
os.unlink(tmpfile)

def test_call_ssu_empty_biom(self):
empty = biom.Table([], [], [])
tre = skbio.TreeNode()
with self.assertRaisesRegex(ValueError, "contain any samples"):
_call_ssu(empty, tre)


if __name__ == "__main__":
unittest.main()

0 comments on commit 8ffb515

Please sign in to comment.