Skip to content

API Reference

I recommend that TACT only be used via its console commands.

While TACT can be used as a Python library, its internal interface is subject to change at any time, even for minor or patch versions. This API documentation is provided merely for the sake of completeness.

Numerical functions

Functions to handle various numerical operations, including optimization.

crown_capture_probability(n, k)

Calculate the probability that a sample of k taxa from a clade of n total taxa includes a root node, under a Yule process.

This equation is taken from:

Sanderson, M. J. 1996. How many taxa must be sampled to identify the root node of a large clade? Systematic Biology 45:168-173

Parameters:

Name Type Description Default
n int

total number of taxa

required
k int

sampled taxa

required

Returns:

Type Description
float

probability

Source code in tact/lib.py
def crown_capture_probability(n, k):
    """
    Calculate the probability that a sample of `k` taxa from a clade
    of `n` total taxa includes a root node, under a Yule process.

    This equation is taken from:

    Sanderson, M. J. 1996. How many taxa must be sampled to identify
    the root node of a large clade? Systematic Biology 45:168-173

    Args:
        n (int): total number of taxa
        k (int): sampled taxa

    Returns:
        float: probability
    """
    if n < k:
        raise Exception(f"n must be greater than or equal to k (n={n}, k={k})")
    if n == 1 and k == 1:
        return 0  # not technically correct but it works for our purposes
    return 1 - 2 * (n - k) / ((n - 1) * (k + 1))

get_bd(r, a)

Converts turnover and relative extinction to birth and death rates.

Parameters:

Name Type Description Default
r float

turnover or net diversification (birth - death)

required
a float

relative extinction (death / birth)

required

Returns:

Type Description
(float, float)

birth, death

Source code in tact/lib.py
def get_bd(r, a):
    """
    Converts turnover and relative extinction to birth and death rates.

    Args:
        r (float): turnover or net diversification (birth - death)
        a (float): relative extinction (death / birth)

    Returns:
        (float, float): birth, death
    """
    return -r / (a - 1), -a * r / (a - 1)

get_new_times(ages, birth, death, missing, told=None, tyoung=None)

Simulates new speciation events in an incomplete phylogeny assuming a constant-rate birth-death process.

Adapted from the R function TreeSim::corsim written by Tanja Stadler.

N. Cusimano, T. Stadler, S. Renner. A new method for handling missing species in diversification analysis applicable to randomly or non-randomly sampled phylogenies. Syst. Biol., 61(5): 785-792, 2012.

Parameters:

Name Type Description Default
ages list

vector of waiting times

required
birth float

birth rate

required
death float

death rate

required
missing int

number of missing taxa to simulate

required
told float

maximum simulated age (default: max(ages))

None
tyoung float

minimum simulated age bound (default: 0)

None

Returns:

Type Description
list

vector of simulated waiting times.

