dendropy.datamodel.treemodel: Trees

The Tree Class

class dendropy.datamodel.treemodel.Tree(*args, **kwargs)[source]

An arborescence, i.e. a fully-connected directed acyclic graph with all edges directing away from the root and toward the tips. The “root” of the tree is represented by the Tree.seed_node attribute. In unrooted trees, this node is an algorithmic artifact. In rooted trees this node is semantically equivalent to the root.

The constructor can optionally construct a Tree object by cloning another Tree object passed as the first positional argument, or out of a data source if stream and schema keyword arguments are passed with a file-like object and a schema-specification string object values respectively.

Parameters:
  • *args (positional argument, optional) – If given, should be exactly one Tree object. The new Tree will then be a structural clone of this argument.

  • **kwargs (keyword arguments, optional) –

    The following optional keyword arguments are recognized and handled by this constructor:

    label

    The label or description of the new Tree object.

    taxon_namespace

    Specifies the TaxonNamespace object to be that the new Tree object will reference.

Examples

Tree objects can be instantiated in the following ways:

# /usr/bin/env python

try:
    from StringIO import StringIO
except ImportError:
    from io import StringIO
from dendropy import Tree, TaxonNamespace

# empty tree
t1 = Tree()

# Tree objects can be instantiated from an external data source
# using the 'get()' factory class method

# From a file-like object
t2 = Tree.get(file=open('treefile.tre', 'r'),
                schema="newick",
                tree_offset=0)

# From a path
t3 = Tree.get(path='sometrees.nexus',
        schema="nexus",
        collection_offset=2,
        tree_offset=1)

# From a string
s = "((A,B),(C,D));((A,C),(B,D));"
# tree will be '((A,B),(C,D))'
t4 = Tree.get(data=s,
        schema="newick")
# tree will be '((A,C),(B,D))'
t5 = Tree.get(data=s,
        schema="newick",
        tree_offset=1)
# passing keywords to underlying tree parser
t7 = dendropy.Tree.get(
        data="((A,B),(C,D));",
        schema="newick",
        taxon_namespace=t3.taxon_namespace,
        suppress_internal_node_taxa=False,
        preserve_underscores=True)

# Tree objects can be written out using the 'write()' method.
t1.write(file=open('treefile.tre', 'r'),
        schema="newick")
t1.write(path='treefile.nex',
        schema="nexus")

# Or returned as a string using the 'as_string()' method.
s = t1.as_string("nexml")

# tree structure deep-copied from another tree
t8 = dendropy.Tree(t7)
assert t8 is not t7                             # Trees are distinct
assert t8.symmetric_difference(t7) == 0         # and structure is identical
assert t8.taxon_namespace is t7.taxon_namespace             # BUT taxa are not cloned.
nds3 = [nd for nd in t7.postorder_node_iter()]  # Nodes in the two trees
nds4 = [nd for nd in t8.postorder_node_iter()]  # are distinct objects,
for i, n in enumerate(nds3):                    # and can be manipulated
    assert nds3[i] is not nds4[i]               # independentally.
egs3 = [eg for eg in t7.postorder_edge_iter()]  # Edges in the two trees
egs4 = [eg for eg in t8.postorder_edge_iter()]  # are also distinct objects,
for i, e in enumerate(egs3):                    # and can also be manipulated
    assert egs3[i] is not egs4[i]               # independentally.
lves7 = t7.leaf_nodes()                         # Leaf nodes in the two trees
lves8 = t8.leaf_nodes()                         # are also distinct objects,
for i, lf in enumerate(lves3):                  # but order is the same,
    assert lves7[i] is not lves8[i]             # and associated Taxon objects
    assert lves7[i].taxon is lves8[i].taxon     # are the same.

# To create deep copy of a tree with a different taxon namespace,
# Use 'copy.deepcopy()'
t9 = copy.deepcopy(t7)

# Or explicitly pass in a new TaxonNamespace instance
taxa = TaxonNamespace()
t9 = dendropy.Tree(t7, taxon_namespace=taxa)
assert t9 is not t7                             # As above, the trees are distinct
assert t9.symmetric_difference(t7) == 0         # and the structures are identical,
assert t9.taxon_namespace is not t7.taxon_namespace         # but this time, the taxa *are* different
assert t9.taxon_namespace is taxa                     # as the given TaxonNamespace is used instead.
lves3 = t7.leaf_nodes()                         # Leaf nodes (and, for that matter other nodes
lves5 = t9.leaf_nodes()                         # as well as edges) are also distinct objects
for i, lf in enumerate(lves3):                  # and the order is the same, as above,
    assert lves7[i] is not lves9[i]             # but this time the associated Taxon
    assert lves7[i].taxon is not lves9[i].taxon # objects are distinct though the taxon
    assert lves7[i].taxon.label == lves9[i].taxon.label # labels are the same.

# to 'switch out' the TaxonNamespace of a tree, replace the reference and
# reindex the taxa:
t11 = Tree.get(data='((A,B),(C,D));', 'newick')
taxa = TaxonNamespace()
t11.taxon_namespace = taxa
t11.reindex_subcomponent_taxa()

# You can also explicitly pass in a seed node:
seed = Node(label="root")
t12 = Tree(seed_node=seed)
assert t12.seed_node is seed
B1()[source]

DEPRECATED: Use dendropy.calculate.treemeasure.B1.

N_bar()[source]

DEPRECATED: Use ‘dendropy.calculate.treemeasure.N_bar()’.

__iter__()[source]

Iterate over nodes on tree in pre-order.

Example

>>> for nd in tree:
...    print(nd.label)
...
Returns:

collections.Iterator [Node] – An iterator yielding the internal nodes of the subtree rooted at this node in post-order sequence.

__len__()[source]

Returns number of tips on tree (could be less than number of taxa in namespace).

__str__()[source]

Dump Newick string.

age_order_node_iter(include_leaves=True, filter_fn=None, descending=False)[source]

Deprecated: use Tree.ageorder_node_iter instead.

ageorder_node_iter(include_leaves=True, filter_fn=None, descending=False)[source]
Iterator over nodes of tree in order of the age of the node (i.e., the

time since the present).

Iterates over nodes in order of age (‘age’ is as given by the age attribute, which is usually the sum of edge lengths from tips to node, i.e., time since present). If include_leaves is True (default), leaves are included in the iteration; if include_leaves is False, leaves will be skipped. If descending is False (default), younger nodes will be returned before older ones; if True, older nodes will be returned before younger ones.

Parameters:
  • include_leaves (boolean, optional) – If True (default), then leaf nodes are included in the iteration. If False, then leaf nodes are skipped.

  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • descending (boolean, optional) – If False (default), then younger nodes are visited before older ones. If True, then older nodes are visited before younger ones.

Returns:

collections.Iterator [Node] – Iterator over age-ordered sequence of nodes of self.

apply(before_fn=None, after_fn=None, leaf_fn=None)[source]

Applies function before_fn and after_fn to all internal nodes and leaf_fn to all terminal nodes in subtree starting with self, with nodes visited in pre-order.

Given a tree with preorder sequence of nodes of [a,b,i,e,j,k,c,g,l,m,f,n,h,o,p,]:

              a
             / \
            /   \
           /     \
          /       \
         /         \
        /           \
       /             c
      b             / \
     / \           /   \
    /   e         /     f
   /   / \       /     / \
  /   /   \     g     /   h
 /   /     \   / \   /   / \
i   j       k l   m n   o   p

the following order of function calls results:

before_fn(a) before_fn(b) leaf_fn(i) before_fn(e) leaf_fn(j) leaf_fn(k) after_fn(e) after_fn(b) before_fn(c) before_fn(g) leaf_fn(l) leaf_fn(m) after_fn(g) before_fn(f) leaf_fn(n) before_fn(h) leaf_fn(o) leaf_fn(p) after_fn(h) after_fn(f) after_fn(c) after_fn(a)

Parameters:
  • before_fn (function object or None) – A function object that takes a Node as its argument.

  • after_fn (function object or None) – A function object that takes a Node as its argument.

  • leaf_fn (function object or None) – A function object that takes a Node as its argument.

Notes

Adapted from work by Mark T. Holder (the peyotl module of the Open Tree of Life Project):

as_ascii_plot(**kwargs)[source]

Returns a string representation a graphic of this tree using ASCII characters. See AsciiTreePlot for details on arguments.

as_python_source(tree_obj_name=None, tree_args=None)[source]

Returns string that will rebuild this tree in Python.

as_string(schema, **kwargs)

Composes and returns string representation of the data.

Mandatory Schema-Specification Keyword Argument:

Optional Schema-Specific Keyword Arguments:

These provide control over how the data is formatted, and supported argument names and values depend on the schema as specified by the value passed as the “schema” argument. See “DendroPy Schemas: Phylogenetic and Evolutionary Biology Data Formats” for more details.

calc_node_ages(ultrametricity_precision=1e-05, is_force_max_age=False, is_force_min_age=False, set_node_age_fn=None, is_return_internal_node_ages_only=False)[source]

Adds an attribute called “age” to each node, with the value equal to the sum of edge lengths from the node to the tips.

NOTE: Consider using the newer and more flexible resolve_node_ages() instead of this.

Parameters:
  • ultrametricity_precision (numeric or bool or None) – If the lengths of different paths to the node differ by more than ultrametricity_precision, then a ValueError exception will be raised indicating deviation from ultrametricity. If ultrametricity_precision is negative or False, then this check will be skipped.

  • is_force_max_age (bool) – If is_force_max_age is True, then each node will be set to the maximum possible age, by being set to the oldest age given its child set and the subtending edge lengths. This option only makes a difference if the tree is not ultrametric, and so the ultrametricity precision check is ignore if this option is set to True.

  • is_force_min_age (bool) – If is_force_min_age is True then each node will be set to the minimum possible age, by being set to the youngest age given its child set and the subtending edge lengths. This option only makes a difference if the tree is not ultrametric, and so the ultrametricity precision check is ignore if this option is set to True.

  • set_node_age_fn (function object) –

    If not None, then this should be a function that takes a node as an argument and returns None or a non-None value. If None, then this indicates that the node’s age should be calculated by this function. If not None, then this is the value that this node’s age should be set to. This can be used to set non-contemporary tip ages by passing something like:

    f = lambda nd: None if not nd.is_leaf else nd.annotations[“height”]

    which returns None if the node is an internal node, but otherwise returns the value in the height annotation.

Returns:

a (iterable[numeric]) – Returns collection of node ages.

calc_node_root_distances(return_leaf_distances_only=True)[source]

Adds attribute “root_distance” to each node, with value set to the sum of edge lengths from the node to the root. Returns list of distances. If return_leaf_distances_only is True, then only leaf distances will be true.

clone(depth=1)

Creates and returns a copy of self.

Parameters:

depth (integer) –

The depth of the copy:

  • 0: shallow-copy: All member objects are references, except for :attr:annotation_set of top-level object and member Annotation objects: these are full, independent instances (though any complex objects in the value field of Annotation objects are also just references).

  • 1: taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon instances: these are references.

  • 2: Exhaustive deep-copy: all objects are cloned.

coalescence_intervals()[source]

Returns list of coalescence intervals of self., i.e., the waiting times between successive coalescence events.

collapse_basal_bifurcation(set_as_unrooted_tree=True)[source]

Converts a degree-2 node at the root to a degree-3 node.

collapse_unweighted_edges(threshold=1e-07, update_bipartitions=False)[source]

Collapse all internal edges with edge lengths less than or equal to threshold (or with None for edge length).

colless_tree_imbalance(normalize='max')[source]

DEPRECATED: Use ‘dendropy.calculate.treemeasure.colless_tree_imbalance()’.

copy_annotations_from(other, attribute_object_mapper=None)

Copies annotations from other, which must be of Annotable type.

