Title: | Applying Landscape Genomic Methods on 'SNP' and 'Silicodart' Data |
---|---|
Description: | Provides landscape genomic functions to analyse 'SNP' (single nuclear polymorphism) data, such as least cost path analysis and isolation by distance. Therefore each sample needs to have coordinate data attached (lat/lon) to be able to run most of the functions. 'dartR.spatial' is a package that belongs to the 'dartRverse' suit of packages and depends on 'dartR.base' and 'dartR.data'. |
Authors: | Bernd Gruber [aut, cre], Arthur Georges [aut], Jose L. Mijangos [aut], Carlo Pacioni [aut], Peter J. Unmack [ctb], Oliver Berry [ctb] |
Maintainer: | Bernd Gruber <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.78 |
Built: | 2025-01-10 04:18:04 UTC |
Source: | https://github.com/cran/dartR.spatial |
Calculates a cost distance matrix, to be used with run.popgensim.
gl.costdistances(landscape, locs, method, NN, verbose = NULL)
gl.costdistances(landscape, locs, method, NN, verbose = NULL)
landscape |
A raster object coding the resistance of the landscape [required]. |
locs |
Coordinates of the subpopulations. If a genlight object is
provided coordinates are taken from @other$latlon and centers for population
(pop(gl)) are calculated. In case you want to calculate costdistances between
individuals redefine pop(gl) via: |
method |
Defines the type of cost distance, types are 'leastcost', 'rSPDistance' or 'commute' (Circuitscape type) [required]. |
NN |
Number of next neighbours recommendation is 8 [required]. |
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, progress log ; 3, progress and results summary; 5, full report [default 2, unless specified using gl.set.verbosity]. |
A costdistance matrix between all pairs of locs.
data(possums.gl) library(raster) #needed for that example landscape.sim <- readRDS(system.file('extdata','landscape.sim.rdata', package='dartR.data')) #calculate mean centers of individuals per population xy <- apply(possums.gl@other$xy, 2, function(x) tapply(x, pop(possums.gl), mean)) cd <- gl.costdistances(landscape.sim, xy, method='leastcost', NN=8) round(cd,3)
data(possums.gl) library(raster) #needed for that example landscape.sim <- readRDS(system.file('extdata','landscape.sim.rdata', package='dartR.data')) #calculate mean centers of individuals per population xy <- apply(possums.gl@other$xy, 2, function(x) tapply(x, pop(possums.gl), mean)) cd <- gl.costdistances(landscape.sim, xy, method='leastcost', NN=8) round(cd,3)
This function calculates the pairwise distances (Euclidean, cost path distances and genetic distances) of populations using a friction matrix and a spatial genind object. The genind object needs to have coordinates in the same projected coordinate system as the friction matrix. The friction matrix can be either a single raster of a stack of several layers. If a stack is provided the specified cost distance is calculated for each layer in the stack. The output of this function can be used with the functions wassermann from package PopGenReport and lgrMMRR from package PopGenReport to test for the significance of a layer on the genetic structure.
gl.genleastcost( x, fric.raster, gen.distance, NN = NULL, pathtype = "leastcost", plotpath = TRUE, theta = 1, verbose = NULL )
gl.genleastcost( x, fric.raster, gen.distance, NN = NULL, pathtype = "leastcost", plotpath = TRUE, theta = 1, verbose = NULL )
x |
A spatial gelight object. [required]. |
fric.raster |
A friction matrix [required]. |
gen.distance |
Specification which genetic distance method should be used to calculate pairwise genetic distances between populations ( 'D', 'Gst.Nei', 'Gst.Hedrick') or individuals ('kosman', 'propShared') [required]. |
NN |
Number of neighbours used when calculating the cost distance (possible values 4, 8 or 16). As the default is NULL a value has to be provided if pathtype='leastcost'. NN=8 is most commonly used. Be aware that linear structures may cause artefacts in the least-cost paths, therefore inspect the actual least-cost paths in the provided output [default NULL]. |
pathtype |
Type of cost distance to be calculated (based on function in
the |
plotpath |
switch if least cost paths should be plotted (works only if pathtype='leastcost'. Be aware this slows down the computation, but it is recommended to do this to check least cost paths visually. |
theta |
value needed for rSPDistance function. See
|
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, progress log; 3, progress and results summary; 5, full report [default 2, unless specified using gl.set.verbosity]. |
Returns a list that consists of four pairwise distance matrices (Euclidean, Cost, length of path and genetic) and the actual paths as spatial line objects.
Bernd Gruber (bugs? Post to https://groups.google.com/d/forum/dartr)
Cushman, S., Wasserman, T., Landguth, E. and Shirk, A. (2013). Re-Evaluating Causal Modeling with Mantel Tests in Landscape Genetics. Diversity, 5(1), 51-72.
Landguth, E. L., Cushman, S. A., Schwartz, M. K., McKelvey, K. S., Murphy, M. and Luikart, G. (2010). Quantifying the lag time to detect barriers in landscape genetics. Molecular ecology, 4179-4191.
Wasserman, T. N., Cushman, S. A., Schwartz, M. K. and Wallin, D. O. (2010). Spatial scaling and multi-model inference in landscape genetics: Martes americana in northern Idaho. Landscape Ecology, 25(10), 1601-1612.
#this example takes about 20 seconds to run... data(possums.gl) library(raster) #needed for that example landscape.sim <- readRDS(system.file('extdata','landscape.sim.rdata', package='dartR.data')) #use only 3 population (first 90 individuals) due to speed glc <- gl.genleastcost(x=possums.gl,fric.raster=landscape.sim , gen.distance = 'D', NN=8, pathtype = 'leastcost',plotpath = TRUE) #### run tests as implemented in PopGenreport (maybe need to install) if (require("PopGenReport", quietly=TRUE)) { PopGenReport::wassermann(eucl.mat = glc$eucl.mat, cost.mat = glc$cost.mats, gen.mat = glc$gen.mat) lgrMMRR(gen.mat = glc$gen.mat, cost.mats = glc$cost.mats, eucl.mat = glc$eucl.mat) }
#this example takes about 20 seconds to run... data(possums.gl) library(raster) #needed for that example landscape.sim <- readRDS(system.file('extdata','landscape.sim.rdata', package='dartR.data')) #use only 3 population (first 90 individuals) due to speed glc <- gl.genleastcost(x=possums.gl,fric.raster=landscape.sim , gen.distance = 'D', NN=8, pathtype = 'leastcost',plotpath = TRUE) #### run tests as implemented in PopGenreport (maybe need to install) if (require("PopGenReport", quietly=TRUE)) { PopGenReport::wassermann(eucl.mat = glc$eucl.mat, cost.mat = glc$cost.mats, gen.mat = glc$gen.mat) lgrMMRR(gen.mat = glc$gen.mat, cost.mats = glc$cost.mats, eucl.mat = glc$eucl.mat) }
This function calculates the mean probability of identity by state (IBS) across loci that would result from all the possible crosses of the individuals analyzed. IBD is calculated by an additive relationship matrix approach developed by Endelman and Jannink (2012) as implemented in the function A.mat (package rrBLUP).
gl.grm2( x, plotheatmap = TRUE, palette_discrete = NULL, palette_convergent = NULL, legendx = 0, legendy = 0.5, plot.file = NULL, plot.dir = NULL, verbose = NULL, ... )
gl.grm2( x, plotheatmap = TRUE, palette_discrete = NULL, palette_convergent = NULL, legendx = 0, legendy = 0.5, plot.file = NULL, plot.dir = NULL, verbose = NULL, ... )
x |
Name of the genlight object containing the SNP data [required]. |
plotheatmap |
A switch if a heatmap should be shown [default TRUE]. |
palette_discrete |
the color of populations [gl.select.colors]. |
palette_convergent |
A convergent palette for the IBD values [default convergent_palette]. |
legendx |
x coordinates for the legend[default 0]. |
legendy |
y coordinates for the legend[default 1]. |
plot.file |
Name for the RDS binary file to save (base name only, exclude extension) [default NULL] |
plot.dir |
Directory in which to save files [default = working directory] |
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, progress log ; 3, progress and results summary; 5, full report [default 2 or as specified using gl.set.verbosity]. |
... |
Parameters passed to function A.mat from package rrBLUP. |
Two or more alleles are identical by descent (IBD) if they are identical copies of the same ancestral allele in a base population. The additive relationship matrix is a theoretical framework for estimating a relationship matrix that is consistent with an approach to estimate the probability that the alleles at a random locus are identical in state (IBS).
This function also plots a heatmap, and a dendrogram, of IBD values where each diagonal element has a mean that equals 1+f, where f is the inbreeding coefficient (i.e. the probability that the two alleles at a randomly chosen locus are IBD from the base population). As this probability lies between 0 and 1, the diagonal elements range from 1 to 2. Because the inbreeding coefficients are expressed relative to the current population, the mean of the off-diagonal elements is -(1+f)/n, where n is the number of loci. Individual names are shown in the margins of the heatmap and colors represent different populations.
An identity by descent matrix
Custodian: Arthur Georges – Post to https://groups.google.com/d/forum/dartr
Endelman, J. B. (2011). Ridge regression and other kernels for genomic selection with r package rrblup. The Plant Genome 4, 250.
Endelman, J. B. , Jannink, J.-L. (2012). Shrinkage estimation of the realized relationship matrix. G3: Genes, Genomics, Genetics 2, 1405.
require(dartR.base) gl.grm2(platypus.gl[1:10,1:100])
require(dartR.base) gl.grm2(platypus.gl[1:10,1:100])
This function performs an isolation by distance analysis based on a Mantel test and also produces an isolation by distance plot. If a genlight object with coordinates is provided, then an Euclidean and genetic distance matrices are calculated.'
gl.ibd( x = NULL, distance = "Fst", coordinates = "latlon", Dgen = NULL, Dgeo = NULL, Dgeo_trans = "Dgeo", Dgen_trans = "Dgen", permutations = 999, plot.out = TRUE, paircols = NULL, plot.theme = theme_dartR(), plot.file = NULL, plot.dir = NULL, verbose = NULL )
gl.ibd( x = NULL, distance = "Fst", coordinates = "latlon", Dgen = NULL, Dgeo = NULL, Dgeo_trans = "Dgeo", Dgen_trans = "Dgen", permutations = 999, plot.out = TRUE, paircols = NULL, plot.theme = theme_dartR(), plot.file = NULL, plot.dir = NULL, verbose = NULL )
x |
Genlight object. If provided a standard analysis on Fst/1-Fst and log(distance) is performed [required]. |
distance |
Type of distance that is calculated and used for the analysis. Can be either population based 'Fst' [stamppFst], 'D' [stamppNeisD] or individual based 'propShared', [gl.propShared], 'euclidean' [gl.dist.ind, method='Euclidean'], 'kosman' [gl.kosman] [default "Fst"]. |
coordinates |
Can be either 'latlon', 'xy' or a two column data.frame
with column names 'lat','lon', 'x', 'y'). Coordinates are provided via
|
Dgen |
Genetic distance matrix if no genlight object is provided [default NULL]. |
Dgeo |
Euclidean distance matrix if no genlight object is provided [default NULL]. |
Dgeo_trans |
Transformation to be used on the Euclidean distances. See Dgen_trans [default "Dgeo"]. |
Dgen_trans |
You can provide a formula to transform the genetic
distance. The transformation can be applied as a formula using Dgen as the
variable to be transformed. For example: |
permutations |
Number of permutations in the Mantel test [default 999]. |
plot.out |
Should an isolation by distance plot be returned [default TRUE]. |
paircols |
Should pairwise dots colored by 'pop'ulation/'ind'ividual pairs [default 'pop']. You can color pairwise individuals by pairwise population colors. |
plot.theme |
Theme for the plot. See details for options [default theme_dartR()]. |
plot.file |
Name for the RDS binary file to save (base name only, exclude extension) [default NULL] |
plot.dir |
Directory in which to save files [default = working directory] |
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, progress log ; 3, progress and results summary; 5, full report [default 2 or as specified using gl.set.verbosity]. |
Currently pairwise Fst and D between populations and 1-propShared and Euclidean distance between individuals are implemented. Coordinates are expected as lat long and converted to Google Earth Mercator projection. If coordinates are already projected, provide them at the x@other$xy slot.
You can provide also your own genetic and Euclidean distance matrices. The function is based on the code provided by the adegenet tutorial (http://adegenet.r-forge.r-project.org/files/tutorial-basics.pdf), using the functions mantel (package vegan), stamppFst, stamppNeisD (package StAMPP) and gl.propShared or gl.dist.ind. For transformation you need to have the dismo package installed. As a new feature you can plot pairwise relationship using double colored points (paircols=TRUE). Pairwise relationship can be visualised via populations or individuals, depending which distance is calculated. Please note: Often a problem arises, if an individual based distance is calculated (e.g. propShared) and some individuals have identical coordinates as this results in distances of zero between those pairs of individuals.
If the standard transformation [log(Dgeo)] is used, this results in an infinite value, because of trying to calculate'log(0)'. To avoid this, the easiest fix is to change the transformation from log(Dgeo) to log(Dgeo+1) or you could add some "noise" to the coordinates of the individuals (e.g. +- 1m, but be aware if you use lat lon then you rather want to add +0.00001 degrees or so).
Returns a list of the following components: Dgen (the genetic distance matrix), Dgeo (the Euclidean distance matrix), Mantel (the statistics of the Mantel test).
Bernd Gruber (bugs? Post to https://groups.google.com/d/forum/dartr)
Rousset, F. (1997). Genetic differentiation and estimation of gene flow from F-statistics under isolation by distance. Genetics, 145(4), 1219-1228.
#because of speed only the first 100 loci ibd <- gl.ibd(bandicoot.gl[,1:100], Dgeo_trans='log(Dgeo)' , Dgen_trans='Dgen/(1-Dgen)') #because of speed only the first 10 individuals) ibd <- gl.ibd(bandicoot.gl[1:10,], distance='euclidean', paircols='pop', Dgeo_trans='Dgeo') #only first 100 loci ibd <- gl.ibd(bandicoot.gl[,1:100], paircols='pop')
#because of speed only the first 100 loci ibd <- gl.ibd(bandicoot.gl[,1:100], Dgeo_trans='log(Dgeo)' , Dgen_trans='Dgen/(1-Dgen)') #because of speed only the first 10 individuals) ibd <- gl.ibd(bandicoot.gl[1:10,], distance='euclidean', paircols='pop', Dgeo_trans='Dgeo') #only first 100 loci ibd <- gl.ibd(bandicoot.gl[,1:100], paircols='pop')
This script calculates an individual based distance matrix.
gl.kosman(x, verbose = NULL)
gl.kosman(x, verbose = NULL)
x |
genlight/dartR object with a unique ploidy |
verbose |
verbosity of the function. |
returns a matrix of [dimensions nInd(x) x nInd(x)] of kosman distances between individuals,
#use only five individuals and seven loci gg <- gl.kosman(possums.gl[1:5,14:21]) gg$kosman gg$nloci
#use only five individuals and seven loci gg <- gl.kosman(possums.gl[1:5,14:21]) gg$kosman gg$nloci
Global spatial autocorrelation is a multivariate approach combining all loci into a single analysis. The autocorrelation coefficient "r" is calculated for each pair of individuals in each specified distance class. For more information see Smouse and Peakall 1999, Peakall et al. 2003 and Smouse et al. 2008.
gl.spatial.autoCorr( x = NULL, Dgeo = NULL, Dgen = NULL, coordinates = "latlon", Dgen_method = "Euclidean", Dgeo_trans = "Dgeo", Dgen_trans = "Dgen", bins = 5, reps = 100, plot.pops.together = FALSE, permutation = TRUE, bootstrap = TRUE, plot.theme = theme_dartR(), plot.colors.pop = NULL, CI.color = "red", plot.out = TRUE, plot.file = NULL, plot.dir = NULL, verbose = NULL )
gl.spatial.autoCorr( x = NULL, Dgeo = NULL, Dgen = NULL, coordinates = "latlon", Dgen_method = "Euclidean", Dgeo_trans = "Dgeo", Dgen_trans = "Dgen", bins = 5, reps = 100, plot.pops.together = FALSE, permutation = TRUE, bootstrap = TRUE, plot.theme = theme_dartR(), plot.colors.pop = NULL, CI.color = "red", plot.out = TRUE, plot.file = NULL, plot.dir = NULL, verbose = NULL )
x |
Genlight object [default NULL]. |
Dgeo |
Geographic distance matrix if no genlight object is provided. This is typically an Euclidean distance but it can be any meaningful (geographical) distance metrics [default NULL]. |
Dgen |
Genetic distance matrix if no genlight object is provided [default NULL]. |
coordinates |
Can be either 'latlon', 'xy' or a two column data.frame
with column names 'lat','lon', 'x', 'y') Coordinates are provided via
|
Dgen_method |
Method to calculate genetic distances. See details [default "Euclidean"]. |
Dgeo_trans |
Transformation to be used on the geographic distances. See Dgen_trans [default "Dgeo"]. |
Dgen_trans |
You can provide a formula to transform the genetic
distance. The transformation can be applied as a formula using Dgen as the
variable to be transformed. For example: |
bins |
The number of bins for the distance classes
(i.e. |
reps |
The number to be used for permutation and bootstrap analyses [default 100]. |
plot.pops.together |
Plot all the populations in one plot. Confidence intervals from permutations are not shown [default FALSE]. |
permutation |
Whether permutation calculations for the null hypothesis of no spatial structure should be carried out [default TRUE]. |
bootstrap |
Whether bootstrap calculations to compute the 95% confidence intervals around r should be carried out [default TRUE]. |
plot.theme |
Theme for the plot. See details [default NULL]. |
plot.colors.pop |
A color palette for populations or a list with as many colors as there are populations in the dataset [default NULL]. |
CI.color |
Color for the shade of the 95% confidence intervals around the r estimates [default "red"]. |
plot.out |
Specify if plot is to be produced [default TRUE]. |
plot.file |
Name for the RDS binary file to save (base name only, exclude extension) [default NULL] temporary directory (tempdir) [default FALSE]. |
plot.dir |
Directory in which to save files [default = working directory] |
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, |
This function executes a modified version
of spautocorr
from the package PopGenReport
. Differently
from PopGenReport
, this function also computes the 95% confidence
intervals around the r via bootstraps, the 95
null hypothesis of no spatial structure and the one-tail test via permutation,
and the correction factor described by Peakall et al 2003.
The input can be i) a genlight object (which has to have the latlon slot
populated), ii) a pair of Dgeo
and Dgen
, which have to be
either
matrix
or dist
objects, or iii) a list
of the
matrix
or dist
objects if the
analysis needs to be carried out for multiple populations (in this case,
all the elements of the list
have to be of the same class (i.e.
matrix
or dist
) and the population order in the two lists has
to be the same.
If the input is a genlight object, the function calculates the linear
distance
for Dgeo
and the relevant Dgen
matrix (see Dgen_method
)
for each population.
When the method selected is a genetic similarity matrix (e.g. "simple"
distance), the matrix is internally transformed with 1 - Dgen
so that
positive values of autocorrelation coefficients indicates more related
individuals similarly as implemented in GenAlEx. If the user provide the
distance matrices, care must be taken in interpreting the results because
similarity matrix will generate negative values for closely related
individuals.
If max(Dgeo)>1000
(e.g. the geographic distances are in thousands of
metres), values are divided by 1000 (in the example before these would then
become km) to facilitate readability of the plots.
If bins
is of length = 1 it is interpreted as the number of (even)
bins to use. In this case the starting point is always the minimum value in
the distance matrix, and the last is the maximum. If it is a numeric vector
of length>1, it is interpreted as the breaking points. In this case, the
first has to be the lowest value, and the last has to be the highest. There
are no internal checks for this and it is user responsibility to ensure that
distance classes are properly set up. If that is not the case, data that fall
outside the range provided will be dropped. The number of bins will be
length(bins) - 1
.
The permutation constructs the 95% confidence intervals around the null hypothesis of no spatial structure (this is a two-tail test). The same data are also used to calculate the probability of the one-tail test (See references below for details).
Bootstrap calculations are skipped and NA
is returned when the number
of possible combinations given the sample size of any given distance class is
< reps
.
Methods available to calculate genetic distances for SNP data:
"propShared" using the function gl.propShared
.
"grm" using the function gl.grm2
.
"Euclidean" using the function gl.dist.ind
.
"Simple" using the function gl.dist.ind
.
"Absolute" using the function gl.dist.ind
.
"Manhattan" using the function gl.dist.ind
.
Methods available to calculate genetic distances for SilicoDArT data:
"Euclidean" using the function gl.dist.ind
.
"Simple" using the function gl.dist.ind
.
"Jaccard" using the function gl.dist.ind
.
"Bray-Curtis" using the function gl.dist.ind
.
Plots and table are saved plot.file in plot.dir if specified.
Examples of other themes that can be used can be consulted in
Returns a data frame with the following columns:
Bin The distance classes
N The number of pairwise comparisons within each distance class
r.uc The uncorrected autocorrelation coefficient
Correction the correction
r The corrected autocorrelation coefficient
L.r The corrected autocorrelation coefficient lower limit
(if bootstap = TRUE
)
U.r The corrected autocorrelation coefficient upper limit
(if bootstap = TRUE
)
L.r.null.uc The uncorrected lower limit for the null hypothesis of no
spatial autocorrelation (if permutation = TRUE
)
U.r.null.uc The uncorrected upper limit for the null hypothesis of no
spatial autocorrelation (if permutation = TRUE
)
L.r.null The corrected lower limit for the null hypothesis of no
spatial autocorrelation (if permutation = TRUE
)
U.r.null The corrected upper limit for the null hypothesis of no
spatial autocorrelation (if permutation = TRUE
)
p.one.tail The p value of the one tail statistical test
Carlo Pacioni, Bernd Gruber & Luis Mijangos (Post to https://groups.google.com/d/forum/dartr)
Smouse PE, Peakall R. 1999. Spatial autocorrelation analysis of individual multiallele and multilocus genetic structure. Heredity 82: 561-573.
Double, MC, et al. 2005. Dispersal, philopatry and infidelity: dissecting local genetic structure in superb fairy-wrens (Malurus cyaneus). Evolution 59, 625-635.
Peakall, R, et al. 2003. Spatial autocorrelation analysis offers new insights into gene flow in the Australian bush rat, Rattus fuscipes. Evolution 57, 1182-1195.
Smouse, PE, et al. 2008. A heterogeneity test for fine-scale genetic structure. Molecular Ecology 17, 3389-3400.
Gonzales, E, et al. 2010. The impact of landscape disturbance on spatial genetic structure in the Guanacaste tree, Enterolobium cyclocarpum (Fabaceae). Journal of Heredity 101, 133-143.
Beck, N, et al. 2008. Social constraint and an absence of sex-biased dispersal drive fine-scale genetic structure in white-winged choughs. Molecular Ecology 17, 4346-4358.
require("dartR.data") res <- gl.spatial.autoCorr(platypus.gl, bins=seq(0,10000,2000)) # using one population, showing sample size test <- gl.keep.pop(platypus.gl,pop.list = "TENTERFIELD") res <- gl.spatial.autoCorr(test, bins=seq(0,10000,2000),CI.color = "green") test <- gl.keep.pop(platypus.gl,pop.list = "TENTERFIELD") res <- gl.spatial.autoCorr(test, bins=seq(0,10000,2000),CI.color = "green")
require("dartR.data") res <- gl.spatial.autoCorr(platypus.gl, bins=seq(0,10000,2000)) # using one population, showing sample size test <- gl.keep.pop(platypus.gl,pop.list = "TENTERFIELD") res <- gl.spatial.autoCorr(test, bins=seq(0,10000,2000),CI.color = "green") test <- gl.keep.pop(platypus.gl,pop.list = "TENTERFIELD") res <- gl.spatial.autoCorr(test, bins=seq(0,10000,2000),CI.color = "green")
This function exports coordinates in a genlight object to a point shape file (including also individual meta data if available). Coordinates are provided under x@other$latlon and assumed to be in WGS84 coordinates, if not proj4 string is provided.
gl2shp( x, type = "shp", proj4 = "+proj=longlat +ellps=WGS84 +datum=WGS84 +no_defs", outfile = "gl", outpath = tempdir(), verbose = NULL )
gl2shp( x, type = "shp", proj4 = "+proj=longlat +ellps=WGS84 +datum=WGS84 +no_defs", outfile = "gl", outpath = tempdir(), verbose = NULL )
x |
Name of the genlight object containing the SNP data and location data, lat longs [required]. |
type |
Type of output 'kml' or 'shp' [default 'shp']. |
proj4 |
Proj4string of data set (see spatialreference.org for projections) [default WGS84]. |
outfile |
Name (path) of the output shape file [default 'gl']. shp extension is added automatically. |
outpath |
Path where to save the output file [default tempdir(), mandated by CRAN]. Use outpath=getwd() or outpath='.' when calling this function to direct output files to your working directory. |
verbose |
Verbosity: 0, silent or fatal errors; 1, begin and end; 2, progress log; 3, progress and results summary; 5, full report [default 2 or as specified using gl.set.verbosity]. |
returns a SpatVector file
Bernd Guber (Post to https://groups.google.com/d/forum/dartr)
out <- gl2shp(testset.gl, outpath=tempdir())
out <- gl2shp(testset.gl, outpath=tempdir())
Carries out calculation for spatial autocorrelation coefficient starting from a genetic and geogaphic distance matrix.
utils.spautocor( GD, GGD, permutation = FALSE, bootstrap = FALSE, bins = 10, reps )
utils.spautocor( GD, GGD, permutation = FALSE, bootstrap = FALSE, bins = 10, reps )
GD |
Genetic distance matrix. |
GGD |
Geographic distance matrix. |
permutation |
Whether permutation calculations for the null hypothesis of no spatial structure should be carried out [default TRUE]. |
bootstrap |
Whether bootstrap calculations to compute the 95% confidence intervals around r should be carried out [default TRUE]. |
bins |
The number of bins for the distance classes
(i.e. |
reps |
The number to be used for permutation and bootstrap analyses [default 100]. |
The code of this function is based one spautocorr
from the package
PopGenReport
, which has been modified to fix a few bugs (as of
PopGenReport v 3.0.4
and allow calculations of bootstraps estimates.
See details from gl.spatial.autoCorr
for a detailed explanation.
Returns a data frame with the following columns:
Bin The distance classes
N The number of pairwise comparisons within each distance class
r.uc The uncorrected autocorrelation coefficient
if both bootstap
and permutation
are FALSE
otherwise only
r
estimates are returned
Carlo Pacioni & Bernd Gruber
Smouse PE, Peakall R. 1999. Spatial autocorrelation analysis of individual multiallele and multilocus genetic structure. Heredity 82: 561-573.
Double, MC, et al. 2005. Dispersal, philopatry and infidelity: dissecting local genetic structure in superb fairy-wrens (Malurus cyaneus). Evolution 59, 625-635.
Peakall, R, et al. 2003. Spatial autocorrelation analysis offers new insights into gene flow in the Australian bush rat, Rattus fuscipes. Evolution 57, 1182-1195.
Smouse, PE, et al. 2008. A heterogeneity test for fine-scale genetic structure. Molecular Ecology 17, 3389-3400.
Gonzales, E, et al. 2010. The impact of landscape disturbance on spatial genetic structure in the Guanacaste tree, Enterolobium cyclocarpum(Fabaceae). Journal of Heredity 101, 133-143.
Beck, N, et al. 2008. Social constraint and an absence of sex-biased dispersal drive fine-scale genetic structure in white-winged choughs. Molecular Ecology 17, 4346-4358.
# See gl.spatial.autoCorr
# See gl.spatial.autoCorr