Source code in tact/lib.py
def get_new_times(ages, birth, death, missing, told=None, tyoung=None):
    """
    Simulates new speciation events in an incomplete phylogeny assuming a
    constant-rate birth-death process.

    Adapted from the R function `TreeSim::corsim` written by Tanja Stadler.

    N. Cusimano, T. Stadler, S. Renner. A new method for handling missing
    species in diversification analysis applicable to randomly or
    non-randomly sampled phylogenies. Syst. Biol., 61(5): 785-792, 2012.

    Args:
        ages (list): vector of waiting times
        birth (float): birth rate
        death (float): death rate
        missing (int): number of missing taxa to simulate
        told (float): maximum simulated age (default: `max(ages)`)
        tyoung (float): minimum simulated age bound (default: `0`)

    Returns:
        list: vector of simulated waiting times.
    """
    if told is None:
        told = max(ages)
    if len(ages) > 0:
        if max(ages) > told and abs(max(ages) - told) > sys.float_info.epsilon:
            raise Exception("Zero or negative branch lengths detected in backbone phylogeny")
    if tyoung is None:
        tyoung = 0

    ages.sort(reverse=True)
    times = [x for x in ages if told >= x >= tyoung]
    times = [told] + times + [tyoung]
    ranks = range(0, len(times))
    only_new = []
    while missing > 0:
        if len(ranks) > 2:
            distrranks = []
            for i in range(1, len(ranks)):
                temp = ranks[i] * (
                    intp1(times[i - 1], birth, death) - intp1(times[i], birth, death)
                )
                distrranks.append(temp)
            try:
                dsum = sum(distrranks)
                distrranks = [x / dsum for x in distrranks]
                for i in range(1, len(distrranks)):
                    distrranks[i] = distrranks[i] + distrranks[i - 1]
                r = random.uniform(0, 1)
                addrank = min([idx for idx, x in enumerate(distrranks) if x > r])
            except ZeroDivisionError:
                addrank = 0
            except ValueError:
                addrank = 0
        else:
            addrank = 0
        r = random.uniform(0, 1)
        const = intp1(times[addrank], birth, death) - intp1(times[addrank + 1], birth, death)
        try:
            temp = intp1(times[addrank + 1], birth, death) / const
        except ZeroDivisionError:
            temp = 0.0
        xnew = (
            1
            / (death - birth)
            * log((1 - (r + temp) * const * birth) / (1 - (r + temp) * const * death))
        )
        only_new.append(xnew)
        missing -= 1
    only_new.sort(reverse=True)
    return only_new

get_ra(b, d)

Converts birth and death to turnover and relative extinction rates.

Parameters:

Name Type Description Default
b float

birth rate

required
d float

extinction rate

required

Returns:

Type Description
(float, float)

turnover, relative extinction

Source code in tact/lib.py
def get_ra(b, d):
    """
    Converts birth and death to turnover and relative extinction rates.

    Args:
        b (float): birth rate
        d (float): extinction rate

    Returns:
        (float, float): turnover, relative extinction
    """
    return (b - d, d / b)

intp1(t, l, m)

Source code in tact/lib.py
def intp1(t, l, m):
    try:
        return (1 - exp(-(l - m) * t)) / (l - m * exp(-(l - m) * t))
    except OverflowError:
        return float(intp1_exact(t, l, m))

intp1_exact(t, l, m)

Exact version of intp1 using Decimal math.

Source code in tact/lib.py
def intp1_exact(t, l, m):
    """Exact version of `intp1` using Decimal math."""
    l = D(l)
    m = D(m)
    t = D(t)
    num = D(1) - (-(l - m) * t).exp()
    denom = l - m * (-(l - m) * t).exp()
    return num / denom

lik_constant(vec, rho, t, root=1, survival=1, p1=<function p1 at 0x7faec1dc0550>)

Calculates the likelihood of a constant-rate birth-death process, conditioned on the waiting times of a phylogenetic tree and degree of incomplete sampling.

Based off of the R function TreePar::LikConstant written by Tanja Stadler.

T. Stadler. On incomplete sampling under birth-death models and connections to the sampling-based coalescent. Jour. Theo. Biol. 261: 58-66, 2009.

Parameters:

Name Type Description Default
vec float, float

two element tuple of birth and death

required
rho float

sampling fraction

required
t list

vector of waiting times

required
root bool

include the root or not? (default: 1)

1
survival bool

assume survival of the process? (default: 1)

1

Returns:

Type Description
float

a likelihood

Source code in tact/lib.py
def lik_constant(vec, rho, t, root=1, survival=1, p1=p1):
    """
    Calculates the likelihood of a constant-rate birth-death process, conditioned
    on the waiting times of a phylogenetic tree and degree of incomplete sampling.

    Based off of the R function `TreePar::LikConstant` written by Tanja Stadler.

    T. Stadler. On incomplete sampling under birth-death models and connections
    to the sampling-based coalescent. Jour. Theo. Biol. 261: 58-66, 2009.

    Args:
        vec (float, float): two element tuple of birth and death
        rho (float): sampling fraction
        t (list): vector of waiting times
        root (bool): include the root or not? (default: 1)
        survival (bool): assume survival of the process? (default: 1)

    Returns:
        float: a likelihood
    """
    l = vec[0]
    m = vec[1]
    t.sort(reverse=True)
    lik = (root + 1) * log(p1(t[0], l, m, rho))
    for tt in t[1:]:
        lik += log(l) + log(p1(tt, l, m, rho))
    if survival == 1:
        lik -= (root + 1) * log(1 - p0(t[0], l, m, rho))
    return -lik