Copies are deep-copies, in that the Annotation objects added to the annotation_set AnnotationSet collection of self are independent copies of those in the annotate_set collection of other. However, dynamic bound-attribute annotations retain references to the original objects as given in other, which may or may not be desirable. This is handled by updated the objects to which attributes are bound via mappings found in attribute_object_mapper. In dynamic bound-attribute annotations, the _value attribute of the annotations object (Annotation._value) is a tuple consisting of “(obj, attr_name)”, which instructs the Annotation object to return “getattr(obj, attr_name)” (via: “getattr(*self._value)”) when returning the value of the Annotation. “obj” is typically the object to which the AnnotationSet belongs (i.e., self). When a copy of Annotation is created, the object reference given in the first element of the _value tuple of dynamic bound-attribute annotations are unchanged, unless the id of the object reference is fo

Parameters:
  • other (Annotable) – Source of annotations to copy.

  • attribute_object_mapper (dict) – Like the memo of __deepcopy__, maps object id’s to objects. The purpose of this is to update the parent or owner objects of dynamic attribute annotations. If a dynamic attribute Annotation gives object x as the parent or owner of the attribute (that is, the first element of the Annotation._value tuple is other) and id(x) is found in attribute_object_mapper, then in the copy the owner of the attribute is changed to attribute_object_mapper[id(x)]. If attribute_object_mapper is None (default), then the following mapping is automatically inserted: id(other): self. That is, any references to other in any Annotation object will be remapped to self. If really no reattribution mappings are desired, then an empty dictionary should be passed instead.

deep_copy_annotations_from(other, memo=None)

Note that all references to other in any annotation value (and sub-annotation, and sub-sub-sub-annotation, etc.) will be replaced with references to self. This may not always make sense (i.e., a reference to a particular entity may be absolute regardless of context).

description(depth=1, indent=0, itemize='', output=None)[source]

Returns description of object, up to level depth.

edges(filter_fn=None)[source]

Returns list of edges on tree.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be included, or False if not. If filter_fn is None (default), then all edges will be included.

Returns:

list [Edge] – List of Edge objects in self.

encode_bipartitions(suppress_unifurcations=True, collapse_unrooted_basal_bifurcation=True, suppress_storage=False, is_bipartitions_mutable=False)[source]

Calculates the bipartitions of this tree.

Parameters:
  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted as they are encountered.

  • collapse_unrooted_basal_bifurcation (bool) – If True, then a basal bifurcation on an unrooted tree will be collapsed to a trifurcation. This mean that an unrooted tree like ‘(A,(B,C))’ will be changed to ‘(A,B,C)’ after this.

  • suppress_storage (bool) – By default, the bipartition encoding is stored as a list (assigned to self.bipartition_encoding) and returned. If suppress_storage is True, then the list is not created.

  • is_bipartitions_mutable (bool) – By default, the Bipartition instances coded will be locked or frozen, allowing their use in hashing containers such as dictionary (keys) and sets. To allow modification of values, the is_mutable attribute must be set to True.

Returns:

list[|Bipartition|] or |None| – A list of Bipartition objects of this Tree representing the structure of this tree, or, if suppress_storage is True, then None.

encode_splits(*args, **kwargs)[source]

Recalculates bipartition hashes for tree.

euclidean_distance(other_tree)[source]

DEPRECATED: Use ‘dendropy.treecompare.euclidean_distance()’.

extract_tree(extraction_source_reference_attr_name='extraction_source', node_filter_fn=None, suppress_unifurcations=True, is_apply_filter_to_leaf_nodes=True, is_apply_filter_to_internal_nodes=False, tree_factory=None, node_factory=None)[source]

Returns a copy of this tree that only includes the basic structure (nodes, edges), and minimal attributes (edge lengths, node labels, and taxon associations). Annotations, comments, and other attributes are not copied.

Parameters:
  • extraction_source_reference_attr_name (str) – Name of attribute to set on cloned nodes that references corresponding original node. If None, then attribute (and reference) will not be created.

  • node_filter_fn (None or function object) – If None, then entire tree structure is cloned. If not None, must be a function object that returns True if a particular Node instance on the original tree should be included in the cloned tree, or False otherwise.

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted. Only will be done if some nodes are excluded from the cloned tree.

  • is_apply_filter_to_leaf_nodes (bool) – If True then the above filter will be applied to leaf nodes. If False then it will not (and all leaf nodes will be automatically included, unless excluded by an ancestral node being filtered out).

  • is_apply_filter_to_internal_nodes (bool) – If True then the above filter will be applied to internal nodes. If False then it will not (internal nodes without children will still be filtered out).

  • tree_factory (function) – If not None, must be a function that optionally takes a TaxonNamespace as an argument and returns a new Tree (or equivalent) instance.

  • node_factory (function) – If not None, must be a function that takes no arguments and returns a new Node (or equivalent) instance.

Examples

A simple clone:

tree0 = dendropy.Tree.get(
            path="mammals.tre",
            schema="newick")
tree1 = tree0.extract_tree()

A clone that only extracts a subtree with taxa in the genus “Rhacophorus”:

tree0 = dendropy.Tree.get(
            path="old_world_frogs.tre",
            schema="newick")
# Include taxa only if label starts with "Rhacophorus"
node_filter_fn = lambda nd: nd.is_internal() or                         nd.taxon.label.startswith("Rhacophorus")
tree1 = tree0.extract_tree(node_filter_fn=node_filter_fn)

# Above is equivalent to, but more efficient than:
#   inclusion_set = [nd.taxon for nd in tree0.leaf_node_iter()
#           if nd.taxon.label.startswith("Rhacophorus)]
#   tree1 = dendropy.Tree(tree0)
#   tree1.retain_taxa(inclusion_set)

A clone that only extracts a subtree with nodes with taxa associated with the habitat “mountain” or “forest”:

tree0 = dendropy.Tree.get(
            path="birds.tre",
            schema="newick")
include_habitats = set(["mountain", "forest"])
node_filter_fn = lambda nd: nd.taxon is None or                         nd.taxon.annotations["habitat"] in include_habitats
tree1 = tree0.extract_tree(node_filter_fn=node_filter_fn)
Returns:

t (|Tree|) – A new tree based on this one, with nodes filtered out if specified.

extract_tree_with_taxa(taxa, extraction_source_reference_attr_name='extraction_source', suppress_unifurcations=True)[source]

Returns a copy of this tree that only includes leaf nodes if they are associated with the taxon objects listed in taxa. Note that this copy will be a “thin” copy, including just the basic structure (nodes, edges) and minimal attributes (edge lengths, node labels, and taxon associations). Annotations, comments, and other attributes are not copied.

Parameters:
  • taxa (iterable of Taxon instances) – List or some other iterable of Taxon objects to include.

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted. Only will be done if some nodes are excluded from the cloned tree.

Examples

A clone that only extracts a subtree with taxa in the genus “Rhacophorus”:

tree0 = dendropy.Tree.get(
            path="old_world_frogs.tre",
            schema="newick")
# Include taxa only if label starts with "Rhacophorus"
taxa_to_retain = set([taxon for taxon in tree0.taxon_namespace
        if taxon.label.startswith("Rhacophorus")])
tree1 = tree0.extract_tree_with_taxa(taxa=taxa_to_retain)

# Above is equivalent to, but more efficient than:
#   inclusion_set = [nd.taxon for nd in tree0.leaf_node_iter()
#           if nd.taxon.label.startswith("Rhacophorus)]
#   tree1 = dendropy.Tree(tree0)
#   tree1.retain_taxa(inclusion_set)
Returns:

t (|Tree|) – A new tree based on this one, with nodes filtered out if specified.

extract_tree_with_taxa_labels(labels, extraction_source_reference_attr_name='extraction_source', suppress_unifurcations=True)[source]

Returns a copy of this tree that only includes leaf nodes if they are associated with taxon objects with labels matching those listed in labels. Note that this copy will be a “thin” copy, including just the basic structure (nodes, edges) and minimal attributes (edge lengths, node labels, and taxon associations). Annotations, comments, and other attributes are not copied.

Parameters:
  • labels (iterable of str instances) – List or some other iterable of strings to match.

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted. Only will be done if some nodes are excluded from the cloned tree.

Examples

A clone that only extracts a subtree with taxa in the genus “Rhacophorus”:

tree0 = dendropy.Tree.get(
            path="old_world_frogs.tre",
            schema="newick")
# Include taxa only if label starts with "Rhacophorus"
labels = set([taxon.label for taxon in tree0.taxon_namespace
        if taxon.label.startswith("Rhacophorus")])
tree1 = tree0.extract_tree_with_taxa_labels(labels=labels)

# Above is equivalent to, but more efficient than:
#   inclusion_set = [nd.taxon for nd in tree0.leaf_node_iter()
#           if nd.taxon.label.startswith("Rhacophorus)]
#   tree1 = dendropy.Tree(tree0)
#   tree1.retain_taxa(inclusion_set)
Returns:

t (|Tree|) – A new tree based on this one, with nodes filtered out if specified.

extract_tree_without_taxa(taxa, extraction_source_reference_attr_name='extraction_source', suppress_unifurcations=True)[source]

Returns a copy of this tree that only includes leaf nodes if they are NOT associated with the taxon objects listed in taxa. Note that this copy will be a “thin” copy, including just the basic structure (nodes, edges) and minimal attributes (edge lengths, node labels, and taxon associations). Annotations, comments, and other attributes are not copied.

Parameters:
  • taxa (iterable of Taxon instances) – List or some other iterable of Taxon objects to exclude.

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted. Only will be done if some nodes are excluded from the cloned tree.

  • is_apply_filter_to_leaf_nodes (bool) – If True then the above filter will be applied to leaf nodes. If False then it will not (and all leaf nodes will be automatically included, unless excluded by an ancestral node being filtered out).

  • is_apply_filter_to_internal_nodes (bool) – If True then the above filter will be applied to internal nodes. If False then it will not (internal nodes without children will still be filtered out).

Examples

A clone that only extracts a subtree with taxa NOT in the genus “Rhacophorus”:

tree0 = dendropy.Tree.get(
            path="old_world_frogs.tre",
            schema="newick")
# Exclude taxa if their name starts with "Rhacophorus"
taxa_to_exclude = set([taxon for taxon in tree0.taxon_namespace
        if taxon.label.startswith("Rhacophorus")])
tree1 = tree0.extract_tree_with_taxa(taxa=taxa_to_retain)

# Above is equivalent to, but more efficient than:
#   inclusion_set = [nd.taxon for nd in tree0.leaf_node_iter()
#           if nd.taxon.label.startswith("Rhacophorus)]
#   tree1 = dendropy.Tree(tree0)
#   tree1.retain_taxa(inclusion_set)
Returns:

t (|Tree|) – A new tree based on this one, with nodes filtered out if specified.

extract_tree_without_taxa_labels(labels, extraction_source_reference_attr_name='extraction_source', suppress_unifurcations=True)[source]

Returns a copy of this tree that only includes leaf nodes if they are NOT associated with the taxon objects listed in taxa. Note that this copy will be a “thin” copy, including just the basic structure (nodes, edges) and minimal attributes (edge lengths, node labels, and taxon associations). Annotations, comments, and other attributes are not copied.

Parameters:
  • labels (iterable of str instances) – List or some other iterable of strings to match.

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted. Only will be done if some nodes are excluded from the cloned tree.

  • is_apply_filter_to_leaf_nodes (bool) – If True then the above filter will be applied to leaf nodes. If False then it will not (and all leaf nodes will be automatically included, unless excluded by an ancestral node being filtered out).

  • is_apply_filter_to_internal_nodes (bool) – If True then the above filter will be applied to internal nodes. If False then it will not (internal nodes without children will still be filtered out).

Examples

A clone that only extracts a subtree with taxa NOT in the genus “Rhacophorus”:

tree0 = dendropy.Tree.get(
            path="old_world_frogs.tre",
            schema="newick")
