| Title: | 3D Morphological Analyses with 'RRphylo' |
|---|---|
| Description: | Combined with 'RRphylo', this package provides a powerful tool to analyse and visualise 3d models (surfaces and meshes) in a phylogenetically explicit context (Melchionna et al., 2024 <doi:10.1038/s42003-024-06710-8>). |
| Authors: | Marina Melchionna [aut], Silvia Castiglione [aut, cre], Carmela Serio [aut], Giorgia Girardi [aut], Alessandro Mondanaro [aut], Pasquale Raia [aut] |
| Maintainer: | Silvia Castiglione <[email protected]> |
| License: | GPL-2 |
| Version: | 0.0.2 |
| Built: | 2026-05-09 07:09:34 UTC |
| Source: | https://github.com/cran/RRmorph |
The function colors a mesh according to a vector of
continuous values related to individual vertices.
col2mesh(mesh,values,pal,from=NULL,to=NULL,NAcol="gray90")col2mesh(mesh,values,pal,from=NULL,to=NULL,NAcol="gray90")
mesh |
a |
values |
a vector of continuous values associated to individual vertices
of the |
pal |
a vector of colors to be passed to
|
from, to
|
lower and upper |
NAcol |
the color associated to NA |
A mesh3d object colored according to values.
Marina Melchionna, Silvia Castiglione
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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) interp<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val,element ="vertices",k = 4) colmesh<-col2mesh(mesh = sur,values = interp,pal = heat.colors(5)) plotLegend(mesh = colmesh,values = interp, main = "Pan troglodytes") open3d() shade3d(colmesh,specular="black")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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) interp<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val,element ="vertices",k = 4) colmesh<-col2mesh(mesh = sur,values = interp,pal = heat.colors(5)) plotLegend(mesh = colmesh,values = interp, main = "Pan troglodytes") open3d() shade3d(colmesh,specular="black")
Given vectors of RW (or PC) scores for some converging species, the function selects the RW (PC) axes which best account for convergence and maps convergent areas on the corresponding 3D surfaces.
conv.map(x1,x2=NULL,scores, pcs, mshape,focal=NULL,mshape_sur=NULL, refmat = NULL,refsur = NULL, k = 4, exclude = NULL, out.rem = TRUE, plot = TRUE, col = "blue", NAcol = "gray", names = TRUE, nsim = 1000)conv.map(x1,x2=NULL,scores, pcs, mshape,focal=NULL,mshape_sur=NULL, refmat = NULL,refsur = NULL, k = 4, exclude = NULL, out.rem = TRUE, plot = TRUE, col = "blue", NAcol = "gray", names = TRUE, nsim = 1000)
x1, x2
|
vectors of convergent species. When convergence within a single
clade was found, |
scores |
data frame (or matrix) with the RW (or PC) scores returned by
RWA/PCA. Species not included in |
pcs |
RW (or PC) vectors (eigenvectors of the covariance matrix) returned by RWA/PCA. |
mshape |
the consensus configuration. |
focal |
vector of species included in |
mshape_sur |
a |
refmat |
a named list of landmark sets corresponding to |
refsur |
a named list of |
k |
the argument |
exclude |
integer: the index numbers of the RWs (or PCs) to be excluded from the comparison. |
out.rem |
logical: if |
plot |
logical: if |
col |
character: the color for plotting. |
NAcol |
the argument |
names |
logical: if |
nsim |
the number of iterations to evaluate significance. |
After selecting the RW (PC) axes which best account for convergence,
conv.map uses such axes (and related scores) within
restoreShapes (Morpho) to reconstruct landmark
matrices for each convergent species (in x1/x2). The reconstruction
of species 3d surfaces is based on mshape_sur, either provided by
the user or generated within the function. Finally, the area differences
between corresponding triangles of reconstructed 3d meshes for each
possible pair of convergent species are calculated. 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).
If the combination of focal species (or species within
refsur/refmat) contains a number equal or lower then 5 items,
conv.map returns a rgl plot mapping the convergence on the 3D
models. If lists of refsur/refmat are not provided, the area
differences are plotted onto reconstructed surfaces. If
refsur/refmat are available, difference values are interpolated by
means of interpolMesh to be plotted onto real surfaces. When
species in either x1 or x2 are missing from focal or
refmat/refsur, conv.map plots the reconstructed surface of
the species having the smallest $selected angle with the focal (see
angle.compare in the description of outputs).
conv.map further gives the opportunity to exclude some RW (or PC)
axes from the analysis because, for example, in most cases the first axes
are mainly related to high-order morphological differences driven by
phylogeny and size variations.
The function returns a list including:
$angle.compare: a data frame including the real angles
between species shape vectors $real.angle, the angles computed between vectors
of the selected RWs (or PCs) $selected, the angles between vectors of the
non-selected RWs (or PCs) $others, the differences selected-others
and its p-values.
$selected.pcs RWs (or PCs) axes selected for convergence.
$average.dist symmetric matrix of pairwise distances between 3D surfaces.
$suface1 list of colored surfaces representing convergence between mesh A and B charted on mesh A.
$suface2 list of colored surfaces representing convergence between mesh A and B charted on mesh B.
$scale the value used to set the color gradient, computed as the maximum of all differences between each surface and the mean shape.
Marina Melchionna, Antonio Profico, Silvia Castiglione, Carmela Serio, Gabriele Sansalone, Pasquale Raia
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.
Melchionna, M., Profico, A., Castiglione, S., Serio, C., Mondanaro, A., Modafferi, M., Tamagnini, D., Maiorano, L. , Raia, P., Witmer, L.M., Wroe, S., & Sansalone, G. (2021). A method for mapping morphological convergence on three-dimensional digital models: the case of the mammalian sabre-tooth. Palaeontology, 64, 573–584. doi:10.1111/pala.12542
search.conv vignette ;
relWarps ; procSym
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) 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"]] # Convergence within group plotted on reconstructed surfaces cm1<-conv.map(x1=c("Pan_troglodytes","Gorilla_gorilla","Pongo_abelii"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, focal=c("Pan_troglodytes","Gorilla_gorilla")) # Convergence between group plotted on reconstructed surfaces cm2<-conv.map(x1=c("Pongo_abelii"),x2=c("Alouatta_caraya"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, focal="Alouatta_caraya") # Convergence within group plotted on real surfaces cm3<-conv.map(x1=c("Homo_sapiens","Gorilla_gorilla","Pongo_abelii"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, refsur=list("Homo_sapiens"=sur_homo), refmat=list("Homo_sapiens"=ldm_homo)) # Convergence between group plotted on real surfaces cm3<-conv.map(x1=c("Homo_sapiens","Pongo_abelii"),x2=c("Macaca_fuscata"), 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))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) 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"]] # Convergence within group plotted on reconstructed surfaces cm1<-conv.map(x1=c("Pan_troglodytes","Gorilla_gorilla","Pongo_abelii"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, focal=c("Pan_troglodytes","Gorilla_gorilla")) # Convergence between group plotted on reconstructed surfaces cm2<-conv.map(x1=c("Pongo_abelii"),x2=c("Alouatta_caraya"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, focal="Alouatta_caraya") # Convergence within group plotted on real surfaces cm3<-conv.map(x1=c("Homo_sapiens","Gorilla_gorilla","Pongo_abelii"), scores=pca$PCscores,pcs=pca$PCs,mshape=pca$mshape, refsur=list("Homo_sapiens"=sur_homo), refmat=list("Homo_sapiens"=ldm_homo)) # Convergence between group plotted on real surfaces cm3<-conv.map(x1=c("Homo_sapiens","Pongo_abelii"),x2=c("Macaca_fuscata"), 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))
The function takes a reconstructed mesh3d object
(sur) with some related values (to either triangles or
vertices of the mesh) and transfers such values to the real mesh
(refsur) from which sur was derived.
interpolMesh(sur,values,refsur,refmat,element=c("triangles","vertices"),k=4)interpolMesh(sur,values,refsur,refmat,element=c("triangles","vertices"),k=4)
sur |
a reconstructed |
values |
the vector of values related to |
refsur |
the reference mesh ( |
refmat |
the landmark set related to |
element |
one of |
k |
the number of nearest neighbor vertices used for interpolation (see details). |
The function starts by locating a set of points (NNps) on
refsur, each being the single nearest neighbor for each vertex of
sur (or barycenter if element="triangles"). Then,
interpolation is performed by identifying the k points among NNps
being the closest to each vertex of refsur and computing the mean of
their values weighted by their distance.
The vector of values related to each vertex of refsur.
Marina Melchionna, Silvia Castiglione
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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val1<-rnorm(ncol(rec$vb)) # Interpolate values associated to vertices val1<-rnorm(ncol(rec$vb)) interp1<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val1,element ="vertices",k = 4) colmesh1<-col2mesh(mesh = sur,values = interp1,pal = heat.colors(5)) open3d() shade3d(colmesh1,specular="black") # Interpolate values associated to triangles val2<-rnorm(ncol(rec$it)) interp2<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val2,element ="triangles",k = 4) colmesh2<-col2mesh(mesh = sur,values = interp2,pal = heat.colors(5)) open3d() shade3d(colmesh2,specular="black")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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val1<-rnorm(ncol(rec$vb)) # Interpolate values associated to vertices val1<-rnorm(ncol(rec$vb)) interp1<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val1,element ="vertices",k = 4) colmesh1<-col2mesh(mesh = sur,values = interp1,pal = heat.colors(5)) open3d() shade3d(colmesh1,specular="black") # Interpolate values associated to triangles val2<-rnorm(ncol(rec$it)) interp2<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val2,element ="triangles",k = 4) colmesh2<-col2mesh(mesh = sur,values = interp2,pal = heat.colors(5)) open3d() shade3d(colmesh2,specular="black")
The function relates PCA loadings of a single PC axis to individual landmarks and plots them on a 3d mesh by means of interpolation.
plotland(pca,sel=1,refsur=NULL,refmat=NULL,k=5,pal=NULL, defo=FALSE,radius=0.001)plotland(pca,sel=1,refsur=NULL,refmat=NULL,k=5,pal=NULL, defo=FALSE,radius=0.001)
pca |
the result of a relative warp analysis. Classes |
sel |
numeric indicating the focal RW/PC axis. |
refsur |
the |
refmat |
the landmark set related to |
k |
the argument |
pal |
a vector of colors to be passed to |
defo |
when |
radius |
argument |
A list including a mesh3d object colored according to landmarks
importance and a matrix of landmarks importance on each RW/PC axis. Additionally,
the function returns a 3d plot of the mesh.
Marina Melchionna, Silvia Castiglione, Carmela Serio, Giorgia Girardi
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(rgl) require(Morpho) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] plotland(pca=pca,sel=1,refsur = sur,refmat = ldm)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(rgl) require(Morpho) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] plotland(pca=pca,sel=1,refsur = sur,refmat = ldm)
Assuming a mesh is colored according to a vector of values,
the function takes the color sequence from the mesh and plots it
associated to values.
plotLegend(mesh,values,main)plotLegend(mesh,values,main)
mesh |
a |
values |
a vector of continuous values associated to individual vertices
of the |
main |
plot title. |
A plot of the color sequence associated to values on the mesh.
Marina Melchionna, Silvia Castiglione
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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) interp<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val,element ="vertices",k = 4) colmesh<-col2mesh(mesh = sur,values = interp,pal = heat.colors(5)) plotLegend(mesh = colmesh,values = interp, main = "Pan troglodytes") open3d() shade3d(colmesh,specular="black")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(rgl) require(Morpho) require(Rvcg) pca<-procSym(endo.set) ldm<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) interp<-interpolMesh(sur = rec,refsur = sur,refmat = ldm, values = val,element ="vertices",k = 4) colmesh<-col2mesh(mesh = sur,values = interp,pal = heat.colors(5)) plotLegend(mesh = colmesh,values = interp, main = "Pan troglodytes") open3d() shade3d(colmesh,specular="black")
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.
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)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)
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 |
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 |
refsur |
a list of |
refmat |
a list of landmark sets to be provided for all species in
|
k |
the argument |
out.rem |
logical: if |
plot |
logical indicating if the 3d plot must be shown. |
pal |
a vector of colors to be passed to
|
NAcol |
the argument |
from, to
|
lower and upper limits to be associated to the ends of
|
show.names |
logical: if |
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).
The function returns a list including:
$selected a list of RWs/PCs axes selected for higher evolutionary rates for each species.
$surfaces a list of reconstructed colored surfaces of the given species/node and the most recent common ancestor.
differences a list area differences between corresponding triangles of species reconstructed 3d mesh and the surface of the ancestor.
lmks if refmat
is not NULL, this is the landmark configuration rotated on the
reconstructed surface .
Marina Melchionna, Antonio Profico, Silvia Castiglione, Gabriele Sansalone, Pasquale Raia
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.
RRphylo
vignette ; relWarps ; procSym
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))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))
The function reconstructs two specimens' meshes (x)
by using their superimposed configurations (from within pca),
calculates the shape difference between them, and plots such differences
on provided meshes (refsur).
shapeDiff(x,pca,refsur,refmat,mshape_sur = NULL, pal=NULL,NAcol="gray90",show.names=TRUE)shapeDiff(x,pca,refsur,refmat,mshape_sur = NULL, pal=NULL,NAcol="gray90",show.names=TRUE)
x |
a vector of specimens pair. |
pca |
the result of a relative warp analysis. Classes |
refsur |
a list of two |
refmat |
a list of two landmark sets related to |
mshape_sur |
a |
pal |
a vector of colors to be passed to |
NAcol |
the color associated to |
show.names |
logical: if |
Two mesh3d objects colored according to shape differences.
Additionally, the function returns 3d plots of the meshes.
Marina Melchionna, Silvia Castiglione
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) 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"]] diffs<-RRmorph::shapeDiff(x=c("Homo_sapiens","Macaca_fuscata"), pca = pca,refsur = list(sur_homo,sur_macaca), refmat = list(ldm_homo,ldm_macaca))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) 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"]] diffs<-RRmorph::shapeDiff(x=c("Homo_sapiens","Macaca_fuscata"), pca = pca,refsur = list(sur_homo,sur_macaca), refmat = list(ldm_homo,ldm_macaca))
The function transfers values associated to triangles of a
mesh3d object to its vertices.
tri2verts(mesh,values)tri2verts(mesh,values)
mesh |
a |
values |
a vector of continuous values associated to individual triangles
of the |
A vector of continuous values associated to individual vertices of
the mesh.
Marina Melchionna, Silvia Castiglione
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<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) vertval<-tri2verts(rec,val)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<-endo.set[,,"Homo_sapiens"] sur<-endo.sur[["Homo_sapiens"]] rec<- vcgBallPivoting(pca$mshape, radius = 0) rec$vb[1:3,]<-t(ldm) val<-rnorm(ncol(rec$vb)) vertval<-tri2verts(rec,val)