optim_bd(ages, sampling, min_bound=1e-09)

Optimizes birth and death parameters given a vector of splitting times and sampling fraction.

Parameters:

Name Type Description Default
ages list

vector of node ages

required
sampling float

sampling fraction (0, 1]

required
min_bound float

minimum birth rate

1e-09

Returns:

Type Description
float, float

birth and death rates

Source code in tact/lib.py
def optim_bd(ages, sampling, min_bound=1e-9):
    """
    Optimizes birth and death parameters given a vector of splitting times and sampling fraction.

    Args:
        ages (list): vector of node ages
        sampling (float): sampling fraction (0, 1]
        min_bound (float): minimum birth rate

    Returns:
        float, float: birth and death rates
    """
    if max(ages) < 0.000001:
        init_r = 1e-3
    else:
        # Magallon-Sanderson crown estimator
        init_r = (log((len(ages) + 1) / sampling) - log(2)) / max(ages)
        init_r = max(1e-3, init_r)
    bounds = ((min_bound, 100), (0, 1 - min_bound))
    result = two_step_optim(wrapped_lik_constant, x0=(init_r, min_bound), bounds=bounds, args=(sampling, ages))
    return get_bd(*result)

optim_yule(ages, sampling, min_bound=1e-09)

Optimizes birth parameter under a Yule model, given a vector of splitting times and sampling fraction.

Parameters:

Name Type Description Default
ages list

vector of node ages

required
sampling float

sampling fraction (0, 1]

required
min_bound float

minimum birth rate

1e-09

Returns:

Type Description
float, float

birth and death rates (where death is always 0)

Source code in tact/lib.py
def optim_yule(ages, sampling, min_bound=1e-9):
    """
    Optimizes birth parameter under a Yule model, given a vector of splitting times and sampling fraction.

    Args:
        ages (list): vector of node ages
        sampling (float): sampling fraction (0, 1]
        min_bound (float): minimum birth rate

    Returns:
        float, float: birth and death rates (where death is always 0)
    """
    bounds = (min_bound, 100)
    result = minimize_scalar(wrapped_lik_constant_yule, bounds=bounds, args=(sampling, ages), method="Bounded")
    if result["success"]:
        return (result["x"], 0.0)

    raise Exception(f"Optimization failed: {result['message']} (code {result['status']})")

p0(t, l, m, rho)

Source code in tact/lib.py
def p0(t, l, m, rho):
    try:
        return 1 - rho * (l - m) / (rho * l + (l * (1 - rho) - m) * exp(-(l - m) * t))
    except FloatingPointError:
        return float(p0_exact(t, l, m, rho))

p0_exact(t, l, m, rho)

Exact version of p0 using Decimal math.

Source code in tact/lib.py
def p0_exact(t, l, m, rho):
    "Exact version of `p0` using Decimal math."
    t = D(t)
    l = D(l)
    m = D(m)
    rho = D(rho)
    return D(1) - rho * (l - m) / (rho * l + (l * (D(1) - rho) - m) * (-(l - m) * t).exp())

p1(t, l, m, rho)

Optimized version of p1_orig using common subexpression elimination and strength reduction from exponentiation to multiplication.

Source code in tact/lib.py
def p1(t, l, m, rho):
    """
    Optimized version of `p1_orig` using common subexpression elimination and strength reduction
    from exponentiation to multiplication.
    """
    try:
        ert = np.exp(-(l - m) * t, dtype=np.float64)
        num = rho * (l - m) ** 2 * ert
        denom = (rho * l + (l * (1 - rho) - m) * ert) ** 2
        res = num / denom
    except (OverflowError, FloatingPointError):
        res = float(p1_exact(t, l, m, rho))
    if res == 0.0:
        return sys.float_info.min
    return res

p1_exact(t, l, m, rho)

Exact version of p1 using Decimal math.