# Exclude taxa if label starts with "Rhacophorus"
labels = set([taxon.label for taxon in tree0.taxon_namespace
        if taxon.label.startswith("Rhacophorus")])
tree1 = tree0.extract_tree_without_taxa_labels(labels=labels)

# Above is equivalent to, but more efficient than:
#   inclusion_set = [nd.taxon for nd in tree0.leaf_node_iter()
#           if nd.taxon.label.startswith("Rhacophorus)]
#   tree1 = dendropy.Tree(tree0)
#   tree1.prune_taxa(inclusion_set)
Returns:

t (|Tree|) – A new tree based on this one, with nodes filtered out if specified.

false_positives_and_negatives(other_tree)[source]

DEPRECATED: Use ‘dendropy.treecompare.false_positives_and_negatives()’.

filter_leaf_nodes(filter_fn, recursive=True, update_bipartitions=False, suppress_unifurcations=True)[source]

Removes all leaves for which filter_fn returns False. If recursive is True, then process is repeated until all leaf nodes in the tree will evaluate to True when passed to filter_fn.

Parameters:
  • filter_fn (function object) – A function that takes a Node object and returns True if the object is to be allowed as a leaf node, and False if otherwise.

  • recursive (bool) – If True, then filter is repeatedly applied until all leaf nodes evaluate to True under filter_fn. If False, then only a single pass is made on the current leaf set. This may result in new leaves for which the filter_fn is False (e.g., the parent node of a cherry in which both children evaluated to False under filter_fn now is a leaf node which may be False under filter_fn).

  • suppress_unifurcations (bool) – If True, nodes of outdegree 1 will be deleted as they are encountered.

  • update_bipartitions (bool) – If True, then bipartitions will be calculated.

Returns:

nds (list[|Node|]) – List of nodes removed.

find_missing_splits(other_tree)[source]

DEPRECATED: Use ‘dendropy.treecompare.find_missing_bipartitions()’.

find_node(filter_fn)[source]

Finds the first node for which filter_fn(node) == True.

For example, if:

filter_fn = lambda n: hasattr(n, 'genes') and n.genes is not None

then:

node = t.find_node(filter_fn=filter_fn)

will return the first node which has the attribute ‘genes’ and this value is not None.

Parameters:

filter_fn (function object) – Takes a single Node object as an argument and returns True if the node should be returned.

Returns:

|Node| or |None| – Returns first Node object for which the filter function filter_fn returns True, or None if no such node exists on this tree.

find_node_for_taxon(taxon)[source]

Returns node associated with Taxon object taxon.

Parameters:

taxon (Taxon object) – Taxon object that should be associated with the node to be returned.

Returns:

|Node| or |None| – Returns first Node object with taxon attribute referencing same object as taxon argument, or None if no such node exists.

find_node_with_label(label)[source]

Returns first node with label attribute matching label argument.

Parameters:

label (string) – Value for label attribute of Node object in this tree.

Returns:

|Node| or |None| – Returns first Node object with label attribute having value given in label, or None if no such node is found.

find_node_with_taxon(taxon_filter_fn=None)[source]

Returns node associated with Taxon object for which taxon_filter_fn returns True.

Parameters:

taxon_filter_fn (function object) – Takes a single Taxon object as an argument and returns True if the node associated with that Taxon should be returned.

Returns:

|Node| or |None| – Returns first Node object with taxon attribute passing filter function taxon_filter_fn, or None if no such node is found.

find_node_with_taxon_label(label)[source]

Returns node associated with Taxon object with the specified label.

Parameters:

label (string) – Label of Taxon object associated with the node to be returned.

Returns:

|Node| or |None| – Returns first Node object with taxon attribute having label label, or|None| if no such node is found.

find_nodes(filter_fn)[source]

Find all nodes for which filter_fn(node) == True.

For example, if:

filter_fn = lambda n: hasattr(n, 'genes') and n.genes is not None

then:

nodes = t.find_node(filter_fn=filter_fn)

will return all nodes which have the attribute ‘genes’ and this value is not None.

Parameters:

filter_fn (function object) – Takes a single Node object as an argument and returns True if the node should be returned.

Returns:

nodes (list of |Node| instances) – Returns list of Node objects for which the filter function filter_fn returns True.

classmethod from_bipartition_encoding(bipartition_encoding, taxon_namespace, is_rooted=False, edge_lengths=None)[source]

Reconstructs a tree from a bipartition encoding.

Parameters:
  • bipartition_encoding (iterable[Bipartition]) – An iterable of Bipartition instances representing a tree. Bipartitions will be added to the tree in the order given by iterating over this. Bipartitions that are incompatible with the tree will be skipped. So if not all bipartitions are compatible with each other, then the sequence of bipartitions given in bipartition_encoding should be in order of their support values or some other preference criteria.

  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to use to manage taxon definitions.

  • is_rooted (bool) – Specifies whether or not the tree is rooted.

  • edge_lengths (iterable or None) – An iterable of edge lengths. This should be in the same order as the bipartitions in the bipartition encoding.

Returns:

|Tree| – The tree reconstructed from the given bipartition encoding.

classmethod from_split_bitmasks(split_bitmasks, taxon_namespace, is_rooted=False, split_edge_lengths=None)[source]

Reconstructs a tree from a collection of splits represented as bitmasks.

Parameters:
  • split_bitmasks (iterable[int]) – An iterable of split bitmasks representing a tree. Splits will be added to the tree in the order given by iterating over this. Splits that are incompatible with the tree will be skipped. So if not all splits are compatible with each other, then the sequence of splits given in bipartition_encoding should be in order of their support values or some other preference criteria.

  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to use to manage taxon definitions.

  • is_rooted (bool) – Specifies whether or not the tree is rooted.

  • edge_lengths (dict or False or None) – If False or None, then no edge lengths will be added. Otherwise, this should be a dictionary mapping splits to edge lengths.

Returns:

|Tree| – The tree reconstructed from the given bipartition encoding.

classmethod get(**kwargs)[source]

Instantiate and return a new Tree object from a data source.

Mandatory Source-Specification Keyword Argument (Exactly One of the Following Required):

  • file (file) – File or file-like object of data opened for reading.

  • path (str) – Path to file of data.

  • url (str) – URL of data.

  • data (str) – Data given directly.

Mandatory Schema-Specification Keyword Argument:

Optional General Keyword Arguments:

  • label (str) – Name or identifier to be assigned to the new object; if not given, will be assigned the one specified in the data source, or None otherwise.

  • taxon_namespace (TaxonNamespace) – The TaxonNamespace instance to use to manage the taxon names. If not specified, a new one will be created.

  • collection_offset (int) – 0-based index of tree block or collection in source to be parsed. If not specified then the first collection (offset = 0) is assumed.

  • tree_offset (int) – 0-based index of tree within the collection specified by collection_offset to be parsed. If not specified, then the first tree (offset = 0) is assumed.

  • ignore_unrecognized_keyword_arguments (bool) – If True, then unsupported or unrecognized keyword arguments will not result in an error. Default is False: unsupported keyword arguments will result in an error.

Optional Schema-Specific Keyword Arguments:

These provide control over how the data is interpreted and processed, and supported argument names and values depend on the schema as specified by the value passed as the “schema” argument. See “DendroPy Schemas: Phylogenetic and Evolutionary Biology Data Formats” for more details.

Examples:

# From a URL
t1 = dendropy.Tree.get(
        url="http://api.opentreeoflife.org/v2/study/pg_1144/tree/tree2324.nex",
        schema="nexus")

# From a file-like object
t2 = Tree.get(file=open('treefile.tre', 'r'),
                schema="newick",
                tree_offset=0)

# From a path
t3 = Tree.get(path='sometrees.nexus',
        schema="nexus",
        collection_offset=2,
        tree_offset=1)

# From a string
s = "((A,B),(C,D));((A,C),(B,D));"
# tree will be '((A,B),(C,D))'
t4 = Tree.get(data=s,
        schema="newick")
# tree will be '((A,C),(B,D))'
t5 = Tree.get(data=s,
        schema="newick",
        tree_offset=1)
# passing keywords to underlying tree parser
t7 = dendropy.Tree.get(
        data="((A,B),(C,D));",
        schema="newick",
        taxon_namespace=t3.taxon_namespace,
        suppress_internal_node_taxa=False,
        preserve_underscores=True)
classmethod get_from_path(src, schema, **kwargs)

Factory method to return new object of this class from file specified by string src.

Parameters:
  • src (string) – Full file path to source of data.

  • schema (string) – Specification of data format (e.g., “nexus”).

  • kwargs (keyword arguments, optional) – Arguments to customize parsing, instantiation, processing, and accession of objects read from the data source, including schema- or format-specific handling. These will be passed to the underlying schema-specific reader for handling.

Returns:

pdo (phylogenetic data object) – New instance of object, constructed and populated from data given in source.

classmethod get_from_stream(src, schema, **kwargs)

Factory method to return new object of this class from file-like object src.

Parameters:
  • src (file or file-like) – Source of data.

  • schema (string) – Specification of data format (e.g., “nexus”).

  • kwargs (keyword arguments, optional) – Arguments to customize parsing, instantiation, processing, and accession of objects read from the data source, including schema- or format-specific handling. These will be passed to the underlying schema-specific reader for handling.

Returns:

pdo (phylogenetic data object) – New instance of object, constructed and populated from data given in source.

classmethod get_from_string(src, schema, **kwargs)

Factory method to return new object of this class from string src.

Parameters:
  • src (string) – Data as a string.

  • schema (string) – Specification of data format (e.g., “nexus”).

  • kwargs (keyword arguments, optional) – Arguments to customize parsing, instantiation, processing, and accession of objects read from the data source, including schema- or format-specific handling. These will be passed to the underlying schema-specific reader for handling.

Returns:

pdo (phylogenetic data object) – New instance of object, constructed and populated from data given in source.

classmethod get_from_url(src, schema, strip_markup=False, **kwargs)

Factory method to return a new object of this class from URL given by src.

Parameters:
  • src (string) – URL of location providing source of data.

  • schema (string) – Specification of data format (e.g., “nexus”).

  • kwargs (keyword arguments, optional) – Arguments to customize parsing, instantiation, processing, and accession of objects read from the data source, including schema- or format-specific handling. These will be passed to the underlying schema-specific reader for handling.

Returns:

pdo (phylogenetic data object) – New instance of object, constructed and populated from data given in source.

infer_taxa()[source]

Creates (and returns) a new TaxonNamespace object for self populated with taxa from this tree.

inorder_edge_iter(filter_fn=None)[source]

In-order iteration over edges of tree.

Visits edges in self, with each edge visited in-between its children. Only valid for strictly-bifurcating trees. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when called with the edge as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

Returns:

collections.Iterator [Edge] – An iterator yielding edges of self in infix or in-order sequence.

inorder_node_iter(filter_fn=None)[source]

In-order iteration over nodes of tree.

Visits nodes in self, with each node visited in-between its children. Only valid for strictly-bifurcating trees. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes of self in infix or in-order sequence.

internal_edges(exclude_seed_edge=False)[source]

Returns list of internal edges on tree.

Parameters:

exclude_seed_node (boolean, optional) – If False (default), then the edge subtending the seed node or root is included. If True, then the seed node is omitted.

Returns:

list [Edge] – List of internal Edge objects in self.

internal_node_ages(ultrametricity_precision=1e-05, is_force_max_age=False, is_force_min_age=False, set_node_age_fn=None)[source]

Returns list of ages of speciation events / coalescence times on tree.

internal_nodes(exclude_seed_node=False)[source]

Returns list of internal nodes in the tree.

Root or seed node is included unless exclude_seed_node is True.

Parameters:

exclude_seed_node (boolean, optional) – If False (default), then the seed node or root is included. If True, then the seed node is omitted.

Returns:

list [Node] – List of internal Node objects in self.

is_compatible_with_bipartition(bipartition, is_bipartitions_updated=False)[source]

