rate.map {RRmorph}R Documentation

Mapping rate and direction of phenotypic change on 3D surfaces.

Description

Given vectors of RW (or PC) scores, the function selects the RW (PC) axes linked to highest (and lowest) evolutionary rate values and reconstructs the morphology weighted on them. In this way, rate.map shows where and how the phenotype changed the most between any pair of taxa.

Usage

rate.map(x, RR, scores, pcs, mshape, mshape_sur=NULL,refsur=NULL,
  refmat=NULL, k=4,out.rem = TRUE, plot=TRUE, pal=NULL,
  NAcol="gray90",from=NULL, to=NULL,show.names=TRUE)

Arguments

x

the species/nodes to be compared; it can be a single species, or a vector containing two species or a species and any of its parental nodes.

RR

an object generated by using the RRphylo function.

scores

RW or PC scores.

pcs

RW (or PC) vectors (eigenvectors of the covariance matrix) returned by RWA/PCA.

mshape

the consensus configuration.

mshape_sur

a mesh3d object used as a reference for mesh reconstruction. The vertices of mshape_sur must be the consensus configuration. If NULL, it is automatically generated by applying vcgBallPivoting on mshape.

refsur

a list of mesh3d to be provided for all species in x.

refmat

a list of landmark sets to be provided for all species in x.

k

the argument k passed to interpolMesh.

out.rem

logical: if TRUE mesh triangles with outlying area difference are removed.

plot

logical indicating if the 3d plot must be shown.

pal

a vector of colors to be passed to colorRampPalette.

NAcol

the argument NAcol passed to col2mesh.

from, to

lower and upper limits to be associated to the ends of pal.

show.names

logical: if TRUE, the names of the species are displayed in the 3d plot.

Details

After selecting the RW (PC) axes, rate.map automatically builds a 3D mesh on the mean shape calculated from the Relative Warp Analysis (RWA) or Principal Component Analysis (PCA) (Schlager 2017) by applying the function vcgBallPivoting (Rvcg). The reconstruction of species 3d surfaces is based on mshape_sur, either provided by the user or generated within the function. Finally, for each species in x, the function computes the area differences between corresponding triangles of its reconstructed 3D mesh and the surface of the ancestor (most recent common ancestor in the case of two species in x). In the calculation of differences, the possibility to find and remove outliers is supplied (out.rem=TRUE, we suggest considering this possibility if the mesh may contain degenerate facets).

Finally, rate.map returns a 3D plot showing such comparisons displayed on shape of the ancestor/mrca used as the reference. The color gradient goes from blue to red, where blue areas represent expansion of the mesh, while the red areas represent contraction of the mesh triangles. If a list refsur (and refmat) is provided, convergence is plotted onto them (see interpolMesh for further details).

Value

The function returns a list including:

Author(s)

Marina Melchionna, Antonio Profico, Silvia Castiglione, Gabriele Sansalone, Pasquale Raia

References

Schlager, S. (2017). Morpho and Rvcg-Shape Analysis in R: R-Packages for geometric morphometrics, shape analysis and surface manipulations. In: Statistical shape and deformation analysis. Academic Press.

Castiglione, S., Melchionna, M., Profico, A., Sansalone, G., Modafferi, M., Mondanaro, A., Wroe, S., Piras, P., & Raia, P. (2021). Human face-off: a new method for mapping evolutionary rates on three-dimensional digital models. Palaeontology, 65, 1. doi:10.1111/pala.12582

Melchionna, M., Castiglione, S., Girardi, G., Serio, C., Esposito, A., Mondanaro, A., Profico, A., Sansalone, G., & Raia, P. (2024). RRmorph—a new R packageto map phenotypic evolutionary rates and patterns on 3D meshes. Communications Biology, 7, 1009.

See Also

RRphylo vignette ; relWarps ; procSym

Examples

  
  da<-"https://github.com/pasraia/RRmorph_example_data/raw/refs/heads/main/RRmorphdata.rda"
  download.file(url=da,destfile = paste0(tempdir(),"/RRmorphdata.rda"))
  load(paste0(tempdir(),"/RRmorphdata.rda"))

  require(Morpho)
  require(Rvcg)

  pca<-procSym(endo.set)
  ldm_homo<-endo.set[,,"Homo_sapiens"]
  sur_homo<-endo.sur[["Homo_sapiens"]]
  ldm_macaca<-endo.set[,,"Macaca_fuscata"]
  sur_macaca<-endo.sur[["Macaca_fuscata"]]

  cc<- 2/parallel::detectCores()
  RR<-RRphylo::RRphylo(tree.prima,pca$PCscores,clus=cc)

  # plotting on reconstructed surfaces
  Rmap1<-rate.map(x=c("Homo_sapiens","Macaca_fuscata"),RR=RR, scores=pca$PCscores,
                  pcs=pca$PCs, mshape=pca$mshape)

  # plotting on real surfaces
  Rmap2<-rate.map(x=c("Homo_sapiens","Macaca_fuscata"),RR=RR, scores=pca$PCscores,
                  pcs=pca$PCs, mshape=pca$mshape,
                  refsur=list("Homo_sapiens"=sur_homo,"Macaca_fuscata"=sur_macaca),
                  refmat=list("Homo_sapiens"=ldm_homo,"Macaca_fuscata"=ldm_macaca))
  

[Package RRmorph version 0.0.1 Index]