Source code in tact/lib.py
def p1_exact(t, l, m, rho):
    """Exact version of `p1` using Decimal math."""
    t = D(t)
    l = D(l)
    m = D(m)
    rho = D(rho)
    num = rho * (l - m) ** D(2) * (-(l - m) * t).exp()
    denom = (rho * l + (l * (1 - rho) - m) * (-(l - m) * t).exp()) ** D(2)
    return num / denom

p1_orig(t, l, m, rho)

Original version of p1, here for testing and comparison purposes.

Source code in tact/lib.py
def p1_orig(t, l, m, rho):
    """Original version of `p1`, here for testing and comparison purposes."""
    try:
        num = rho * (l - m) ** 2 * np.exp(-(l - m) * t)
        denom = (rho * l + (l * (1 - rho) - m) * np.exp(-(l - m) * t)) ** 2
        res = num / denom
    except (OverflowError, FloatingPointError):
        res = float(p1_exact(t, l, m, rho))
    if res == 0.0:
        return sys.float_info.min
    return res

two_step_optim(func, x0, bounds, args)

Conduct a two-step function optimization, first by using the fast L-BFGS-B method, and if that fails, use simulated annealing.

Parameters:

Name Type Description Default
func callable

function to optimize

required
x0 tuple

initial conditions

required
bounds tuple

boundary conditions

required
args lsit

additional argumnets to pass to func

required

Returns:

Type Description
tuple

optimized parameter values

Source code in tact/lib.py
def two_step_optim(func, x0, bounds, args):
    """
    Conduct a two-step function optimization, first by using the fast L-BFGS-B method,
    and if that fails, use simulated annealing.

    Args:
        func (callable): function to optimize
        x0 (tuple): initial conditions
        bounds (tuple): boundary conditions
        args (lsit): additional argumnets to pass to `func`

    Returns:
        tuple: optimized parameter values
    """
    try:
        result = minimize(func, x0=x0, bounds=bounds, args=args, method="L-BFGS-B")
        if result["success"]:
            return result["x"].tolist()
    except FloatingPointError:
        pass

    result = dual_annealing(func, x0=x0, bounds=bounds, args=args)
    if result["success"]:
        return result["x"].tolist()

    raise Exception(f"Optimization failed: {result['message']} (code {result['status']})")

wrapped_lik_constant(x, sampling, ages)

Wrapper for birth-death likelihood to make optimizing more convenient.

Parameters:

Name Type Description Default
x float, float

turnover, relative extinction

required
sampling float

sampling fraction (0, 1]

required
ages list

vector of node ages

required

Returns:

Type Description
float

a likelihood

Source code in tact/lib.py
def wrapped_lik_constant(x, sampling, ages):
    """
    Wrapper for birth-death likelihood to make optimizing more convenient.

    Args:
        x (float, float): turnover, relative extinction
        sampling (float): sampling fraction (0, 1]
        ages (list): vector of node ages

    Returns:
        float: a likelihood
    """
    return lik_constant(get_bd(*x), sampling, ages)

wrapped_lik_constant_yule(x, sampling, ages)

Wrapper for Yule likelihood to make optimizing more convenient.

Parameters:

Name Type Description Default
x float

birth rate

required
sampling float

sampling fraction (0, 1]

required
ages list

vector of node ages

required

Returns:

Type Description
float

a likelihood

Source code in tact/lib.py
def wrapped_lik_constant_yule(x, sampling, ages):
    """
    Wrapper for Yule likelihood to make optimizing more convenient.

    Args:
        x (float): birth rate
        sampling (float): sampling fraction (0, 1]
        ages (list): vector of node ages

    Returns:
        float: a likelihood
    """
    return lik_constant((x, 0.0), sampling, ages)

Tree functions

Functions specifically to handle DendroPy tree objects.

compute_node_depths(tree)

Returns a dictionary of node depths for each node with a label.

Source code in tact/tree_util.py
def compute_node_depths(tree):
    """Returns a dictionary of node depths for each node with a label."""
    res = {}
    for leaf in tree.leaf_node_iter():
        cnt = 0
        for anc in leaf.ancestor_iter():
            if anc.label:
                cnt += 1
        res[leaf.taxon.label] = cnt
    return res