Returns true if the Bipartition bipartition is compatible with all bipartitions of this tree.

ladderize(ascending=True)[source]

Sorts child nodes in ascending (if ascending is True) or descending (if ascending is False) order in terms of the number of children each child node has.

Ladderize sort is stable. To control order between nodes with same child count, call reorder prior to ladderization.

leaf_edge_iter(filter_fn=None)[source]

Iterate over all tips or leaves of tree.

Visits all leaf or tip in self. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when called with the edge as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

Returns:

collections.Iterator [Edge] – An iterator yielding leaf edges in self.

leaf_edges()[source]

Returns list of leaf edges on the tree.

Returns:

list [Edge] – List of leaf Edge objects in self.

leaf_iter(filter_fn=None)[source]

Deprecated: use Tree.leaf_node_iter instead.

leaf_node_iter(filter_fn=None)[source]

Iterate over all tips or leaves of tree.

Visits all leaf or tip in self. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding leaf nodes in self.

leaf_nodes()[source]

Returns list of leaf nodes on the tree.

Returns:

list [Node] – List of leaf Node objects in self.

length()[source]

Returns sum of edge lengths of self. Edges with no lengths defined (None) will be considered to have a length of 0. Note that we do not overrride __len__ as this requires an integer return value.

level_order_edge_iter(filter_fn=None)[source]

Deprecated: use Tree.levelorder_edge_iter instead.

level_order_node_iter(filter_fn=None)[source]

Deprecated: use Tree.levelorder_node_iter instead.

levelorder_edge_iter(filter_fn=None)[source]

Level-order iteration over edges of tree.

Visits edges in self, with each edge and other edges at the same level (distance from root) visited before their children. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when called with the edge as an argument are visited.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

Returns:

collections.Iterator [Edge] – An iterator yielding edges of self in level-order sequence.

levelorder_node_iter(filter_fn=None)[source]

Level-order iteration over nodes of tree.

Visits nodes in self, with each node and other nodes at the same level (distance from root) visited before their children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are visited.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes of self in level-order sequence.

max_distance_from_root()[source]

Returns distance of node furthest from root.

migrate_taxon_namespace(taxon_namespace, unify_taxa_by_label=True, taxon_mapping_memo=None)

Move this object and all members to a new operational taxonomic unit concept namespace scope.

Current self.taxon_namespace value will be replaced with value given in taxon_namespace if this is not None, or a new TaxonNamespace object. Following this, reconstruct_taxon_namespace() will be called: each distinct Taxon object associated with self or members of self that is not alread in taxon_namespace will be replaced with a new Taxon object that will be created with the same label and added to self.taxon_namespace. Calling this method results in the object (and all its member objects) being associated with a new, independent taxon namespace.

Label mapping case sensitivity follows the self.taxon_namespace.is_case_sensitive setting. If False and unify_taxa_by_label is also True, then the establishment of correspondence between Taxon objects in the old and new namespaces with be based on case-insensitive matching of labels. E.g., if there are four Taxon objects with labels ‘Foo’, ‘Foo’, ‘FOO’, and ‘FoO’ in the old namespace, then all objects that reference these will reference a single new Taxon object in the new namespace (with a label some existing casing variant of ‘foo’). If True: if unify_taxa_by_label is True, Taxon objects with labels identical except in case will be considered distinct.

Parameters:
  • taxon_namespace (TaxonNamespace) – The TaxonNamespace into the scope of which this object will be moved.

  • unify_taxa_by_label (boolean, optional) – If True, then references to distinct Taxon objects with identical labels in the current namespace will be replaced with a reference to a single Taxon object in the new namespace. If False: references to distinct Taxon objects will remain distinct, even if the labels are the same.

  • taxon_mapping_memo (dictionary) – Similar to memo of deepcopy, this is a dictionary that maps Taxon objects in the old namespace to corresponding Taxon objects in the new namespace. Mostly for interal use when migrating complex data to a new namespace. Note that any mappings here take precedence over all other options: if a Taxon object in the old namespace is found in this dictionary, the counterpart in the new namespace will be whatever value is mapped, regardless of, e.g. label values.

Examples

Use this method to move an object from one taxon namespace to another.

For example, to get a copy of an object associated with another taxon namespace and associate it with a different namespace:

# Get handle to the new TaxonNamespace
other_taxon_namespace = some_other_data.taxon_namespace

# Get a taxon-namespace scoped copy of a tree
# in another namespace
t2 = Tree(t1)

# Replace taxon namespace of copy
t2.migrate_taxon_namespace(other_taxon_namespace)

You can also use this method to get a copy of a structure and then move it to a new namespace:

t2 = Tree(t1) t2.migrate_taxon_namespace(TaxonNamespace())

# Note: the same effect can be achived by: t3 = copy.deepcopy(t1)

minmax_leaf_distance_from_root()[source]

Returns pair of values, representing the distance of the leaf closest to a furthest from the root.

mrca(**kwargs)[source]

Returns most-recent common ancestor node of a set of taxa on the tree.

Returns the shallowest node in the tree (the node nearest the tips) that has all of the taxa that:

  • are specified by the leafset bitmask given by the keyword argument leafset_bitmask

  • are in the list of Taxon objects given by the keyword argument taxa

  • have the labels specified by the list of strings given by the keyword argument taxon_labels

Returns None if no appropriate node is found. Assumes that bipartitions have been encoded on the tree. It is possible that the leafset bitmask is not compatible with the subtree that is returned! (compatibility tests are not fully performed). This function is used to find the “insertion point” for a new bipartition via a root to tip search.

Parameters:

**kwargs (keyword arguments) –

Exactly one of the following must be specified:

leafset_bitmaskinteger

Node object subtended by the first edge compatible with this leafset bitmask will be returned.

taxacollections.Iterable [Taxon]

Shallowest node object with descendent nodes associated with all the Taxon objects specified will be returned.

taxon_labelscollections.Iterable [string]

Shallowest node object with descendent nodes associated with the minimal set of Taxon objects that collectively have all the labels specified in taxon_labels will be returned.

In addition, the following optional keywords are supported:

start_nodeNode, optional

If given, specifies the node at which to start searching. If not, defaults to the root or seed_node.

Returns:

|Node| or |None| – The most-recent common ancestor of the nodes specified, or None if no such node exists.

node_ages(ultrametricity_precision=1e-05, is_force_max_age=False, is_force_min_age=False, set_node_age_fn=None, internal_only=False)[source]

Returns list of ages of all nodes on tree. NOTE: Changed from DendroPy3: this function now returns the ages of ALL nodes. To get only internal node ages, use Tree.internal_node_ages.

classmethod node_factory(**kwargs)[source]

Creates and returns a Node object.

Derived classes can override this method to provide support for specialized or different types of nodes on the tree.

Parameters:

**kwargs (keyword arguments) – Passed directly to constructor of Node.

Returns:

|Node| – A new Node object.

nodes(filter_fn=None)[source]

Returns list of nodes on tree.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be included in the list, or False if not. If filter_fn is None (default), then all nodes visited will be included.

Returns:

list [Node] – List of Node objects in the tree.

num_lineages_at(distance_from_root)[source]

Returns the number of lineages on the tree at a particular distance from the root.

phylogenetic_distance_matrix(*args, **kwargs)[source]

Returns a PhylogeneticDistanceMatrix instance based on the tree (in its current state).

Returns:

pdc (a |PhylogeneticDistanceMatrix| instance) – A PhylogeneticDistanceMatrix instance corresponding to the tree in its current state.

poll_taxa(taxa=None)[source]

Returns a set populated with all of Taxon instances associated with self.

Parameters:

taxa (set()) – Set to populate. If not specified, a new one will be created.

Returns:

set[|Taxon|] – Set of taxa associated with self.

postorder_edge_iter(filter_fn=None)[source]

Post-order iterator over edges of tree.

Visits edges in self, with each edge visited after its children. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when called with the edge as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

Returns:

collections.Iterator [Edge] – An iterator yielding the edges in self in post-order sequence.

postorder_internal_edge_iter(filter_fn=None, exclude_seed_edge=False)[source]

Pre-order iterator over internal edges tree.

Visits internal edges in self, with each edge visited after its children. In DendroPy, “internal edges” are edges that have at least one child edge, and thus the root or seed edge is typically included unless exclude_seed_edge is True. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when passed the edge as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

  • exclude_seed_edge (boolean, optional) – If False (default), then the seed edge or root is visited. If True, then the seed edge is skipped.

Returns:

collections.Iterator [Edge] – An iterator yielding the internal edges of self in post-order sequence.

postorder_internal_node_iter(filter_fn=None, exclude_seed_node=False)[source]

Pre-order iterator over internal nodes tree.

Visits internal nodes in self, with each node visited after its children. In DendroPy, “internal nodes” are nodes that have at least one child node, and thus the root or seed node is typically included unless exclude_seed_node is True. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when passed the node as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • exclude_seed_node (boolean, optional) – If False (default), then the seed node or root is visited. If True, then the seed node is skipped.

Returns:

collections.Iterator [Node] – An iterator yielding the internal nodes of self in post-order sequence.

postorder_node_iter(filter_fn=None)[source]

Post-order iterator over nodes of tree.

Visits self and all descendant nodes, with each node visited after its children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding the nodes in self in post-order sequence.

preorder_edge_iter(filter_fn=None)[source]

Pre-order iterator over nodes in tree.

Visits nodes in self, with each node visited before its children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes in self in pre-order sequence.

preorder_internal_edge_iter(filter_fn=None, exclude_seed_edge=False)[source]

Pre-order iterator over internal edges in tree.

Visits internal edges in self, with each edge visited before its children. In DendroPy, “internal edges” are edges that have at least one child edge, and thus the root or seed edge is typically included unless exclude_seed_edge is True. Edges can optionally be filtered by filter_fn: only edges for which filter_fn returns True when passed the edge as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

  • exclude_seed_edge (boolean, optional) – If False (default), then the edge subtending the seed node or root is visited. If True, then this edge is skipped.

Returns:

collections.Iterator [Edge] – An iterator yielding the internal edges of self.

preorder_internal_node_iter(filter_fn=None, exclude_seed_node=False)[source]

Pre-order iterator over internal nodes in tree.

Visits internal nodes in self, with each node visited before its children. In DendroPy, “internal nodes” are nodes that have at least one child node, and thus the root or seed node is typically included unless exclude_seed_node is True. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when passed the node as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • exclude_seed_node (boolean, optional) – If False (default), then the seed node or root is visited. If True, then the seed node is skipped.

Returns:

collections.Iterator [Node] – An iterator yielding the internal nodes of self.

preorder_node_iter(filter_fn=None)[source]

Pre-order iterator over nodes in tree.

Visits nodes in self, with each node visited before its children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes in self in pre-order sequence.

print_plot(**kwargs)[source]

Writes an ASCII text graphic of this tree to standard output.

prune_leaves_without_taxa(recursive=True, update_bipartitions=False, suppress_unifurcations=True)[source]

Removes all terminal nodes that have their taxon attribute set to None.

prune_subtree(node, update_bipartitions=False, suppress_unifurcations=True)[source]

Removes subtree starting at node from tree.

prune_taxa(taxa, update_bipartitions=False, suppress_unifurcations=True, is_apply_filter_to_leaf_nodes=True, is_apply_filter_to_internal_nodes=False)[source]

Removes terminal nodes associated with Taxon objects given by the container taxa (which can be any iterable, including a TaxonNamespace object) from self.

prune_taxa_with_labels(labels, update_bipartitions=False, suppress_unifurcations=True, is_apply_filter_to_leaf_nodes=True, is_apply_filter_to_internal_nodes=False)[source]

Removes terminal nodes that are associated with Taxon objects with labels given by labels.

purge_taxon_namespace()

Remove all Taxon instances in self.taxon_namespace that are not associated with self or any item in self.

