Skip to content

simulate_alignment

piqtree.simulate_alignment(tree, model, length=1000, rand_seed=None, insertion_rate=0.0, deletion_rate=0.0, insertion_size_distribution='POW{1.7/100}', deletion_size_distribution='POW{1.7/100}', root_seq=None, num_threads=None)

Uses AliSim to simulate an Alignment through IQ-TREE.

PARAMETER DESCRIPTION
tree

A tree to simulate an alignment over.

TYPE: PhyloNode

model

The substitution model's specification.

TYPE: Model | str

length

The length of the alignment (by default 1000). Alignment may be longer when indel model is used due to insertion events.

TYPE: int DEFAULT: 1000

rand_seed

The random seed - None means no seed is used, by default None.

TYPE: int | None DEFAULT: None

insertion_rate

The insertion rate relative to substitution rate (by default 0.0).

TYPE: float DEFAULT: 0.0

deletion_rate

The deletion rate relative to substitution rate (by default 0.0).

TYPE: float DEFAULT: 0.0

insertion_size_distribution

The insertion size distribution (by default the Zipfian distribution with a=1.7 and maximum size 100).

TYPE: IndelDistribution | str DEFAULT: 'POW{1.7/100}'

deletion_size_distribution

The deletion size distribution (by default the Zipfian distribution with a=1.7 and maximum size 100).

TYPE: IndelDistribution | str DEFAULT: 'POW{1.7/100}'

root_seq

The root sequence (by default None).

TYPE: str | None DEFAULT: None

num_threads

Number of threads for IQ-TREE to use, by default None (single-threaded).

TYPE: int | None DEFAULT: None

RETURNS DESCRIPTION
AlignedSeqsType

The simulated alignment.

Source code in src/piqtree/iqtree/_alignment.py
def simulate_alignment(
    tree: PhyloNode,
    model: Model | str,
    length: int = 1000,
    rand_seed: int | None = None,
    insertion_rate: float = 0.0,
    deletion_rate: float = 0.0,
    insertion_size_distribution: IndelDistribution | str = "POW{1.7/100}",
    deletion_size_distribution: IndelDistribution | str = "POW{1.7/100}",
    root_seq: str | None = None,
    num_threads: int | None = None,
) -> Alignment:
    """Uses AliSim to simulate an Alignment through IQ-TREE.

    Parameters
    ----------
    tree: list[cogent3.PhyloNode]
        A tree to simulate an alignment over.
    model: str
        The substitution model's specification.
    length: int
        The length of the alignment (by default 1000).
        Alignment may be longer when indel model is used due to insertion events.
    rand_seed : int | None, optional
        The random seed - None means no seed is used, by default None.
    insertion_rate: float | None, optional
        The insertion rate relative to substitution rate (by default 0.0).
    deletion_rate: float | None, optional
        The deletion rate relative to substitution rate (by default 0.0).
    insertion_size_distribution: str | None, optional
        The insertion size distribution (by default the Zipfian
        distribution with a=1.7 and maximum size 100).
    deletion_size_distribution: str | None, optional
        The deletion size distribution (by default the Zipfian
        distribution with a=1.7 and maximum size 100).
    root_seq: str | None, optional
        The root sequence (by default None).
    num_threads: int | None, optional
        Number of threads for IQ-TREE to use, by default None (single-threaded).

    Returns
    -------
    c3_types.AlignedSeqsType
        The simulated alignment.

    """
    if rand_seed is None:
        rand_seed = make_rand_seed()

    if root_seq is None:
        root_seq = ""

    if num_threads is None:
        num_threads = 1

    if isinstance(model, str):
        model = make_model(model)

    if model.submod_type.base_model in UNSUPPORTED_MODELS:
        msg = f"Lie Model {cast('LieModel', model.submod_type.base_model).value} is unsupported."
        raise ValueError(msg)

    newick_tree = get_newick(tree)

    yaml_result = yaml.safe_load(
        iq_simulate_alignment(
            newick_tree,
            str(model),
            rand_seed,
            "",
            "",
            length,
            insertion_rate,
            deletion_rate,
            root_seq,
            num_threads,
            str(insertion_size_distribution),
            str(deletion_size_distribution),
            -1,
        ),
    )

    seqs = _parse_yaml_alignment(yaml_result["alignment"])

    return make_aligned_seqs(seqs, moltype=model.submod_type.get_moltype())

Usage

For usage, see "Simulate alignments with AliSim".