count_locked(node)

How many edges under node are locked?

Source code in tact/tree_util.py
def count_locked(node):
    """How many edges under `node` are locked?"""
    sum([x.label == "locked" for x in edge_iter(node)])

edge_iter(node, filter_fn=None)

Iterates over the child edge of node and all its descendants. Can optionally be filtered by filter_fn.

Source code in tact/tree_util.py
def edge_iter(node, filter_fn=None):
    """
    Iterates over the child edge of `node` and all its descendants.
    Can optionally be filtered by `filter_fn`.
    """
    stack = list(node.child_edge_iter())
    while stack:
        edge = stack.pop()
        if filter_fn is None or filter_fn(edge):
            yield edge
        stack.extend(edge.head_node.child_edge_iter())

ensure_tree_node_depths(tree)

Source code in tact/tree_util.py
def ensure_tree_node_depths(tree):
    node_depths = compute_node_depths(tree)
    stats = collections.defaultdict(int)
    for v in node_depths.values():
        stats[v] += 1
    msg = ""
    if len(stats) > 1:
        msg += "The tips of your taxonomy tree do not have equal numbers of ranked clades in their ancestor chain:\n"
        for k in sorted(stats.keys()):
            msg += f"* {stats[k]} tips have {k} ranked ancestors\n"
    return msg

get_age_intervals(node)

Gets the (possibly disjoint) interval that could be generated in the clade under node, assuming that grafts to locked edges are restricted.

Source code in tact/tree_util.py
def get_age_intervals(node):
    """
    Gets the (possibly disjoint) interval that could be generated in the
    clade under `node`, assuming that grafts to locked edges are restricted.
    """
    acc = portion.empty()
    for edge in edge_iter(node, lambda x: x.label != "locked"):
        acc = acc | portion.closed(edge.head_node.age, edge.tail_node.age)
    return acc

get_ages(node, include_root=False)

Returns the list of ages of the children of a given node, optionally including the node's age if include_root is True.

Source code in tact/tree_util.py
def get_ages(node, include_root=False):
    """
    Returns the list of ages of the children of a given `node`,
    optionally including the `node`'s age if `include_root` is True.
    """
    ages = [x.age for x in node.ageorder_iter(include_leaves=False, descending=True)]
    if include_root:
        ages += [node.age]
    return ages

get_birth_death_rates(node, sampfrac, yule=False, include_root=False)

Estimates the birth and death rates for the subtree descending from node with sampling fraction sampfrac. Optionally restrict to a Yule pure-birth model.

Source code in tact/tree_util.py
def get_birth_death_rates(node, sampfrac, yule=False, include_root=False):
    """
    Estimates the birth and death rates for the subtree descending from
    `node` with sampling fraction `sampfrac`. Optionally restrict to a
    Yule pure-birth model.
    """
    if yule:
        return optim_yule(get_ages(node, include_root), sampfrac)

    return optim_bd(get_ages(node, include_root), sampfrac)

get_min_age(node)

Gets the minimum possible age that could be generated in a clade under node, assuming that grafts to locked edges are restricted.

Source code in tact/tree_util.py
def get_min_age(node):
    """
    Gets the minimum possible age that could be generated in a clade under `node`,
    assuming that grafts to locked edges are restricted.
    """
    interval = get_age_intervals(node)

    if not interval.atomic:
        raise DisjointConstraintError(f"Constraint on {node} implies disjoint interval {interval}")

    if interval.empty:
        return 0.0

    return interval.lower

get_monophyletic_node(tree, species)

Returns the node or None that is the MRCA of the species in tree.

Source code in tact/tree_util.py
def get_monophyletic_node(tree, species):
    """Returns the node or None that is the MRCA of the `species` in `tree`."""
    mrca = tree.mrca(taxon_labels=species)
    if mrca and species.issuperset(get_tip_labels(mrca)):
        return mrca

    return None

get_short_branches(node)

Yields an iterator of especially short edges under node.