pybus_harvey_gamma(prec=1e-05)[source]

DEPRECATED: Use ‘dendropy.calculate.treemeasure.pybus_harvey_gamma()’.

randomly_assign_taxa(create_required_taxa=True, rng=None)[source]

Randomly assigns taxa to leaf nodes. If the number of taxa defined in the taxon set of the tree is more than the number of tips, then a random subset of taxa in taxon_namespace will be assigned to the tips of tree. If the number of tips is more than the number of taxa in the taxon_namespace, and add_extra_taxa is not True [default], then new Taxon objects will be created and added to the taxon_namespace; if create_required_taxa is False, then an exception is raised.

In addition, a Random() object or equivalent can be passed using rng; otherwise GLOBAL_RNG is used.

randomly_reorient(rng=None, update_bipartitions=False)[source]

Randomly picks a new rooting position and rotates the branches around all internal nodes in the self. If update_bipartitions is True, the the bipartition_bitmask and bipartition_edge_map attributes kept valid.

randomly_rotate(rng=None)[source]

Randomly rotates the branches around all internal nodes in self

reconstruct_taxon_namespace(unify_taxa_by_label=True, taxon_mapping_memo=None)[source]

Repopulates the current taxon namespace with new taxon objects, preserving labels. Each distinct Taxon object associated with self or members of self that is not already in self.taxon_namespace will be replaced with a new Taxon object that will be created with the same label and added to self.taxon_namespace.

Label mapping case sensitivity follows the self.taxon_namespace.is_case_sensitive setting. If False and unify_taxa_by_label is also True, then the establishment of correspondence between Taxon objects in the old and new namespaces with be based on case-insensitive matching of labels. E.g., if there are four Taxon objects with labels ‘Foo’, ‘Foo’, ‘FOO’, and ‘FoO’ in the old namespace, then all objects that reference these will reference a single new Taxon object in the new namespace (with a label some existing casing variant of ‘foo’). If True: if unify_taxa_by_label is True, Taxon objects with labels identical except in case will be considered distinct.

Note

Existing Taxon objects in self.taxon_namespace are not removed. This method should thus only be called only when self.taxon_namespace has been changed. In fact, typical usage would not involve calling this method directly, but rather through

Parameters:
  • unify_taxa_by_label (boolean, optional) – If True, then references to distinct Taxon objects with identical labels in the current namespace will be replaced with a reference to a single Taxon object in the new namespace. If False: references to distinct Taxon objects will remain distinct, even if the labels are the same.

  • taxon_mapping_memo (dictionary) – Similar to memo of deepcopy, this is a dictionary that maps Taxon objects in the old namespace to corresponding Taxon objects in the new namespace. Mostly for interal use when migrating complex data to a new namespace.

reindex_subcomponent_taxa()[source]

Remaps node taxon objects

reindex_taxa(taxon_namespace=None, clear=False)

DEPRECATED: Use migrate_taxon_namespace() instead. Rebuilds taxon_namespace from scratch, or assigns Taxon objects from given TaxonNamespace object taxon_namespace based on label values.

reorder(ascending=True, key=<function Tree.<lambda>>)[source]

Reorder the children of each node in the tree, by default in ascending order by Taxon with missing taxa treated as labeled as empty string. Does not alter tree topology.

Specify key to sort by a different attribute or function of nodes.

reroot_at_edge(edge, length1=None, length2=None, update_bipartitions=False, suppress_unifurcations=True)[source]

Takes an internal edge, edge, adds a new node to it, and then roots the tree on the new node. length1 will be assigned to the new (sub-)edge leading to the original parent of the original edge. length2 will be assigned to the new (sub-)edge leading to the original child of the original edge. If update_bipartitions is True, then the edges’ bipartition and the tree’s bipartition_encoding attributes will be updated. If the old root of the tree had an outdegree of 2, then after this operation, it will have an outdegree of one. In this case, unless suppress_unifurcations is False, then it will be removed from the tree.

reroot_at_midpoint(update_bipartitions=False, suppress_unifurcations=True, collapse_unrooted_basal_bifurcation=True)[source]

Reroots the tree at the the mid-point of the longest distance between two taxa in a tree. Sets the rooted flag on the tree to True. If update_bipartitions is True, then the edges’ bipartition and the tree’s bipartition_encoding attributes will be updated. If the old root of the tree had an outdegree of 2, then after this operation, it will have an outdegree of one. In this case, unless suppress_unifurcations is False, then it will be removed from the tree.

reroot_at_node(new_root_node, update_bipartitions=False, suppress_unifurcations=True, collapse_unrooted_basal_bifurcation=True)[source]

Takes an internal node, new_seed_node that must already be in the tree and roots the tree at that node. This is a ‘hard’ rerooting – i.e., changes the tree representation so tree traversal behaves as if the tree is rooted at ‘new_seed_node’, and changes the tree’s rooting state. If update_bipartitions is True, then the edges’ bipartition and the tree’s bipartition_encoding attributes will be updated. If the old root of the tree had an outdegree of 2, then after this operation, it will have an outdegree of one. In this case, unless suppress_unifurcations is False, then it will be removed from the tree.

reseed_at(new_seed_node, update_bipartitions=False, collapse_unrooted_basal_bifurcation=True, suppress_unifurcations=True)[source]

Reseeds the tree at a different (existing) node.

Takes an internal node, new_seed_node that must already be in the tree and rotates the tree such that new_seed_node is the seed_node of the tree. This is a ‘soft’ rerooting – i.e., changes the tree representation so tree traversal behaves as if the tree is rooted at ‘new_seed_node’, but it does not actually change the tree’s rooting state. If update_bipartitions is True, then the edges’ bipartition_bitmask and the tree’s bipartition_edge_map attributes will be updated. If the old root of the tree had an outdegree of 2, then after this operation, it will have an outdegree of one. In this case, unless suppress_unifurcations is False, then it will be removed from the tree.

resolve_node_ages(node_callback_fn=None, node_edge_length_fn=None, attr_name='age')[source]

Adds an attribute called “age” to each node, with the value equal to the time elapsed since the present.

This is calculated by: (1) setting the age of the root node to the sum of path lengths to the most distant tip (2) setting the age of each other node as the sum of path lengths from the root.

Unlike the (legacy) calc_node_ages() there is no ultrametricity requirement or check.

resolve_node_depths(node_callback_fn=None, node_edge_length_fn=None, attr_name='depth')[source]

Adds an attribute given by attr_name to each node, with the value equal to the sum of edge lengths from the root.

resolve_polytomies(limit=2, update_bipartitions=False, rng=None)[source]

Arbitrarily resolve polytomies using 0-length edges.

Parameters:
  • limit (int) – The maximum number of children a node can have before being resolved.

  • update_bipartitions (bool) – If True, then bipartitions will be calculated.

  • rng (random.Random object or None) – If rng is an object with a sample() method then the polytomy will be resolved by sequentially adding, generating all tree topologies equiprobably. rng.sample() should behave like random.sample() If rng is None, then polytomy is broken deterministically by repeatedly joining pairs of children.

retain_taxa(taxa, update_bipartitions=False, suppress_unifurcations=True)[source]

Removes terminal nodes that are not associated with any of the Taxon objects given by taxa (which can be any iterable, including a TaxonNamespace object) from the self.

retain_taxa_with_labels(labels, update_bipartitions=False, suppress_unifurcations=True)[source]

Removes terminal nodes that are not associated with Taxon objects with labels given by labels.

robinson_foulds_distance(other_tree)[source]

DEPRECATED: Use ‘dendropy.treecompare.weighted_robinson_foulds_distance()’.

sackin_index(normalize=True)[source]

DEPRECATED: Use ‘dendropy.calculate.treemeasure.sackin_index()’.

scale_edges(edge_len_multiplier)[source]

Multiplies every edge length in self by edge_len_multiplier

set_edge_lengths_from_node_ages(minimum_edge_length=0.0, error_on_negative_edge_lengths=False)[source]

Sets the edge lengths of the tree so that the path lengths from the tips equal the value of the age attribute of the nodes.

Parameters:
  • minimum_edge_length (numeric) – All edge lengths calculated to have a value less than this will be set to this.

  • error_on_negative_edge_lengths (bool) – If True, an inferred edge length that is less than 0 will result in a ValueError.

shuffle_taxa(include_internal_nodes=False, rng=None)[source]

Randomly re-assigns taxa associated with nodes. Note that in the case of not all nodes being associated with taxa, this will NOT assign taxa to nodes that currently do not have them, nor will nodes currently associated with taxa end up not being associated with taxa. Returns a dictionary mapping the old taxa to their new counterparts.

strip_comments()[source]

Remove comments from tree/nodes.

suppress_unifurcations(update_bipartitions=False)[source]

Delete all nodes of outdegree-one from this tree.

Parameters:

update_bipartitions (bool) – If True then the bipartitions encoding will be calculated.

symmetric_difference(other_tree)[source]

DEPRECATED: Use ‘dendropy.treecompare.symmetric_difference()’.

taxon_namespace_scoped_copy(memo=None)[source]

Cloning level: 1. Taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon objects: these are preserved as references.

to_outgroup_position(outgroup_node, update_bipartitions=False, suppress_unifurcations=True)[source]

Reroots the tree at the parent of outgroup_node and makes outgroup_node the first child of the new root. This is just a convenience function to make it easy to place a clade as the first child under the root. Assumes that outgroup_node and outgroup_node._parent_node and are in the tree/ If update_bipartitions is True, then the edges’ bipartition and the tree’s bipartition_encoding attributes will be updated. If the old root of the tree had an outdegree of 2, then after this operation, it will have an outdegree of one. In this case, unless suppress_unifurcations is False, then it will be removed from the tree.

treeness()[source]

DEPRECATED: Use ‘dendropy.calculate.treemeasure.treeness()’.

unassign_taxa(exclude_leaves=False, exclude_internal=False)[source]

Strips taxon assignments from tree. If exclude_leaves is True, then taxa on leaves will be retained. If exclude_internal is True, then taxa on internal nodes will be retained. The taxon_namespace is not affected by this operation.

update_bipartitions(*args, **kwargs)[source]

Recalculates bipartition hashes for tree.

update_splits(*args, **kwargs)[source]

Recalculates bipartition hashes for tree.

update_taxon_namespace()[source]

All Taxon objects in self that are not in self.taxon_namespace will be added.

write(**kwargs)

Writes out self in schema format.

Mandatory Destination-Specification Keyword Argument (Exactly One of the Following Required):

  • file (file) – File or file-like object opened for writing.

  • path (str) – Path to file to which to write.

Mandatory Schema-Specification Keyword Argument:

Optional Schema-Specific Keyword Arguments:

These provide control over how the data is formatted, and supported argument names and values depend on the schema as specified by the value passed as the “schema” argument. See “DendroPy Schemas: Phylogenetic and Evolutionary Biology Data Formats” for more details.

Examples

# Using a file path:
d.write(path="path/to/file.dat", schema="nexus")

# Using an open file:
with open("path/to/file.dat", "w") as f:
    d.write(file=f, schema="nexus")
write_as_dot(out, **kwargs)[source]

Writes the tree to out as a DOT formatted digraph

write_ascii_plot(stream, **kwargs)[source]

Writes an ASCII text graphic of this tree to stream.

write_to_path(dest, schema, **kwargs)

Writes to file specified by dest.

write_to_stream(dest, schema, **kwargs)

Writes to file-like object dest.

classmethod yield_from_files(files, schema, taxon_namespace=None, **kwargs)[source]

Iterates over trees from files, returning them one-by-one instead of instantiating all of them in memory at once.

For operations where it is sufficient to process each tree individually (e.g., performing a calculation or set of calculations on a tree and storing the result, after the which the entire tree itself is not needed), this approach is far more performant that reading in the trees using a TreeList. This is because a full tree structure requires significant memory overhead, and as memory gets used up and the OS starts page faulting, performance starts taking some serious hits.

