| Title: | 'R' Access to the 'tskit C' API |
|---|---|
| Description: | 'Tskit' enables efficient storage, manipulation, and analysis of ancestral recombination graphs (ARGs) using succinct tree sequence encoding. The tree sequence encoding of an ARG is described in Wong et al. (2024) <doi:10.1093/genetics/iyae100>, while 'tskit' project is described in Jeffrey et al. (2026) <doi:10.48550/arXiv.2602.09649>. See also <https://tskit.dev> for project news, documentation, and tutorials. 'Tskit' provides 'Python', 'C', and 'Rust' application programming interfaces (APIs). The 'Python' API can be called from 'R' via the 'reticulate' package to load and analyse tree sequences as described at <https://tskit.dev/tutorials/tskitr.html>. 'RcppTskit' provides 'R' access to the 'tskit C' API for cases where the 'reticulate' option is not optimal; for example, high-performance or low-level work with tree sequences. Currently, 'RcppTskit' provides a limited set of functions because the 'Python' API and 'reticulate' already cover most needs. The provided 'RcppTskit R' API mirrors the 'tskit Python' API, while the 'RcppTskit C++' API mirrors the 'tskit C' API. Users should explore the 'RcppTskit' help pages of 'R' functions, while developers should explore the 'RcppTskit:::rtsk_*' low-level 'R' and 'C++'' functions. |
| Authors: | Gregor Gorjanc [aut, cre, cph] (ORCID: <https://orcid.org/0000-0001-8008-2787>), Tskit Developers [cph] (Authors of included tskit C library) |
| Maintainer: | Gregor Gorjanc <[email protected]> |
| License: | MIT + file LICENSE |
| Version: | 0.3.0 |
| Built: | 2026-05-10 06:41:03 UTC |
| Source: | https://github.com/highlanderlab/RcppTskit |
tskit Python moduleThis function imports the reticulate Python tskit module.
If it is not yet installed, it attempts to install it first.
get_tskit_py(object_name = "tskit", force = FALSE) check_tskit_py(object, stop = FALSE)get_tskit_py(object_name = "tskit", force = FALSE) check_tskit_py(object, stop = FALSE)
object_name |
character name of the object holding the reticulate
|
force |
logical; force installation and/or import before returning the reticulate Python module. |
object |
reticulate Python module. |
stop |
logical; whether to throw an error in |
This function is meant for users running tskit <- get_tskit_py()
or similar code, and for other functions in this package that need the
tskit reticulate Python module. The point of get_tskit_py is
to avoid importing the module repeatedly; if it has been imported already,
we reuse that instance. This process can be sensitive to the reticulate
Python setup, module availability, and internet access.
get_tskit_py returns the reticulate Python module tskit
if successful. Otherwise it throws an error (when object_name exists
but is not a reticulate Python module) or returns simpleError
(when installation or import failed). check_tskit_py returns
TRUE if object is a reticulate Python module or FALSE
otherwise.
check_tskit_py(): Test whether get_tskit_py returned a reticulate Python module object
## Not run: tskit <- get_tskit_py() is(tskit) if (check_tskit_py(tskit)) { tskit$ALLELES_01 } ## End(Not run)## Not run: tskit <- get_tskit_py() is(tskit) if (check_tskit_py(tskit)) { tskit$ALLELES_01 } ## End(Not run)
Report the version of installed kastore C API
kastore_version()kastore_version()
The version is stored in the installed header kastore.h.
A named vector with three elements major, minor, and
patch.
kastore_version()kastore_version()
An R6 class holding an external pointer to
a table collection object. As an R6 class, method-calling looks Pythonic
and hence resembles the tskit Python API. Since the class only
holds the pointer, it is lightweight. Currently there is a limited set of
R methods for working with the table collection object.
xptrexternal pointer to the table collection
new()
Create a TableCollection from a file or an external pointer.
TableCollection$new( file, skip_tables = FALSE, skip_reference_sequence = FALSE, xptr = NULL )
filea string specifying the full path of the tree sequence file.
skip_tableslogical; if TRUE, load only non-table information.
skip_reference_sequencelogical; if TRUE, skip loading
reference genome sequence information.
xptran external pointer (externalptr) to a table collection.
See the tskit Python equivalent at
https://github.com/tskit-dev/tskit/blob/dc394d72d121c99c6dcad88f7a4873880924dd72/python/tskit/tables.py#L3463.
TODO: Update URL to TableCollection.load() method #104
https://github.com/HighlanderLab/RcppTskit/issues/104
A TableCollection object.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- TableCollection$new(file = ts_file)
is(tc)
tc
dump()
Write a table collection to a file.
TableCollection$dump(file)
filea string specifying the full path of the tree sequence file.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TableCollection.dump.
No return value; called for side effects.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- TableCollection$new(file = ts_file)
dump_file <- tempfile()
tc$dump(dump_file)
tc$write(dump_file) # alias
\dontshow{file.remove(dump_file)}
write()
Alias for TableCollection$dump.
TableCollection$write(file)
filesee TableCollection$dump.
tree_sequence()
Create a TreeSequence from this table collection.
TableCollection$tree_sequence()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TableCollection.tree_sequence.
A TreeSequence object.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- TableCollection$new(file = ts_file)
ts <- tc$tree_sequence()
is(ts)
num_provenances()
Get the number of provenances in a table collection.
TableCollection$num_provenances()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_provenances()
num_populations()
Get the number of populations in a table collection.
TableCollection$num_populations()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_populations()
num_migrations()
Get the number of migrations in a table collection.
TableCollection$num_migrations()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_migrations()
num_individuals()
Get the number of individuals in a table collection.
TableCollection$num_individuals()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_individuals()
individual_table_add_row()
Add a row to the individuals table.
TableCollection$individual_table_add_row( flags = 0L, location = NULL, parents = NULL, metadata = NULL )
flagsinteger scalar flags for the new individual.
locationnumeric vector with the location of the new individual;
can be NULL if unknown.
parentsinteger vector with parent individual IDs (0-based);
can be NULL if unknown
metadatafor the new individual; accepts NULL,
a raw vector, or a character of length 1.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/stable/python-api.html#tskit.IndividualTable.add_row.
Integer row ID (0-based) of the newly added individual.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(ts_file)
(n_before <- tc$num_individuals())
new_id <- tc$individual_table_add_row()
new_id <- tc$individual_table_add_row(location = c(5, 8))
new_id <- tc$individual_table_add_row(flags = 0L)
new_id <- tc$individual_table_add_row(parents = c(0L, 2L))
new_id <- tc$individual_table_add_row(metadata = "abc")
new_id <- tc$individual_table_add_row(metadata = charToRaw("cba"))
(n_after <- tc$num_individuals())
num_nodes()
Get the number of nodes in a table collection.
TableCollection$num_nodes()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_nodes()
node_table_add_row()
Add a row to the nodes table.
TableCollection$node_table_add_row( flags = 0L, time = 0, population = -1L, individual = -1L, metadata = NULL )
flagsinteger scalar flags for the new node.
timenumeric scalar time value for the new node.
populationinteger scalar population row ID (0-based);
use -1 if not known - NULL maps to -1 (TSK_NULL).
individualinteger scalar individual row ID (0-based);
use -1 if not known - NULL maps to -1 (TSK_NULL).
metadatafor the new node; accepts NULL,
a raw vector, or a character of length 1.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/stable/python-api.html#tskit.NodeTable.add_row.
Integer row ID (0-based) of the newly added node.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(ts_file)
(n_before <- tc$num_nodes())
new_id <- tc$node_table_add_row()
new_id <- tc$node_table_add_row(time = 2.5)
new_id <- tc$node_table_add_row(flags = 1L, time = 3.5, population = 0L)
new_id <- tc$node_table_add_row(flags = 1L, time = 4.5, individual = 0L)
new_id <- tc$node_table_add_row(metadata = "abc")
new_id <- tc$node_table_add_row(metadata = charToRaw("cba"))
(n_after <- tc$num_nodes())
num_edges()
Get the number of edges in a table collection.
TableCollection$num_edges()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_edges()
edge_table_add_row()
Add a row to the edges table.
TableCollection$edge_table_add_row(left, right, parent, child, metadata = NULL)
leftnumeric scalar left coordinate for the new edge.
rightnumeric scalar right coordinate for the new edge.
parentinteger scalar parent node row ID (0-based).
childinteger scalar child node row ID (0-based).
metadatafor the new edge; accepts NULL,
a raw vector, or a character of length 1.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/stable/python-api.html#tskit.EdgeTable.add_row.
Integer row ID (0-based) of the newly added edge.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(ts_file)
child <- tc$node_table_add_row(time = 0.0)
(n_before <- tc$num_edges())
new_id <- tc$edge_table_add_row(
left = 0, right = 50, parent = 16L, child = child
)
new_id <- tc$edge_table_add_row(
left = 50, right = 75, parent = 17L, child = child, metadata = "abc"
)
new_id <- tc$edge_table_add_row(
left = 75, right = 100, parent = 18L, child = child, metadata = charToRaw("cba")
)
(n_after <- tc$num_edges())
num_sites()
Get the number of sites in a table collection.
TableCollection$num_sites()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_sites()
site_table_add_row()
Add a row to the sites table.
TableCollection$site_table_add_row(position, ancestral_state, metadata = NULL)
positionnumeric scalar site position.
ancestral_statecharacter string for the new site.
metadatafor the new site; accepts NULL,
a raw vector, or a character of length 1.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/stable/python-api.html#tskit.SiteTable.add_row.
Integer row ID (0-based) of the newly added site.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(ts_file)
(n_before <- tc$num_sites())
new_id <- tc$site_table_add_row(position = 0.5, ancestral_state = "A")
new_id <- tc$site_table_add_row(position = 2.5, ancestral_state = "T", metadata = "abc")
(n_after <- tc$num_sites())
num_mutations()
Get the number of mutations in a table collection.
TableCollection$num_mutations()
A signed 64 bit integer bit64::integer64.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$num_mutations()
mutation_table_add_row()
Add a row to the mutations table.
TableCollection$mutation_table_add_row( site, node, derived_state, parent = -1L, metadata = NULL, time = NaN )
siteinteger scalar site row ID (0-based).
nodeinteger scalar node row ID (0-based).
derived_statecharacter string for the new mutation.
parentinteger scalar parent mutation row ID (0-based);
use -1 if not known - NULL maps to -1 (TSK_NULL).
metadatafor the new mutation; accepts NULL,
a raw vector, or a character of length 1.
timenumeric scalar mutation time;
use NaN if not known - NULL maps to NaN (TSK_UNKNOWN_TIME).
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/stable/python-api.html#tskit.MutationTable.add_row.
Integer row ID (0-based) of the newly added mutation.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(ts_file)
(n_before <- tc$num_mutations())
# From inspection of tc we have:
# node13(time=0) <- node16(time=0.02...) <- node20(time=0.08...)
# Add mutation above 16L
m0 <- tc$mutation_table_add_row(site = 0L, node = 16L, derived_state = "T", time = 0.03)
# Add mutation above 13L
m1 <- tc$mutation_table_add_row(
site = 0L,
node = 13L,
parent = m0,
time = 0.01,
derived_state = "C",
metadata = "abc"
)
(n_after <- tc$num_mutations())
sequence_length()
Get the sequence length.
TableCollection$sequence_length()
A numeric.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$sequence_length()
time_units()
Get the time units string.
TableCollection$time_units()
A character.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$time_units()
has_index()
Get whether the table collection has edge indexes.
TableCollection$has_index()
A logical.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$has_index()
build_index()
Build edge indexes for this table collection.
TableCollection$build_index()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TableCollection.build_index.
No return value; called for side effects.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$has_index()
tc$drop_index()
tc$has_index()
tc$build_index()
tc$has_index()
drop_index()
Drop edge indexes for this table collection.
TableCollection$drop_index()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TableCollection.drop_index.
No return value; called for side effects.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$has_index()
tc$drop_index()
tc$has_index()
has_reference_sequence()
Get whether the table collection has a reference genome sequence.
TableCollection$has_reference_sequence()
A logical.
tc_file1 <- system.file("examples/test.trees", package = "RcppTskit")
tc_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit")
tc1 <- tc_load(tc_file1)
tc1$has_reference_sequence()
tc2 <- tc_load(tc_file2)
tc2$has_reference_sequence()
file_uuid()
Get the UUID string of the file the table collection was loaded from.
TableCollection$file_uuid()
A character; NA_character_ when file is information is
unavailable.
tc_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(tc_file)
tc$file_uuid()
r_to_py()
This function saves a table collection from R to disk
and loads it into reticulate Python for use with the
tskit Python API.
TableCollection$r_to_py(tskit_module = get_tskit_py(), cleanup = TRUE)
tskit_modulereticulate Python module of tskit.
By default, it calls get_tskit_py to obtain the module.
cleanuplogical; delete the temporary file at the end of the function?
See https://tskit.dev/tutorials/tables_and_editing.html#tables-and-editing on what you can do with the tables.
TableCollection object in reticulate Python.
\dontrun{
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc_r <- tc_load(ts_file)
is(tc_r)
tc_r$print()
# Transfer the table collection to reticulate Python and use tskit Python API
tskit <- get_tskit_py()
if (check_tskit_py(tskit)) {
tc_py <- tc_r$r_to_py()
is(tc_py)
tmp <- tc_py$simplify(samples = c(0L, 1L, 2L, 3L))
tmp
tc_py$individuals$num_rows # 2
tc_py$nodes$num_rows # 8
tc_py$nodes$time # 0.0 ... 5.0093910
}
}
print()
Print a summary of a table collection and its contents.
TableCollection$print()
A list with two data.frames; the first contains table collection properties and their values; the second contains the number of rows in each table and the length of their metadata. All columns are characters since output types differ across the entries. Use individual getters to obtain raw values before they are converted to character.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
tc <- tc_load(file = ts_file)
tc$print()
tc
clone()
The objects of this class are cloneable with this method.
TableCollection$clone(deep = FALSE)
deepWhether to make a deep clone.
tc_py_to_r, tc_load, and
TableCollection$dump.
## ------------------------------------------------ ## Method `TableCollection$new` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) is(tc) tc ## ------------------------------------------------ ## Method `TableCollection$dump` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) dump_file <- tempfile() tc$dump(dump_file) tc$write(dump_file) # alias ## ------------------------------------------------ ## Method `TableCollection$tree_sequence` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) ts <- tc$tree_sequence() is(ts) ## ------------------------------------------------ ## Method `TableCollection$num_provenances` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_provenances() ## ------------------------------------------------ ## Method `TableCollection$num_populations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_populations() ## ------------------------------------------------ ## Method `TableCollection$num_migrations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_migrations() ## ------------------------------------------------ ## Method `TableCollection$num_individuals` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_individuals() ## ------------------------------------------------ ## Method `TableCollection$individual_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_individuals()) new_id <- tc$individual_table_add_row() new_id <- tc$individual_table_add_row(location = c(5, 8)) new_id <- tc$individual_table_add_row(flags = 0L) new_id <- tc$individual_table_add_row(parents = c(0L, 2L)) new_id <- tc$individual_table_add_row(metadata = "abc") new_id <- tc$individual_table_add_row(metadata = charToRaw("cba")) (n_after <- tc$num_individuals()) ## ------------------------------------------------ ## Method `TableCollection$num_nodes` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_nodes() ## ------------------------------------------------ ## Method `TableCollection$node_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_nodes()) new_id <- tc$node_table_add_row() new_id <- tc$node_table_add_row(time = 2.5) new_id <- tc$node_table_add_row(flags = 1L, time = 3.5, population = 0L) new_id <- tc$node_table_add_row(flags = 1L, time = 4.5, individual = 0L) new_id <- tc$node_table_add_row(metadata = "abc") new_id <- tc$node_table_add_row(metadata = charToRaw("cba")) (n_after <- tc$num_nodes()) ## ------------------------------------------------ ## Method `TableCollection$num_edges` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_edges() ## ------------------------------------------------ ## Method `TableCollection$edge_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) child <- tc$node_table_add_row(time = 0.0) (n_before <- tc$num_edges()) new_id <- tc$edge_table_add_row( left = 0, right = 50, parent = 16L, child = child ) new_id <- tc$edge_table_add_row( left = 50, right = 75, parent = 17L, child = child, metadata = "abc" ) new_id <- tc$edge_table_add_row( left = 75, right = 100, parent = 18L, child = child, metadata = charToRaw("cba") ) (n_after <- tc$num_edges()) ## ------------------------------------------------ ## Method `TableCollection$num_sites` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_sites() ## ------------------------------------------------ ## Method `TableCollection$site_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_sites()) new_id <- tc$site_table_add_row(position = 0.5, ancestral_state = "A") new_id <- tc$site_table_add_row(position = 2.5, ancestral_state = "T", metadata = "abc") (n_after <- tc$num_sites()) ## ------------------------------------------------ ## Method `TableCollection$num_mutations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_mutations() ## ------------------------------------------------ ## Method `TableCollection$mutation_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_mutations()) # From inspection of tc we have: # node13(time=0) <- node16(time=0.02...) <- node20(time=0.08...) # Add mutation above 16L m0 <- tc$mutation_table_add_row(site = 0L, node = 16L, derived_state = "T", time = 0.03) # Add mutation above 13L m1 <- tc$mutation_table_add_row( site = 0L, node = 13L, parent = m0, time = 0.01, derived_state = "C", metadata = "abc" ) (n_after <- tc$num_mutations()) ## ------------------------------------------------ ## Method `TableCollection$sequence_length` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$sequence_length() ## ------------------------------------------------ ## Method `TableCollection$time_units` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$time_units() ## ------------------------------------------------ ## Method `TableCollection$has_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$build_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() tc$drop_index() tc$has_index() tc$build_index() tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$drop_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() tc$drop_index() tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$has_reference_sequence` ## ------------------------------------------------ tc_file1 <- system.file("examples/test.trees", package = "RcppTskit") tc_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit") tc1 <- tc_load(tc_file1) tc1$has_reference_sequence() tc2 <- tc_load(tc_file2) tc2$has_reference_sequence() ## ------------------------------------------------ ## Method `TableCollection$file_uuid` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$file_uuid() ## ------------------------------------------------ ## Method `TableCollection$r_to_py` ## ------------------------------------------------ ## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc_r <- tc_load(ts_file) is(tc_r) tc_r$print() # Transfer the table collection to reticulate Python and use tskit Python API tskit <- get_tskit_py() if (check_tskit_py(tskit)) { tc_py <- tc_r$r_to_py() is(tc_py) tmp <- tc_py$simplify(samples = c(0L, 1L, 2L, 3L)) tmp tc_py$individuals$num_rows # 2 tc_py$nodes$num_rows # 8 tc_py$nodes$time # 0.0 ... 5.0093910 } ## End(Not run) ## ------------------------------------------------ ## Method `TableCollection$print` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(file = ts_file) tc$print() tc## ------------------------------------------------ ## Method `TableCollection$new` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) is(tc) tc ## ------------------------------------------------ ## Method `TableCollection$dump` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) dump_file <- tempfile() tc$dump(dump_file) tc$write(dump_file) # alias ## ------------------------------------------------ ## Method `TableCollection$tree_sequence` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- TableCollection$new(file = ts_file) ts <- tc$tree_sequence() is(ts) ## ------------------------------------------------ ## Method `TableCollection$num_provenances` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_provenances() ## ------------------------------------------------ ## Method `TableCollection$num_populations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_populations() ## ------------------------------------------------ ## Method `TableCollection$num_migrations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_migrations() ## ------------------------------------------------ ## Method `TableCollection$num_individuals` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_individuals() ## ------------------------------------------------ ## Method `TableCollection$individual_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_individuals()) new_id <- tc$individual_table_add_row() new_id <- tc$individual_table_add_row(location = c(5, 8)) new_id <- tc$individual_table_add_row(flags = 0L) new_id <- tc$individual_table_add_row(parents = c(0L, 2L)) new_id <- tc$individual_table_add_row(metadata = "abc") new_id <- tc$individual_table_add_row(metadata = charToRaw("cba")) (n_after <- tc$num_individuals()) ## ------------------------------------------------ ## Method `TableCollection$num_nodes` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_nodes() ## ------------------------------------------------ ## Method `TableCollection$node_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_nodes()) new_id <- tc$node_table_add_row() new_id <- tc$node_table_add_row(time = 2.5) new_id <- tc$node_table_add_row(flags = 1L, time = 3.5, population = 0L) new_id <- tc$node_table_add_row(flags = 1L, time = 4.5, individual = 0L) new_id <- tc$node_table_add_row(metadata = "abc") new_id <- tc$node_table_add_row(metadata = charToRaw("cba")) (n_after <- tc$num_nodes()) ## ------------------------------------------------ ## Method `TableCollection$num_edges` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_edges() ## ------------------------------------------------ ## Method `TableCollection$edge_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) child <- tc$node_table_add_row(time = 0.0) (n_before <- tc$num_edges()) new_id <- tc$edge_table_add_row( left = 0, right = 50, parent = 16L, child = child ) new_id <- tc$edge_table_add_row( left = 50, right = 75, parent = 17L, child = child, metadata = "abc" ) new_id <- tc$edge_table_add_row( left = 75, right = 100, parent = 18L, child = child, metadata = charToRaw("cba") ) (n_after <- tc$num_edges()) ## ------------------------------------------------ ## Method `TableCollection$num_sites` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_sites() ## ------------------------------------------------ ## Method `TableCollection$site_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_sites()) new_id <- tc$site_table_add_row(position = 0.5, ancestral_state = "A") new_id <- tc$site_table_add_row(position = 2.5, ancestral_state = "T", metadata = "abc") (n_after <- tc$num_sites()) ## ------------------------------------------------ ## Method `TableCollection$num_mutations` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$num_mutations() ## ------------------------------------------------ ## Method `TableCollection$mutation_table_add_row` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) (n_before <- tc$num_mutations()) # From inspection of tc we have: # node13(time=0) <- node16(time=0.02...) <- node20(time=0.08...) # Add mutation above 16L m0 <- tc$mutation_table_add_row(site = 0L, node = 16L, derived_state = "T", time = 0.03) # Add mutation above 13L m1 <- tc$mutation_table_add_row( site = 0L, node = 13L, parent = m0, time = 0.01, derived_state = "C", metadata = "abc" ) (n_after <- tc$num_mutations()) ## ------------------------------------------------ ## Method `TableCollection$sequence_length` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$sequence_length() ## ------------------------------------------------ ## Method `TableCollection$time_units` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$time_units() ## ------------------------------------------------ ## Method `TableCollection$has_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$build_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() tc$drop_index() tc$has_index() tc$build_index() tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$drop_index` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$has_index() tc$drop_index() tc$has_index() ## ------------------------------------------------ ## Method `TableCollection$has_reference_sequence` ## ------------------------------------------------ tc_file1 <- system.file("examples/test.trees", package = "RcppTskit") tc_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit") tc1 <- tc_load(tc_file1) tc1$has_reference_sequence() tc2 <- tc_load(tc_file2) tc2$has_reference_sequence() ## ------------------------------------------------ ## Method `TableCollection$file_uuid` ## ------------------------------------------------ tc_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(tc_file) tc$file_uuid() ## ------------------------------------------------ ## Method `TableCollection$r_to_py` ## ------------------------------------------------ ## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc_r <- tc_load(ts_file) is(tc_r) tc_r$print() # Transfer the table collection to reticulate Python and use tskit Python API tskit <- get_tskit_py() if (check_tskit_py(tskit)) { tc_py <- tc_r$r_to_py() is(tc_py) tmp <- tc_py$simplify(samples = c(0L, 1L, 2L, 3L)) tmp tc_py$individuals$num_rows # 2 tc_py$nodes$num_rows # 8 tc_py$nodes$time # 0.0 ... 5.0093910 } ## End(Not run) ## ------------------------------------------------ ## Method `TableCollection$print` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(file = ts_file) tc$print() tc
Load a table collection from a file
tc_load(file, skip_tables = FALSE, skip_reference_sequence = FALSE) tc_read(file, skip_tables = FALSE, skip_reference_sequence = FALSE)tc_load(file, skip_tables = FALSE, skip_reference_sequence = FALSE) tc_read(file, skip_tables = FALSE, skip_reference_sequence = FALSE)
file |
a string specifying the full path to a tree sequence file. |
skip_tables |
logical; if |
skip_reference_sequence |
logical; if |
See the tskit Python equivalent at
https://github.com/tskit-dev/tskit/blob/dc394d72d121c99c6dcad88f7a4873880924dd72/python/tskit/tables.py#L3463.
TODO: Update URL to TableCollection.load() method #104
https://github.com/HighlanderLab/RcppTskit/issues/104
A TableCollection object.
tc_read(): Alias for tc_load()
ts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) is(tc) tcts_file <- system.file("examples/test.trees", package = "RcppTskit") tc <- tc_load(ts_file) is(tc) tc
This function saves a table collection from reticulate Python
to temporary file on disk and reads it into R for use with RcppTskit.
tc_py_to_r(tc, cleanup = TRUE)tc_py_to_r(tc, cleanup = TRUE)
tc |
table collection in reticulate Python. |
cleanup |
logical; delete the temporary file at the end of the function? |
Because this transfer is via a temporary file, the file UUID property changes.
A TableCollection object.
TableCollection$r_to_py
tc_load, and TableCollection$dump.
## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") # Use the tskit Python API to work with a table collection (via reticulate) tskit <- get_tskit_py() if (check_tskit_py(tskit)) { tc_py <- tskit$TableCollection$load(ts_file) is(tc_py) tc_py$individuals$num_rows # 8 tmp <- tc_py$simplify(samples = c(0L, 1L, 2L, 3L)) tmp tc_py$individuals$num_rows # 2 tc_py$nodes$num_rows # 8 tc_py$nodes$time # 0.0 ... 5.0093910 # Transfer the table collection to R and use RcppTskit tc_r <- tc_py_to_r(tc_py) is(tc_r) tc_r$print() } ## End(Not run)## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") # Use the tskit Python API to work with a table collection (via reticulate) tskit <- get_tskit_py() if (check_tskit_py(tskit)) { tc_py <- tskit$TableCollection$load(ts_file) is(tc_py) tc_py$individuals$num_rows # 8 tmp <- tc_py$simplify(samples = c(0L, 1L, 2L, 3L)) tmp tc_py$individuals$num_rows # 2 tc_py$nodes$num_rows # 8 tc_py$nodes$time # 0.0 ... 5.0093910 # Transfer the table collection to R and use RcppTskit tc_r <- tc_py_to_r(tc_py) is(tc_r) tc_r$print() } ## End(Not run)
An R6 class holding an external pointer to
a tree sequence object. As an R6 class, method-calling looks Pythonic
and hence resembles the tskit Python API. Since the class only
holds the pointer, it is lightweight. Currently there is a limited set of
R methods for working with the tree sequence.
xptrexternal pointer to the tree sequence
new()
Create a TreeSequence from a file or an external pointer.
See ts_load for details and examples.
TreeSequence$new( file, skip_tables = FALSE, skip_reference_sequence = FALSE, xptr = NULL )
filea string specifying the full path of the tree sequence file.
skip_tableslogical; if TRUE, load only non-table information.
skip_reference_sequencelogical; if TRUE, skip loading
reference genome sequence information.
xptran external pointer (externalptr) to a tree sequence.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.load.
A TreeSequence object.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- TreeSequence$new(file = ts_file)
is(ts)
ts
ts$num_nodes()
# Also
ts <- ts_load(ts_file)
is(ts)
dump()
Write a tree sequence to a file.
TreeSequence$dump(file)
filea string specifying the full path of the tree sequence file.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.dump.
No return value; called for side effects.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
dump_file <- tempfile()
ts$dump(dump_file)
ts$write(dump_file) # alias
\dontshow{file.remove(dump_file)}
write()
Alias for TreeSequence$dump.
TreeSequence$write(file)
filesee TreeSequence$dump.
dump_tables()
Copy the tables into a TableCollection.
TreeSequence$dump_tables()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.dump_tables.
A TableCollection object.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
tc <- ts$dump_tables()
is(tc)
print()
Print a summary of a tree sequence and its contents.
TreeSequence$print()
A list with two data.frames; the first contains tree sequence properties and their values; the second contains the number of rows in each table and the length of their metadata. All columns are characters since output types differ across the entries. Use individual getters to obtain raw values before they are converted to character.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$print()
ts
r_to_py()
This function saves a tree sequence from R to disk
and loads it into reticulate Python for use with the
tskit Python API.
TreeSequence$r_to_py(tskit_module = get_tskit_py(), cleanup = TRUE)
tskit_modulereticulate Python module of tskit.
By default, it calls get_tskit_py to obtain the module.
cleanuplogical; delete the temporary file at the end of the function?
TreeSequence object in reticulate Python.
\dontrun{
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts_r <- ts_load(ts_file)
is(ts_r)
ts_r$num_individuals() # 8
# Transfer the tree sequence to reticulate Python and use tskit Python API
tskit <- get_tskit_py()
if (check_tskit_py(tskit)) {
ts_py <- ts_r$r_to_py()
is(ts_py)
ts_py$num_individuals # 8
ts2_py <- ts_py$simplify(samples = c(0L, 1L, 2L, 3L))
ts_py$num_individuals # 8
ts2_py$num_individuals # 2
ts2_py$num_nodes # 8
ts2_py$tables$nodes$time # 0.0 ... 5.0093910
}
}
num_provenances()
Get the number of provenances in a tree sequence.
TreeSequence$num_provenances()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_provenances.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_provenances()
num_populations()
Get the number of populations in a tree sequence.
TreeSequence$num_populations()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_populations.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_populations()
num_migrations()
Get the number of migrations in a tree sequence.
TreeSequence$num_migrations()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_migrations.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_migrations()
num_individuals()
Get the number of individuals in a tree sequence.
TreeSequence$num_individuals()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_individuals.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_individuals()
num_samples()
Get the number of samples (of nodes) in a tree sequence.
TreeSequence$num_samples()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_samples.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_samples()
num_nodes()
Get the number of nodes in a tree sequence.
TreeSequence$num_nodes()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_nodes.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_nodes()
num_edges()
Get the number of edges in a tree sequence.
TreeSequence$num_edges()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_nodes.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_edges()
num_trees()
Get the number of trees in a tree sequence.
TreeSequence$num_trees()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_trees.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_trees()
num_sites()
Get the number of sites in a tree sequence.
TreeSequence$num_sites()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_sites.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_sites()
num_mutations()
Get the number of mutations in a tree sequence.
TreeSequence$num_mutations()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.num_mutations.
A signed 64 bit integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$num_mutations()
sequence_length()
Get the sequence length.
TreeSequence$sequence_length()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.sequence_length.
A numeric.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$sequence_length()
discrete_genome()
Get the discrete genome status.
TreeSequence$discrete_genome()
Returns TRUE if all genomic coordinates in the tree
sequence are discrete integer values.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.discrete_genome.
A logical.
ts_file1 <- system.file("examples/test.trees", package = "RcppTskit")
ts_file2 <- system.file("examples/test_non_discrete_genome.trees", package = "RcppTskit")
ts1 <- ts_load(ts_file1)
ts1$discrete_genome()
ts2 <- ts_load(ts_file2)
ts2$discrete_genome()
has_reference_sequence()
Get whether the tree sequence has a reference genome sequence.
TreeSequence$has_reference_sequence()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.has_reference_sequence.
A logical.
ts_file1 <- system.file("examples/test.trees", package = "RcppTskit")
ts_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit")
ts1 <- ts_load(ts_file1)
ts1$has_reference_sequence()
ts2 <- ts_load(ts_file2)
ts2$has_reference_sequence()
time_units()
Get the time units string.
TreeSequence$time_units()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.time_units.
A character.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$time_units()
discrete_time()
Get the discrete time status.
TreeSequence$discrete_time()
Returns TRUE if all time values in the tree sequence are
discrete integer values.
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.discrete_time.
A logical.
ts_file1 <- system.file("examples/test.trees", package = "RcppTskit")
ts_file2 <- system.file("examples/test_discrete_time.trees", package = "RcppTskit")
ts1 <- ts_load(ts_file1)
ts1$discrete_time()
ts2 <- ts_load(ts_file2)
ts2$discrete_time()
min_time()
Get the min time in node table and mutation table.
TreeSequence$min_time()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.min_time.
A numeric.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$min_time()
max_time()
Get the max time in node table and mutation table.
TreeSequence$max_time()
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.TreeSequence.max_time.
A numeric.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$max_time()
metadata_length()
Get the length of metadata in a tree sequence and its tables.
TreeSequence$metadata_length()
A named list with the length of metadata, each as a signed 64 bit
integer bit64::integer64.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$metadata_length()
file_uuid()
Get the UUID string of the file the tree sequence was loaded from.
TreeSequence$file_uuid()
A character; NA_character_ when file is information is
unavailable.
ts_file <- system.file("examples/test.trees", package = "RcppTskit")
ts <- ts_load(ts_file)
ts$file_uuid()
clone()
The objects of this class are cloneable with this method.
TreeSequence$clone(deep = FALSE)
deepWhether to make a deep clone.
ts_py_to_r, ts_load, and
TreeSequence$dump.
## ------------------------------------------------ ## Method `TreeSequence$new` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- TreeSequence$new(file = ts_file) is(ts) ts ts$num_nodes() # Also ts <- ts_load(ts_file) is(ts) ## ------------------------------------------------ ## Method `TreeSequence$dump` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) dump_file <- tempfile() ts$dump(dump_file) ts$write(dump_file) # alias ## ------------------------------------------------ ## Method `TreeSequence$dump_tables` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) tc <- ts$dump_tables() is(tc) ## ------------------------------------------------ ## Method `TreeSequence$print` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$print() ts ## ------------------------------------------------ ## Method `TreeSequence$r_to_py` ## ------------------------------------------------ ## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts_r <- ts_load(ts_file) is(ts_r) ts_r$num_individuals() # 8 # Transfer the tree sequence to reticulate Python and use tskit Python API tskit <- get_tskit_py() if (check_tskit_py(tskit)) { ts_py <- ts_r$r_to_py() is(ts_py) ts_py$num_individuals # 8 ts2_py <- ts_py$simplify(samples = c(0L, 1L, 2L, 3L)) ts_py$num_individuals # 8 ts2_py$num_individuals # 2 ts2_py$num_nodes # 8 ts2_py$tables$nodes$time # 0.0 ... 5.0093910 } ## End(Not run) ## ------------------------------------------------ ## Method `TreeSequence$num_provenances` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_provenances() ## ------------------------------------------------ ## Method `TreeSequence$num_populations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_populations() ## ------------------------------------------------ ## Method `TreeSequence$num_migrations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_migrations() ## ------------------------------------------------ ## Method `TreeSequence$num_individuals` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_individuals() ## ------------------------------------------------ ## Method `TreeSequence$num_samples` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_samples() ## ------------------------------------------------ ## Method `TreeSequence$num_nodes` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_nodes() ## ------------------------------------------------ ## Method `TreeSequence$num_edges` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_edges() ## ------------------------------------------------ ## Method `TreeSequence$num_trees` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_trees() ## ------------------------------------------------ ## Method `TreeSequence$num_sites` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_sites() ## ------------------------------------------------ ## Method `TreeSequence$num_mutations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_mutations() ## ------------------------------------------------ ## Method `TreeSequence$sequence_length` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$sequence_length() ## ------------------------------------------------ ## Method `TreeSequence$discrete_genome` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_non_discrete_genome.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$discrete_genome() ts2 <- ts_load(ts_file2) ts2$discrete_genome() ## ------------------------------------------------ ## Method `TreeSequence$has_reference_sequence` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$has_reference_sequence() ts2 <- ts_load(ts_file2) ts2$has_reference_sequence() ## ------------------------------------------------ ## Method `TreeSequence$time_units` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$time_units() ## ------------------------------------------------ ## Method `TreeSequence$discrete_time` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_discrete_time.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$discrete_time() ts2 <- ts_load(ts_file2) ts2$discrete_time() ## ------------------------------------------------ ## Method `TreeSequence$min_time` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$min_time() ## ------------------------------------------------ ## Method `TreeSequence$max_time` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$max_time() ## ------------------------------------------------ ## Method `TreeSequence$metadata_length` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$metadata_length() ## ------------------------------------------------ ## Method `TreeSequence$file_uuid` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$file_uuid()## ------------------------------------------------ ## Method `TreeSequence$new` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- TreeSequence$new(file = ts_file) is(ts) ts ts$num_nodes() # Also ts <- ts_load(ts_file) is(ts) ## ------------------------------------------------ ## Method `TreeSequence$dump` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) dump_file <- tempfile() ts$dump(dump_file) ts$write(dump_file) # alias ## ------------------------------------------------ ## Method `TreeSequence$dump_tables` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) tc <- ts$dump_tables() is(tc) ## ------------------------------------------------ ## Method `TreeSequence$print` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$print() ts ## ------------------------------------------------ ## Method `TreeSequence$r_to_py` ## ------------------------------------------------ ## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts_r <- ts_load(ts_file) is(ts_r) ts_r$num_individuals() # 8 # Transfer the tree sequence to reticulate Python and use tskit Python API tskit <- get_tskit_py() if (check_tskit_py(tskit)) { ts_py <- ts_r$r_to_py() is(ts_py) ts_py$num_individuals # 8 ts2_py <- ts_py$simplify(samples = c(0L, 1L, 2L, 3L)) ts_py$num_individuals # 8 ts2_py$num_individuals # 2 ts2_py$num_nodes # 8 ts2_py$tables$nodes$time # 0.0 ... 5.0093910 } ## End(Not run) ## ------------------------------------------------ ## Method `TreeSequence$num_provenances` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_provenances() ## ------------------------------------------------ ## Method `TreeSequence$num_populations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_populations() ## ------------------------------------------------ ## Method `TreeSequence$num_migrations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_migrations() ## ------------------------------------------------ ## Method `TreeSequence$num_individuals` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_individuals() ## ------------------------------------------------ ## Method `TreeSequence$num_samples` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_samples() ## ------------------------------------------------ ## Method `TreeSequence$num_nodes` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_nodes() ## ------------------------------------------------ ## Method `TreeSequence$num_edges` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_edges() ## ------------------------------------------------ ## Method `TreeSequence$num_trees` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_trees() ## ------------------------------------------------ ## Method `TreeSequence$num_sites` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_sites() ## ------------------------------------------------ ## Method `TreeSequence$num_mutations` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$num_mutations() ## ------------------------------------------------ ## Method `TreeSequence$sequence_length` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$sequence_length() ## ------------------------------------------------ ## Method `TreeSequence$discrete_genome` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_non_discrete_genome.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$discrete_genome() ts2 <- ts_load(ts_file2) ts2$discrete_genome() ## ------------------------------------------------ ## Method `TreeSequence$has_reference_sequence` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_with_ref_seq.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$has_reference_sequence() ts2 <- ts_load(ts_file2) ts2$has_reference_sequence() ## ------------------------------------------------ ## Method `TreeSequence$time_units` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$time_units() ## ------------------------------------------------ ## Method `TreeSequence$discrete_time` ## ------------------------------------------------ ts_file1 <- system.file("examples/test.trees", package = "RcppTskit") ts_file2 <- system.file("examples/test_discrete_time.trees", package = "RcppTskit") ts1 <- ts_load(ts_file1) ts1$discrete_time() ts2 <- ts_load(ts_file2) ts2$discrete_time() ## ------------------------------------------------ ## Method `TreeSequence$min_time` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$min_time() ## ------------------------------------------------ ## Method `TreeSequence$max_time` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$max_time() ## ------------------------------------------------ ## Method `TreeSequence$metadata_length` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$metadata_length() ## ------------------------------------------------ ## Method `TreeSequence$file_uuid` ## ------------------------------------------------ ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) ts$file_uuid()
Load a tree sequence from a file
ts_load(file, skip_tables = FALSE, skip_reference_sequence = FALSE) ts_read(file, skip_tables = FALSE, skip_reference_sequence = FALSE)ts_load(file, skip_tables = FALSE, skip_reference_sequence = FALSE) ts_read(file, skip_tables = FALSE, skip_reference_sequence = FALSE)
file |
a string specifying the full path to a tree sequence file. |
skip_tables |
logical; if |
skip_reference_sequence |
logical; if |
See the tskit Python equivalent at
https://tskit.dev/tskit/docs/latest/python-api.html#tskit.load.
A TreeSequence object.
ts_read(): Alias for ts_load()
ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) is(ts) ts ts$num_nodes() # Also ts <- TreeSequence$new(file = ts_file) is(ts)ts_file <- system.file("examples/test.trees", package = "RcppTskit") ts <- ts_load(ts_file) is(ts) ts ts$num_nodes() # Also ts <- TreeSequence$new(file = ts_file) is(ts)
This function saves a tree sequence from reticulate Python to
temporary file on disk and reads it into R for use with RcppTskit.
ts_py_to_r(ts, cleanup = TRUE)ts_py_to_r(ts, cleanup = TRUE)
ts |
tree sequence in reticulate Python. |
cleanup |
logical; delete the temporary file at the end of the function? |
Because this transfer is via a temporary file, the file UUID property changes.
A TreeSequence object.
TreeSequence$r_to_py
ts_load, and TreeSequence$dump.
## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") # Use the tskit Python API to work with a tree sequence (via reticulate) tskit <- get_tskit_py() if (check_tskit_py(tskit)) { ts_py <- tskit$load(ts_file) is(ts_py) ts_py$num_individuals # 8 ts2_py <- ts_py$simplify(samples = c(0L, 1L, 2L, 3L)) ts_py$num_individuals # 8 ts2_py$num_individuals # 2 ts2_py$num_nodes # 8 ts2_py$tables$nodes$time # 0.0 ... 5.0093910 # Transfer the tree sequence to R and use RcppTskit ts2_r <- ts_py_to_r(ts2_py) is(ts2_r) ts2_r$num_individuals() # 2 } ## End(Not run)## Not run: ts_file <- system.file("examples/test.trees", package = "RcppTskit") # Use the tskit Python API to work with a tree sequence (via reticulate) tskit <- get_tskit_py() if (check_tskit_py(tskit)) { ts_py <- tskit$load(ts_file) is(ts_py) ts_py$num_individuals # 8 ts2_py <- ts_py$simplify(samples = c(0L, 1L, 2L, 3L)) ts_py$num_individuals # 8 ts2_py$num_individuals # 2 ts2_py$num_nodes # 8 ts2_py$tables$nodes$time # 0.0 ... 5.0093910 # Transfer the tree sequence to R and use RcppTskit ts2_r <- ts_py_to_r(ts2_py) is(ts2_r) ts2_r$num_individuals() # 2 } ## End(Not run)
tskit C APIReport the version of installed tskit C API
tskit_version()tskit_version()
The version is defined in the installed header tskit/core.h.
A named vector with three elements major, minor, and
patch.
tskit_version()tskit_version()