Source code in tact/tree_util.py
def get_short_branches(node):
    """Yields an iterator of especially short edges under `node`."""
    for edge in edge_iter(node):
        if edge.length <= 0.001:
            yield edge

get_tip_labels(tree_or_node)

Returns a set of tip labels for a node or tree.

Source code in tact/tree_util.py
def get_tip_labels(tree_or_node):
    """Returns a `set` of tip labels for a node or tree."""
    try:
        return {x.taxon.label for x in tree_or_node.leaf_node_iter()}
    except AttributeError:
        return {x.taxon.label for x in tree_or_node.leaf_iter()}

get_tree(path, namespace=None)

Gets a DendroPy tree from a path and precalculate its node ages and bipartition bitmask.

Source code in tact/tree_util.py
def get_tree(path, namespace=None):
    """
    Gets a DendroPy tree from a path and precalculate its node ages and bipartition bitmask.
    """
    tree = dendropy.Tree.get_from_path(
        path, schema="newick", taxon_namespace=namespace, rooting="default-rooted"
    )
    update_tree_view(tree)
    return tree

graft_node(graft_recipient, graft, stem=False)

Grafts a node graft randomly in the subtree below node graft_recipient. The attribute graft.age must be set so we know where is the best place to graft the node. The node graft can optionally have child nodes, in this case the edge.length attribute should be set on all child nodes if the tree is to remain ultrametric.

We graft things "below" a node by picking one of the children of that node and forcing it to be sister to the grafted node and adjusting the edge lengths accordingly. Therefore, the node above which the graft lives (i.e., the one that will be the child of the new graft) must fulfill the following requirements:

  1. Must not be the crown node (cannot graft things above crown node)
  2. Must be younger than the graft node (no negative branches)
  3. Seed node must be older than graft node (no negative branches)
  4. Must not be locked (intruding on monophyly)
Source code in tact/tree_util.py
def graft_node(graft_recipient, graft, stem=False):
    """
    Grafts a node `graft` randomly in the subtree below node
    `graft_recipient`. The attribute `graft.age` must be set so
    we know where is the best place to graft the node. The node
    `graft` can optionally have child nodes, in this case the
    `edge.length` attribute should be set on all child nodes if
    the tree is to remain ultrametric.

    We graft things "below" a node by picking one of the children
    of that node and forcing it to be sister to the grafted node
    and adjusting the edge lengths accordingly. Therefore, the node
    *above* which the graft lives (i.e., the one that will be the child
    of the new graft) must fulfill the following requirements:

    1. Must not be the crown node (cannot graft things above crown node)
    2. Must be younger than the graft node (no negative branches)
    3. Seed node must be older than graft node (no negative branches)
    4. Must not be locked (intruding on monophyly)
    """
    def filter_fn(x):
        return x.head_node.age <= graft.age and x.head_node.parent_node.age >= graft.age and x.label != "locked"

    all_edges = list(edge_iter(graft_recipient))
    if stem:
        # also include the crown node's subtending edge
        all_edges.append(graft_recipient.edge)
    eligible_edges = [x for x in all_edges if filter_fn(x)]

    if not eligible_edges:
        raise Exception(f"could not place node {graft} in clade {graft_recipient}")

    focal_node = random.choice([x.head_node for x in eligible_edges])
    seed_node = focal_node.parent_node
    sisters = focal_node.sibling_nodes()

    # pick a child edge and detach its corresponding node
    #
    # DendroPy's Node.remove_child() messes with the edge lengths.
    # But, Node.clear_child_nodes() simply cuts that bit of the tree out.
    seed_node.clear_child_nodes()

    # set the correct edge length on the grafted node and make the grafted
    # node a child of the seed node
    graft.edge.length = seed_node.age - graft.age
    if graft.edge.length < 0:
        raise Exception("negative branch length")
    sisters.append(graft)
    seed_node.set_child_nodes(sisters)

    # make the focal node a child of the grafted node and set edge length
    focal_node.edge.length = graft.age - focal_node.age
    if focal_node.edge.length < 0:
        raise Exception("negative branch length")
    graft.add_child(focal_node)

    # return the (potentially new) crown of the clade
    if graft_recipient.parent_node == graft:
        return graft
    return graft_recipient