Parameters:
  • files (iterable of file paths or file-like objects.) – Iterable of sources, which can either be strings specifying file paths or file-like objects open for reading. If a source element is a string (isinstance(i,str) == True), then it is assumed to be a path to a file. Otherwise, the source is assumed to be a file-like object.

  • schema (string) – The name of the data format (e.g., “newick” or “nexus”).

  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to use to manage taxon definitions.

  • **kwargs (keyword arguments) – These will be passed directly to the schema-parser implementation.

Yields:

t (|Tree|) – Trees as read from the file.

Examples

taxon_namespace = dendropy.TaxonNamespace()
f1 = open("path/to/trees1.nex", "r")
f2 = open("path/to/trees2.nex", "r")
tree_yielder = dendropy.Tree.yield_from_files(
        files=[f1, f2, "path/to/trees3.nex", "path/to/trees4.nex"],
        schema="nexus",
        taxon_namespace=taxon_namespace,
        store_tree_weights=True,
        preserve_underscores=True,
        rooting="default-unrooted",
        ignore_unrecognized_keyword_arguments=True,
        )
lengths = []
root_ages = []
for tree in tree_yielder:
    length = 0.0
    for edge in tree:
        length += edge.length
    lengths.append(length)
    tree.calc_node_ages()
    root_ages.append(tree.seed_node.age)

The Node Class

class dendropy.datamodel.treemodel.Node(**kwargs)[source]

A :term:Node on a :term:Tree.

Keyword Arguments:
  • taxon (Taxon, optional) – The Taxon instance representing the operational taxonomic unit concept associated with this Node.

  • label (string, optional) – A label for this node.

  • edge_length (numeric, optional) – Length or weight of the edge subtending this node.

add_child(node)[source]

Adds a child node to this node if it is not already a child.

Results in the parent_node attribute of node as well as the tail_node attribute of node.edge being assigned to self.

Parameters:

node (Node) – The node to be added as a child of this node.

Returns:

|Node| – The node that was added.

adjacent_nodes()[source]

Return parent and child nodes.

Returns:

list [Node] – A list with all child nodes and parent node of this node.

age_order_iter(include_leaves=True, filter_fn=None, descending=False)[source]

Deprecated: use Node.ageorder_iter instead.

ageorder_iter(filter_fn=None, include_leaves=True, descending=False)[source]

Iterator over nodes of subtree rooted at this node in order of the age of the node (i.e., the time since the present).

Iterates over nodes in order of age (‘age’ is as given by the age attribute, which is usually the sum of edge lengths from tips to node, i.e., time since present). If include_leaves is True (default), leaves are included in the iteration; if include_leaves is False, leaves will be skipped. If descending is False (default), younger nodes will be returned before older ones; if True, older nodes will be returned before younger ones.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (defau

  • include_leaves (boolean, optional) – If True (default), then leaf nodes are included in the iteration. If False, then leaf nodes are skipped.lt), then all nodes visited will be yielded.

  • descending (boolean, optional) – If False (default), then younger nodes are visited before older ones. If True, then older nodes are visited before younger ones.

Returns:

collections.Iterator [Node] – Iterator over age-ordered sequence of nodes in subtree rooted at this node.

ancestor_iter(filter_fn=None, inclusive=False)[source]

Iterator over all ancestors of this node.

Visits all nodes that are the ancestors of this node. If inclusive is True, self is returned as the first item of the sequence; otherwise self is skipped. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when passed the node as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • inclusive (boolean, optional) – If True, includes this node in the sequence. If False, this is skipped.

Returns:

collections.Iterator [Node] – Iterator over all predecessor/ancestor nodes of this node.

apply(before_fn=None, after_fn=None, leaf_fn=None)[source]

Applies function before_fn and after_fn to all internal nodes and leaf_fn to all terminal nodes in subtree starting with self, with nodes visited in pre-order.

Given a tree with preorder sequence of nodes of [a,b,i,e,j,k,c,g,l,m,f,n,h,o,p,]:

              a
             / \
            /   \
           /     \
          /       \
         /         \
        /           \
       /             c
      b             / \
     / \           /   \
    /   e         /     f
   /   / \       /     / \
  /   /   \     g     /   h
 /   /     \   / \   /   / \
i   j       k l   m n   o   p

the following order of function calls results:

before_fn(a) before_fn(b) leaf_fn(i) before_fn(e) leaf_fn(j) leaf_fn(k) after_fn(e) after_fn(b) before_fn(c) before_fn(g) leaf_fn(l) leaf_fn(m) after_fn(g) before_fn(f) leaf_fn(n) before_fn(h) leaf_fn(o) leaf_fn(p) after_fn(h) after_fn(f) after_fn(c) after_fn(a)

Parameters:
  • before_fn (function object or None) – A function object that takes a Node as its argument.

  • after_fn (function object or None) – A function object that takes a Node as its argument.

  • leaf_fn (function object or None) – A function object that takes a Node as its argument.

Notes

Adapted from work by Mark T. Holder (the peyotl module of the Open Tree of Life Project):

property bipartition

Returns the bipartition for the edge subtending this node.

child_edge_iter(filter_fn=None)[source]

Iterator over all edges that are the (immediate) children of this edge.

Parameters:

filter_fn (function object, optional) – A function object that takes a Edge object as an argument and returns True if the Edge object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all edges visited will be yielded.

Returns:

collections.Iterator [Edge] – An iterator yielding edges that have this edge as a parent.

child_edges()[source]

Returns a shallow-copy list of all child edges of this node.

Note

Unless an actual list is needed, iterating over the child edges using Node.child_edge_iter is preferable to avoid the overhead of list construction.

Returns:

list [Edge] – A list of Edge objects that have self as a tail node.

child_node_iter(filter_fn=None)[source]

Iterator over all nodes that are the (immediate) children of this node.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes that have this node as a parent.

child_nodes()[source]

Returns a shallow-copy list of all child nodes of this node.

Note

Unless an actual list is needed, iterating over the child nodes using Node.child_node_iter is preferable to avoid the overhead of list construction.

Returns:

list [Node] – A list of Node objects that have self as a parent.

clear_child_nodes()[source]

Removes all child nodes.

clone(depth=1)

Creates and returns a copy of self.

Parameters:

depth (integer) –

The depth of the copy:

  • 0: shallow-copy: All member objects are references, except for :attr:annotation_set of top-level object and member Annotation objects: these are full, independent instances (though any complex objects in the value field of Annotation objects are also just references).

  • 1: taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon instances: these are references.

  • 2: Exhaustive deep-copy: all objects are cloned.

collapse_clade()[source]

Collapses all internal edges that are descendants of self.

collapse_conflicting(bipartition)[source]

Collapses every edge in the subtree that conflicts with the given bipartition. This can include the edge subtending subtree_root.

copy_annotations_from(other, attribute_object_mapper=None)

Copies annotations from other, which must be of Annotable type.

Copies are deep-copies, in that the Annotation objects added to the annotation_set AnnotationSet collection of self are independent copies of those in the annotate_set collection of other. However, dynamic bound-attribute annotations retain references to the original objects as given in other, which may or may not be desirable. This is handled by updated the objects to which attributes are bound via mappings found in attribute_object_mapper. In dynamic bound-attribute annotations, the _value attribute of the annotations object (Annotation._value) is a tuple consisting of “(obj, attr_name)”, which instructs the Annotation object to return “getattr(obj, attr_name)” (via: “getattr(*self._value)”) when returning the value of the Annotation. “obj” is typically the object to which the AnnotationSet belongs (i.e., self). When a copy of Annotation is created, the object reference given in the first element of the _value tuple of dynamic bound-attribute annotations are unchanged, unless the id of the object reference is fo

Parameters:
  • other (Annotable) – Source of annotations to copy.

  • attribute_object_mapper (dict) – Like the memo of __deepcopy__, maps object id’s to objects. The purpose of this is to update the parent or owner objects of dynamic attribute annotations. If a dynamic attribute Annotation gives object x as the parent or owner of the attribute (that is, the first element of the Annotation._value tuple is other) and id(x) is found in attribute_object_mapper, then in the copy the owner of the attribute is changed to attribute_object_mapper[id(x)]. If attribute_object_mapper is None (default), then the following mapping is automatically inserted: id(other): self. That is, any references to other in any Annotation object will be remapped to self. If really no reattribution mappings are desired, then an empty dictionary should be passed instead.

deep_copy_annotations_from(other, memo=None)

Note that all references to other in any annotation value (and sub-annotation, and sub-sub-sub-annotation, etc.) will be replaced with references to self. This may not always make sense (i.e., a reference to a particular entity may be absolute regardless of context).

description(depth=1, indent=0, itemize='', output=None, taxon_namespace=None)[source]

Returns description of object, up to level depth.

distance_from_root()[source]

Weighted path length of self from root.

Returns:

numeric – Total weight of all edges connecting self with the root of the tree.

distance_from_tip()[source]

Maximum weighted length of path of self to tip.

If tree is not ultrametric (i.e., descendent edges have different lengths), then count the maximum of edge lengths. Note that Tree.calc_node_ages is a more efficient way of doing this over the whole tree if this value is need for many or all the nodes on the tree.

Returns:

numeric – Maximum weight of edges connecting self to tip.

property edge

Returns the edge subtending this node.

classmethod edge_factory(**kwargs)[source]

Creates and returns a Edge object.

Derived classes can override this method to provide support for specialized or different types of edges on the tree.

Parameters:

kwargs (keyword arguments) – Passed directly to constructor of Edge.

Returns:

|Edge| – A new Edge object.

property edge_length

Returns the length of the edge subtending this node.

extract_subtree(extraction_source_reference_attr_name='extraction_source', node_filter_fn=None, suppress_unifurcations=True, is_apply_filter_to_leaf_nodes=True, is_apply_filter_to_internal_nodes=False, node_factory=None)[source]

Returns a clone of the structure descending from this node.

Parameters:
  • extraction_source_reference_attr_name (str) – Name of attribute to set on cloned nodes that references corresponding original node. If None, then attribute (and reference) will not be created.

  • node_filter_fn (None or function object) – If None, then entire tree structure is cloned. If not None, must be a function object that returns True if a particular Node instance on the original tree should be included in the cloned tree, or False otherwise.

  • is_apply_filter_to_leaf_nodes (bool) – If True then the above filter will be applied to leaf nodes. If False then it will not (and all leaf nodes will be automatically included, unless excluded by an ancestral node being filtered out).

  • is_apply_filter_to_internal_nodes (bool) – If True then the above filter will be applied to internal nodes. If False then it will not (internal nodes without children will still be filtered out).

  • node_factory (function) – If not None, must be a function that takes no arguments and returns a new Node (or equivalent) instance.

Returns:

nd (|Node|) – A node with descending subtree mirroring this one.

get_adjacent_nodes()[source]

Legacy synonym for Node.adjacent_edges

get_incident_edges()[source]

Legacy synonym for Node.incident_edges.

incident_edges()[source]

Return parent and child edges.

Returns:

list [Edge] – A list of edges linking to this node, with outgoing edges (edges connecting to child nodes) followed by the edge connecting this node to its parent.

inorder_iter(filter_fn=None)[source]

In-order iteration over nodes of subtree rooted at this node.

Visits self and all descendant nodes, with each node visited in-between its children. Only valid for strictly-bifurcating trees. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes of the subtree rooted at this node in infix or in-order sequence.

insert_child(index, node)[source]

Adds a child node to this node.

If the node is already a child of this node, then it is moved to the specified position. Results in the parent_node attribute of node as well as the tail_node attribute of node.edge being assigned to self.

Parameters:
  • index (integer) – The index before which to insert the new node.

  • node (Node) – The node to be added as a child of this node.

Returns:

|Node| – The node that was added.

insert_new_child(index, **kwargs)[source]

Create and add a new child to this node at a particular position.

