Title: | Tools for Ecological Niche Evolution Assessment Considering Uncertainty |
---|---|
Description: | A collection of tools that allow users to perform critical steps in the process of assessing ecological niche evolution over phylogenies, with uncertainty incorporated explicitly in reconstructions. The method proposed here for ancestral reconstruction of ecological niches characterizes species' niches using a bin-based approach that incorporates uncertainty in estimations. Compared to other existing methods, the approaches presented here reduce risk of overestimation of amounts and rates of ecological niche evolution. The main analyses include: initial exploration of environmental data in occurrence records and accessible areas, preparation of data for phylogenetic analyses, executing comparative phylogenetic analyses of ecological niches, and plotting for interpretations. Details on the theoretical background and methods used can be found in: Owens et al. (2020) <doi:10.1002/ece3.6359>, Peterson et al. (1999) <doi:10.1126/science.285.5431.1265>, Soberón and Peterson (2005) <doi:10.17161/bi.v2i0.4>, Peterson (2011) <doi:10.1111/j.1365-2699.2010.02456.x>, Barve et al. (2011) <doi:10.1111/ecog.02671>, Machado-Stredel et al. (2021) <doi:10.21425/F5FBG48814>, Owens et al. (2013) <doi:10.1016/j.ecolmodel.2013.04.011>, Saupe et al. (2018) <doi:10.1093/sysbio/syx084>, and Cobos et al. (2021) <doi:10.1111/jav.02868>. |
Authors: | Marlon E. Cobos [aut, cre], Hannah L. Owens [aut], A. Townsend Peterson [aut] |
Maintainer: | Marlon E. Cobos <[email protected]> |
License: | GPL-3 |
Version: | 0.1.20 |
Built: | 2024-10-30 03:15:06 UTC |
Source: | https://github.com/marlonecobos/nichevol |
Helper function to prepare bin tables
bin_env(overall_range, M_range, sp_range, bin_size)
bin_env(overall_range, M_range, sp_range, bin_size)
overall_range |
(numeric) minimum and maximum values of all species and Ms to be analyzed. |
M_range |
matrix of ranges of environmental values in M for all species. Columns must be minimum and maximum, and rows correspond to species. |
sp_range |
matrix of ranges of environmental values in occurrences for all species. Columns must be minimum and maximum, and rows correspond to species. |
bin_size |
(numeric) size of bins. Range of environmental values to be considered when creating each character in bin tables. See details. |
The argument bin_size
helps to create characters that represent not
only one value of an environmental variable, but a range of environmental
conditions. For instance, if a variable of precipitation in mm is used, a
value of 10 for bin_size
indicates that each character will represent
a class that correspond to 10 continuous values of precipitation (e.g., from
100 to 110 mm).
A character matrix (table of characters) containing bins for a given variable
and for all species considered. See more details in bin_tables
.
# example o_range <- c(1, 25) m_range <- rbind(c(5, 15), c(10, 23), c(4, 20)) s_range <- rbind(c(7, 15), c(12, 21), c(3, 18)) # bin preparation bins <- bin_env(overall_range = o_range, M_range = m_range, sp_range = s_range, bin_size = 1)
# example o_range <- c(1, 25) m_range <- rbind(c(5, 15), c(10, 23), c(4, 20)) s_range <- rbind(c(7, 15), c(12, 21), c(3, 18)) # bin preparation bins <- bin_env(overall_range = o_range, M_range = m_range, sp_range = s_range, bin_size = 1)
Maximum likelihood reconstruction of ancestral character states
bin_ml_rec(tree_data, ...)
bin_ml_rec(tree_data, ...)
tree_data |
a list of two elements (phy and data) resulting from using the
function |
... |
other arguments from |
Reconstructions are done using the function ace
from the
ape
package. The argument method is set as "ML" and the type
of variable is "discrete".
A table with columns representing bins, rows representing first tip states and then reconstructed nodes.
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum likelihood reconstruction ml_rec <- bin_ml_rec(treeWdata)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum likelihood reconstruction ml_rec <- bin_ml_rec(treeWdata)
Maximum parsimony reconstruction of ancestral character states
bin_par_rec(tree_data, ...)
bin_par_rec(tree_data, ...)
tree_data |
a list of two elements (phy and data) resulting from using
the function |
... |
other arguments from |
Reconstructions are done using the asr_max_parsimony
function from the castor
package.
A table with columns representing bins, rows representing first tip states and then reconstructed nodes.
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction par_rec <- bin_par_rec(treeWdata)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction par_rec <- bin_par_rec(treeWdata)
bin_table helps in creating a bin table of environmental conditions in accessible areas (M) and for species occurrence records (i.e., table of characters).
bin_table(Ms, occurrences, species, longitude, latitude, variable, percentage_out = 5, n_bins = 20, bin_size, verbose = TRUE)
bin_table(Ms, occurrences, species, longitude, latitude, variable, percentage_out = 5, n_bins = 20, bin_size, verbose = TRUE)
Ms |
a list of SpatVector objects representing the accessible area
(M) for all species to be analyzed. The order of species represented by each
object here must coincide with the one in |
occurrences |
a list of data.frames of occurrence records for all species.
The order of species represented by each data.frame must coincide with the one
in |
species |
(character) name of the column in occurrence data.frames that contains the name of the species. |
longitude |
(character) name of the column in occurrence files containing values of longitude. |
latitude |
(character) name of the column in occurrence files containing values of latitude. |
variable |
a single SpatRaster layer representing an environmental variable of interest. See details. |
percentage_out |
(numeric) percentage of extreme environmental data in M to be excluded in bin creation for further analyses. See details. Default = 5. |
n_bins |
(numeric) number of bins to be created from the range of environmental values considered when creating each character in bin tables. Default = 20. See details. |
bin_size |
(numeric) argument deprecated, use n_bins instead. |
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Coordinates in occurrences
, SpatVector objects in Ms
, and
SpatRaster in variable
must coincide in the geographic projection in
which they are represented. WGS84 with no planar projection is recommended.
Accessible area (M) is understood as the geographic area that has been accessible for a species for relevant periods of time. Defining M is usually a hard task, but also a very important one, because it allows identifying uncertainties about the ability of a species to maintain populations in certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
The percentage to be defined in percentage_out
excludes a percentage
of extreme environmental values to prevent from considering extremely rare
environmental values in the accessible area for the species (M). Being too
rare, these values may have never been explored by the species; therefore,
including them in the process of preparation of the table of characters
(bin table) is risky.
The argument n_bins
helps to define how many characters (bins) will be
considered for the range of values in each variable. This is, a value of 20
determines that a range of temperature (5-25) will be split approximately
every 1 degree. The argument bin_size
has been deprecated.
A list containing a table of characters to represent ecological niches of the species of interest.
Potential values for characters are:
"1" = the species is present in those environmental conditions.
"0" = the species is not present in those environmental conditions. This is, those environmental conditions inside the accessible area (M) are more extreme than the ones used for the species.
"?" = there is no certainty about the species presence in those environmental conditions. This happens in environmental combinations more extreme than the ones found in the accessible area (M), when environmental conditions in species records are as extreme as the most extreme ones in M.
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # preparing bins char_table <- bin_table(Ms = m_list, occurrences = occ_list, species = "species", longitude = "x", latitude = "y", variable = temp, percentage_out = 5, n_bins = 20)
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # preparing bins char_table <- bin_table(Ms = m_list, occurrences = occ_list, species = "species", longitude = "x", latitude = "y", variable = temp, percentage_out = 5, n_bins = 20)
bin_tables helps in creating bin tables of environmental conditions in accessible areas (M) and species occurrence records (i.e., table of characters). This is done using results from previous analyses, and can be applied to various species and multiple variables.
bin_tables(ranges, percentage_out = 5, n_bins = 20, bin_size, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
bin_tables(ranges, percentage_out = 5, n_bins = 20, bin_size, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
ranges |
list of ranges of environmental values in M and in species
occurrences derived from using the function |
percentage_out |
(numeric) percentage of extreme environmental data in M to be excluded in bin creation for further analyses. See details. Default = 5. |
n_bins |
(numeric) number of bins to be created from the range of environmental values considered when creating each character in bin tables. Default = 20. See details. |
bin_size |
(numeric) argument deprecated, use n_bins instead. |
save |
(logical) whether or not to save the results in working directory. Default = FALSE. |
output_directory |
(character) name of the folder in which results will be written. |
overwrite |
(logical) whether or not to overwrite existing results
in |
verbose |
(logical) whether messages should be printed. Default = TRUE. |
The percentage to be defined in percentage_out
must correspond with
one of the confidence limits defined in histograms_env
(argument CL_lines
). For instance, if CL_lines
= 95, then
percentage_out
can only be either 5 (keeping data inside the 95 CL) or
0 (to avoid exclusion of extreme values in M).
Excluding a certain percentage of extreme environmental values prevents the algorithm from considering extremely rare environmental values in the accessible area for the species (M). Being too rare, these values may have never been explored by the species; therefore, including them in the process of preparation of the table of characters (bin table) is risky.
The argument n_bins
helps to define how many characters (bins) will be
considered for the range of values in each variable. This is, a value of 20
determines that a range of temperature (5-25) will be split approximately
every 1 degree. The argument bin_size
has been deprecated.
A list named as in ranges
containing the table(s) of characters.
A folder named as in output_directory
containing all resulting csv
files with the tables of characters will be created if save
is set as
TRUE.
Potential values for characters are:
"1" = the species is present in those environmental conditions.
"0" = the species is not present in those environmental conditions. This is, those environmental conditions inside the accessible area (M) are more extreme than the ones used for the species.
"?" = there is no certainty about the species presence in those environmental conditions. This happens if environmental combinations are more extreme than the ones found in the accessible area (M), when environmental conditions in species records are as extreme as the most extreme ones in M.
# simple list of ranges ranges <- list(temp = data.frame(Species = c("sp1", "sp2", "sp3"), Species_lower = c(120, 56, 59.75), Species_upper = c(265, 333, 333), M_lower = c(93, 39, 56), M_upper = c(302, 333, 333), M_95_lowerCL = c(158, 91, 143), M_95_upperCL = c(292, 290, 326)), prec = data.frame(Species = c("sp1", "sp2", "sp3"), Species_lower = c(597, 3, 3), Species_upper = c(3492, 2673, 6171), M_lower = c(228, 3, 3), M_upper = c(6369, 7290, 6606), M_95_lowerCL = c(228, 3, 3), M_95_upperCL = c(3114, 2376, 2568))) # bin preparation bins <- bin_tables(ranges, percentage_out = 5, n_bins = 20) # see arguments save and output_directory to write results in local directory
# simple list of ranges ranges <- list(temp = data.frame(Species = c("sp1", "sp2", "sp3"), Species_lower = c(120, 56, 59.75), Species_upper = c(265, 333, 333), M_lower = c(93, 39, 56), M_upper = c(302, 333, 333), M_95_lowerCL = c(158, 91, 143), M_95_upperCL = c(292, 290, 326)), prec = data.frame(Species = c("sp1", "sp2", "sp3"), Species_lower = c(597, 3, 3), Species_upper = c(3492, 2673, 6171), M_lower = c(228, 3, 3), M_upper = c(6369, 7290, 6606), M_95_lowerCL = c(228, 3, 3), M_95_upperCL = c(3114, 2376, 2568))) # bin preparation bins <- bin_tables(ranges, percentage_out = 5, n_bins = 20) # see arguments save and output_directory to write results in local directory
bin_tables0 helps in creating bin tables of environmental conditions in accessible areas (M) and species occurrence records (i.e., table of characters). This is done using data read directly from a local directory, and can be applied to various species and multiple variables.
bin_tables0(M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, round = FALSE, round_names, multiplication_factor = 1, percentage_out = 5, n_bins = 20, bin_size, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
bin_tables0(M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, round = FALSE, round_names, multiplication_factor = 1, percentage_out = 5, n_bins = 20, bin_size, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
M_folder |
(character) name of the folder containing files representing the accessible area (M) for all species to be analyzed. See details. |
M_format |
format of files representing the accessible area (M) for the
species. Names of M files must match the ones for occurrence files in
|
occ_folder |
(character) name of the folder containing csv files of
occurrence data for all species. Names of csv files must match the ones of M
files in |
longitude |
(character) name of the column in occurrence files containing values of longitude. |
latitude |
(character) name of the column in occurrence files containing values of latitude. |
var_folder |
(character) name of the folder containing layers to represent environmental variables. |
var_format |
format of layers to represent environmental variables.
Format options are all the ones supported by |
round |
(logical) whether or not to round the values of one or more
variables after multiplying them times the value in |
round_names |
(character) names of the variables to be rounded.
Default = NULL. If |
multiplication_factor |
(numeric) value to be used to multiply the
variables defined in |
percentage_out |
(numeric) percentage of extreme environmental data in M to be excluded in bin creation for further analyses. See details. Default = 5. |
n_bins |
(numeric) number of bins to be created from the range of environmental values considered when creating each character in bin tables. Default = 20. See details. |
bin_size |
(numeric) argument deprecated, use n_bins instead. |
save |
(logical) whether or not to save the results in working directory. Default = FALSE. |
output_directory |
(character) name of the folder in which results will be written. |
overwrite |
(logical) whether or not to overwrite existing results in
|
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Coordinates in csv files in occ_folder
, SpatVector-like files in
M_folder
, and raster layers in var_folder
must coincide in the
geographic projection in which they are represented. WGS84 with no planar
projection is recommended.
Accessible area (M) is understood as the geographic area that has been accessible for a species for relevant periods of time. Defining M is usually a hard task, but also a very important one, because it allows identifying uncertainties about the ability of a species to maintain populations in certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
Rounding variables may be useful when multiple variables are considered and
the values of some or all of them are too small (e.g., when using principal
components). To round specific variables arguments round
,
round_names
, and multiplication_factor
, must be used accordingly.
The percentage to be defined in percentage_out
excludes a percentage
of extreme environmental values to prevent from considering extremely rare
environmental values in the accessible area for the species (M). Being too
rare, these values may have never been explored by the species; therefore,
including them in the process of preparation of the table of characters
(bin table) is risky.
The argument n_bins
helps to define how many characters (bins) will be
considered for the range of values in each variable. This is, a value of 20
determines that a range of temperature (5-25) will be split approximately
every 1 degree. The argument bin_size
has been deprecated.
A list named as the variables present in var_folder
, containing all
tables of characters. A folder named as in output_directory
containing
all resultant csv files with the tables of characters will be created if
save
is set as TRUE.
Potential values for characters are:
"1" = the species is present in those environmental conditions.
"0" = the species is not present in those environmental conditions. This is, those environmental conditions inside the accessible area (M) are more extreme than the ones used for the species.
"?" = there is no certainty about the species presence in those environmental conditions. This happens in environmental combinations more extreme than the ones found in the accessible area (M), when environmental conditions in species records are as extreme as the most extreme ones in M.
# preparing data and directories for example ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) # preparing tables tabs <- bin_tables0(M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif")
# preparing data and directories for example ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) # preparing tables tabs <- bin_tables0(M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif")
A character table representing species ecological niches derived from previous preparation processes. Each row represents a species and each column a binary character in which one or more values of the environmental variable are categorized as used "1", non used "0", or uncertain "?".
character_table
character_table
A character matrix with 6 rows and 28 columns.
data("character_table", package = "nichevol") head(character_table)
data("character_table", package = "nichevol") head(character_table)
hist_evalues helps in creating histograms to explore environmental conditions in M, lines for the confidence limits of values in M, and the location of values in occurrence records, for one species at the time.
hist_evalues(M, occurrences, species, longitude, latitude, variable, CL_lines = c(95, 99), col = NULL)
hist_evalues(M, occurrences, species, longitude, latitude, variable, CL_lines = c(95, 99), col = NULL)
M |
a SpatVector object representing the accessible area (M) for one species. See details. |
occurrences |
a data.frame of occurrence records for one species. See details. |
species |
(character) name of the column in |
longitude |
(character) name of the column in |
latitude |
(character) name of the column in |
variable |
a single SpatRaster layer representing an environmental variable of interest. See details. |
CL_lines |
(numeric) confidence limits of environmental values in M to be plotted as lines in the histograms. See details. Default = c(95, 99). |
col |
colors for lines representing confidence limits. If NULL, colors are selected from a gray palette. Default = NULL. |
Coordinates in occurrences
, SpatVector object in M
, and
SpatRaster in variable
must coincide in the geographic projection in
which they are represented. WGS84 with no planar projection is recommended.
The accessible area (M) is understood as the geographic area that has been accessible to a species over relevant periods of time. Defining M is usually a hard task, but also a very important one because it allows identifying uncertainties about the ability of a species to maintain populations under certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # running stats hist_evalues(M = m_list[[1]], occurrences = occ_list[[1]], species = "species", longitude = "x", latitude = "y", variable = temp, CL_lines = c(95, 99), col = c("blue", "red"))
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # running stats hist_evalues(M = m_list[[1]], occurrences = occ_list[[1]], species = "species", longitude = "x", latitude = "y", variable = temp, CL_lines = c(95, 99), col = c("blue", "red"))
histograms_env creates PDF files with histogram plots of environmental conditions in M, lines for the confidence limits of values in M, and the location of values in occurrence records. This is done using data read directly from a local directory, and can be applied to various species and multiple variables.
histograms_env(M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, CL_lines = c(95, 99), col = NULL, round = FALSE, round_names = NULL, multiplication_factor = 1, save_ranges = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
histograms_env(M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, CL_lines = c(95, 99), col = NULL, round = FALSE, round_names = NULL, multiplication_factor = 1, save_ranges = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
M_folder |
(character) name of the folder containing files representing the accessible area (M) for all species to be analyzed. See details. |
M_format |
format of files representing the accessible area (M) for the
species. Names of M files must match the ones for occurrence files in
|
occ_folder |
(character) name of the folder containing csv files of
occurrence data for all species. Names of csv files must match the ones of M
files in |
longitude |
(character) name of the column in occurrence files containing values of longitude. |
latitude |
(character) name of the column in occurrence files containing values of latitude. |
var_folder |
(character) name of the folder containing layers to represent environmental variables. |
var_format |
format of layers to represent environmental variables.
Format options are all the ones supported by |
CL_lines |
(numeric) confidence limits of environmental values in M to be plotted as lines in the histograms. See details. Default = c(95, 99). |
col |
colors for lines representing confidence limits. If NULL, colors are selected from a gray palette. Default = NULL. |
round |
(logical) whether or not to round values of one or more
variables after multiplying them times the value in |
round_names |
(character) names of the variables to be rounded.
Default = NULL. If |
multiplication_factor |
(numeric) value to be used to multiply the
variables defined in |
save_ranges |
(logical) whether or not to save the values identified as
ranges considering the whole set of values and confidence limits defined in
|
output_directory |
(character) name of the folder in which results will be written. |
overwrite |
(logical) whether or not to overwrite existing results in
|
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Coordinates in csv files in occ_folder
, SpatVector-like files in
M_folder
, and raster layers in var_folder
must coincide in the
geographic projection in which they are represented. WGS84 with no planar
projection is recommended.
Accessible area (M) is understood as the geographic area that has been accessible for a species for relevant periods of time. Defining M is usually a hard task, but also a very important one, because it allows identifying uncertainties about the ability of a species to maintain populations under certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
Rounding variables may be useful when multiple variables are considered and
the values of some or all of them are too small (e.g., when using principal
components). To round specific variables arguments round
,
round_names
, and multiplication_factor
, must be used accordingly.
A list of data.frames containing intervals of environmental values in species
occurrences and accessible areas (M), as well as values corresponding to the
confidence limits defined in CL_lines
. A folder named as
in output_directory
containing all resulting PDF files (one per
variable) with histograms for all species. Files (csv) of ranges found during
the analyses will be also written in output_directory
if
save_ranges
is set as TRUE.
# preparing data and directories for examples ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) histdir <- paste0(tempdir, "/Hists") ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) # running analysis to produce plots hists <- histograms_env(M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif", output_directory = histdir)
# preparing data and directories for examples ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) histdir <- paste0(tempdir, "/Hists") ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) # running analysis to produce plots hists <- histograms_env(M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif", output_directory = histdir)
A SpatVector object representing the accessible area for a species.
A SpatVector object.
No return value, used with function vect
to
bring an example of an accessible area for a species.
m1 <- terra::vect(system.file("extdata", "m1.gpkg", package = "nichevol")) terra::plot(m1)
m1 <- terra::vect(system.file("extdata", "m1.gpkg", package = "nichevol")) terra::plot(m1)
map_nichevol produces a SpatRaster layer representing geographic areas corresponding to environmental bins of niche or events of niche evolution detected in reconstructions.
map_nichevol(whole_rec_table, variable, return = "niche", from, to = NULL, id_unknown = TRUE, verbose = TRUE)
map_nichevol(whole_rec_table, variable, return = "niche", from, to = NULL, id_unknown = TRUE, verbose = TRUE)
whole_rec_table |
matrix of environmental bins for all tips and nodes
derived from functions |
variable |
a SpatRaster layer corresponding to the variable for which
the reconstruction was performed (represented in |
return |
(character) type of result to return. Options are: "niche",
"evolution", or "nichevol" (a combination of both). Default = "niche". If
"niche", values correspond to that defined in |
from |
(character) if |
to |
(character) valid if |
id_unknown |
(logical) whether to identify areas of unknown or uncertain change. Default = TRUE. See details. |
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Mapping is done following Cobos et al. (2021) doi:10.1111/jav.02868. This allows to represent geographic areas with environments where niche expanded, retracted, or stayed stable (evolution). Niche is represented as presence, absence, or unknown.
Defining id_unknown
= TRUE allows to map areas where niche or niche
change are uncertain. id_unknown
= FALSE returns NA in areas with these
characteristics, hence they will not be visible when plotting the resulting
map.
A SpatRaster object classified according to values of niche in
whole_rec_table
, and/or according to niche changes detected in
comparisons between an ancestor and a tip, or another more recent ancestor.
Options of values resulting from classifications are as follow:
If return
= "niche":
ID | category |
0 | Absent |
10 | Unknown |
100 | Present |
If return
= "evolution":
ID | category |
0 | Stable |
1 | Expansion low |
3 | Expansion high |
2 | Retraction high |
4 | Retraction low |
10 | Unknown |
If return
= "nichevol":
ID | category |
0 | Stable |
1 | Expansion low |
3 | Expansion high |
10 | Unknown |
100 | Present |
102 | Retraction high |
104 | Retraction low |
# a tree data("tree", package = "nichevol") # raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # results from reconstruction data("par_rec_table", package = "nichevol") # rename tree tips tree$tip.label <- rownames(par_rec_table)[1:6] # check in plot plot.phylo(tree, label.offset = 0.02) nodelabels() nichevol_labels(tree, par_rec_table) # mapping nichevol nevol_map <- map_nichevol(whole_rec_table = par_rec_table, variable = temp, return = "nichevol", from = "9", to = "RD 6933") terra::plot(nevol_map)
# a tree data("tree", package = "nichevol") # raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # results from reconstruction data("par_rec_table", package = "nichevol") # rename tree tips tree$tip.label <- rownames(par_rec_table)[1:6] # check in plot plot.phylo(tree, label.offset = 0.02) nodelabels() nichevol_labels(tree, par_rec_table) # mapping nichevol nevol_map <- map_nichevol(whole_rec_table = par_rec_table, variable = temp, return = "nichevol", from = "9", to = "RD 6933") terra::plot(nevol_map)
niche_bars produces bar plots that represent species ecological niches in one environmental variable. Bars are exported as png figures to an output directory for posterior use.
niche_bars(tree, whole_rec_table, present = "1", unknown = "?", present_col = "#e41a1c", unknown_col = "#969696", absent_col = "#377eb8", width = 50, height = 5, res = 300, output_directory, overwrite = FALSE)
niche_bars(tree, whole_rec_table, present = "1", unknown = "?", present_col = "#e41a1c", unknown_col = "#969696", absent_col = "#377eb8", width = 50, height = 5, res = 300, output_directory, overwrite = FALSE)
tree |
an object of class "phylo". |
whole_rec_table |
matrix of environmental bins for all tips and nodes
derived from functions |
present |
(character) code indicating environmental bins in which the species is present. Default = "1". |
unknown |
(character) code indicating environmental bins in which the species presence is unknown (uncertain). Default = "?". |
present_col |
color for area of the bar representing environments where the species is present. Default = "#e41a1c". |
unknown_col |
color for area of the bar representing environments where the species presence is unknown (uncertain). Default = "#969696". |
absent_col |
color for area of the bar representing environments where no change has been detected. Default = "#377eb8". |
width |
(numeric) width of the device in mm to be passed to the
|
height |
(numeric) height of the device in mm to be passed to the
|
res |
(numeric) nominal resolution in ppi to be passed to the
|
output_directory |
(character) name of the folder in which results will be written. The directory will be created as part of the process. |
overwrite |
(logical) whether or not to overwrite existing results in
|
Ecological niches are represented in one environmental dimension with vertical bars that indicate if the species is present, absent, or if its presence is uncertain in the range of environmental conditions. Lower values of environmental variables are represented in the left part of the bar, and the opposite part of the bar represents higher values.
A folder named as in output_directory
containing all bar figures
produced, as well as a legend to describe what is plotted.
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # the running (before running, define a working directory) niche_bars(tree5, rec_tab, output_directory = file.path(tempdir(), "nichebars"))
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # the running (before running, define a working directory) niche_bars(tree5, rec_tab, output_directory = file.path(tempdir(), "nichebars"))
niche_labels helps in adding bar-type labels that represent species ecological niches in one environmental variable.
niche_labels(tree, whole_rec_table, label_type = "tip_node", tip_offset = 0.015, present = "1", unknown = "?", present_col = "#e41a1c", unknown_col = "#969696", absent_col = "#377eb8", width = 1, height = 1)
niche_labels(tree, whole_rec_table, label_type = "tip_node", tip_offset = 0.015, present = "1", unknown = "?", present_col = "#e41a1c", unknown_col = "#969696", absent_col = "#377eb8", width = 1, height = 1)
tree |
an object of class "phylo". |
whole_rec_table |
matrix of environmental bins for all tips and nodes
derived from functions |
label_type |
(character) type of label; options are: "tip", "node", and "tip_node". Default = "tip_node". |
tip_offset |
(numeric) space between tips and the labels. Default = 0.015. |
present |
(character) code indicating environmental bins in which the species is present. Default = "1". |
unknown |
(character) code indicating environmental bins in which the species presence is unknown (uncertain). Default = "?". |
present_col |
color for area of the bar representing environments where the species is present. Default = "#e41a1c". |
unknown_col |
color for area of the bar representing environments where the species presence is unknown (uncertain). Default = "#969696". |
absent_col |
color for area of the bar representing environments where no change has been detected. Default = "#377eb8". |
width |
value defining the width of niche bars; default = 1. |
height |
value defining the height of niche bars; default = 1. |
For the moment, only plots of type "phylogram" with "rightwards" or "leftwards"
directions, created with the function plot.phylo
from the
package ape
are supported.
Ecological niches are represented in one environmental dimension with vertical bars that indicate if the species is present, absent, or if its presence is uncertain in the range of environmental conditions. Lower values of environmental variables are represented in the lower part of the bar, and the opposite part of the bar represents higher values.
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels ape::plot.phylo(tree5, label.offset = 0.04) niche_labels(tree5, rec_tab, height = 0.6)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels ape::plot.phylo(tree5, label.offset = 0.04) niche_labels(tree5, rec_tab, height = 0.6)
Legends for niche labels in phylogenetic trees
niche_legend(position, legend = c("Uncertain", "Present", "Not present"), pch = 22, pt.bg = c("#969696", "#e41a1c", "#377eb8"), col = "transparent", pt.cex = 2.2, bty = "n", ...)
niche_legend(position, legend = c("Uncertain", "Present", "Not present"), pch = 22, pt.bg = c("#969696", "#e41a1c", "#377eb8"), col = "transparent", pt.cex = 2.2, bty = "n", ...)
position |
(character or numeric) position of legend. If character,
part of the plot (e.g., "topleft"), see |
legend |
(character) vector of length = three indicating the text to identify environments with uncertain presence, presence, and absence of the species. Default = c("Uncertain", "Present", "Not present"). |
pch |
point type as in |
pt.bg |
colors to represent what is in |
col |
border of symbol (points). Default = "transparent". |
pt.cex |
size of symbol (points). Default = 2.2. |
bty |
legend border type. Default = "n". |
... |
Other arguments from function |
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels and legend ape::plot.phylo(tree5, label.offset = 0.04) niche_labels(tree5, rec_tab, height = 0.6) niche_legend(position = "topleft", cex = 0.7)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels and legend ape::plot.phylo(tree5, label.offset = 0.04) niche_labels(tree5, rec_tab, height = 0.6) niche_legend(position = "topleft", cex = 0.7)
nichevol is a collection of tools that allow users to perform critical steps in the process of assessing ecological niche evolution over phylogenies, with uncertainty incorporated explicitly in reconstructions. The method proposed here for ancestral reconstruction of ecological niches characterizes species' niches using a bin-based approach that incorporates uncertainty in estimations. Compared to other existing methods, the approaches presented here reduce risk of overestimation of amounts and rates of ecological niche evolution. The main analyses include: initial exploration of environmental data in occurrence records and accessible areas, preparation of data for phylogenetic analyses, executing comparative phylogenetic analyses of ecological niches, and plotting for interpretations.
bin_ml_rec
, bin_par_rec
, bin_table
,
bin_tables
, bin_tables0
, hist_evalues
,
histograms_env
, map_nichevol
,
niche_bars
, nichevol_bars
,
niche_labels
, nichevol_labels
,
niche_legend
, nichevol_legend
,
set_uncertainty
, smooth_rec
,
stats_eval
, stats_evalues
Other functions (important helpers)
bin_env
, pdf_histograms
,
rename_tips
, score_tip
,
score_tree
, sig_sq
nichevol_bars produces bar plots that represent how species' niches (considering one environmental variable at a time) have evolved. Bars are exported as png figures to an output directory for posterior use.
nichevol_bars(tree, whole_rec_table, ancestor_line = FALSE, present = "1", absent = "0", unknown = "?", present_col = "#252525", unknown_col = "#d9d9d9", no_change_col = "#b2df8a", retraction_col = "#984ea3", expansion_col = "#4daf4a", width = 50, height = 5, res = 300, output_directory, overwrite = FALSE)
nichevol_bars(tree, whole_rec_table, ancestor_line = FALSE, present = "1", absent = "0", unknown = "?", present_col = "#252525", unknown_col = "#d9d9d9", no_change_col = "#b2df8a", retraction_col = "#984ea3", expansion_col = "#4daf4a", width = 50, height = 5, res = 300, output_directory, overwrite = FALSE)
tree |
an object of class "phylo". |
whole_rec_table |
matrix of reconstructed bins for nodes and species derived from a process of maximum parsimony reconstruction. |
ancestor_line |
controls whether ancestor line is plotted. Default = FALSE. |
present |
(character) code indicating environmental bins in which the species is present. Default = "1". |
absent |
(character) code indicating environmental bins in which the species is absent. Default = "0". |
unknown |
(character) code indicating environmental bins in which the species presence is unknown (uncertain). Default = "?". |
present_col |
color for line representing environments where the species is present. Default = "#252525". |
unknown_col |
color for line representing environments where the species presence is unknown (uncertain). Default = "#d9d9d9". |
no_change_col |
color for area of the bar representing environments where no change has been detected. Default = "#b2df8a". |
retraction_col |
color for area of the bar representing environments where niche retraction has been detected. Default = "#984ea3". |
expansion_col |
color for area of the bar representing environments where niche expansion has been detected. Default = "#4daf4a". |
width |
(numeric) width of the device in mm to be passed to the
|
height |
(numeric) height of the device in mm to be passed to the
|
res |
(numeric) nominal resolution in ppi to be passed to the
|
output_directory |
(character) name of the folder in which results will be written. The directory will be created as part of the process. |
overwrite |
(logical) whether or not to overwrite existing results in
|
Evolution of ecological niches is represented in one environmental dimension with horizontal bars indicating if the niche of the descendant has expanded, retracted, or has not changed compared to its ancestor. Lower values of environmental variables are represented in the left part of the bar, higher values at the right.
Changes in niches (evolution) are defined as follows:
if (ancestor == present & descendant == absent) change <- "retraction"
if (ancestor == present & descendant == present) change <- "no_change"
if (ancestor == present & descendant == unknown) change <- "no_change"
if (ancestor == absent & descendant == present) change <- "expansion"
if (ancestor == absent & descendant == absent) change <- "no_change"
if (ancestor == absent & descendant == unknown) change <- "no_change"
if (ancestor == unknown & descendant == absent) change <- "no_change"
if (ancestor == unknown & descendant == present) change <- "no_change"
if (ancestor == unknown & descendant == unknown) change <- "no_change"
If ancestor_line
is TRUE, the ancestor line will be plotted on the bar
representing niche evolution. The line will represent where, in the range of
environmental conditions, the ancestor was present, and where its presence is
uncertain (unknown).
A folder named as in output_directory
containing all bar figures
produced, as well as a legend to describe what is plotted.
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # the running (before running, define a working directory) nichevol_bars(tree5, rec_tab, output_directory = file.path(tempdir(), "evolbars"))
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # the running (before running, define a working directory) nichevol_bars(tree5, rec_tab, output_directory = file.path(tempdir(), "evolbars"))
nichevol_labels helps in adding bar-type labels that represent how species' niches changed from ancestors to descendants.
nichevol_labels(tree, whole_rec_table, ancestor_line = FALSE, present = "1", absent = "0", unknown = "?", present_col = "#252525", unknown_col = "#d9d9d9", no_change_col = "#b2df8a", retraction_col = "#984ea3", expansion_col = "#4daf4a", width = 1, height = 1)
nichevol_labels(tree, whole_rec_table, ancestor_line = FALSE, present = "1", absent = "0", unknown = "?", present_col = "#252525", unknown_col = "#d9d9d9", no_change_col = "#b2df8a", retraction_col = "#984ea3", expansion_col = "#4daf4a", width = 1, height = 1)
tree |
an object of class "phylo". |
whole_rec_table |
matrix of reconstructed bins for nodes and species
derived from a process of maximum parsimony or maximum likelihood reconstruction.
See functions |
ancestor_line |
controls whether ancestor line is plotted. Default = FALSE. |
present |
(character) code indicating environmental bins in which the species is present. Default = "1". |
absent |
(character) code indicating environmental bins in which the species is absent. Default = "0". |
unknown |
(character) code indicating environmental bins in which the species presence is unknown (uncertain). Default = "?". |
present_col |
color for line representing environments where the species is present. Default = "#252525". |
unknown_col |
color for line representing environments where the species presence is unknown (uncertain). Default = "#d9d9d9". |
no_change_col |
color for area of the bar representing environments where no change has been detected. Default = "#b2df8a". |
retraction_col |
color for area of the bar representing environments where niche retraction has been detected. Default = "#984ea3". |
expansion_col |
color for area of the bar representing environments where niche expansion has been detected. Default = "#4daf4a". |
width |
value defining the width of bars representing changes in niches; default = 1. |
height |
value defining the height of bars representing changes in niches; default = 1. |
For the moment, only plots of type "phylogram" with "rightwards" or "leftwards"
directions, created with the function plot.phylo
from the
package ape
are supported.
Evolution of ecological niches is represented in one environmental dimension, with vertical bars indicating if the niche of the descendant has expanded, retracted, or has not changed compared to its ancestor's niche. Lower values of environmental variables are represented in the lower part of the bar, and the opposite part of the bar represents higher values.
Changes in niches (evolution) are defined as follows:
if (ancestor == present & descendant == absent) change <- "retraction"
if (ancestor == present & descendant == present) change <- "no_change"
if (ancestor == present & descendant == unknown) change <- "no_change"
if (ancestor == absent & descendant == present) change <- "expansion"
if (ancestor == absent & descendant == absent) change <- "no_change"
if (ancestor == absent & descendant == unknown) change <- "no_change"
if (ancestor == unknown & descendant == absent) change <- "no_change"
if (ancestor == unknown & descendant == present) change <- "no_change"
if (ancestor == unknown & descendant == unknown) change <- "no_change"
If ancestor_line
is TRUE, the ancestor line will be plotted on the bar
representing niche evolution. The line will represent where, in the range of
environmental conditions, the ancestor was present, and where its presence is
uncertain (unknown).
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels ape::plot.phylo(tree5, label.offset = 0.04) nichevol_labels(tree5, rec_tab, height = 0.6)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels ape::plot.phylo(tree5, label.offset = 0.04) nichevol_labels(tree5, rec_tab, height = 0.6)
Legends for niche evolution labels in phylogenetic trees
nichevol_legend(position, ancestor_line = FALSE, ancestor_legend = c("Uncertain", "Present"), evol_legend = c("No change", "Retraction", "Expansion"), ancestor_col = c("#d9d9d9", "#252525"), evol_col = c("#b2df8a", "#984ea3", "#4daf4a"), pch = 22, pt.cex = 2.2, lty = 1, lwd = 1, cex = 1, bty = "n", ...)
nichevol_legend(position, ancestor_line = FALSE, ancestor_legend = c("Uncertain", "Present"), evol_legend = c("No change", "Retraction", "Expansion"), ancestor_col = c("#d9d9d9", "#252525"), evol_col = c("#b2df8a", "#984ea3", "#4daf4a"), pch = 22, pt.cex = 2.2, lty = 1, lwd = 1, cex = 1, bty = "n", ...)
position |
(character or numeric) position of legend. If character,
part of the plot (e.g., "topleft"), see |
ancestor_line |
whether or not ancestor line was plotted. Default = FALSE. |
ancestor_legend |
(character) vector of length = two indicating the text to identify environments with uncertain presence and true presence of the species. Default = c("Uncertain", "Present"). |
evol_legend |
(character) vector of length = three indicating the text to identify environments where niches have not changed, have retracted or expanded. Default = c("No change", "Retraction", "Expansion"). |
ancestor_col |
vector of two colors to represent what is indicated in
|
evol_col |
vector of three colors to represent what is indicated in
|
pch |
point type as in |
pt.cex |
size of symbol (points). Default = 2.2. |
lty |
line type see |
lwd |
line width see |
cex |
size of all elements in legend see |
bty |
legend border type. Default = "n". |
... |
Other arguments from function |
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels and legend ape::plot.phylo(tree5, label.offset = 0.04) nichevol_labels(tree5, rec_tab, height = 0.6) nichevol_legend(position = "bottomleft", cex = 0.7)
# a simple tree data("tree5", package = "nichevol") # a matrix of niche charactes (1 = present, 0 = absent, ? = unknown) dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label # list with two objects (tree and character table) treeWdata <- geiger::treedata(tree5, dataTable) # Maximum parsimony reconstruction rec_tab <- smooth_rec(bin_par_rec(treeWdata)) # plotting and adding labels and legend ape::plot.phylo(tree5, label.offset = 0.04) nichevol_labels(tree5, rec_tab, height = 0.6) nichevol_legend(position = "bottomleft", cex = 0.7)
A list of 6 data.frames containing name and geographic coordinates for 6 species.
occ_list
occ_list
A list of 6 data.frames:
species name, a code in this example
longitude, longitude value
latitude, latitude value
data("occ_list", package = "nichevol") str(occ_list)
data("occ_list", package = "nichevol") str(occ_list)
A character table representing species ecological niches derived from previous preparation processes and reconstructed niches for ancestors. Each row represents a species or a node and each column a binary character in which one or more values of the environmental variable are categorized as used "1", non used "0", or uncertain "?".
par_rec_table
par_rec_table
A character matrix with 11 rows and 20 columns.
data("par_rec_table", package = "nichevol") head(par_rec_table)
data("par_rec_table", package = "nichevol") head(par_rec_table)
Helper function to create PDF files with histograms
pdf_histograms(env_data, occ_data, y_values, sp_names, variable_name, CL_lines, limits, col, output_directory)
pdf_histograms(env_data, occ_data, y_values, sp_names, variable_name, CL_lines, limits, col, output_directory)
env_data |
list of environmental values in M for all species. |
occ_data |
list of environmental values in occurrences for all species. |
y_values |
list of values for the y axis to be used to represent where occurrences are distributed across the environmental values in M. |
sp_names |
(character) names of the species for which the process will be performed. |
variable_name |
(character) name of the variable to be plotted. |
CL_lines |
(numeric) confidence limits to be plotted in the histograms. |
limits |
numeric matrix containing the actual values for the confidence limits of M. |
col |
color for lines representing the confidence limits of M. |
output_directory |
(character) name of the folder in which results will be written. |
A PDF file written in the output directory containing all resulting figures.
# example data e_data <- list(rnorm(1000, 15, 7), rnorm(800, 20, 6), rnorm(1000, 12, 3)) o_data <- list(sample(seq(5, 29, 0.1), 45), sample(seq(10, 33, 0.1), 40), sample(seq(1, 16, 0.1), 50)) for (i in 1:3) { names(e_data[[i]]) <- e_data[[i]] names(o_data[[i]]) <- o_data[[i]] } y_val <- list(rep(3, length(o_data)), rep(4, length(o_data)), rep(2, length(o_data))) s_names <- c("sp1", "sp2", "sp3") lims <- rbind(c(3.5, 26.47), c(10.83, 29.66), c(6.92, 16.91)) tmpd <- file.path(tempdir(), "Hist_to_check") # temporal directory dir.create(tmpd) # the running (before running, create output_directory in current directory) bins <- pdf_histograms(env_data = e_data, occ_data = o_data, y_values = y_val, sp_names = s_names, variable_name = "Temperature", CL_lines = 95, limits = lims, col = "green", output_directory = tmpd)
# example data e_data <- list(rnorm(1000, 15, 7), rnorm(800, 20, 6), rnorm(1000, 12, 3)) o_data <- list(sample(seq(5, 29, 0.1), 45), sample(seq(10, 33, 0.1), 40), sample(seq(1, 16, 0.1), 50)) for (i in 1:3) { names(e_data[[i]]) <- e_data[[i]] names(o_data[[i]]) <- o_data[[i]] } y_val <- list(rep(3, length(o_data)), rep(4, length(o_data)), rep(2, length(o_data))) s_names <- c("sp1", "sp2", "sp3") lims <- rbind(c(3.5, 26.47), c(10.83, 29.66), c(6.92, 16.91)) tmpd <- file.path(tempdir(), "Hist_to_check") # temporal directory dir.create(tmpd) # the running (before running, create output_directory in current directory) bins <- pdf_histograms(env_data = e_data, occ_data = o_data, y_values = y_val, sp_names = s_names, variable_name = "Temperature", CL_lines = 95, limits = lims, col = "green", output_directory = tmpd)
Read one or multiple tables binary niche characters from directory.
read_bin_table(file) read_bin_tables(directory)
read_bin_table(file) read_bin_tables(directory)
file |
(character) name of CSV file containing a table of binary niche characters. |
directory |
(character) name of directory where tables of binary niche characters were written as CSV files. |
A matrix if read_bin_table
is used.
A list of matrices if read_bin_tables
is used.
Helper function to rename tips of trees for simulations
rename_tips(tree, names)
rename_tips(tree, names)
tree |
an object of class "phylo". |
names |
(character) vector of new names. Length must be equal to number of tips. They will be assigned in the order given. |
Tree of class "phylo" with specified names
# a simple tree data("tree5", package = "nichevol") # renaming tips renamedTree <- rename_tips(tree5, c("a", "b", "c", "d", "e"))
# a simple tree data("tree5", package = "nichevol") # renaming tips renamedTree <- rename_tips(tree5, c("a", "b", "c", "d", "e"))
Helper function to calculate the median bin score for a given species
score_tip(character_table, species_name, include_unknown = FALSE)
score_tip(character_table, species_name, include_unknown = FALSE)
character_table |
data.frame containing bin scores for all species. NOTE: row names must be species' names. |
species_name |
(character) name of the species to be analyzed. |
include_unknown |
(logical) whether or not unknown bin status should be included. |
Median bin value for a given species (for inferring sigma squared or other comparative phylogenetic analyses requiring a single continuous variable).
# Simulate data for single number bin labels dataTable <- cbind("241" = rep("1", 5), "242" = rep("1", 5), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Simulate data for bin labels as strings dataTableStringLabel <- cbind("241 to 244" = rep("1", 5), "244 to 246" = c("1", "1", "0", "0", "0"), "246 to 248" = c("1", "?", "0", "0", "0")) rownames(dataTableStringLabel) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Use function score_tip(character_table = dataTable, species_name = "GadusMorhua", include_unknown = TRUE) score_tip(character_table = dataTableStringLabel, species_name = "GadusMorhua", include_unknown = FALSE)
# Simulate data for single number bin labels dataTable <- cbind("241" = rep("1", 5), "242" = rep("1", 5), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Simulate data for bin labels as strings dataTableStringLabel <- cbind("241 to 244" = rep("1", 5), "244 to 246" = c("1", "1", "0", "0", "0"), "246 to 248" = c("1", "?", "0", "0", "0")) rownames(dataTableStringLabel) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Use function score_tip(character_table = dataTable, species_name = "GadusMorhua", include_unknown = TRUE) score_tip(character_table = dataTableStringLabel, species_name = "GadusMorhua", include_unknown = FALSE)
Helper function to assign bin scores to every tip in a given tree
score_tree(tree_data, include_unknown = FALSE)
score_tree(tree_data, include_unknown = FALSE)
tree_data |
a list of two elements (phy and data) resulting from using
the function |
include_unknown |
(logical) whether or not there are unknown tips. |
a list of two elements (phy and data). Data is the median bin scored as present or present + unknown.
# Simulate data table dataTable <- cbind("241" = rep("1", 5), "242" = rep("1", 5), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # a simple tree data("tree5", package = "nichevol") tree5$tip.label <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Unite data treeWithData <- geiger::treedata(tree5, dataTable) # Get a new tree with tips scored from median bin scores score_tree(treeWithData, include_unknown = TRUE)
# Simulate data table dataTable <- cbind("241" = rep("1", 5), "242" = rep("1", 5), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # a simple tree data("tree5", package = "nichevol") tree5$tip.label <- c("GadusMorhua", "GadusMacrocephalus", "GadusChalcogrammus", "ArctogadusGlacials", "BoreogadusSaida") # Unite data treeWithData <- geiger::treedata(tree5, dataTable) # Get a new tree with tips scored from median bin scores score_tree(treeWithData, include_unknown = TRUE)
set_uncertainty allows to define uncertainty ("?") values around values denoting presence ("1") towards one or both ends of the variable in a table of binary characters.
set_uncertainty(character_table, species, end)
set_uncertainty(character_table, species, end)
character_table |
a matrix of characters to represent ecological niches
of the species of interest. A matrix containing values "1" = presence,
"0" = absence, and "?" = uncertain. See |
species |
(character) name of the species in the table for which values of uncertainty will be set. |
end |
(character) end towards which uncertainty values ("?") will be set. Options are: "high", "low", or "both". |
Values of characters around those denoting presence ("1") are manually transformed to uncertain ("?") to help producing more conservative reconstructions of ancestral ecological niches. This increases uncertainty in reconstructions and further niche comparisons, which reduces the events of niche change that can be detected. This may be especially useful when dealing with species with one or just a few known records.
A modified matrix of characters to represent ecological niches of the species of interest.
Potential values for characters are:
"1" = the species is present in those environmental conditions.
"0" = the species is not present in those environmental conditions. This is, those environmental conditions inside the accessible area (M) are more extreme than the ones used for the species.
"?" = there is no certainty about the species presence in those environmental conditions.
# a character table data("character_table", package = "nichevol") character_table[, 20:28] # set values of uncertainty towards the lower end of the variable for species t3 char_tableu <- set_uncertainty(character_table, species = "t2", end = "low") char_tableu[, 20:28]
# a character table data("character_table", package = "nichevol") character_table[, 20:28] # set values of uncertainty towards the lower end of the variable for species t3 char_tableu <- set_uncertainty(character_table, species = "t2", end = "low") char_tableu[, 20:28]
Sigma squared values for a single niche summary statistic
are calculated using fitContinuous
.
sig_sq(tree_data, model = "BM")
sig_sq(tree_data, model = "BM")
tree_data |
a list of two elements (phy and data) resulted from using
the function |
model |
model to fit to comparative data; see
|
the sigma squared value (evolutionary rate) for the data, given the tree.
# a simple tree data("tree5", package = "nichevol") # simple data data <- rnorm(n = length(tree5$tip.label)) names(data) <- tree5$tip.label # tree with data treeWdata <- geiger::treedata(tree5, data) # Estimating sigma squared for the dataset sig_sq(treeWdata)
# a simple tree data("tree5", package = "nichevol") # simple data data <- rnorm(n = length(tree5$tip.label)) names(data) <- tree5$tip.label # tree with data treeWdata <- geiger::treedata(tree5, data) # Estimating sigma squared for the dataset sig_sq(treeWdata)
Smooth character table values resulted from ancestral character state reconstructions
smooth_rec(whole_rec_table)
smooth_rec(whole_rec_table)
whole_rec_table |
matrix containing all reconstructed characters for all
tips and nodes. It results from using the functions |
The matrix of reconstructed characters with smoothed values.
# a simple tree data("tree5", package = "nichevol") # simple matrix of data dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label treeWdata <- geiger::treedata(tree5, dataTable) # ancestral reconstruction parsimonyReconstruction <- bin_par_rec(treeWdata) # smoothing reconstructions smooth_rec(parsimonyReconstruction)
# a simple tree data("tree5", package = "nichevol") # simple matrix of data dataTable <- cbind("241" = rep("1", length(tree5$tip.label)), "242" = rep("1", length(tree5$tip.label)), "243" = c("1", "1", "0", "0", "0"), "244" = c("1", "1", "0", "0", "0"), "245" = c("1", "?", "0", "0", "0")) rownames(dataTable) <- tree5$tip.label treeWdata <- geiger::treedata(tree5, dataTable) # ancestral reconstruction parsimonyReconstruction <- bin_par_rec(treeWdata) # smoothing reconstructions smooth_rec(parsimonyReconstruction)
stats_eval helps in creating tables of descriptive statistics of environmental conditions in accessible areas (M) and occurrence records for one environmental variable at a time.
stats_eval(stats = c("median", "range"), Ms, occurrences, species, longitude, latitude, variable, percentage_out = 0, verbose = TRUE)
stats_eval(stats = c("median", "range"), Ms, occurrences, species, longitude, latitude, variable, percentage_out = 0, verbose = TRUE)
stats |
(character) name or vector of names of functions to be applied to get basic statistics of environmental values. |
Ms |
a list of SpatVector objects representing the accessible area
(M) for each species to be analyzed. The order of species represented by each
object here must coincide with the one in |
occurrences |
a list of data.frames of occurrence records for all species.
The order of species represented by each data.frame must coincide with the one
in |
species |
(character) name of the column in occurrence data.frames that contains the name of the species. |
longitude |
(character) name of the column in occurrence files containing values of longitude. |
latitude |
(character) name of the column in occurrence files containing values of latitude. |
variable |
a single SpatRaster layer of an environmental variable of interest. See details. |
percentage_out |
(numeric) percentage of extreme environmental data in M to be excluded in bin creation for further analyses. See details. Default = 0. |
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Coordinates in occurrences
, SpatVector objects in Ms
, and
SpatRaster in variable
must coincide in the geographic projection in
which they are represented. WGS84 with no planar projection is recommended.
Accessible area (M) is understood as the geographic area that has been accessible for a species for relevant periods of time. Defining M is usually a hard task, but also a very important one, because it allows identifying uncertainties about the ability of a species to maintain populations in certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
The percentage to be defined in percentage_out
excludes a percentage
of extreme environmental values to prevent from considering extremely rare
environmental values in the accessible area for the species (M). Being too
rare, these values may have never been explored by the species; therefore,
including them in the process of preparation of the table of characters
(bin table) is risky.
A list containing tables with statistics of the values in variable
,
for the species M and occurrences.
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # running stats stat <- stats_eval(stats = c("mean", "sd", "median", "range", "quantile"), Ms = m_list, occurrences = occ_list, species = "species", longitude = "x", latitude = "y", variable = temp, percentage_out = 0)
# example data ## list of species records data("occ_list", package = "nichevol") ## list of species accessible areas m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) m_list <- lapply(m_files, terra::vect) ## raster variable temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) # running stats stat <- stats_eval(stats = c("mean", "sd", "median", "range", "quantile"), Ms = m_list, occurrences = occ_list, species = "species", longitude = "x", latitude = "y", variable = temp, percentage_out = 0)
stats_evalues helps in creating csv files with statistics of environmental conditions in accessible areas (M) and species occurrence records. This is done using data read directly from a local directory, and can be applied to various species and multiple variables.
stats_evalues(stats = c("median", "range"), M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, round = FALSE, round_names, multiplication_factor = 1, percentage_out = 0, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
stats_evalues(stats = c("median", "range"), M_folder, M_format, occ_folder, longitude, latitude, var_folder, var_format, round = FALSE, round_names, multiplication_factor = 1, percentage_out = 0, save = FALSE, output_directory, overwrite = FALSE, verbose = TRUE)
stats |
(character) name or vector of names of functions to be applied to get basic statistics of environmental values. |
M_folder |
(character) name of the folder containing files representing the accessible area (M) for each species to be analyzed. See details. |
M_format |
format of files representing the accessible area (M) for the
species. Names of M files must match the ones for occurrence files in
|
occ_folder |
(character) name of the folder containing csv files of
occurrence data for all species. Names of csv files must match the ones of M
files in |
longitude |
(character) name of the column in occurrence files containing values of longitude. |
latitude |
(character) name of the column in occurrence files containing values of latitude. |
var_folder |
(character) name of the folder containing layers to represent environmental variables. |
var_format |
format of layers to represent environmental variables.
Format options are all the ones supported by |
round |
(logical) whether or not to round the values of one or more
variables after multiplying them times the value in |
round_names |
(character) names of the variables to be rounded.
Default = NULL. If |
multiplication_factor |
(numeric) value to be used to multiply the
variables defined in |
percentage_out |
(numeric) percentage of extreme environmental data in M to be excluded in bin creation for further analyses. See details. Default = 0. |
save |
(logical) whether or not to save the results in working directory. Default = FALSE. |
output_directory |
(character) name of the folder in which results will be written. |
overwrite |
(logical) whether or not to overwrite existing results in
|
verbose |
(logical) whether messages should be printed. Default = TRUE. |
Coordinates in csv files in occ_folder
, SpatVector-like files in
M_folder
, and raster layers in var_folder
must coincide in the
geographic projection in which they are represented. WGS84 with no planar
projection is recommended.
Accessible area (M) is understood as the geographic area that has been accessible for a species for relevant periods of time. Defining M is usually a hard task, but also a very important one, because it allows identifying uncertainties about the ability of a species to maintain populations in certain environmental conditions. For further details on this topic, see Barve et al. (2011) doi:10.1016/j.ecolmodel.2011.02.011 and Machado-Stredel et al. (2021) doi:10.21425/F5FBG48814.
Rounding variables may be useful when multiple variables are considered and
the values of some or all of them are too small (e.g., when using principal
components). To round specific variables arguments round
,
round_names
, and multiplication_factor
, must be used accordingly.
The percentage to be defined in percentage_out
excludes a percentage
of extreme environmental values to prevent the algorithm from considering
extremely rare environmental values in the accessible area for the species (M).
Being too rare, these values may have never been explored by the species;
therefore, including them in the process of preparation of the table of
characters (bin table) is risky.
A list named as the variables present in var_folder
, containing all
tables with statistics of environmental values in M and in species records.
A folder named as in output_directory
containing all resultant csv
files with the tables of statistics will be created if save
is set as
TRUE.
# preparing data and directories for examples ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) stats <- stats_evalues(stats = c("median", "range"), M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif", percentage_out = 5)
# preparing data and directories for examples ## directories tempdir <- file.path(tempdir(), "nevol_test") dir.create(tempdir) cvariables <- paste0(tempdir, "/variables") dir.create(cvariables) records <- paste0(tempdir, "/records") dir.create(records) m_areas <- paste0(tempdir, "/M_areas") dir.create(m_areas) ## data data("occ_list", package = "nichevol") temp <- system.file("extdata", "temp.tif", package = "nichevol") m_files <- list.files(system.file("extdata", package = "nichevol"), pattern = "m\\d.gpkg", full.names = TRUE) ## writing data in temporal directories spnames <- sapply(occ_list, function (x) as.character(x[1, 1])) ocnames <- paste0(records, "/", spnames, ".csv") occs <- lapply(1:length(spnames), function (x) { write.csv(occ_list[[x]], ocnames[x], row.names = FALSE) }) to_replace <- paste0(system.file("extdata", package = "nichevol"), "/") otemp <- gsub(to_replace, "", temp) file.copy(from = temp, to = paste0(cvariables, "/", otemp)) file.copy(from = m_files, to = paste0(m_areas, "/", spnames, ".gpkg")) stats <- stats_evalues(stats = c("median", "range"), M_folder = m_areas, M_format = "gpkg", occ_folder = records, longitude = "x", latitude = "y", var_folder = cvariables, var_format = "tif", percentage_out = 5)
A SpatRaster object representing the variable temperature.
A SpatRaster object.
No return value, used with function rast
to
bring an example of an environmental variable used in analysis.
temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) terra::plot(temp)
temp <- terra::rast(system.file("extdata", "temp.tif", package = "nichevol")) terra::plot(temp)
A phylogenetic tree with 6 species and their relationships.
tree
tree
An object of class phylo for 6 species.
data("tree", package = "nichevol") str(tree)
data("tree", package = "nichevol") str(tree)
A list of 2 elements (phy and data) resulting from using the function
treedata
.
tree_data
tree_data
A list of 2 elements:
object of class phylo for 6 species
matrix of 6 rows and 28 columns
data("tree_data", package = "nichevol") str(tree_data)
data("tree_data", package = "nichevol") str(tree_data)
A phylogenetic tree with 5 species and their relationships.
tree5
tree5
An object of class phylo for 5 species.
data("tree5", package = "nichevol") str(tree5)
data("tree5", package = "nichevol") str(tree5)