is_binary(node)

Is the subtree under node a fully bifurcating tree?

Source code in tact/tree_util.py
def is_binary(node):
    """Is the subtree under `node` a fully bifurcating tree?"""
    for internal_node in node.preorder_internal_node_iter():
        if len(internal_node.child_nodes()) != 2:
            return False
    return True

is_fully_locked(node)

Are all the edges below node locked?

Source code in tact/tree_util.py
def is_fully_locked(node):
    """
    Are all the edges below `node` locked?
    """
    return all(x.label == "locked" for x in edge_iter(node))

is_ultrametric(tree, tolerance=1e-06)

Is the tree ultrametric, within a specified tolerance?

Uses the relative difference between minimum and maximum root-to-tip distances.

Source code in tact/tree_util.py
def is_ultrametric(tree, tolerance=1e-6):
    """Is the `tree` ultrametric, within a specified `tolerance`?

    Uses the relative difference between minimum and maximum root-to-tip distances.
    """
    tree.calc_node_root_distances()
    lengths = {}
    for leaf in tree.leaf_node_iter():
        lengths[leaf.taxon.label] = leaf.root_distance
    t_min = min(lengths.items(), key=lambda x: x[1])
    t_max = max(lengths.items(), key=lambda x: x[1])
    return (math.isclose(t_min[1], t_max[1], rel_tol=tolerance), (t_min, t_max))

lock_clade(node)

Locks a clade descending from node so future grafts will avoid locked edges.

Source code in tact/tree_util.py
def lock_clade(node):
    """
    Locks a clade descending from `node` so future grafts will avoid locked edges.
    """
    for edge in edge_iter(node):
        edge.label = "locked"

unlock_clade(node)

Unlocks a clade descending from node so new tips can be grafted to its edges.

Source code in tact/tree_util.py
def unlock_clade(node):
    """
    Unlocks a clade descending from `node` so new tips can be grafted to its edges.
    """
    for edge in edge_iter(node):
        edge.label = ""

update_tree_view(tree)

Mutates a DendroPy tree object with updated node ages and bipartition bitmask. We also correct for minor ultrametricity errors.

Returns a list of tip labels.

Source code in tact/tree_util.py
def update_tree_view(tree):
    """
    Mutates a DendroPy tree object with updated node ages and bipartition bitmask. We also
    correct for minor ultrametricity errors.

    Returns a list of tip labels.
    """
    tree.calc_node_ages(is_force_max_age=True)
    tree.update_bipartitions()
    return get_tip_labels(tree)

FastMRCA

Singleton object that helps speed up MRCA lookups.

bitmask(labels)

Gets a bitmask for the taxa in labels, potentially in parallel.

Source code in tact/fastmrca.py
def bitmask(labels):
    """
    Gets a bitmask for the taxa in `labels`, potentially in parallel.
    """
    global tree
    tn = tree.taxon_namespace
    return tn.taxa_bitmask(labels=labels)

fastmrca_getter(tn, x)

Helper function for submitting stuff.

Source code in tact/fastmrca.py
def fastmrca_getter(tn, x):
    """Helper function for submitting stuff."""
    taxa = tn.get_taxa(labels=x)
    mask = 0
    for taxon in taxa:
        mask |= tn.taxon_bitmask(taxon)
    return mask

get(labels)

Pulls a MRCA node out for the taxa in labels.

Source code in tact/fastmrca.py
def get(labels):
    """Pulls a MRCA node out for the taxa in `labels`."""
    global tree
    labels = set(labels)
    mrca = tree.mrca(leafset_bitmask=bitmask(labels))
    if mrca and labels.issuperset(get_tip_labels(mrca)):
        return mrca
    return None

initialize(phy)

Initialize the fastmrca singleton with a tree.

Source code in tact/fastmrca.py
def initialize(phy):
    """
    Initialize the fastmrca singleton with a tree.
    """
    global tree
    tree = phy

Last update: October 11, 2021