Results in the parent_node attribute of node as well as the tail_node attribute of node.edge being assigned to self.

Parameters:
  • index (integer) – The index before which to insert the new node.

  • kwargs (keyword arguments, optional) – Keyword arguments will be passed directly to the Node constructor (Node.__init()__).

Returns:

|Node| – The new child node that was created and added.

is_internal()[source]

Returns True if the node is not a tip or a leaf node.

Returns:

booleanTrue if the node is not a leaf. False otherwise.

is_leaf()[source]

Returns True if the node is a tip or a leaf node, i.e. has no child nodes.

Returns:

booleanTrue if the node is a leaf, i.e., has no child nodes. False otherwise.

leaf_iter(filter_fn=None)[source]

Iterate over all tips or leaves that ultimately descend from this node.

Visits all leaf or tip nodes descended from this node. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding leaf nodes of the subtree rooted at this node.

leaf_nodes()[source]

Returns list of all leaf_nodes descended from this node (or just list with self as the only member if self is a leaf).

Note

Usage of leaf_iter() is preferable for efficiency reasons unless actual list is required.

Returns:

list [Node] – A list of Node objects descended from this node (inclusive of self) that are the leaves.

level()[source]

Returns the number of nodes between self and the seed node of the tree.

Returns:

integer – The number of nodes between self and the seed node of the tree, or 0 if self has no parent.

level_order_iter(filter_fn=None)[source]

DEPRECATED: Use Node.levelorder_iter instead.

levelorder_iter(filter_fn=None)[source]

Level-order iteration over nodes of subtree rooted at this node.

Visits self and all descendant nodes, with each node and other nodes at the same level (distance from root) visited before their children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are visited.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes of the subtree rooted at this node in level-order sequence.

new_child(**kwargs)[source]

Create and add a new child to this node.

Parameters:

kwargs (keyword arguments) – Keyword arguments will be passed directly to the Node constructor (Node.__init()__).

Returns:

|Node| – The new child node that was created and added.

num_child_nodes()[source]

Returns number of child nodes.

Returns:

int – Number of children in self.

property parent_node

Returns the parent node of this node.

postorder_internal_node_iter(filter_fn=None, exclude_seed_node=False)[source]

Pre-order iterator over internal nodes of subtree rooted at this node.

Visits self and all internal descendant nodes, with each node visited after its children. In DendroPy, “internal nodes” are nodes that have at least one child node, and thus the root or seed node is typically included unless exclude_seed_node is True. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when passed the node as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • exclude_seed_node (boolean, optional) – If False (default), then the seed node or root is visited. If True, then the seed node is skipped.

Returns:

collections.Iterator [Node] – An iterator yielding the internal nodes of the subtree rooted at this node in post-order sequence.

postorder_iter(filter_fn=None)[source]

Post-order iterator over nodes of subtree rooted at this node.

Visits self and all descendant nodes, with each node visited after its children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding the nodes of the subtree rooted at this node in post-order sequence.

preorder_internal_node_iter(filter_fn=None, exclude_seed_node=False)[source]

Pre-order iterator over internal nodes of subtree rooted at this node.

Visits self and all internal descendant nodes, with each node visited before its children. In DendroPy, “internal nodes” are nodes that have at least one child node, and thus the root or seed node is typically included unless exclude_seed_node is True. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when passed the node as an argument are yielded.

Parameters:
  • filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

  • exclude_seed_node (boolean, optional) – If False (default), then the seed node or root is visited. If True, then the seed node is skipped.

Returns:

collections.Iterator [Node] – An iterator yielding the internal nodes of the subtree rooted at this node in pre-order sequence.

preorder_iter(filter_fn=None)[source]

Pre-order iterator over nodes of subtree rooted at this node.

Visits self and all descendant nodes, with each node visited before its children. Nodes can optionally be filtered by filter_fn: only nodes for which filter_fn returns True when called with the node as an argument are yielded.

Parameters:

filter_fn (function object, optional) – A function object that takes a Node object as an argument and returns True if the Node object is to be yielded by the iterator, or False if not. If filter_fn is None (default), then all nodes visited will be yielded.

Returns:

collections.Iterator [Node] – An iterator yielding nodes of the subtree rooted at this node in pre-order sequence.

reinsert_nodes(nd_connection_list)[source]

