View source: R/input_creation.R
create_input_landscape | R Documentation |
create an landscape input from a named list of rasters or raster files
create_input_landscape(
landscapes,
cost_function,
directions,
output_directory,
timesteps = NULL,
calculate_full_distance_matrices = FALSE,
crs = NULL,
overwrite_output = FALSE,
verbose = FALSE
)
landscapes |
list of named list(s) of raster(s) or raster file(s) name(s). Starting from the present towards the past. NOTE: the list names are important since these are the environmental names |
cost_function |
function that returns a cost value between a pair of sites (neighbors) that should have the following signature:
|
directions |
4, 8 or 16 neighbors, dictates the connection of cell neighbors on adjacency matrix (see gistance package) |
output_directory |
path for storing the gen3sis ready landscape (i.e. landscape.rds, metadata.txt and full- and/or local_distance folders) |
timesteps |
vector of names for every time-step to represent the time-step at gen3sis ready landscape. If timesteps=NULL (default), time-steps are sequentially numbered from 0 to the latest time-step. |
calculate_full_distance_matrices |
should a full distance matrix be calculated? TRUE or FALSE? If TRUE calculates the entire distance matrix for every time-step and between all habitable cells (faster CPU time, higher storage required). If FALSE (default), only local distances are calculated (slower CPU time when simulating but smaller gen3sis landscape size) |
crs |
the coordinate reference system in crs format (see raster::crs) |
overwrite_output |
TRUE or FALSE |
verbose |
print distance calculation progress (default: FALSE) |
This function creates the input landscapes files needed by the run_simulation function. It uses as input the dynamic landscape rasters and user defined geodesimal corrections as well as rules to define the connection costs between sites
no return object. This function saves the landscape input files for gen3sis at the output_directory
run_simulation
# load needed library
library(raster)
# get path containing example rasters
datapath <- system.file(file.path("extdata", "WorldCenter"), package="gen3sis")
# create raster bricks
temperature_brick <- brick(file.path(datapath, "input_rasters/temp_rasters.grd"))
aridity_brick <- brick(file.path(datapath, "input_rasters/arid_rasters.grd"))
area_brick <- brick(file.path(datapath, "input_rasters/area_rasters.grd"))
# create sub-list of environmental variables for fast example
# (i.e. 4 time-steps)
landscapes_sub_list <- list(temp=NULL, arid=NULL, area=NULL)
for(i in 1:4){
landscapes_sub_list$temp <- c(landscapes_sub_list$temp, temperature_brick[[i]])
landscapes_sub_list$arid <- c(landscapes_sub_list$arid, aridity_brick[[i]])
landscapes_sub_list$area <- c(landscapes_sub_list$area, area_brick[[i]])
}
# define cost function, crossing water as double as land sites
cost_function_water <- function(source, habitable_src, dest, habitable_dest) {
if(!all(habitable_src, habitable_dest)) {
return(2/1000)
} else {
return(1/1000)
}
}
## Not run:
# create input landscape ready for gen3sis from sub-list
# (i.e. 10 time-steps) and only local-distances.
create_input_landscape(
landscapes = landscapes_sub_list,
cost_function = cost_function_water,
output_directory = file.path(tempdir(), "landscape_sub"),
directions = 8, # surrounding sites for each site
timesteps = paste0(round(150:147,2), "Ma"),
calculate_full_distance_matrices = FALSE) # full distance matrix
# create list of all environmental variables available
landscapes_list <- list(temp=NULL, arid=NULL, area=NULL)
for(i in 1:nlayers(temperature_brick)){
landscapes_list$temp <- c(landscapes_list$temp, temperature_brick[[i]])
landscapes_list$arid <- c(landscapes_list$arid, aridity_brick[[i]])
landscapes_list$area <- c(landscapes_list$area, area_brick[[i]])
}
# create input landscape ready for gen3sis (~ 3min run-time)
# and full distance matrix
create_input_landscape(
landscapes = landscapes_list,
cost_function = cost_function_water,
output_directory = file.path(tempdir(), "landscape_WorldCenter_5"),
directions = 8, # surrounding sites for each site
timesteps = paste0(round(150:100,2), "Ma"),
crs="+proj=longlat +datum=WGS84 +no_defs +ellps=WGS84 +towgs84=0,0,0",
calculate_full_distance_matrices = FALSE) # full distance matrix
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.