This function should be used to “undo” the effects of Node.reversible_remove_child NOTE: the behavior is only guaranteed if the tree has not been modified between the remove_child and reinsert_nodes calls! (or the tree has been restored such that the node/edge identities are identical to the state before the remove_child call.

The order of info in each tuple is:

0 - node removed 1 - parent of node removed 2 - pos in parent matrix 3 - children of node removed that were “stolen” 4 - edge that was lengthened by “stealing” length from node’s edge

remove_child(node, suppress_unifurcations=False)[source]

Removes a node from the child set of this node.

Results in the parent of the node being removed set to None. If suppress_unifurcations is True, if this node ends up having only one child after removal of the specified node, then this node will be removed from the tree, with its single child added to the child node set of its parent and the edge length adjusted accordingly. suppress_unifurcations should only be True for unrooted trees.

Parameters:
  • node (Node) – The node to be removed.

  • suppress_unifurcations (boolean, optional) – If False (default), no action is taken. If True, then if the node removal results in a node with degree of two (i.e., a single parent and a single child), then it will be removed from the tree and its (sole) child will be added as a child of its parent (with edge lengths adjusted accordingly).

Returns:

|Node| – The node removed.

reversible_remove_child(node, suppress_unifurcations=False)[source]

This function is a (less-efficient) version of remove_child that also returns the data needed by reinsert_nodes to “undo” the removal.

Returns a list of tuples. The first element of each tuple is the node removed, the other elements are the information needed by reinsert_nodes in order to restore the tree to the same topology as it was before the call to remove_child. If suppress_unifurcations is False then the returned list will contain only one item.

suppress_unifurcations should only be called on unrooted trees.

set_child_nodes(child_nodes)[source]

Assigns the set of child nodes for this node.

Results in the parent_node attribute of each Node in nodes as well as the tail_node attribute of corresponding Edge objects being assigned to self.

Parameters:

child_nodes (collections.Iterable[Node]) – The (iterable) collection of child nodes to be assigned this node as a parent.

set_children(child_nodes)[source]

Deprecated: use Node.set_child_nodes instead.

sibling_nodes()[source]

Return all other children of parent, excluding self.

Returns:

list [Node] – A list of all nodes descended from the same parent as self, excluding self.

sister_nodes()[source]

Legacy synonym for Node.sister_nodes

taxon_namespace_scoped_copy(memo=None)[source]

Cloning level: 1. Taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon objects: these are preserved as references.

The Edge Class

class dendropy.datamodel.treemodel.Edge(**kwargs)[source]

An :term:edge on a :term:tree.

Keyword Arguments:
  • head_node (Node, optional) – Node from to which this edge links, i.e., the child node of this node tail_node.

  • length (numerical, optional) – A value representing the weight of the edge.

  • rootedge (boolean, optional) – Is the child node of this edge the root or seed node of the tree?

  • label (string, optional) – Label for this edge.

property adjacent_edges

Returns a list of all edges that “share” a node with self.

clone(depth=1)

Creates and returns a copy of self.

Parameters:

depth (integer) –

The depth of the copy:

  • 0: shallow-copy: All member objects are references, except for :attr:annotation_set of top-level object and member Annotation objects: these are full, independent instances (though any complex objects in the value field of Annotation objects are also just references).

  • 1: taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon instances: these are references.

  • 2: Exhaustive deep-copy: all objects are cloned.

collapse(adjust_collapsed_head_children_edge_lengths=False)[source]

Inserts all children of the head_node of self as children of the tail_node of self in the same place in the child_node list that head_node had occupied. The edge length and head_node will no longer be part of the tree unless adjust_collapsed_head_children_edge_lengths. is True.

copy_annotations_from(other, attribute_object_mapper=None)

Copies annotations from other, which must be of Annotable type.

Copies are deep-copies, in that the Annotation objects added to the annotation_set AnnotationSet collection of self are independent copies of those in the annotate_set collection of other. However, dynamic bound-attribute annotations retain references to the original objects as given in other, which may or may not be desirable. This is handled by updated the objects to which attributes are bound via mappings found in attribute_object_mapper. In dynamic bound-attribute annotations, the _value attribute of the annotations object (Annotation._value) is a tuple consisting of “(obj, attr_name)”, which instructs the Annotation object to return “getattr(obj, attr_name)” (via: “getattr(*self._value)”) when returning the value of the Annotation. “obj” is typically the object to which the AnnotationSet belongs (i.e., self). When a copy of Annotation is created, the object reference given in the first element of the _value tuple of dynamic bound-attribute annotations are unchanged, unless the id of the object reference is fo

Parameters:
  • other (Annotable) – Source of annotations to copy.

  • attribute_object_mapper (dict) – Like the memo of __deepcopy__, maps object id’s to objects. The purpose of this is to update the parent or owner objects of dynamic attribute annotations. If a dynamic attribute Annotation gives object x as the parent or owner of the attribute (that is, the first element of the Annotation._value tuple is other) and id(x) is found in attribute_object_mapper, then in the copy the owner of the attribute is changed to attribute_object_mapper[id(x)]. If attribute_object_mapper is None (default), then the following mapping is automatically inserted: id(other): self. That is, any references to other in any Annotation object will be remapped to self. If really no reattribution mappings are desired, then an empty dictionary should be passed instead.

deep_copy_annotations_from(other, memo=None)

Note that all references to other in any annotation value (and sub-annotation, and sub-sub-sub-annotation, etc.) will be replaced with references to self. This may not always make sense (i.e., a reference to a particular entity may be absolute regardless of context).

description(depth=1, indent=0, itemize='', output=None, taxon_namespace=None)[source]

Returns description of object, up to level depth.

get_adjacent_edges()[source]

Returns a list of all edges that “share” a node with self.

invert(update_bipartitions=False)[source]

Changes polarity of edge.

is_internal()[source]

Returns True if the head node has children

is_leaf()[source]

Returns True if the head node has no children

taxon_namespace_scoped_copy(memo=None)[source]

Cloning level: 1. Taxon-namespace-scoped copy: All member objects are full independent instances, except for TaxonNamespace and Taxon objects: these are preserved as references.

The Bipartition Class

class dendropy.datamodel.treemodel.Bipartition(**kwargs)[source]

A bipartition on a tree.

A bipartition of a tree is a division or sorting of the leaves/tips of a tree into two mutually-exclusive and collectively-comprehensive subsets, obtained by bisecting the tree at a particular edge. There is thus a one-to-one correspondence with an edge of a tree and a bipartition. The term “split” is often also used to refer to the same concept, though this is typically applied to unrooted trees.

A bipartition is modeled using a bitmask. This is a a bit array representing the membership of taxa, with the least-significant bit corresponding to the first taxon, the next least-signficant bit corresponding to the second taxon, and so on, till the last taxon corresponding to the most-significant bit. Taxon membership in one of two arbitrary groups, ‘0’ or ‘1’, is indicated by its corresponding bit being unset or set, respectively.

To allow comparisons and correct identification of the same bipartition across different rotational and orientiational representations of unrooted trees, we normalize the bipartition such that the first taxon is always assigned to group ‘0’ for bipartition representations of unrooted trees.

The normalization of the bitmask loses information about the actual descendents of a particular edge. Thus in addition to the Bipartition.bitmask attribute, each Bipartition object also maintains a Bipartition.leafset_bitmask attribute which is unnormalized. This is a bit array representing the presence or absence of taxa in the subtree descending from the child node of the edge of which this bipartition is associated. The least-significant bit corresponds to the first taxon, the next least-signficant bit corresponds to the second taxon, and so on, with the last taxon corresponding to the most-significant bit. For rooted trees, the value of Bipartition.bitmask and Bipartition.leafset_bitmask are identical. For unrooted trees, they may or may not be equal.

In general, we use Bipartition.bitmask data to establish the identity of a split or bipartition across different trees: for example, when computing the Robinson-Foulds distances between trees, or in assessing the support for different bipartitions given an MCMC or bootstrap sample of trees. Here the normalization of the bitmask in unrooted trees allows for the (arbitrarily-labeled) group ‘0’ to be consistent across different representations, rotations, and orientations of trees.

On the other hand, we use Bipartition.leafset_bitmask data to work with various ancestor-descendent relationships within the same tree: for example, to quickly assess if a taxon descends from a particular node in a given tree, or if a particular node is a common ancestor of two taxa in a given tree.

The Bipartition object might be used in keys in dictionaries and look-up tables implemented as sets to allow for, e.g., calculation of support in terms of the number times a particular bipartition is observed. The Bipartition.bitmask is used as hash value for this purpose. As such, it is crucial that this value does not change once a particular Bipartition object is stored in a dictionary or set. To this end, we impose the constraint that Bipartition objects are immutable unless the is_mutable attribute is explicitly set to True as a sort of waiver signed by the client code. Client code does this at its risk, with the warning that anything up to and including the implosion of the universe may occur if the Bipartition object is a member of an set of dictionary at the time (or, at the very least, the modified Bipartition object may not be accessible from dictionaries and sets in which it is stored, or may occlude other Bipartition objects in the container).

Note

There are two possible ways of mapping taxa to bits in a bitarray or bitstring.

In the “Least-Signficiant-Bit” (LSB) scheme, the first taxon corresponds to the least-significant, or left-most bit. So, given four taxa, indexed from 1 to 4, taxon 1 would map to 0b0001, taxon 2 would map to 0b0010, taxon 3 would map to 0b0100, and taxon 4 would map to 0b1000.

In the “Most-Significant-Bit” (MSB) scheme, on the other hand, the first taxon corresponds to the most-significant, or right-most bit. So, given four taxa, indexed from 1 to 4, taxon 1 would map to 0b1000, taxon 2 would map to 0b0100, taxon 3 would map to 0b0010, and taxon 4 would map to 0b0001.

We selected the Least Significant Bit (LSB) approach because the MSB scheme requires the size of the taxon namespace to fixed before the index can be assigned to any taxa. For example, under the MSB scheme, if there are 4 taxa, the bitmask for taxon 1 is 0b1000 == 8, but if another taxon is added, then the bitmask for taxon 1 will become 0b10000 == 16. On the other hand, under the LSB scheme, the bitmask for taxon 1 will be 0b0001 == 1 if there are 4 taxa, and 0b00001 == 1 if there 5 taxa, and so on. This stability of taxon indexes even as the taxon namespace grows is a strongly desirable property, and thus the adoption of the LSB scheme.

Constraining the first taxon to be in group 0 (LSB-0) rather than group 1 (LSB-1) is motivated by the fact that, in the former, we would combine the bitmasks of child nodes using OR (logical addition) operations when calculating the bitmask for a parent node, whereas, with the latter, we would need to use AND operations. The former strikes us as more intuitive.

Keyword Arguments:
  • bitmask (integer) – A bit array representing the membership of taxa, with the least-significant bit corresponding to the first taxon, the next least-significant bit corresponding to the second taxon, and so on, till the last taxon corresponding to the most-significant bit. Taxon membership in one of two arbitrary groups, ‘0’ or ‘1’, is indicated by its corresponding bit being unset or set, respectively.

  • leafset_bitmask (integer) – A bit array representing the presence or absence of taxa in the subtree descending from the child node of the edge of which this bipartition is associated. The least-significant bit corresponds to the first taxon, the next least-significant bit corresponds to the second taxon, and so on, with the last taxon corresponding to the most-significant bit.

  • tree_leafset_bitmask (integer) – The leafset_bitmask of the root edge of the tree with which this bipartition is associated. In, general, this will be $0b1111…n$, where $n$ is the number of taxa, except in cases of trees with incomplete leaf-sets, where the positions corresponding to the missing taxa will have the bits unset.

  • is_rooted (bool) – Specifies whether or not the tree with which this bipartition is associated is rooted.

  • is_mutable (bool) – Specifies whether or not the tree is mutable.

compile_bipartition(is_mutable=None)[source]

Updates the values of the various masks specified and calculates the normalized bipartition bitmask.

If a rooted bipartition, then this is set to the value of the leafset bitmask. If an unrooted bipartition, then the leafset bitmask is normalized such that the lowest-significant bit (i.e., the group to which the first taxon belongs) is set to ‘0’.

Also makes this bipartition immutable (unless is_mutable is False), which facilitates it being used in dictionaries and sets.

Note that this requires full population of the following fields:
  • self._leafset_bitmask

  • self._tree_leafset_bitmask

compile_split_bitmask(leafset_bitmask=None, tree_leafset_bitmask=None, is_rooted=None, is_mutable=True)[source]

Updates the values of the various masks specified and calculates the normalized bipartition bitmask.

If a rooted bipartition, then this is set to the value of the leafset bitmask. If an unrooted bipartition, then the leafset bitmask is normalized such that the lowest-significant bit (i.e., the group to which the first taxon belongs) is set to ‘0’.

Also makes this bipartition immutable (unless is_mutable is False), which facilitates it being used in dictionaries and sets.

Parameters:
  • leafset_bitmask (integer) – A bit array representing the presence or absence of taxa in the subtree descending from the child node of the edge of which this bipartition is associated. The least-significant bit corresponds to the first taxon, the next least-signficant bit corresponds to the second taxon, and so on, with the last taxon corresponding to the most-significant bit. If not specified or None, the current value of self.leafset_bitmask is used.

  • tree_leafset_bitmask (integer) – The leafset_bitmask of the root edge of the tree with which this bipartition is associated. In, general, this will be $0b1111…n$, where $n$ is the number of taxa, except in cases of trees with incomplete leaf-sets, where the positions corresponding to the missing taxa will have the bits unset. If not specified or None, the current value of self.tree_leafset_bitmask is used.

  • is_rooted (bool) – Specifies whether or not the tree with which this bipartition is associated is rooted. If not specified or None, the current value of self.is_rooted is used.

Returns:

integer – The bipartition bitmask.

compile_tree_leafset_bitmask(tree_leafset_bitmask, lowest_relevant_bit=None)[source]

Avoids recalculation of lowest_relevant_bit if specified.

static is_compatible_bitmasks(m1, m2, fill_bitmask)[source]

Returns True if m1 is compatible with m2

Parameters:
  • m1 (int) – A bitmask representing a split.

  • m2 (int) – A bitmask representing a split.

Returns:

boolTrue if m1 is compatible with m2. False otherwise.

is_compatible_with(other)[source]

Returns True if other is compatible with self.

Parameters:

other (Bipartition) – The bipartition to check for compatibility.

Returns:

boolTrue if other is compatible with self; False otherwise.

is_incompatible_with(other)[source]

Returns True if other conflicts with self.

Parameters:

other (Bipartition) – The bipartition to check for conflicts.

Returns:

boolTrue if other conflicts with self; False otherwise.

is_leafset_nested_within(other)[source]

Returns True if the leafset of self is a subset of the leafset of other.

Parameters:

other (Bipartition) – The bipartition to check for compatibility.

Returns:

boolTrue if the leafset of self is contained in other.

is_nested_within(other, is_other_masked_for_tree_leafset=False)[source]

Returns True if the current bipartition is contained within other.

Parameters:

other (Bipartition) – The bipartition to check.

Returns:

boolTrue if the the bipartition is “contained” within other

is_trivial()[source]
Returns:

boolTrue if this bipartition divides a leaf and the rest of the tree.

static is_trivial_bitmask(bitmask, fill_bitmask)[source]

Returns True if the bitmask occurs in any tree of the taxa mask – if there is only fewer than two 1’s or fewer than two 0’s in bitmask (among all of the that are 1 in mask).

leafset_as_bitstring(symbol0='0', symbol1='1', reverse=False)[source]

Composes and returns and representation of the bipartition leafset as a bitstring.

Parameters:
  • symbol1 (str) – The symbol to represent group ‘0’ in the bitmask.

  • symbol1 – The symbol to represent group ‘1’ in the bitmask.

  • reverse (bool) – If True, then the first taxon will correspond to the most-significant bit, instead of the least-significant bit, as is the default.

Returns:

str – The bitstring representing the bipartition.

Example

To represent a bipartition in the same scheme used by, e.g. PAUP* or Mr. Bayes:

print(bipartition.leafset_as_bitstring('.', '*', reverse=True))
leafset_as_newick_string(taxon_namespace, preserve_spaces=False, quote_underscores=True)[source]

Represents this bipartition leafset as a newick string.

Parameters:
  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to reference.

  • preserve_spaces (boolean, optional) – If False (default), then spaces in taxon labels will be replaced by underscores. If True, then taxon labels with spaces will be wrapped in quotes.

  • quote_underscores (boolean, optional) – If True (default), then taxon labels with underscores will be wrapped in quotes. If False, then the labels will not be wrapped in quotes.

Returns:

string – NEWICK representation of split specified by bitmask.

leafset_taxa(taxon_namespace, index=0)[source]

Returns list of Taxon objects in the leafset of this bipartition.

Parameters:
  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to reference.

  • index (integer, optional) – Start from this Taxon object instead of the first Taxon object in the collection.

Returns:

list [Taxon] – List of Taxon objects specified or spanned by bitmask.

normalize(bitmask, convention='lsb0')[source]

Return bitmask ensuring that the bit corresponding to the first taxon is 1.

split_as_bitstring(symbol0='0', symbol1='1', reverse=False)[source]

Composes and returns and representation of the bipartition as a bitstring.

Parameters:
  • symbol1 (str) – The symbol to represent group ‘0’ in the bitmask.

  • symbol1 – The symbol to represent group ‘1’ in the bitmask.

  • reverse (bool) – If True, then the first taxon will correspond to the most-significant bit, instead of the least-significant bit, as is the default.

Returns:

str – The bitstring representing the bipartition.

Example

To represent a bipartition in the same scheme used by, e.g. PAUP* or Mr. Bayes:

print(bipartition.split_as_bitstring('.', '*', reverse=True))
split_as_newick_string(taxon_namespace, preserve_spaces=False, quote_underscores=True)[source]

Represents this bipartition split as a newick string.

Parameters:
  • taxon_namespace (TaxonNamespace instance) – The operational taxonomic unit concept namespace to reference.

  • preserve_spaces (boolean, optional) – If False (default), then spaces in taxon labels will be replaced by underscores. If True, then taxon labels with spaces will be wrapped in quotes.

  • quote_underscores (boolean, optional) – If True (default), then taxon labels with underscores will be wrapped in quotes. If False, then the labels will not be wrapped in quotes.

Returns:

string – NEWICK representation of split specified by bitmask.

The AsciiTreePlot Class

class dendropy.datamodel.treemodel.AsciiTreePlot(**kwargs)[source]
Keyword Arguments:
  • plot_metric (str) – A string which specifies how branches should be scaled, one of: ‘age’ (distance from tips), ‘depth’ (distance from root), ‘level’ (number of branches from root) or ‘length’ (edge length/weights).

  • show_internal_node_labels (bool) – Whether or not to write out internal node labels.

  • leaf_spacing_factor (int) – Positive integer: number of rows between each leaf.

  • width (int) – Force a particular display width, in terms of number of columns.

  • node_label_compose_fn (function object) – A function that takes a Node object as an argument and returns the string to be used to display it.

exception NullEdgeLengthError(*args, **kwargs)[source]
with_traceback()

Exception.with_traceback(tb) – set self.__traceback__ to tb and return self.

calc_plot(node, edge_scale_factor)[source]

First pass through tree, post-order traversal to calculate coordinates of each node.

draw_node(node)[source]

Second pass through tree, plotting nodes onto given self.grid.